Loading...
Searching...
No Matches
ReducedUnsteadyBB.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 "ReducedUnsteadyBB.H"
37
38
39// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
40
41// Constructor
45
47 :
48 problem(&FOMproblem)
49{
50 N_BC_t = problem->inletIndexT.rows();
51 N_BC = problem->inletIndex.rows();
52 Nphi_u = problem->B_matrix.rows();
53 Nphi_prgh = 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 LUmodes.append((problem->liftfield[k]).clone());
60 }
61
62 for (int k = 0; k < problem->NUmodes; k++)
63 {
64 LUmodes.append((problem->Umodes[k]).clone());
65 }
66
67 for (int k = 0; k < problem->NSUPmodes; k++)
68 {
69 LUmodes.append((problem->supmodes[k]).clone());
70 }
71
72 // Create locally the temperature modes including BC with liftfield
73 for (int k = 0; k < problem->liftfieldT.size(); k++)
74 {
75 LTmodes.append((problem->liftfieldT[k]).clone());
76 }
77
78 for (int k = 0; k < problem->NTmodes; k++)
79 {
80 LTmodes.append((problem->Tmodes[k]).clone());
81 }
82
85 FOMproblem);
88 FOMproblem);
89}
90
91
92// * * * * * * * * * * * * * * * Operators supremizer * * * * * * * * * * * * * //
93//Operator to evaluate the residual for the supremizer approach
94int newton_unsteadyBB_sup::operator()(const Eigen::VectorXd& x,
95 Eigen::VectorXd& fvec) const
96{
97 Eigen::VectorXd a_dot(Nphi_u);
98 Eigen::VectorXd a_tmp(Nphi_u);
99 Eigen::VectorXd b_tmp(Nphi_prgh);
100 Eigen::VectorXd c_dot(Nphi_t);
101 Eigen::VectorXd c_tmp(Nphi_t);
102 a_tmp = x.head(Nphi_u);
103 a_dot = (x.head(Nphi_u) - y_old.head(Nphi_u)) / dt;
104 b_tmp = x.segment(Nphi_u, Nphi_prgh);
105 c_tmp = x.tail(Nphi_t);
106 c_dot = (x.tail(Nphi_t) - y_old.tail(Nphi_t)) / dt;
107 // Convective term
108 Eigen::MatrixXd cc(1, 1);
109 // Diffusive Term
110 Eigen::VectorXd M1 = problem->B_matrix * a_tmp * nu;
111 // Mass Term Velocity
112 Eigen::VectorXd M5 = problem->M_matrix * a_dot;
113 // Gradient of pressure
114 Eigen::VectorXd M2 = problem->K_matrix * b_tmp;
115 // Continuity
116 Eigen::VectorXd M3 = problem->P_matrix * a_tmp;
117 // Buoyancy Term
118 Eigen::VectorXd M10 = problem->H_matrix * c_tmp;
119 // Convective term temperature
120 Eigen::MatrixXd qq(1, 1);
121 // diffusive term temperature
122 Eigen::VectorXd M6 = problem->Y_matrix * c_tmp * (nu / Pr);
123 // Mass Term Temperature
124 Eigen::VectorXd M8 = problem->W_matrix * c_dot;
125
126 for (int i = 0; i < Nphi_u; i++)
127 {
128 cc = a_tmp.transpose() * Eigen::SliceFromTensor(problem->C_tensor, 0,
129 i) * a_tmp;
130 fvec(i) = - M5(i) + M1(i) - cc(0, 0) - M10(i) - M2(i);
131 }
132
133 for (int j = 0; j < Nphi_prgh; j++)
134 {
135 int k = j + Nphi_u;
136 fvec(k) = M3(j);
137 }
138
139 for (int j = 0; j < Nphi_t; j++)
140 {
141 int k = j + Nphi_u + Nphi_prgh;
142 qq = a_tmp.transpose() * problem->Q_matrix[j] * c_tmp;
143 fvec(k) = -M8(j) + M6(j) - qq(0, 0);
144 }
145
146 for (int j = 0; j < N_BC; j++)
147 {
148 fvec(j) = x(j) - BC(j);
149 }
150
151 for (int j = 0; j < N_BC_t; j++)
152 {
153 int k = j + Nphi_u + Nphi_prgh;
154 fvec(k) = x(k) - BC_t(j);
155 }
156
157 return 0;
158}
159
160// Operator to evaluate the Jacobian for the supremizer approach
161int newton_unsteadyBB_sup::df(const Eigen::VectorXd& x,
162 Eigen::MatrixXd& fjac) const
163{
164 Eigen::NumericalDiff<newton_unsteadyBB_sup> numDiff(*this);
165 numDiff.df(x, fjac);
166 return 0;
167}
168
169// * * * * * * * * * * * * * * * Operators PPE * * * * * * * * * * * * * * * //
170
171// Operator to evaluate the residual for the supremizer approach
172int newton_unsteadyBB_PPE::operator()(const Eigen::VectorXd& x,
173 Eigen::VectorXd& fvec) const
174{
175 Eigen::VectorXd a_dot(Nphi_u);
176 Eigen::VectorXd a_tmp(Nphi_u);
177 Eigen::VectorXd b_tmp(Nphi_prgh);
178 Eigen::VectorXd c_dot(Nphi_t);
179 Eigen::VectorXd c_tmp(Nphi_t);
180 a_tmp = x.head(Nphi_u);
181 b_tmp = x.segment(Nphi_u, Nphi_prgh);
182 a_dot = (x.head(Nphi_u) - y_old.head(Nphi_u)) / dt;
183 c_tmp = x.tail(Nphi_t);
184 c_dot = (x.tail(Nphi_t) - y_old.tail(Nphi_t)) / dt;
185 // Convective terms
186 Eigen::MatrixXd cc(1, 1);
187 Eigen::MatrixXd gg(1, 1);
188 Eigen::MatrixXd bb(1, 1);
189 // Convective term temperature
190 Eigen::MatrixXd qq(1, 1);
191 // Eigen::MatrixXd st(1, 1);
192 // Mom Term
193 Eigen::VectorXd M1 = problem->B_matrix * a_tmp * nu;
194 // Gradient of pressure
195 Eigen::VectorXd M2 = problem->K_matrix * b_tmp;
196 // Mass Term
197 Eigen::VectorXd M5 = problem->M_matrix * a_dot;
198 // Pressure Term
199 Eigen::VectorXd M3 = problem->D_matrix * b_tmp;
200 // BC PPE
201 Eigen::VectorXd M6 = problem->BC1_matrix * a_tmp * nu;
202 // BC PPE
203 // Buoyancy Term
204 Eigen::VectorXd M10 = problem->H_matrix * c_tmp;
205 Eigen::VectorXd M11 = problem->HP_matrix * c_tmp;
206 Eigen::VectorXd M7 = problem->BC3_matrix * a_tmp * nu;
207 // diffusive term temperature
208 Eigen::VectorXd M9 = problem->Y_matrix * c_tmp * (nu / Pr);
209 // Mass Term Temperature
210 Eigen::VectorXd M8 = problem->W_matrix * c_dot;
211
212 for (int i = 0; i < Nphi_u; i++)
213 {
214 cc = a_tmp.transpose() * Eigen::SliceFromTensor(problem->C_tensor, 0,
215 i) * a_tmp;
216 fvec(i) = - M5(i) + M1(i) - cc(0, 0) - M10(i) - M2(i);
217 }
218
219 for (int j = 0; j < Nphi_prgh; j++)
220 {
221 int k = j + Nphi_u;
222 gg = a_tmp.transpose() * problem->G_matrix[j] * a_tmp;
223 bb = a_tmp.transpose() * problem->BC2_matrix[j] * a_tmp;
224 //fvec(k) = M3(j, 0) - gg(0, 0) - M6(j, 0) + bb(0, 0);
225 fvec(k) = M3(j, 0) + gg(0, 0) + M11(j, 0) - M7(j, 0);
226 }
227
228 for (int j = 0; j < Nphi_t; j++)
229 {
230 int k = j + Nphi_u + Nphi_prgh;
231 qq = a_tmp.transpose() * problem->Q_matrix[j] * c_tmp;
232 fvec(k) = -M8(j) + M9(j) - qq(0, 0);
233 }
234
235 // for (int j = 0; j < N_BC; j++)
236 //{
237 // fvec(j) = x(j) - BC(j);
238 // }
239 return 0;
240}
241
242// Operator to evaluate the Jacobian for the supremizer approach
243int newton_unsteadyBB_PPE::df(const Eigen::VectorXd& x,
244 Eigen::MatrixXd& fjac) const
245{
246 Eigen::NumericalDiff<newton_unsteadyBB_PPE> numDiff(*this);
247 numDiff.df(x, fjac);
248 return 0;
249}
250
251// * * * * * * * * * * * * * * * Solve Functions * * * * * * * * * * * * * //
252Eigen::MatrixXd ReducedUnsteadyBB::solveOnline_sup(Eigen::MatrixXd& temp_now_BC,
253 Eigen::MatrixXd& vel_now_BC, int NParaSet, int startSnap)
254{
255 std::cout << "################## Online solve N° " << NParaSet <<
256 " ##################" << std::endl;
257 std::cout << "Solving for the parameter: " << temp_now_BC << std::endl;
258 // Count number of time steps
259 int counter = 0;
260 time = tstart;
261
262 while (time < finalTime - 0.5 * dt)
263 {
264 time = time + dt;
265 counter ++;
266 }
267
268 std::cerr << "File: ReducedUnsteadyBB.C, Line: 264" << std::endl;
269 // Set size of online solution
270 online_solutiont.resize(Nphi_u + Nphi_prgh + Nphi_t + 1, counter + 1);
271 // Set initial condition for online solve
272 volScalarField T_IC("T_IC", problem->Tfield[0].clone());
273 Info << T_IC.boundaryField() << endl;
274
275 for (int j = 0; j < T_IC.boundaryField().size(); j++)
276 {
277 for (int i = 0; i < N_BC_t; i++)
278 {
279 if (j == problem->inletIndexT(i, 0))
280 {
281 T_IC.boundaryFieldRef()[problem->inletIndexT(i, 0)][j] = temp_now_BC(i, 0);
282 std::cerr << "File: ReducedUnsteadyBB.C, Line: 277" << std::endl;
283 }
284 else
285 {
286 }
287 }
288 }
289
290 std::cerr << "File: ReducedUnsteadyBB.C, Line: 284" << std::endl;
291 // Create and resize the solution vector
292 y.resize(Nphi_u + Nphi_prgh + Nphi_t, 1);
293 y.setZero();
294 // Calculate the time-dependent coefficients
296 LUmodes);
297
298 if (Nphi_prgh != 0)
299 {
301 problem->Prghfield[2],
303 }
304
305 std::cerr << "File: ReducedUnsteadyBB.C, Line: 299" << std::endl;
307 std::cerr << "File: ReducedUnsteadyBB.C, Line: 302" << std::endl;
308 // Set some properties of the newton object
314 newton_object_sup.BC.resize(N_BC);
315
316 // Change initial condition for the lifting function
317 for (int j = 0; j < N_BC_t; j++)
318 {
319 newton_object_sup.BC_t(j) = temp_now_BC(j, 0);
320 }
321
322 for (int j = 0; j < N_BC; j++)
323 {
324 newton_object_sup.BC(j) = vel_now_BC(j, 0);
325 }
326
327 // Set the initial time
328 time = tstart;
329 // Create vector to store temporal solution and save initial condition as first solution
330 Eigen::MatrixXd tmp_sol(Nphi_u + Nphi_prgh + Nphi_t + 1, 1);
331 tmp_sol(0) = time;
332 tmp_sol.col(0).tail(y.rows()) = y;
333 online_solutiont.col(0) = tmp_sol;
334 // Create nonlinear solver object
335 Eigen::HybridNonLinearSolver<newton_unsteadyBB_sup> hnls(newton_object_sup);
336 // Set output colors for fancy output
340
341 // Start the time loop
342 for (int i = 1; i < online_solutiont.cols(); i++)
343 {
344 time = time + dt;
345 Eigen::VectorXd res(y);
346 res.setZero();
347 hnls.solve(y);
348
349 for (int j = 0; j < N_BC; j++)
350 {
351 y(j) = vel_now_BC(j, 0);
352 }
353
354 for (int j = 0; j < N_BC_t; j++)
355 {
356 int k = j + Nphi_prgh + Nphi_u;
357 y(k) = temp_now_BC(j, 0);
358 }
359
360 newton_object_sup.operator()(y, res);
362 Info << "Time = " << time << endl;
363
364 if (res.norm() < 1e-5)
365 {
366 std::cout << green << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
367 hnls.iter << " iterations " << def << std::endl << std::endl;
368 }
369 else
370 {
371 std::cout << red << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
372 hnls.iter << " iterations " << def << std::endl << std::endl;
373 }
374
375 tmp_sol(0) = time;
376 tmp_sol.col(0).tail(y.rows()) = y;
377 online_solutiont.col(i) = tmp_sol;
378 }
379
380 // Save the current solution
381 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "python",
382 "./ITHACAoutput/red_coeff/" + name(NParaSet) + "/");
383 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "matlab",
384 "./ITHACAoutput/red_coeff/" + name(NParaSet) + "/");
385 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "eigen",
386 "./ITHACAoutput/red_coeff/" + name(NParaSet) + "/");
387 return online_solutiont;
388}
389
390// * * * * * * * * * * * * * * * Solve Functions * * * * * * * * * * * * * //
391Eigen::MatrixXd ReducedUnsteadyBB::solveOnline_PPE(Eigen::MatrixXd&
392 temp_now_BC,
393 Eigen::MatrixXd& vel_now_BC, int NParaSet, int startSnap)
394{
395 std::cout << "################## Online solve N° " << NParaSet <<
396 " ##################" << std::endl;
397 std::cout << "Solving for the parameter: " << temp_now_BC << std::endl;
398 // Count number of time steps
399 int counter = 0;
400 time = tstart;
401
402 while (time < finalTime - 0.5 * dt)
403 {
404 time = time + dt;
405 counter ++;
406 }
407
408 // Set size of online solution
409 online_solutiont.resize(Nphi_u + Nphi_prgh + Nphi_t + 1, counter + 1);
410 //Average Method
411 //if (problem->AveMethod == "mean")
412 //{
413 // for (int j = 0; j < temp_now_BC.cols(); j++)
414 //{
415 // for (int i = 0; i < N_BC_t; i++)
416 // {
417 //int patche = problem->inletIndexT(i, 0);
418 // temp_now_BC(i, j)= temp_now_BC(i, j)-problem->Tsub[0].boundaryField()[problem->inletIndexT(i, 0)][0];
419 // }
420 // }
421 //volScalarField T_IC("T_IC", problem->Tsub[0]);
422 //for (int j = 0; j < T_IC.boundaryField().size(); j++)
423 // {
424 // for (int i = 0; i < N_BC_t; i++)
425 // {
426 // T_IC.boundaryFieldRef()[problem->inletIndexT(i, 0)][j] =
427 // temp_now_BC(i, 0)-problem->Tsub[0].boundaryField()[problem->inletIndexT(i, 0)][0];
428 //
429 // }
430 // }
431 volScalarField T_IC("T_IC", problem->Tfield[0]);
432
433 for (int j = 0; j < T_IC.boundaryField().size(); j++)
434 {
435 for (int i = 0; i < N_BC_t; i++)
436 {
437 if (j == problem->inletIndexT(i, 0))
438 {
439 T_IC.boundaryFieldRef()[problem->inletIndexT(i, 0)][j] = temp_now_BC(i, 0);
440 }
441 else
442 {
443 }
444 }
445 }
446
447 // Create and resize the solution vector
448 y.resize(Nphi_u + Nphi_prgh + Nphi_t, 1);
449 y.setZero();
450 // Calculate the time-dependent coefficients
452 LUmodes);
453
454 if (Nphi_prgh != 0)
455 {
458 }
459
461 // Set some properties of the newton object
467 newton_object_PPE.BC.resize(N_BC);
468
469 //Eigen::MatrixXd Ncoeff = ITHACAutilities::getCoeffs(problem->nutFields
470 // , problem->nuTmodes);
471 // std::vector<double> tv;
472 // tv.resize(1);
473 // tv[0] = time;
474
475 //for (int l = 0; l < Nphi_nut; l++)
476 //{
477 // newton_object_PPE.nu_c(l) = problem->rbfsplines[l]->eval(tv);
478 //}
479
480 //Change initial condition for the lifting function
481 for (int j = 0; j < N_BC_t; j++)
482 {
483 newton_object_PPE.BC_t(j) = temp_now_BC(j, 0);
484 }
485
486 for (int j = 0; j < N_BC; j++)
487 {
488 newton_object_PPE.BC(j) = vel_now_BC(j, 0);
489 }
490
491 // Set the initial time
492 time = tstart;
493 // Create vector to store temporal solution and save initial condition as first solution
494 Eigen::MatrixXd tmp_sol(Nphi_u + Nphi_prgh + Nphi_t + 1, 1);
495 tmp_sol(0) = time;
496 tmp_sol.col(0).tail(y.rows()) = y;
497 online_solutiont.col(0) = tmp_sol;
498 // Create nonlinear solver object
499 Eigen::HybridNonLinearSolver<newton_unsteadyBB_PPE> hnls(newton_object_PPE);
500 // Set output colors for fancy output
504
505 // Start the time loop
506 for (int i = 1; i < online_solutiont.cols(); i++)
507 {
508 time = time + dt;
509 Eigen::VectorXd res(y);
510 res.setZero();
511 hnls.solve(y);
512
513 //if (problem->BCmethod == "lift")
514 // {
515
516 //}
517
518 for (int j = 0; j < N_BC; j++)
519 {
520 y(j) = vel_now_BC(j, 0);
521 }
522
523 for (int j = 0; j < N_BC_t; j++)
524 {
525 int k = j + Nphi_prgh + Nphi_u;
526 y(k) = temp_now_BC(j, 0);
527 }
528
529 newton_object_PPE.operator()(y, res);
531 Info << "Time = " << time << endl;
532
533 if (res.norm() < 1e-5)
534 {
535 std::cout << green << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
536 hnls.iter << " iterations " << def << std::endl << std::endl;
537 }
538 else
539 {
540 std::cout << red << "|F(x)| = " << res.norm() << " - Minimun reached in " <<
541 hnls.iter << " iterations " << def << std::endl << std::endl;
542 }
543
544 tmp_sol(0) = time;
545 tmp_sol.col(0).tail(y.rows()) = y;
546 online_solutiont.col(i) = tmp_sol;
547 }
548
549 // Save the current solution
550 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "python",
551 "./ITHACAoutput/red_coeff/" + name(NParaSet) + "/");
552 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "python",
553 ".");
554 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "matlab",
555 "./ITHACAoutput/red_coeff/" + name(NParaSet) + "/");
556 ITHACAstream::exportMatrix(online_solutiont, "red_coeff", "eigen",
557 "./ITHACAoutput/red_coeff/" + name(NParaSet) + "/");
558 return online_solutiont;
559}
560
561void ReducedUnsteadyBB::reconstruct_sup(fileName folder, int printevery)
562{
564 {
565 }
566 else
567 {
568 mkDir(folder);
570 }
571
572 int counter = 0;
573 int nextwrite = 0;
574 int counter2 = 1 + TREC.size();
575
576 for (int i = 0; i < online_solutiont.cols(); i++)
577 {
578 if (counter == nextwrite)
579 {
580 volVectorField U_rec("U_rec", LUmodes[0] * 0);
581
582 for (int j = 0; j < Nphi_u; j++)
583 {
584 U_rec += LUmodes[j] * online_solutiont(j + 1, i);
585 }
586
587 ITHACAstream::exportSolution(U_rec, name(counter2), folder);
588
589 if (Nphi_prgh != 0)
590 {
591 volScalarField P_rec("P_rec", problem->Prghmodes[0] * 0);
592
593 for (int j = 0; j < Nphi_prgh; j++)
594 {
595 P_rec += problem->Prghmodes[j] * online_solutiont(j + Nphi_u + 1, i);
596 }
597
598 ITHACAstream::exportSolution(P_rec, name(counter2), folder);
599 PREC.append((P_rec).clone());
600 }
601
602 volScalarField T_rec("T_rec", LTmodes[0] * 0);
603
604 for (int j = 0; j < Nphi_t; j++)
605 {
606 T_rec += LTmodes[j] * online_solutiont(j + Nphi_prgh + Nphi_u + 1, i);
607 }
608
609 ITHACAstream::exportSolution(T_rec, name(counter2), folder);
610 nextwrite += printevery;
611 double timenow = online_solutiont(0, i);
612 std::ofstream of(folder + name(counter2) + "/" + name(timenow));
613 counter2 ++;
614 UREC.append((U_rec).clone());
615 TREC.append((T_rec).clone());
616 }
617
618 counter++;
619 }
620}
621
622
623// ************************************************************************* //
Header file of the ReducedUnsteadyBB class.
Class to change color to the output stream.
Definition colormod.H:24
int Nphi_t
Number of temperature modes.
int N_BC_t
Number of parametrized boundary conditions related to temperature field.
PtrList< volScalarField > PREC
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.
newton_unsteadyBB_PPE newton_object_PPE
Function object to call the non linear solver PPE approach.
scalar Pr
DimensionedScalar Pr;.
PtrList< volVectorField > LUmodes
List of pointers to store the modes for velocity including lift field modes.
Eigen::MatrixXd online_solutiont
List of Eigen matrices to store current online solution for temperature equation.
Eigen::MatrixXd solveOnline_sup(Eigen::MatrixXd &temp_now_BC, Eigen::MatrixXd &vel_now_BC, int NParaSet=0, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
PtrList< volScalarField > LTmodes
List of pointers to store the modes for temperature including lift field modes.
int Nphi_prgh
Number of shifted pressure modes.
PtrList< volScalarField > TREC
Reconstructed temperature field.
Eigen::MatrixXd solveOnline_PPE(Eigen::MatrixXd &temp_now_BC, Eigen::MatrixXd &vel_now_BC, int NParaSet=0, int startSnap=0)
Method to perform an online solve using a PPE stabilisation method.
newton_unsteadyBB_sup newton_object_sup
Function object to call the non linear solver sup approach.
ReducedUnsteadyBB()
Construct Null.
UnsteadyBB * problem
Pointer to the FOM problem.
Implementation of a parametrized full order unsteady Boussinesq problem and preparation of the the ...
Definition UnsteadyBB.H:60
Eigen::MatrixXd HP_matrix
Buoyancy term - PPE equation.
Definition UnsteadyBB.H:115
Eigen::MatrixXd Y_matrix
Diffusive term - energy equation.
Definition UnsteadyBB.H:118
PtrList< volScalarField > Prghmodes
List of pointers used to form the shifted pressure modes.
Definition UnsteadyBB.H:79
Eigen::MatrixXd W_matrix
Mass Matrix - energy equation.
Definition UnsteadyBB.H:125
PtrList< volScalarField > Prghfield
List of pointers used to form the shifted pressure snapshots matrix.
Definition UnsteadyBB.H:73
Eigen::MatrixXd H_matrix
Buoyancy term - momentum equation.
Definition UnsteadyBB.H:112
List< Eigen::MatrixXd > Q_matrix
Non linear convective term - energy equation.
Definition UnsteadyBB.H:121
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
Definition UnsteadyBB.H:76
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of lifting functions.
Definition UnsteadyBB.H:94
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
Definition UnsteadyBB.H:82
label NTmodes
Number of temperature modes used for the projection.
Definition UnsteadyBB.H:106
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
scalar nu
Reduced viscosity in case of parametrized viscosity.
int N_BC
Number of parametrized boundary conditions.
PtrList< volVectorField > UREC
Recontructed velocity field.
int Nphi_u
Number of velocity modes.
scalar finalTime
Scalar to store the final time if the online simulation.
scalar tstart
Scalar to store the initial time if the online simulation.
scalar time
Scalar to store the current time.
double dt
Scalar to store the time increment.
Eigen::MatrixXi inletIndexT
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
Eigen::MatrixXd BC1_matrix
PPE BC1.
Definition steadyNS.H:182
Eigen::MatrixXd BC3_matrix
PPE BC3.
Definition steadyNS.H:191
Eigen::Tensor< double, 3 > C_tensor
Diffusion term.
Definition steadyNS.H:167
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
List< Eigen::MatrixXd > G_matrix
Divergence of momentum PPE.
Definition steadyNS.H:176
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
Eigen::MatrixXd D_matrix
Laplacian term PPE.
Definition steadyNS.H:173
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
Eigen::MatrixXd P_matrix
Div of velocity.
Definition steadyNS.H:170
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition steadyNS.H:160
List< Eigen::MatrixXd > BC2_matrix
PPE BC2.
Definition steadyNS.H:185
@ FG_GREEN
Definition colormod.H:14
@ FG_DEFAULT
Definition colormod.H:16
@ FG_RED
Definition colormod.H:13
Matrix< VectorType, Dynamic, Dynamic > SliceFromTensor(Eigen::Tensor< VectorType, 3 > &tensor, label dim, label index1)
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.
bool check_folder(word folder)
Checks if a folder exists.
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
label i
Definition pEqn.H:46
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) 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