Loading...
Searching...
No Matches
ReducedSimpleSteadyNS.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
33
35
36// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
37
38// Constructor
42
44 :
45 problem(&FOMproblem)
46{
47 // Create a new Umodes set where the first ones are the lift functions
48 for (int i = 0; i < problem->inletIndex.rows(); i++)
49 {
50 ULmodes.append((problem->liftfield[i]).clone());
51 }
52
53 for (int i = 0; i < problem->Umodes.size(); i++)
54 {
55 ULmodes.append((problem->Umodes.toPtrList()[i]).clone());
56 }
57}
58
59// * * * * * * * * * * * * * * * Solve Functions * * * * * * * * * * * * * //
60
61
63 int NmodesUproj, int NmodesPproj, int NmodesNut, int NmodesSup,
64 word Folder)
65{
66 ULmodes.resize(0);
67
68 for (int i = 0; i < problem->inletIndex.rows(); i++)
69 {
70 ULmodes.append((problem->liftfield[i]).clone());
71 }
72
73 for (int i = 0; i < NmodesUproj; i++)
74 {
75 ULmodes.append((problem->Umodes.toPtrList()[i]).clone());
76 }
77
78 for (int i = 0; i < NmodesSup; i++)
79 {
80 ULmodes.append((problem->supmodes.toPtrList()[i]).clone());
81 }
82
83 counter++;
84
85 if (NmodesUproj == 0)
86 {
87 UprojN = ULmodes.size();
88 }
89 else
90 {
91 UprojN = NmodesUproj + NmodesSup;
92 }
93
94 if (NmodesPproj == 0)
95 {
96 PprojN = problem->Pmodes.size();
97 }
98 else
99 {
100 PprojN = NmodesPproj;
101 }
102
103 if (NmodesNut == 0)
104 {
105 NmodesNut = problem->nutModes.size();
106 }
107
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);
116 a(0) = vel_now(0, 0);
117 float residualJumpLim =
118 problem->para->ITHACAdict->lookupOrDefault<float>("residualJumpLim", 1e-5);
119 float normalizedResidualLim =
120 problem->para->ITHACAdict->lookupOrDefault<float>("normalizedResidualLim",
121 1e-5);
122 maxIterOn = problem->para->ITHACAdict->lookupOrDefault<int>("maxIterOn",
123 1000);
124 scalar residual_jump(1 + residualJumpLim);
125 //problem->restart();
126 volScalarField& P = problem->_p();
127 volVectorField& U = problem->_U();
128 fvMesh& mesh = problem->_mesh();
129 Time& runTime = problem->_runTime();
130 P.rename("p");
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();
135 int iter = 0;
136 simpleControl& simple = problem->_simple();
137
139 {
140 Eigen::MatrixXd nutCoeff;
141 nutCoeff.resize(NmodesNut, 1);
142
143 for (int i = 0; i < NmodesNut; i++)
144 {
145 Eigen::MatrixXd muEval;
146 muEval.resize(1, 1);
147 muEval(0, 0) = mu_now;
148 nutCoeff(i, 0) = problem->rbfSplines[i]->eval(muEval);
149 }
150
151 volScalarField& nut = const_cast<volScalarField&>
152 (problem->_mesh().lookupObject<volScalarField>("nut"));
153 problem->nutModes.reconstruct(nut, nutCoeff, "nut");
155 }
156
157 PtrList<volVectorField> gradModP;
158
159 for (int i = 0; i < NmodesPproj; i++)
160 {
161 gradModP.append(fvc::grad(problem->Pmodes[i]));
162 }
163
164 projGradModP = ULmodes.project(gradModP, NmodesUproj);
165
166 while ((residual_jump > residualJumpLim
167 || std::max(U_norm_res, P_norm_res) > normalizedResidualLim)
168 && iter < maxIterOn)
169 {
170 iter++;
171 std::cout << "Iteration " << iter << std::endl;
172#if defined(OFVER) && (OFVER == 6)
173 simple.loop(runTime);
174#else
175 simple.loop();
176#endif
177 volScalarField nueff = problem->turbulence->nuEff().ref();
178 fvVectorMatrix UEqn
179 (
180 fvm::div(phi, U)
181 - fvm::laplacian(nueff, U)
182 - fvc::div(nueff * dev2(T(fvc::grad(U))))
183 );
184 UEqn.relax();
185 List<Eigen::MatrixXd> RedLinSysU = ULmodes.project(UEqn, UprojN);
186 RedLinSysU[1] = RedLinSysU[1] - projGradModP * b;
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);
194
195 if (simple.consistent())
196 {
197 rAtU = 1.0 / (1.0 / rAU - UEqn.H1());
198 phiHbyA +=
199 fvc::interpolate(rAtU() - rAU) * fvc::snGrad(P) * mesh.magSf();
200 HbyA -= (rAU - rAtU()) * fvc::grad(P);
201 }
202
203 List<Eigen::MatrixXd> RedLinSysP;
204
205 while (simple.correctNonOrthogonal())
206 {
207 fvScalarMatrix pEqn
208 (
209 fvm::laplacian(rAtU(), P) == fvc::div(phiHbyA)
210 );
211 RedLinSysP = problem->Pmodes.project(pEqn, PprojN);
212 b = reducedProblem::solveLinearSys(RedLinSysP, b, presidual);
213 problem->Pmodes.reconstruct(P, b, "p");
214
215 if (simple.finalNonOrthogonalIter())
216 {
217 phi = phiHbyA - pEqn.flux();
218 }
219 }
220
221 P.relax();
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();
235
236 if (problem->para->debug)
237 {
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;
241 }
242 }
243
244 std::cout << "Solution " << counter << " converged in " << iter <<
245 " iterations." << std::endl;
246 std::cout << "Final normalized residual for velocity: " << U_norm_res <<
247 std::endl;
248 std::cout << "Final normalized residual for pressure: " << P_norm_res <<
249 std::endl;
250 ULmodes.reconstruct(U, a, "Uaux");
251 P.rename("Paux");
252 problem->Pmodes.reconstruct(P, b, "Paux");
253 ITHACAstream::exportSolution(U, name(counter), Folder);
254 ITHACAstream::exportSolution(P, name(counter), Folder);
255 runTime.setTime(runTime.startTime(), 0);
256}
257
259{
260 M_Assert(problem->inletIndex.rows() == vel.size(),
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());
264
265 for (int k = 0; k < problem->inletIndex.rows(); k++)
266 {
267 int p = problem->inletIndex(k, 0);
268 int l = problem->inletIndex(k, 1);
269 scalar area = gSum(problem->liftfield[0].mesh().magSf().boundaryField()[p]);
270 scalar u_lf = gSum(problem->liftfield[k].mesh().magSf().boundaryField()[p] *
271 problem->liftfield[k].boundaryField()[p]).component(l) / area;
272 vel_scal(k, 0) = vel(k, 0) / u_lf;
273 }
274
275 vel_now = vel_scal;
276}
fvScalarMatrix pEqn
Definition CFM.H:31
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
#define M_Assert(Expr, Msg)
scalar uresidual
scalar presidual
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...
Definition Modes.C:275
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.
Definition Modes.C:63
PtrList< GeometricField< Type, PatchField, GeoMesh > > & toPtrList()
Function that returns the Modes object as a standard PtrList.
Definition Modes.H:86
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.
Definition steadyNS.H:293
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:284
volVectorModes supmodes
List of pointers used to form the supremizer modes.
Definition steadyNS.H:104
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:290
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition steadyNS.H:101
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:281
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
ITHACAparameters * para
Definition steadyNS.H:82
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition steadyNS.H:299
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
volScalarField & T
Definition createT.H:46
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)
surfaceScalarField & phi
volVectorField & U
volScalarField & nut
fvVectorMatrix & UEqn
Definition UEqn.H:37
volScalarField & p
adjustPhi(phiHbyA, U, p)
tmp< volScalarField > rAtU(rAU)
surfaceScalarField phiHbyA("phiHbyA", fvc::flux(HbyA))
volVectorField HbyA(constrainHbyA(rAU *UEqn.H(), U, p))
label i
Definition pEqn.H:46
volScalarField rAU(1.0/UEqn.A())