Loading...
Searching...
No Matches
ReducedUnsteadyNST.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2v/*---------------------------------------------------------------------------*\
3 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
4 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
5 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
6 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
7 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
8 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
9
10 * In real Time Highly Advanced Computational Applications for Finite Volumes
11 * Copyright (C) 2017 by the ITHACA-FV authors
12-------------------------------------------------------------------------------
13
14License
15 This file is part of ITHACA-FV
16
17 ITHACA-FV is free software: you can redistribute it and/or modify
18 it under the terms of the GNU Lesser General Public License as published by
19 the Free Software Foundation, either version 3 of the License, or
20 (at your option) any later version.
21
22 ITHACA-FV is distributed in the hope that it will be useful,
23 but WITHOUT ANY WARRANTY; without even the implied warranty of
24 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25 GNU Lesser General Public License for more details.
26
27 You should have received a copy of the GNU Lesser General Public License
28 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
29
30\*---------------------------------------------------------------------------*/
31
34
35
36#include "ReducedUnsteadyNST.H"
37
38
39// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
40
41// Constructor
45
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
56 // Create locally the velocity modes
57 for (int k = 0; k < problem->liftfield.size(); k++)
58 {
59 Umodes.append((problem->liftfield[k]).clone());
60 }
61
62 for (int k = 0; k < problem->NUmodes; k++)
63 {
64 Umodes.append((problem->Umodes[k]).clone());
65 }
66
67 for (int k = 0; k < problem->NSUPmodes; k++)
68 {
69 Umodes.append((problem->supmodes[k]).clone());
70 }
71
72 // Create locally the pressure modes
73 for (int k = 0; k < problem->NPmodes; k++)
74 {
75 Pmodes.append((problem->Pmodes[k]).clone());
76 }
77
78 for (int k = 0; k < problem->liftfieldT.size(); k++)
79 {
80 Tmodes.append((problem->liftfieldT[k]).clone());
81 }
82
83 // Create locally the temperature modes
84 for (int k = 0; k < problem->NTmodes; k++)
85 {
86 Tmodes.append((problem->Tmodes[k]).clone());
87 }
88
89 // Store locally the snapshots for projections
90 for (int k = 0; k < problem->Ufield.size(); k++)
91 {
92 Usnapshots.append((problem->Ufield[k]).clone());
93 Psnapshots.append((problem->Pfield[k]).clone());
94 }
95
96 for (int k = 0; k < problem->Tfield.size(); k++)
97 {
98 Tsnapshots.append((problem->Tfield[k]).clone());
99 }
100
102 FOMproblem);
104}
105
106
107// * * * * * * * * * * * * * * * Operators supremizer * * * * * * * * * * * * * //
108//Operator to evaluate the residual for the supremizer approach
109int newton_unsteadyNST_sup::operator()(const Eigen::VectorXd& x,
110 Eigen::VectorXd& fvec) const
111{
112 Eigen::VectorXd a_dot(Nphi_u);
113 Eigen::VectorXd a_tmp(Nphi_u);
114 Eigen::VectorXd b_tmp(Nphi_p);
115 a_tmp = x.head(Nphi_u);
116 b_tmp = x.tail(Nphi_p);
117 a_dot = (x.head(Nphi_u) - y_old.head(Nphi_u)) / dt;
118 // Convective term
119 Eigen::MatrixXd cc(1, 1);
120 // Momentum Term
121 Eigen::VectorXd M1 = problem->B_matrix * a_tmp * nu;
122 // Gradient of pressure
123 Eigen::VectorXd M2 = problem->K_matrix * b_tmp;
124 // Mass Term Velocity
125 Eigen::VectorXd M5 = problem->M_matrix * a_dot;
126 // Pressure Term
127 Eigen::VectorXd M3 = problem->P_matrix * a_tmp;
128
129 for (int i = 0; i < Nphi_u; i++)
130 {
131 cc = a_tmp.transpose() * problem->C_matrix[i] * a_tmp;
132 fvec(i) = - M5(i) + M1(i) - cc(0, 0) - M2(i);
133 }
134
135 for (int j = 0; j < Nphi_p; j++)
136 {
137 int k = j + Nphi_u;
138 fvec(k) = M3(j);
139 }
140
141 for (int j = 0; j < N_BC; j++)
142 {
143 fvec(j) = x(j) - BC(j);
144 }
145
146 return 0;
147}
148
149// Operator to evaluate the Jacobian for the supremizer approach
150int newton_unsteadyNST_sup::df(const Eigen::VectorXd& x,
151 Eigen::MatrixXd& fjac) const
152{
153 Eigen::NumericalDiff<newton_unsteadyNST_sup> numDiff(*this);
154 numDiff.df(x, fjac);
155 return 0;
156}
157
158int newton_unsteadyNST_sup_t::operator()(const Eigen::VectorXd& t,
159 Eigen::VectorXd& fvect) const
160{
161 Eigen::VectorXd c_dot(Nphi_t);
162 Eigen::VectorXd c_tmp(Nphi_t);
163 c_tmp = t.head(Nphi_t);
164 c_dot = (t.head(Nphi_t) - z_old.head(Nphi_t)) / dt;
165 // Convective term temperature
166 Eigen::MatrixXd qq(1, 1);
167 // diffusive term temperature
168 Eigen::VectorXd M6 = problem->Y_matrix * c_tmp * DT;
169 // Mass Term Temperature
170 Eigen::VectorXd M8 = problem->MT_matrix * c_dot;
171
172 for (int i = 0; i < Nphi_t; i++)
173 {
174 qq = a_tmp.transpose() * problem-> Q_matrix[i] * c_tmp;
175 fvect(i) = -M8(i) + M6(i) - qq(0, 0);
176 }
177
178 for (int j = 0; j < N_BC_t; j++)
179 {
180 fvect(j) = t(j) - BC_t(j);
181 }
182
183 return 0;
184}
185int newton_unsteadyNST_sup_t::df(const Eigen::VectorXd& t,
186 Eigen::MatrixXd& fjact) const
187{
188 Eigen::NumericalDiff<newton_unsteadyNST_sup_t> numDiff(*this);
189 numDiff.df(t, fjact);
190 return 0;
191}
192
193// * * * * * * * * * * * * * * * Solve Functions * * * * * * * * * * * * * //
194void reducedUnsteadyNST::solveOnline_sup(Eigen::MatrixXd& vel_now,
195 Eigen::MatrixXd& temp_now, int startSnap)
196{
197 // Create and resize the solution vector
198 y.resize(Nphi_u + Nphi_p, 1);
199 y.setZero();
200 z.resize(Nphi_t, 1);
201 z.setZero();
202 // Set Initial Conditions
206
207 // Change initial condition for the lifting function
208 for (int j = 0; j < N_BC; j++)
209 {
210 y(j) = vel_now(j, 0);
211 }
212
213 for (int j = 0; j < N_BC_t; j++)
214 {
215 z(j) = temp_now(j, 0);
216 }
217
218 // Set some properties of the newton object
225 newton_object_sup.BC.resize(N_BC);
227
228 for (int j = 0; j < N_BC_t; j++)
229 {
230 newton_object_sup_t.BC_t(j) = temp_now(j, 0);
231 }
232
233 for (int j = 0; j < N_BC; j++)
234 {
235 newton_object_sup.BC(j) = vel_now(j, 0);
236 }
237
238 // Set number of online solutions
239 int Ntsteps = static_cast<int>((finalTime - tstart) / dt);
240 online_solution.resize(Ntsteps);
241 online_solutiont.resize(Ntsteps);
242 // Set the initial time
243 time = tstart;
244 // Counting variable
245 int counter = 0;
246 // Create vector to store temporal solution and save initial condition as first solution
247 Eigen::MatrixXd tmp_sol(Nphi_u + Nphi_p + 1, 1);
248 tmp_sol(0) = time;
249 tmp_sol.col(0).tail(y.rows()) = y;
250 online_solution[counter] = tmp_sol;
251 Eigen::MatrixXd tmp_solt(Nphi_t + 1, 1);
252 tmp_solt(0) = time;
253 tmp_solt.col(0).tail(z.rows()) = z;
254 online_solutiont[counter] = tmp_solt;
255 counter ++;
256 // Create nonlinear solver object
257 Eigen::HybridNonLinearSolver<newton_unsteadyNST_sup> hnls(newton_object_sup);
258 Eigen::HybridNonLinearSolver<newton_unsteadyNST_sup_t> hnlst(
260 // Set output colors for fancy output
264
265 // Start the time loop
266 while (time < finalTime)
267 {
268 time = time + dt;
269 Eigen::VectorXd res(y);
270 Eigen::VectorXd rest(z);
271 res.setZero();
272 rest.setZero();
273 hnls.solve(y);
274
275 for (int j = 0; j < N_BC; j++)
276 {
277 y(j) = vel_now(j, 0);
278 }
279
280 // set the a_temp
281 // solve for temperature
283 hnlst.solve(z);
284
285 for (int j = 0; j < N_BC_t; j++)
286 {
287 z(j) = temp_now(j, 0);
288 }
289
290 newton_object_sup.operator()(y, res);
292 newton_object_sup_t.operator()(z, rest);
294 std::cout << "################## Online solve N° " << count_online_solve <<
295 " ##################" << std::endl;
296 Info << "Time = " << time << endl;
297 std::cout << "Solving for the parameter: " << vel_now << std::endl;
298 std::cout << "Solving for the parameter: " << temp_now << std::endl;
299
300 if (res.norm() < 1e-5)
301 {
302 std::cout << green << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
303 hnls.iter << " iterations " << def << std::endl << std::endl;
304 }
305 else
306 {
307 std::cout << red << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
308 hnls.iter << " iterations " << def << std::endl << std::endl;
309 }
310
311 if (rest.norm() < 1e-5)
312 {
313 std::cout << green << "|F(x)| = " << rest.norm() << " - Minimun reached in " <<
314 hnlst.iter << " iterations " << def << std::endl << std::endl;
315 }
316 else
317 {
318 std::cout << red << "|F(x)| = " << rest.norm() << " - Minimun reached in " <<
319 hnlst.iter << " iterations " << def << std::endl << std::endl;
320 }
321
323 tmp_sol(0) = time;
324 tmp_sol.col(0).tail(y.rows()) = y;
325
326 if (counter >= online_solution.size())
327 {
328 online_solution.append(tmp_sol);
329 }
330 else
331 {
332 online_solution[counter] = tmp_sol;
333 }
334
335 tmp_solt(0) = time;
336 tmp_solt.col(0).tail(z.rows()) = z;
337
338 if (counter >= online_solutiont.size())
339 {
340 online_solutiont.append(tmp_solt);
341 }
342 else
343 {
344 online_solutiont[counter] = tmp_solt;
345 }
346
347 counter ++;
348 }
349
350 // Save the solution
351 ITHACAstream::exportMatrix(online_solution, "red_coeff", "python",
352 "./ITHACAoutput/red_coeff");
353 ITHACAstream::exportMatrix(online_solution, "red_coeff", "matlab",
354 "./ITHACAoutput/red_coeff");
355 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "python",
356 "./ITHACAoutput/red_coeff_t");
357 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "matlab",
358 "./ITHACAoutput/red_coeff_t");
360}
361
362
363void reducedUnsteadyNST::reconstruct_sup(fileName folder, int printevery)
364{
365 mkDir(folder);
367 int counter = 0;
368 int nextwrite = 0;
369 int counter2 = 1;
370
371 for (int i = 0; i < online_solution.size(); i++)
372 {
373 if (counter == nextwrite)
374 {
375 volVectorField U_rec("U_rec", Umodes[0] * 0);
376
377 for (int j = 0; j < Nphi_u; j++)
378 {
379 U_rec += Umodes[j] * online_solution[i](j + 1, 0);
380 }
381
382 //problem.exportSolution(U_rec, name(online_solution[i](0, 0)), folder);
383 ITHACAstream::exportSolution(U_rec, name(counter2), folder);
384 volScalarField P_rec("P_rec", Pmodes[0] * 0);
385
386 for (int j = 0; j < Nphi_p; j++)
387 {
388 P_rec += Pmodes[j] * online_solution[i](j + Nphi_u + 1, 0);
389 }
390
391 //problem.exportSolution(P_rec, name(online_solution[i](0, 0)), folder);
392 ITHACAstream::exportSolution(P_rec, name(counter2), folder);
393 nextwrite += printevery;
394 counter2 ++;
395 UREC.append((U_rec).clone());
396 PREC.append((P_rec).clone());
397 }
398
399 counter++;
400 }
401}
402
403
404
405void reducedUnsteadyNST::reconstruct_supt(fileName folder, int printevery)
406{
407 mkDir(folder);
409 int counter = 0;
410 int nextwrite = 0;
411 int counter2 = 1;
412
413 for (int i = 0; i < online_solutiont.size(); i++)
414 {
415 if (counter == nextwrite)
416 {
417 volScalarField T_rec("T_rec", Tmodes[0] * 0);
418
419 for (int j = 0; j < Nphi_t; j++)
420 {
421 T_rec += Tmodes[j] * online_solutiont[i](j + 1, 0);
422 }
423
424 ITHACAstream::exportSolution(T_rec, name(counter2), folder);
425 nextwrite += printevery;
426 counter2 ++;
427 TREC.append((T_rec).clone());
428 }
429
430 counter++;
431 }
432}
433
434// ************************************************************************* //
Header file of the reducedUnsteadyNST class.
Class to change color to the output stream.
Definition colormod.H:24
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.
scalar finalTime
Scalar to store the final time if the online simulation.
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.
void reconstruct_supt(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique for t...
PtrList< volScalarField > T_rec
Reconstructed temperature field.
reducedUnsteadyNST()
Construct Null.
List< Eigen::MatrixXd > online_solutiont
List of Eigen matrices to store the online solution for temperature equation.
newton_unsteadyNST_sup_t newton_object_sup_t
Functor object to call the non linear solver sup. approach for temperature case.
scalar tstart
Scalar to store the final time if the online simulation.
unsteadyNST * problem
pointer to the FOM problem
int Nphi_t
Number of temperature modes.
void reconstruct_sup(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique.
scalar dt
Scalar to store the time increment.
newton_unsteadyNST_sup newton_object_sup
Functor object to call the non linear solver sup. approach for velocity case.
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.
scalar time
Scalar to store the current time.
void solveOnline_sup(Eigen::MatrixXd &vel_now, Eigen::MatrixXd &temp_now, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
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
Implementation of a parametrized full order unsteady NS problem weakly coupled with the energy equat...
Definition unsteadyNST.H:61
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
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 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.
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
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
Newton object for the resolution of the reduced problem using a supremizer approach.
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const