36#include "viscosityModel.H"
49 NUmodesOut = para->ITHACAdict->lookupOrDefault<label>(
"NmodesUout", 15);
51 NPmodesOut = para->ITHACAdict->lookupOrDefault<label>(
"NmodesPout", 15);
53 NNutModesOut = para->ITHACAdict->lookupOrDefault<label>(
"NmodesNutOut", 15);
55 NUmodes = para->ITHACAdict->lookupOrDefault<label>(
"NmodesUproj", 10);
57 NPmodes = para->ITHACAdict->lookupOrDefault<label>(
"NmodesPproj", 10);
59 NNutModes = para->ITHACAdict->lookupOrDefault<label>(
"NmodesNutProj", 0);
72 Eigen::MatrixXd weights;
76 word weightName =
"wRBF_M" + name(i + 1);
80 samples[i] =
new SPLINTER::DataTable(1, 1);
82 for (label j = 0; j <
coeffL2.cols();
90 SPLINTER::RadialBasisFunctionType::GAUSSIAN, weights);
92 Info <<
"Constructing RadialBasisFunction for mode " << i + 1 << endl;
96 samples[i] =
new SPLINTER::DataTable(1, 1);
98 for (label j = 0; j <
coeffL2.cols();
105 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
107 "./ITHACAoutput/weights/", weightName);
109 Info <<
"Constructing RadialBasisFunction for mode " << i + 1 << endl;
117 volScalarField& p =
_p();
118 volVectorField& U =
_U();
119 fvMesh& mesh =
_mesh();
120 surfaceScalarField& phi =
_phi();
121 simpleControl& simple =
_simple();
124 scalar uresidual = 1;
125 Vector<double> uresidual_v(0, 0, 0);
126 scalar presidual = 1;
129 std::ofstream res_os;
130 std::ofstream snaps_os;
134 res_os.open(Folder +
"/residuals", std::ios_base::app);
135 snaps_os.open(Folder +
"/snaps", std::ios_base::app);
136 iters.open(Folder +
"/iters", std::ios_base::app);
137 res_U.open(Folder + name(
counter) +
"/res_U", std::ios_base::app);
138 res_P.open(Folder + name(
counter) +
"/res_P", std::ios_base::app);
141 middleStep = para->ITHACAdict->lookupOrDefault<label>(
"middleStep", 20);
142 middleExport = para->ITHACAdict->lookupOrDefault<
bool>(
"middleExport",
true);
143#if defined(OFVER) && (OFVER == 6)
152 Info <<
"Time = " << runTime.timeName() << nl << endl;
153 volScalarField nueff =
turbulence->nuEff().ref();
157 - fvm::laplacian(nueff, U)
158 - fvc::div(nueff * dev2(T(fvc::grad(U))))
162 if (simple.momentumPredictor())
164 uresidual_v = solve(UEqn == - fvc::grad(p)).initialResidual();
167 volVectorField U2(U);
170 for (label i = 0; i < 3; i++)
172 if (C < uresidual_v[i])
179 volScalarField rAU(1.0 / UEqn.A());
180 volVectorField HbyA(constrainHbyA(1.0 / UEqn.A() * UEqn.H(), U, p));
181 surfaceScalarField phiHbyA(
"phiHbyA", fvc::flux(HbyA));
182 adjustPhi(phiHbyA, U, p);
183 tmp<volScalarField> rAtU(rAU);
185 if (simple.consistent())
187 rAtU = 1.0 / (1.0 / rAU - UEqn.H1());
189 fvc::interpolate(rAtU() - rAU) * fvc::snGrad(p) * mesh.magSf();
190 HbyA -= (rAU - rAtU()) * fvc::grad(p);
195 while (simple.correctNonOrthogonal())
199 fvm::laplacian(rAtU(), p) == fvc::div(phiHbyA)
205 presidual = pEqn.solve().initialResidual();
209 pEqn.solve().initialResidual();
212 if (simple.finalNonOrthogonalIter())
214 phi = phiHbyA - pEqn.flux();
220#include "continuityErrs.H"
223 U = HbyA - rAtU() * fvc::grad(p);
224 U.correctBoundaryConditions();
225 residual = max(presidual, uresidual);
226 Info <<
"Time = " << runTime.timeName() << nl << endl;
227 laminarTransport.correct();
245 res_U << uresidual << endl;
246 res_P << presidual << endl;
250 auto nut = mesh.lookupObject<volScalarField>(
"nut");
257 snaps_os <<
folderN + 1 << endl;
258 iters << csolve << endl;
259 res_os << residual << endl;
265 runTime.setTime(runTime.startTime(), 0);
280 auto nut = mesh.lookupObject<volScalarField>(
"nut");
292 for (label i = 0; i < mu_now.size(); i++)
Header file of the steadyNS class.
std::vector< SPLINTER::DataTable * > samples
Create a samples for interpolation.
std::vector< SPLINTER::RBFSpline * > rbfSplines
Create a RBF splines for interpolation.
SteadyNSSimple()
Null constructor.
label folderN
Counter to save intermediate steps in the correct folder (for turbulent case only).
void getTurbRBF(label NNutModes)
Function to calculate RBF weights for turbulence.
void truthSolve2(List< scalar > mu_now, word Folder="./ITHACAoutput/Offline/")
Offline solver for the whole Navier-Stokes problem.
volScalarModes nutModes
List of POD modes for eddy viscosity.
bool middleExport
Export also intermediate fields.
Eigen::MatrixXd coeffL2
The matrix of L2 projection coefficients for the eddy viscosity.
label saver
Counter to check if the middleStep has been reached or not (for turbulent case only).
label middleStep
Distancing between intermediate steps (for turbulent case only).
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
void writeMu(List< scalar > mu_now)
Write out a list of scalar corresponding to the parameters used in the truthSolve.
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
label counter
Counter used for the output of the full order solutions.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
Eigen::MatrixXd mu
Row matrix of parameters.
scalar maxIter
Number of maximum iterations to be done for the computation of the truth solution.
label NPmodes
Number of pressure modes used for the projection.
autoPtr< surfaceScalarField > _phi
Flux.
autoPtr< simpleControl > _simple
simpleControl
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
autoPtr< Time > _runTime
Time.
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
steadyNS()
Null constructor.
scalar tolerance
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and p...
autoPtr< fvMesh > _mesh
Mesh.
autoPtr< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model).
label NUmodes
Number of velocity modes used for the projection.
scalar pRefValue
Reference pressure value.
label pRefCell
Reference pressure cell.
label NNutModesOut
Number of nut modes to be calculated.
label NNutModes
Number of nut modes used for the projection.
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
label NUmodesOut
Number of velocity modes to be calculated.
label NPmodesOut
Number of pressure modes to be calculated.
autoPtr< volVectorField > _U
Velocity field.
autoPtr< volScalarField > _p
Pressure field.
void ReadDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Read a dense matrix from a binary format file.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
void SaveDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Save a dense matrix to a binary format file.
Eigen::VectorXd getCoeffs(GeometricField< Type, PatchField, GeoMesh > &snapshot, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Projects a snapshot on a basis function and gets the coefficients of the projection.
bool isTurbulent()
This function checks if the case is turbulent.
bool check_file(std::string fileName)
Function that returns true if a file exists.