63 int NmodesUproj,
int NmodesPproj,
int NmodesNut,
int NmodesSup,
68 for (
int i = 0; i <
problem->inletIndex.rows(); i++)
73 for (
int i = 0; i < NmodesUproj; i++)
78 for (
int i = 0; i < NmodesSup; i++)
91 UprojN = NmodesUproj + NmodesSup;
96 PprojN =
problem->Pmodes.size();
100 PprojN = NmodesPproj;
105 NmodesNut =
problem->nutModes.size();
108 Eigen::VectorXd uresidualOld = Eigen::VectorXd::Zero(UprojN);
109 Eigen::VectorXd presidualOld = Eigen::VectorXd::Zero(PprojN);
110 Eigen::VectorXd uresidual = Eigen::VectorXd::Zero(UprojN);
111 Eigen::VectorXd presidual = Eigen::VectorXd::Zero(PprojN);
112 scalar U_norm_res(1);
113 scalar P_norm_res(1);
114 Eigen::MatrixXd a = Eigen::VectorXd::Zero(UprojN);
115 Eigen::MatrixXd b = Eigen::VectorXd::Zero(PprojN);
117 float residualJumpLim =
118 problem->para->ITHACAdict->lookupOrDefault<
float>(
"residualJumpLim", 1e-5);
119 float normalizedResidualLim =
120 problem->para->ITHACAdict->lookupOrDefault<
float>(
"normalizedResidualLim",
124 scalar residual_jump(1 + residualJumpLim);
126 volScalarField& P =
problem->_p();
127 volVectorField& U =
problem->_U();
128 fvMesh& mesh =
problem->_mesh();
129 Time& runTime =
problem->_runTime();
131 surfaceScalarField& phi(
problem->_phi());
132 ULmodes.reconstruct(U, a,
"U");
133 problem->Pmodes.reconstruct(P, b,
"p");
134 phi = fvc::interpolate(U) & U.mesh().Sf();
136 simpleControl& simple =
problem->_simple();
140 Eigen::MatrixXd nutCoeff;
141 nutCoeff.resize(NmodesNut, 1);
143 for (
int i = 0; i < NmodesNut; i++)
145 Eigen::MatrixXd muEval;
147 muEval(0, 0) = mu_now;
148 nutCoeff(i, 0) =
problem->rbfSplines[i]->eval(muEval);
151 volScalarField& nut =
const_cast<volScalarField&
>
152 (
problem->_mesh().lookupObject<volScalarField>(
"nut"));
153 problem->nutModes.reconstruct(nut, nutCoeff,
"nut");
157 PtrList<volVectorField> gradModP;
159 for (
int i = 0; i < NmodesPproj; i++)
161 gradModP.append(fvc::grad(
problem->Pmodes[i]));
166 while ((residual_jump > residualJumpLim
167 || std::max(U_norm_res, P_norm_res) > normalizedResidualLim)
171 Info <<
"Iteration " << iter << endl;
172#if defined(OFVER) && (OFVER == 6)
173 simple.loop(runTime);
177 volScalarField nueff =
problem->turbulence->nuEff().ref();
181 - fvm::laplacian(nueff, U)
182 - fvc::div(nueff * dev2(T(fvc::grad(U))))
185 List<Eigen::MatrixXd> RedLinSysU =
ULmodes.project(UEqn, UprojN);
188 ULmodes.reconstruct(U, a,
"U");
189 volScalarField rAU(1.0 / UEqn.A());
190 volVectorField HbyA(constrainHbyA(1.0 / UEqn.A() * UEqn.H(), U, P));
191 surfaceScalarField phiHbyA(
"phiHbyA", fvc::flux(HbyA));
192 adjustPhi(phiHbyA, U, P);
193 tmp<volScalarField> rAtU(rAU);
195 if (simple.consistent())
197 rAtU = 1.0 / (1.0 / rAU - UEqn.H1());
199 fvc::interpolate(rAtU() - rAU) * fvc::snGrad(P) * mesh.magSf();
200 HbyA -= (rAU - rAtU()) * fvc::grad(P);
203 List<Eigen::MatrixXd> RedLinSysP;
205 while (simple.correctNonOrthogonal())
209 fvm::laplacian(rAtU(), P) == fvc::div(phiHbyA)
211 RedLinSysP =
problem->Pmodes.project(pEqn, PprojN);
213 problem->Pmodes.reconstruct(P, b,
"p");
215 if (simple.finalNonOrthogonalIter())
217 phi = phiHbyA - pEqn.flux();
222 U = HbyA - rAtU() * fvc::grad(P);
223 U.correctBoundaryConditions();
224 uresidualOld = uresidualOld - uresidual;
225 presidualOld = presidualOld - presidual;
226 uresidualOld = uresidualOld.cwiseAbs();
227 presidualOld = presidualOld.cwiseAbs();
228 residual_jump = std::max(uresidualOld.sum(), presidualOld.sum());
229 uresidualOld = uresidual;
230 presidualOld = presidual;
231 uresidual = uresidual.cwiseAbs();
232 presidual = presidual.cwiseAbs();
233 U_norm_res = uresidual.sum() / (RedLinSysU[1].cwiseAbs()).sum();
234 P_norm_res = presidual.sum() / (RedLinSysP[1].cwiseAbs()).sum();
238 Info <<
"Residual jump = " << residual_jump << endl;
239 Info <<
"Normalized residual = " << std::max(U_norm_res,
244 Info <<
"Solution " <<
counter <<
" converged in " << iter <<
245 " iterations." << endl;
246 Info <<
"Final normalized residual for velocity: " << U_norm_res <<
248 Info <<
"Final normalized residual for pressure: " << P_norm_res <<
250 ULmodes.reconstruct(U, a,
"Uaux");
252 problem->Pmodes.reconstruct(P, b,
"Paux");
255 runTime.setTime(runTime.startTime(), 0);