63 int NmodesUproj,
int NmodesPproj,
int NmodesNut,
int NmodesSup,
73 for (
int i = 0;
i < NmodesUproj;
i++)
78 for (
int i = 0;
i < NmodesSup;
i++)
91 UprojN = NmodesUproj + NmodesSup;
108 Eigen::VectorXd uresidualOld = Eigen::VectorXd::Zero(
UprojN);
109 Eigen::VectorXd presidualOld = 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 =
119 float normalizedResidualLim =
124 scalar residual_jump(1 + residualJumpLim);
134 phi = fvc::interpolate(
U) &
U.mesh().Sf();
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;
151 volScalarField&
nut =
const_cast<volScalarField&
>
157 PtrList<volVectorField> gradModP;
159 for (
int i = 0;
i < NmodesPproj;
i++)
166 while ((residual_jump > residualJumpLim
167 || std::max(U_norm_res, P_norm_res) > normalizedResidualLim)
171 std::cout <<
"Iteration " << iter << std::endl;
172#if defined(OFVER) && (OFVER == 6)
181 - fvm::laplacian(nueff,
U)
182 - fvc::div(nueff * dev2(
T(fvc::grad(
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));
199 fvc::interpolate(
rAtU() -
rAU) * fvc::snGrad(P) *
mesh.magSf();
203 List<Eigen::MatrixXd> RedLinSysP;
205 while (
simple.correctNonOrthogonal())
215 if (
simple.finalNonOrthogonalIter())
223 U.correctBoundaryConditions();
226 uresidualOld = uresidualOld.cwiseAbs();
227 presidualOld = presidualOld.cwiseAbs();
228 residual_jump = std::max(uresidualOld.sum(), presidualOld.sum());
233 U_norm_res =
uresidual.sum() / (RedLinSysU[1].cwiseAbs()).sum();
234 P_norm_res =
presidual.sum() / (RedLinSysP[1].cwiseAbs()).sum();
238 std::cout <<
"Residual jump = " << residual_jump << std::endl;
239 std::cout <<
"Normalized residual = " << std::max(U_norm_res,
240 P_norm_res) << std::endl;
244 std::cout <<
"Solution " <<
counter <<
" converged in " << iter <<
245 " iterations." << std::endl;
246 std::cout <<
"Final normalized residual for velocity: " << U_norm_res <<
248 std::cout <<
"Final normalized residual for pressure: " << P_norm_res <<
261 "Imposed boundary conditions dimensions do not match given values matrix dimensions");
262 Eigen::MatrixXd vel_scal;
263 vel_scal.resize(vel.rows(), vel.cols());
272 vel_scal(k, 0) = vel(k, 0) / u_lf;
#define M_Assert(Expr, Msg)
Header file of the reducedSteadyNS class.
IOdictionary * ITHACAdict
Dictionary for input objects from file.
GeometricField< Type, PatchField, GeoMesh > reconstruct(GeometricField< Type, PatchField, GeoMesh > &inputField, Eigen::MatrixXd Coeff, word Name)
Function to reconstruct the solution starting from the coefficients, in this case the field is passed...
List< Eigen::MatrixXd > project(fvMatrix< Type > &Af, label numberOfModes=0, word projType="G")
A function that project an FvMatrix (OpenFoam linear System) on the modes.
PtrList< GeometricField< Type, PatchField, GeoMesh > > & toPtrList()
Function that returns the Modes object as a standard PtrList.
Implementation of a parametrized full order steady NS problem and preparation of the the reduced ma...
std::vector< SPLINTER::RBFSpline * > rbfSplines
Create a RBF splines for interpolation.
volScalarModes nutModes
List of POD modes for eddy viscosity.
static Eigen::MatrixXd solveLinearSys(List< Eigen::MatrixXd > LinSys, Eigen::MatrixXd x, Eigen::VectorXd &residual, const Eigen::MatrixXd &bc=Eigen::MatrixXd::Zero(0, 0), const std::string solverType="fullPivLu")
Linear system solver for the online problem.
void solveOnline_Simple(scalar mu_now, int NmodesUproj, int NmodesPproj, int NmodesNut=0, int NmodesSup=0, word Folder="./ITHACAoutput/Reconstruct/")
Method to perform an online solve using a PPE stabilisation method.
int maxIterOn
Maximum iterations number for the online step.
SteadyNSSimple * problem
Full problem.
Eigen::MatrixXd vel_now
Imposed boundary conditions.
reducedSimpleSteadyNS()
Construct Null.
Eigen::MatrixXd projGradModP
Projected gradient of the pressure modes.
volVectorModes ULmodes
Lifted velocity modes.
void setOnlineVelocity(Eigen::MatrixXd vel)
It checks if the number of imposed boundary conditions is correct and set the inlet velocity equal to...
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
autoPtr< surfaceScalarField > _phi
Flux.
autoPtr< simpleControl > _simple
simpleControl
volVectorModes supmodes
List of pointers used to form the supremizer modes.
autoPtr< Time > _runTime
Time.
volVectorModes Umodes
List of pointers used to form the velocity modes.
autoPtr< fvMesh > _mesh
Mesh.
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
autoPtr< volVectorField > _U
Velocity field.
volScalarModes Pmodes
List of pointers used to form the pressure modes.
autoPtr< volScalarField > _p
Pressure field.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
bool isTurbulent()
This function checks if the case is turbulent.
simpleControl simple(mesh)
tmp< volScalarField > rAtU(rAU)
surfaceScalarField phiHbyA("phiHbyA", fvc::flux(HbyA))
volVectorField HbyA(constrainHbyA(rAU *UEqn.H(), U, p))
volScalarField rAU(1.0/UEqn.A())