Loading...
Searching...
No Matches
ReducedUnsteadyNSTTurb.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
34
36
37
38// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
39
40// Constructor
44
46
47//problem(&FOMproblem)
48{
49 problem = &FOMproblem;
50 N_BC = problem->inletIndex.rows();
51 N_BC_t = problem->inletIndexT.rows();
52 Nphi_u = problem->B_matrix.rows();
53 Nphi_p = problem->K_matrix.cols();
54 Nphi_t = problem->Y_matrix.rows();
55 Nphi_nut = problem->CT2_matrix[0].rows();
56
57 // Create locally the velocity modes
58 for (int k = 0; k < problem->liftfield.size(); k++)
59 {
60 Umodes.append((problem->liftfield[k]).clone());
61 }
62
63 for (int k = 0; k < problem->NUmodes; k++)
64 {
65 Umodes.append((problem->Umodes[k]).clone());
66 }
67
68 for (int k = 0; k < problem->NSUPmodes; k++)
69 {
70 Umodes.append((problem->supmodes[k]).clone());
71 }
72
73 // Create locally the pressure modes
74 for (int k = 0; k < problem->NPmodes; k++)
75 {
76 Pmodes.append((problem->Pmodes[k]).clone());
77 }
78
79 for (int k = 0; k < problem->liftfieldT.size(); k++)
80 {
81 Tmodes.append((problem->liftfieldT[k]).clone());
82 }
83
84 // Create locally the temperature modes
85 for (int k = 0; k < problem->NTmodes; k++)
86 {
87 Tmodes.append((problem->Tmodes[k]).clone());
88 }
89
90 // Store locally the snapshots for projections
91 for (int k = 0; k < problem->Ufield.size(); k++)
92 {
93 Usnapshots.append((problem->Ufield[k]).clone());
94 Psnapshots.append((problem->Pfield[k]).clone());
95 }
96
97 for (int k = 0; k < problem->Tfield.size(); k++)
98 {
99 Tsnapshots.append((problem->Tfield[k]).clone());
100 }
101
103 FOMproblem);
105}
106
107// * * * * * * * * * * * * * * * Operators supremizer * * * * * * * * * * * * * //
108
109// Operator to evaluate the residual for the supremizer approach
110int newton_unsteadyNSTTurb_sup::operator()(const Eigen::VectorXd& x,
111 Eigen::VectorXd& fvec) const
112{
113 Eigen::VectorXd a_dot(Nphi_u);
114 Eigen::VectorXd a_tmp(Nphi_u);
115 Eigen::VectorXd b_tmp(Nphi_p);
116 a_tmp = x.head(Nphi_u);
117 b_tmp = x.tail(Nphi_p);
118 a_dot = (x.head(Nphi_u) - y_old.head(Nphi_u)) / dt;
119 // Convective term
120 Eigen::MatrixXd cc(1, 1);
121 // Mom Term
122 Eigen::VectorXd M1 = problem ->B_total_matrix * a_tmp * nu;
123 // Gradient of pressure
124 Eigen::VectorXd M2 = problem->K_matrix * b_tmp;
125 // Mass Term
126 Eigen::VectorXd M5 = problem->M_matrix * a_dot;
127 // Pressure Term
128 Eigen::VectorXd M3 = problem->P_matrix * a_tmp;
129
130 //std::cerr << "I am here 5" << std::endl;
131 for (int i = 0; i < Nphi_u; i++)
132 {
133 cc = a_tmp.transpose() * problem->C_matrix[i] * a_tmp - nu_c.transpose() *
134 problem->C_total_matrix[i] * a_tmp;
135 fvec(i) = - M5(i) + M1(i) - cc(0, 0) - M2(i);
136 //Info << "Non-turb part is " << a_tmp.transpose() * C_matrix[i] * a_tmp << endl; //Info << "Turb part is " << nu_c.transpose() * C_total_matrix[i] * a_tmp << endl
137 }
138
139 for (int j = 0; j < Nphi_p; j++)
140 {
141 int k = j + Nphi_u;
142 fvec(k) = M3(j);
143 }
144
145 for (int j = 0; j < N_BC; j++)
146 {
147 fvec(j) = x(j) - BC(j);
148 }
149
150 return 0;
151}
152
153// Operator to evaluate the Jacobian for the supremizer approach
154int newton_unsteadyNSTTurb_sup::df(const Eigen::VectorXd& x,
155 Eigen::MatrixXd& fjac) const
156{
157 Eigen::NumericalDiff<newton_unsteadyNSTTurb_sup> numDiff(*this);
158 numDiff.df(x, fjac);
159 return 0;
160}
161
162
163int newton_unsteadyNSTTurb_sup_t::operator()(const Eigen::VectorXd& t,
164 Eigen::VectorXd& fvect) const
165{
166 // Eigen::VectorXd a_tmp(Nphi_u);
167 Eigen::VectorXd c_dot(Nphi_t);
168 Eigen::VectorXd c_tmp(Nphi_t);
169 c_tmp = t.head(Nphi_t);
170 c_dot = (t.head(Nphi_t) - z_old.head(Nphi_t)) / dt;
171 // Convective term temperature
172 Eigen::MatrixXd qq(1, 1);
173 Eigen::MatrixXd st(1, 1);
174 // diffusive term temperature
175 Eigen::VectorXd M6 = problem->Y_matrix * c_tmp * nu / Pr;
176 // Mass Term Temperature
177 Eigen::VectorXd M8 = problem->MT_matrix * c_dot;
178
179 for (int i = 0; i < Nphi_t; i++)
180 {
181 qq = a_tmp.transpose() * problem->Q_matrix[i] * c_tmp;
182 st = nu_c.transpose() * problem->S_matrix[i] * c_tmp;
183 fvect(i) = -M8(i) + M6(i) - qq(0, 0) + st(0, 0) / Prt;
184 }
185
186 for (int j = 0; j < N_BC_t; j++)
187 {
188 fvect(j) = t(j) - BC_t(j);
189 }
190
191 return 0;
192}
193int newton_unsteadyNSTTurb_sup_t::df(const Eigen::VectorXd& t,
194 Eigen::MatrixXd& fjact) const
195{
196 Eigen::NumericalDiff<newton_unsteadyNSTTurb_sup_t> numDiff(*this);
197 numDiff.df(t, fjact);
198 return 0;
199}
200
201
202// * * * * * * * * * * * * * * * Solve Functions * * * * * * * * * * * * * //
203void ReducedUnsteadyNSTTurb::solveOnlineSup(Eigen::MatrixXd& vel_now,
204 Eigen::MatrixXd& temp_now, int startSnap)
205{
206 // Create and resize the solution vector
207 y.resize(Nphi_u + Nphi_p, 1);
208 y.setZero();
209 z.resize(Nphi_t, 1);
210 z.setZero();
211 volScalarField T_IC("T_IC", problem->Tfield[0]);
212
213 for (int j = 0; j < T_IC.boundaryField().size(); j++)
214 {
215 for (int i = 0; i < N_BC_t; i++)
216 {
217 if (j == problem->inletIndexT(i, 0))
218 {
219 T_IC.boundaryFieldRef()[problem->inletIndexT(i, 0)][j] = temp_now(i, 0);
220 }
221 else
222 {
223 }
224 }
225 }
226
230
231 // Change initial condition for the lifting function
232 for (int j = 0; j < N_BC; j++)
233 {
234 y(j) = vel_now(j, 0);
235 }
236
237 for (int j = 0; j < N_BC_t; j++)
238 {
239 z(j) = temp_now(j, 0);
240 }
241
242 // Set some properties of the newton object
251 newton_object_sup.BC.resize(N_BC);
253
254 for (int j = 0; j < N_BC; j++)
255 {
256 newton_object_sup.BC(j) = vel_now(j, 0);
257 }
258
259 for (int j = 0; j < N_BC_t; j++)
260 {
261 newton_object_sup_t.BC_t(j) = temp_now(j, 0);
262 }
263
264 // Set number of online solutions
265 int Ntsteps = static_cast<int>((finalTime - tstart) / dt);
266 online_solution.resize(Ntsteps);
267 online_solutiont.resize(Ntsteps);
268 // Set the initial time
269 time = tstart;
270 // Counting variable
271 int counter = 0;
272 // Create vector to store temporal solution and save initial condition as first solution
273 Eigen::MatrixXd tmp_sol(Nphi_u + Nphi_p + 1, 1);
274 tmp_sol(0) = time;
275 tmp_sol.col(0).tail(y.rows()) = y;
276 online_solution[counter] = tmp_sol;
277 Eigen::MatrixXd tmp_solt(Nphi_t + 1, 1);
278 tmp_solt(0) = time;
279 tmp_solt.col(0).tail(z.rows()) = z;
280 online_solutiont[counter] = tmp_solt;
281 counter ++;
282 Eigen::HybridNonLinearSolver<newton_unsteadyNSTTurb_sup> hnls(
284 Eigen::HybridNonLinearSolver<newton_unsteadyNSTTurb_sup_t> hnlst(
286 // Set output colors for fancy output
290
291 // Start the time loop
292 while (time < finalTime)
293 {
294 time = time + dt;
295 std::vector<double> tv;
296 tv.resize(1 + vel_now.size());
297 tv[0] = time;
298
299 for (int i = 1; i < tv.size(); i++)
300 {
301 tv[i] = vel_now(i - 1);
302 }
303
304 for (int i = 0; i < Nphi_nut; i++)
305 {
307 }
308
309 volScalarField nut_rec("nut_rec", problem->nuTmodes[0] * 0);
310
311 for (int j = 0; j < Nphi_nut; j++)
312 {
313 nut_rec += problem->nuTmodes[j] * newton_object_sup.nu_c(j);
314 }
315
316 nutREC.append((nut_rec).clone());
318 Eigen::VectorXd res(y);
319 Eigen::VectorXd rest(z);
320 res.setZero();
321 rest.setZero();
322 hnls.solve(y);
323
324 for (int j = 0; j < N_BC; j++)
325 {
326 y(j) = vel_now(j, 0);
327 }
328
330 hnlst.solve(z);
331
332 for (int j = 0; j < N_BC_t; j++)
333 {
334 z(j) = temp_now(j, 0);
335 }
336
337 newton_object_sup.operator()(y, res);
339 newton_object_sup_t.operator()(z, rest);
341 std::cout << "################## Online solve N° " << count_online_solve <<
342 " ##################" << std::endl;
343 Info << "Time = " << time << endl;
344 std::cout << "Solving for the parameter: " << vel_now << std::endl;
345
346 if (res.norm() < 1e-5)
347 {
348 std::cout << green << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
349 hnls.iter << " iterations " << def << std::endl << std::endl;
350 }
351 else
352 {
353 std::cout << red << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
354 hnls.iter << " iterations " << def << std::endl << std::endl;
355 }
356
358 tmp_sol(0) = time;
359 tmp_sol.col(0).tail(y.rows()) = y;
360
361 if (counter >= online_solution.size())
362 {
363 online_solution.append(tmp_sol);
364 }
365 else
366 {
367 online_solution[counter] = tmp_sol;
368 }
369
370 tmp_solt(0) = time;
371 tmp_solt.col(0).tail(z.rows()) = z;
372
373 if (counter >= online_solutiont.size())
374 {
375 online_solutiont.append(tmp_solt);
376 }
377 else
378 {
379 online_solutiont[counter] = tmp_solt;
380 }
381
382 counter ++;
383 }
384
385 // Save the solution
386 ITHACAstream::exportMatrix(online_solution, "red_coeff", "python",
387 "./ITHACAoutput/red_coeff");
388 ITHACAstream::exportMatrix(online_solution, "red_coeff", "matlab",
389 "./ITHACAoutput/red_coeff");
390 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "python",
391 "./ITHACAoutput/red_coeff_t");
392 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "matlab",
393 "./ITHACAoutput/red_coeff_t");
394 ITHACAstream::exportFields(nutREC, "nutfield", "./ITHACAoutput/nutfield/");
396}
397
398// * * * * * * * * * * * * * * * Solve Functions * * * * * * * * * * * * * //
399
400void ReducedUnsteadyNSTTurb::reconstructSup(fileName folder, int printevery)
401{
402 mkDir(folder);
403 system("ln -s ../../constant " + folder + "/constant");
404 system("ln -s ../../0 " + folder + "/0");
405 system("ln -s ../../system " + folder + "/system");
406 int counter = 0;
407 int nextwrite = 0;
408 int counter2 = 1;
409
410 for (int i = 0; i < online_solution.size(); i++)
411 {
412 if (counter == nextwrite)
413 {
414 volVectorField U_rec("U_rec", Umodes[0] * 0);
415
416 for (int j = 0; j < Nphi_u; j++)
417 {
418 U_rec += Umodes[j] * online_solution[i](j + 1, 0);
419 }
420
421 ITHACAstream::exportSolution(U_rec, name(counter2), folder);
422 volScalarField P_rec("P_rec", Pmodes[0] * 0);
423
424 for (int j = 0; j < Nphi_p; j++)
425 {
426 P_rec += Pmodes[j] * online_solution[i](j + Nphi_u + 1, 0);
427 }
428
429 ITHACAstream::exportSolution(P_rec, name(counter2), folder);
430 ITHACAstream::exportSolution(nutREC[nextwrite], name(counter2), folder);
431 nextwrite += printevery;
432 counter2 ++;
433 UREC.append((U_rec).clone());
434 PREC.append((P_rec).clone());
435 }
436
437 counter++;
438 }
439}
440
441void ReducedUnsteadyNSTTurb::reconstructSupt(fileName folder, int printevery)
442{
443 mkDir(folder);
444 system("ln -s ../../constant " + folder + "/constant");
445 system("ln -s ../../0 " + folder + "/0");
446 system("ln -s ../../system " + folder + "/system");
447 int counter = 0;
448 int nextwrite = 0;
449 int counter2 = 1;
450
451 for (int i = 0; i < online_solutiont.size(); i++)
452 {
453 if (counter == nextwrite)
454 {
455 volScalarField T_rec("T_rec", Tmodes[0] * 0);
456
457 for (int j = 0; j < Nphi_t; j++)
458 {
459 T_rec += Tmodes[j] * online_solutiont[i](j + 1, 0);
460 }
461
462 ITHACAstream::exportSolution(T_rec, name(counter2), folder);
463 nextwrite += printevery;
464 counter2 ++;
465 TREC.append((T_rec).clone());
466 }
467
468 counter++;
469 }
470}
471// ************************************************************************* //
Header file of the reducedUnsteadyNSTTurb class.
Class to change color to the output stream.
Definition colormod.H:24
ReducedUnsteadyNSTTurb()
Construct Null.
void solveOnlineSup(Eigen::MatrixXd &vel_now, Eigen::MatrixXd &temp_now, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
double time
Scalar to store the current time.
newton_unsteadyNSTTurb_sup_t newton_object_sup_t
Functor object to call the non linear solver sup. approach for temperature field.
double dt
Scalar to store the time increment.
newton_unsteadyNSTTurb_sup newton_object_sup
Functor object to call the non linear solver sup. approach.
scalar Prt
Scalar to store the turbulent Prandtl number.
int Nphi_nut
Number of nut field modes.
UnsteadyNSTTurb * problem
Pointer to the FOM problem.
void reconstructSup(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution for velocity and pressure from an online solve with a supremizer sta...
PtrList< volScalarField > nutREC
Reconstructed eddy viscosity field.
double finalTime
Scalar to store the final time if the online simulation.
void reconstructSupt(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution for temperature from an online solve with a supremizer stabilisation...
scalar Pr
Scalar to store the Prandtl number.
double tstart
Scalar to store the final time if the online simulation.
Implementation of a parametrized full order unsteady NST problem weakly coupled with the energy equa...
List< Eigen::MatrixXd > S_matrix
Turbulent diffusivity term.
List< Eigen::MatrixXd > C_total_matrix
Total C Matrix.
List< Eigen::MatrixXd > CT2_matrix
Turbulent viscosity term.
Eigen::MatrixXd B_total_matrix
Total B Matrix.
std::vector< SPLINTER::RBFSpline * > rbfsplines
Create a SAMPLES for interpolation.
PtrList< volScalarField > nuTmodes
List of POD modes for eddy viscosity.
PtrList< volScalarField > PREC
Reconstructed pressure field.
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
PtrList< volScalarField > Pmodes
List of pointers to store the modes for pressure.
List< Eigen::MatrixXd > online_solution
List of Eigen matrices to store the online solution.
int Nphi_p
Number of pressure modes.
scalar nu
Reduced viscosity in case of parametrized viscosity.
int count_online_solve
Counter to count the online solutions.
Eigen::MatrixXd vel_now
Online inlet velocity vector.
PtrList< volScalarField > Psnapshots
List of pointers to store the snapshots for pressure.
PtrList< volVectorField > Usnapshots
List of pointers to store the snapshots for velocity.
PtrList< volVectorField > Umodes
List of pointers to store the modes for velocity.
int N_BC
Number of parametrized boundary conditions.
PtrList< volVectorField > UREC
Recontructed velocity field.
int Nphi_u
Number of velocity modes.
PtrList< volScalarField > TREC
Reconstructed temperature field.
int N_BC_t
Number of parametrized boundary conditions related to temperature field.
PtrList< volScalarField > Tmodes
List of pointers to store the modes for temperature.
PtrList< volScalarField > T_rec
Reconstructed temperature field.
List< Eigen::MatrixXd > online_solutiont
List of Eigen matrices to store the online solution for temperature equation.
int Nphi_t
Number of temperature modes.
Eigen::VectorXd z
Vector to store the temperature solution during the Newton procedure.
PtrList< volScalarField > Tsnapshots
List of pointers to store the snapshots for temperature.
Eigen::MatrixXi inletIndexT
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
label NPmodes
Number of pressure modes used for the projection.
Definition steadyNS.H:143
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
volVectorModes supmodes
List of pointers used to form the supremizer modes.
Definition steadyNS.H:104
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition steadyNS.H:101
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition steadyNS.H:89
label NUmodes
Number of velocity modes used for the projection.
Definition steadyNS.H:140
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
Eigen::MatrixXd B_matrix
Diffusion term.
Definition steadyNS.H:157
label NSUPmodes
Number of supremizer modes used for the projection.
Definition steadyNS.H:146
Eigen::MatrixXd K_matrix
Gradient of pressure matrix.
Definition steadyNS.H:163
List< Eigen::MatrixXd > C_matrix
Non linear term.
Definition steadyNS.H:166
Eigen::MatrixXd P_matrix
Div of velocity.
Definition steadyNS.H:170
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition steadyNS.H:160
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
Eigen::MatrixXd Y_matrix
Gradient of pressure matrix.
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
Definition unsteadyNST.H:99
List< Eigen::MatrixXd > Q_matrix
Non linear term.
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of the temperature lifting functions.
label NTmodes
Number of temperature modes used for the projection.
Eigen::MatrixXd MT_matrix
Mass Matrix T.
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
@ FG_GREEN
Definition colormod.H:14
@ FG_DEFAULT
Definition colormod.H:16
@ FG_RED
Definition colormod.H:13
void exportFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &field, word folder, word fieldname)
Function to export a scalar of vector 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.
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 ...
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.
label i
Definition pEqn.H:46
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const