Loading...
Searching...
No Matches
UnsteadyBB.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
13 License
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
31
34
35#include "UnsteadyBB.H"
36#include <cmath>
37
38// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
39
40// Constructor
42
43UnsteadyBB::UnsteadyBB(int argc, char* argv[])
44{
45 _args = autoPtr<argList>
46 (
47 new argList(argc, argv)
48 );
49
50 if (!_args->checkRootCase())
51 {
52 Foam::FatalError.exit();
53 }
54
55 argList& args = _args();
56#include "createTime.H"
57#include "createMesh.H"
58 _pimple = autoPtr<pimpleControl>
59 (
60 new pimpleControl
61 (
62 mesh
63 )
64 );
65#pragma GCC diagnostic push
66#pragma GCC diagnostic ignored "-Wunused-variable"
67#include "createFields.H"
68#pragma GCC diagnostic pop
69#include "createFvOptions.H"
73}
74
75
76// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
77#include "fvCFD.H"
78
79// Method to performa a truthSolve
80void UnsteadyBB::truthSolve(List<scalar> mu_now)
81{
82 Time& runTime = _runTime();
83 fvMesh& mesh = _mesh();
84#include "initContinuityErrs.H"
85 volScalarField& p = _p();
86 volVectorField& U = _U();
87 volScalarField& p_rgh = _p_rgh();
88 volScalarField& T = _T();
89 volScalarField& alphat = _alphat();
90 volScalarField& rhok = _rhok();
91 volScalarField& gh = _gh();
92 surfaceScalarField& ghf = _ghf();
93 surfaceScalarField& phi = _phi();
94 fv::options& fvOptions = _fvOptions();
95 pimpleControl& pimple = _pimple();
96 IOMRFZoneList& MRF = _MRF();
97 singlePhaseTransportModel& laminarTransport = _laminarTransport();
98 dimensionedScalar& beta = _beta();
99 dimensionedScalar& TRef = _TRef();
100 dimensionedScalar& Pr = _Pr();
101 dimensionedScalar& Prt = _Prt();
102 instantList Times = runTime.times();
103 runTime.setEndTime(finalTime);
104 runTime.setTime(Times[1], 1);
105 runTime.setDeltaT(timeStep);
107 // save initial condition in folder 0
108 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
109 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
110 ITHACAstream::exportSolution(p_rgh, name(counter), "./ITHACAoutput/Offline/");
111 ITHACAstream::exportSolution(T, name(counter), "./ITHACAoutput/Offline/");
112 std::ofstream of("./ITHACAoutput/Offline/" + name(counter) + "/" +
113 runTime.timeName());
114 Ufield.append(U.clone());
115 Pfield.append(p.clone());
116 Prghfield.append(p_rgh.clone());
117 Tfield.append(T.clone());
118 counter++;
120
121 // Start the time loop
122 while (runTime.run())
123 {
124#include "readTimeControls.H"
125#include "CourantNo.H"
126#include "setDeltaT.H"
127 runTime.setEndTime(finalTime + timeStep);
128 Info << "Time = " << runTime.timeName() << nl << endl;
129
130 // --- Pressure-velocity PIMPLE corrector loop
131 while (pimple.loop())
132 {
133#include "UEqn.H"
134#include "TEqn.H"
135
136 // --- Pressure corrector loop
137 while (pimple.correct())
138 {
139#include "pEqn.H"
140 }
141
142 if (pimple.turbCorr())
143 {
144 laminarTransport.correct();
145 turbulence->correct();
146 }
147 }
148
149 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
150 << " ClockTime = " << runTime.elapsedClockTime() << " s"
151 << nl << endl;
152
153 if (checkWrite(runTime))
154 {
155 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
156 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
157 ITHACAstream::exportSolution(p_rgh, name(counter), "./ITHACAoutput/Offline/");
158 ITHACAstream::exportSolution(T, name(counter), "./ITHACAoutput/Offline/");
159 std::ofstream of("./ITHACAoutput/Offline/" + name(counter) + "/" +
160 runTime.timeName());
161 Ufield.append(U.clone());
162 Pfield.append(p.clone());
163 Prghfield.append(p_rgh.clone());
164 Tfield.append(T.clone());
165 counter++;
167 writeMu(mu_now);
168 }
169
170 runTime++;
171 }
172}
173
174void UnsteadyBB::truthSolve(fileName folder)
175{
176#include "initContinuityErrs.H"
177 Time& runTime = _runTime();
178 fvMesh& mesh = _mesh();
179 volScalarField& p = _p();
180 volVectorField& U = _U();
181 volScalarField& p_rgh = _p_rgh();
182 volScalarField& T = _T();
183 // volScalarField& nut = _nut();
184 volScalarField& alphat = _alphat();
185 volScalarField& rhok = _rhok();
186 // dimensionedVector& g = _g();
187 volScalarField& gh = _gh();
188 surfaceScalarField& ghf = _ghf();
189 surfaceScalarField& phi = _phi();
190 fv::options& fvOptions = _fvOptions();
191 pimpleControl& pimple = _pimple();
192 IOMRFZoneList& MRF = _MRF();
193 singlePhaseTransportModel& laminarTransport = _laminarTransport();
194 // dimensionedScalar& nu = _nu();
195 dimensionedScalar& beta = _beta();
196 // dimensionedScalar& hRef = _hRef();
197 // dimensionedScalar& ghRef = _ghRef();
198 dimensionedScalar& TRef = _TRef();
199 dimensionedScalar& Pr = _Pr();
200 dimensionedScalar& Prt = _Prt();
201 instantList Times = runTime.times();
202 runTime.setEndTime(finalTime);
203 runTime.setTime(Times[1], 1);
204 runTime.setDeltaT(timeStep);
206 // Save initial condition
207 ITHACAstream::exportSolution(U, name(counter), folder);
208 ITHACAstream::exportSolution(p, name(counter), folder);
209 ITHACAstream::exportSolution(T, name(counter), folder);
211 std::ofstream of(folder + name(counter) + "/" + runTime.timeName());
212 counter++;
214
215 // Start the time loop
216 while (runTime.run())
217 {
218#include "readTimeControls.H"
219#include "CourantNo.H"
220#include "setDeltaT.H"
221 runTime.setEndTime(finalTime + timeStep);
222 Info << "Time = " << runTime.timeName() << nl << endl;
223
224 // --- Pressure-velocity PIMPLE corrector loop
225 while (pimple.loop())
226 {
227#include "UEqn.H"
228#include "TEqn.H"
229
230 // --- Pressure corrector loop
231 while (pimple.correct())
232 {
233#include "pEqn.H"
234 }
235
236 if (pimple.turbCorr())
237 {
238 laminarTransport.correct();
239 turbulence->correct();
240 }
241 }
242
243 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
244 << " ClockTime = " << runTime.elapsedClockTime() << " s"
245 << nl << endl;
246
247 if (checkWrite(runTime))
248 {
249 ITHACAstream::exportSolution(U, name(counter), folder);
250 ITHACAstream::exportSolution(p, name(counter), folder);
251 ITHACAstream::exportSolution(T, name(counter), folder);
253 std::ofstream of(folder + name(counter) + "/" + runTime.timeName());
254 counter++;
256 }
257
258 runTime++;
259 }
260}
261
262
263
264
265bool UnsteadyBB::checkWrite(Time& timeObject)
266{
267 scalar diffnow = mag(nextWrite - atof(timeObject.timeName().c_str()));
268 scalar diffnext = mag(nextWrite - atof(timeObject.timeName().c_str()) -
269 timeObject.deltaTValue());
270
271 if ( diffnow < diffnext)
272 {
273 return true;
274 }
275 else
276 {
277 return false;
278 }
279}
280
281// Method to solve the supremizer problem
283{
284 PtrList<volScalarField> P_sup;
285
286 if (type == "modes")
287 {
288 P_sup = Prghmodes;
289 }
290 else if (type == "snapshots")
291 {
292 P_sup = Prghfield;
293 }
294 else
295 {
296 std::cout << "You must specify the variable type with either snapshots or modes"
297 << std::endl;
298 exit(0);
299 }
300
301 if (supex == 1)
302 {
303 volVectorField U = _U();
304 volVectorField Usup
305 (
306 IOobject
307 (
308 "Usup",
309 U.time().timeName(),
310 U.mesh(),
311 IOobject::NO_READ,
312 IOobject::AUTO_WRITE
313 ),
314 U.mesh(),
315 dimensionedVector("zero", U.dimensions(), vector::zero)
316 );
317
318 if (type == "snapshots")
319 {
320 ITHACAstream::read_fields(supfield, Usup, "./ITHACAoutput/supfield/");
321 }
322 else if (type == "modes")
323 {
324 ITHACAstream::read_fields(supmodes, Usup, "./ITHACAoutput/supremizer/");
325 }
326 else
327 {
328 }
329 }
330 else
331 {
332 volVectorField U = _U();
333 volVectorField Usup
334 (
335 IOobject
336 (
337 "Usup",
338 U.time().timeName(),
339 U.mesh(),
340 IOobject::NO_READ,
341 IOobject::AUTO_WRITE
342 ),
343 U.mesh(),
344 dimensionedVector("zero", U.dimensions(), vector::zero)
345 );
346 dimensionedScalar nu_fake
347 (
348 "nu_fake",
349 dimensionSet(0, 2, -1, 0, 0, 0, 0),
350 scalar(1)
351 );
352 Vector<double> v(0, 0, 0);
353
354 for (label i = 0; i < Usup.boundaryField().size(); i++)
355 {
356 ITHACAutilities::changeBCtype(Usup, "fixedValue", i);
357 assignBC(Usup, i, v);
358 assignIF(Usup, v);
359 }
360
361 if (type == "snapshots")
362 {
363 for (label i = 0; i < P_sup.size(); i++)
364 {
365 fvVectorMatrix u_sup_eqn
366 (
367 - fvm::laplacian(nu_fake, Usup)
368 );
369 solve
370 (
371 u_sup_eqn == fvc::grad(P_sup[i])
372 );
373 supfield.append(Usup.clone());
374 ITHACAstream::exportSolution(Usup, name(i + 1), "./ITHACAoutput/supfield/");
375 }
376
377 label systemRet =
378 system("ln -s ../../constant ./ITHACAoutput/supfield/constant");
379 systemRet += system("ln -s ../../0 ./ITHACAoutput/supfield/0");
380 systemRet += system("ln -s ../../system ./ITHACAoutput/supfield/system");
381
382 if (systemRet < 0)
383 {
384 Info << "System Command Failed in steadyNS.C" << endl;
385 exit(0);
386 }
387 }
388 else
389 {
390 for (label i = 0; i < Prghmodes.size(); i++)
391 {
392 fvVectorMatrix u_sup_eqn
393 (
394 - fvm::laplacian(nu_fake, Usup)
395 );
396 solve
397 (
398 u_sup_eqn == fvc::grad(Prghmodes[i])
399 );
400 supmodes.append(Usup.clone());
401 ITHACAstream::exportSolution(Usup, name(i + 1), "./ITHACAoutput/supremizer/");
402 }
403
404 label systemRet =
405 system("ln -s ../../constant ./ITHACAoutput/supremizer/constant");
406 systemRet += system("ln -s ../../0 ./ITHACAoutput/supremizer/0");
407 systemRet += system("ln -s ../../system ./ITHACAoutput/supremizer/system");
408
409 if (systemRet < 0)
410 {
411 Info << "System Command Failed in steadyNS.C" << endl;
412 exit(0);
413 }
414 }
415 }
416}
417
418
419// * * * * * * * * * * * * * * lifting function temperature * * * * * * * * * * * * * //
420
421
422// Method to compute the lifting function for temperature
424{
425 for (label k = 0; k < inletIndexT.rows(); k++)
426 {
427 Time& runTime = _runTime();
428 fvMesh& mesh = _mesh();
429 volScalarField& T = _T();
430 volVectorField& U = _U();
431 surfaceScalarField& phi = _phi();
432 phi = linearInterpolate(U) & mesh.Sf();
433 simpleControl simple(mesh);
434 // IOMRFZoneList& MRF = _MRF();
435 // singlePhaseTransportModel& laminarTransport = _laminarTransport();
436 // volScalarField& nut = _nut();
437 volScalarField& alphat = _alphat();
438 // dimensionedScalar& nu = _nu();
439 dimensionedScalar& Pr = _Pr();
440 dimensionedScalar& Prt = _Prt();
441 label BCind = inletIndexT(k, 0);
442 volScalarField Tlift("Tlift" + name(k), T);
443 instantList Times = runTime.times();
444 runTime.setTime(Times[1], 1);
445 Info << "Solving a lifting Problem" << endl;
446 scalar t1 = 1;
447 scalar t0 = 0;
448 alphat = turbulence->nut() / Prt;
449 alphat.correctBoundaryConditions();
450 volScalarField alphaEff("alphaEff", turbulence->nu() / Pr + alphat);
451
452 for (label j = 0; j < T.boundaryField().size(); j++)
453 {
454 if (j == BCind)
455 {
456 assignBC(Tlift, j, t1);
457 assignIF(Tlift, t0);
458 }
459 else if (T.boundaryField()[BCind].type() == "fixedValue")
460 {
461 assignBC(Tlift, j, t0);
462 assignIF(Tlift, t0);
463 }
464 else
465 {
466 }
467 }
468
469 while (simple.correctNonOrthogonal())
470 {
471 fvScalarMatrix TEqn
472 (
473 fvm::div(phi, Tlift)
474 - fvm::laplacian(alphaEff, Tlift)
475 );
476 TEqn.solve();
477 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
478 << " ClockTime = " << runTime.elapsedClockTime() << " s"
479 << nl << endl;
480 }
481
482 Tlift.write();
483 liftfieldT.append(Tlift.clone());
484 }
485}
486
487// * * * * * * * * * * * * * * Projection Methods * * * * * * * * * * * * * * //
488
489void UnsteadyBB::projectSUP(fileName folder, label NU, label NP, label NT,
490 label NSUP)
491{
492 NUmodes = NU;
493 NTmodes = NT;
494 NSUPmodes = NSUP;
495 NPrghmodes = NP;
496
497 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
498 {
499 word M_str = "M_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
500 NSUPmodes);
501
502 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + M_str))
503 {
504 ITHACAstream::ReadDenseMatrix(M_matrix, "./ITHACAoutput/Matrices/", M_str);
505 }
506 else
507 {
509 }
510
511 word B_str = "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
512 NSUPmodes);
513
514 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + B_str))
515 {
516 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", B_str);
517 }
518 else
519 {
521 }
522
523 word C_str = "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
524 NSUPmodes) + "_t";
525
526 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + C_str))
527 {
528 ITHACAstream::ReadDenseTensor(C_tensor, "./ITHACAoutput/Matrices/", C_str);
529 }
530 else
531 {
533 }
534
535 word K_str = "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
536 NSUPmodes) + "_" + name(NPrghmodes);
537
538 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + K_str))
539 {
540 ITHACAstream::ReadDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/", K_str);
541 }
542 else
543 {
545 }
546
547 word H_str = "H_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
548 NSUPmodes) + "_" + name(liftfieldT.size()) + "_" + name(NTmodes);
549
550 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + H_str))
551 {
552 ITHACAstream::ReadDenseMatrix(H_matrix, "./ITHACAoutput/Matrices/", H_str);
553 }
554 else
555 {
557 }
558
559 word W_str = "W_" + name(liftfieldT.size()) + "_" + name(NTmodes);
560
561 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + W_str))
562 {
563 ITHACAstream::ReadDenseMatrix(W_matrix, "./ITHACAoutput/Matrices/", W_str);
564 }
565 else
566 {
568 }
569
571 word Y_str = "Y_" + name(liftfieldT.size()) + "_" + name(NTmodes);
572
573 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + Y_str))
574 {
575 ITHACAstream::ReadDenseMatrix(Y_matrix, "./ITHACAoutput/Matrices/", Y_str);
576 }
577 else
578 {
580 }
581
582 word P_str = "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
583 NSUPmodes) + "_" + name(NPrghmodes);
584
585 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + P_str))
586 {
587 ITHACAstream::ReadDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/", P_str);
588 }
589 else
590 {
592 }
593 }
594 else
595 {
596 L_U_SUPmodes.resize(0);
597
598 if (liftfield.size() != 0)
599 {
600 for (label k = 0; k < liftfield.size(); k++)
601 {
602 L_U_SUPmodes.append(liftfield[k].clone());
603 }
604 }
605
606 if (NUmodes != 0)
607 {
608 for (label k = 0; k < NUmodes; k++)
609 {
610 L_U_SUPmodes.append(Umodes[k].clone());
611 }
612 }
613
614 if (NSUPmodes != 0)
615 {
616 for (label k = 0; k < NSUPmodes; k++)
617 {
618 L_U_SUPmodes.append(supmodes[k].clone());
619 }
620 }
621
622 L_T_modes.resize(0);
623
624 if (liftfieldT.size() != 0)
625 {
626 for (label k = 0; k < liftfieldT.size(); k++)
627 {
628 L_T_modes.append(liftfieldT[k].clone());
629 }
630 }
631
632 if (NTmodes != 0)
633 {
634 for (label k = 0; k < NTmodes; k++)
635 {
636 L_T_modes.append(Tmodes[k].clone());
637 }
638 }
639
649 }
650}
651
652
653void UnsteadyBB::projectPPE(fileName folder, label NU, label NP, label NT,
654 label NSUP)
655{
656 NUmodes = NU;
657 NTmodes = NT;
658 NSUPmodes = 0;
659 NPrghmodes = NP;
660
661 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
662 {
663 word M_str = "M_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
664 NSUPmodes);
665
666 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + M_str))
667 {
668 ITHACAstream::ReadDenseMatrix(M_matrix, "./ITHACAoutput/Matrices/", M_str);
669 }
670 else
671 {
673 }
674
675 word B_str = "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
676 NSUPmodes);
677
678 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + B_str))
679 {
680 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", B_str);
681 }
682 else
683 {
685 }
686
687 ITHACAstream::exportMatrix(B_matrix, "B_matrix", "eigen",
688 "./ITHACAoutput/Matrices/");
689 word C_str = "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
690 NSUPmodes) + "_t";
691
692 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + C_str))
693 {
694 ITHACAstream::ReadDenseTensor(C_tensor, "./ITHACAoutput/Matrices/", C_str);
695 }
696 else
697 {
699 }
700
701 word K_str = "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
702 NSUPmodes) + "_" + name(NPrghmodes);
703
704 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + K_str))
705 {
706 ITHACAstream::ReadDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/", K_str);
707 }
708 else
709 {
711 }
712
713 word H_str = "H_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
714 NSUPmodes) + "_" + name(liftfieldT.size()) + "_" + name(NTmodes);
715
716 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + H_str))
717 {
718 ITHACAstream::ReadDenseMatrix(H_matrix, "./ITHACAoutput/Matrices/", H_str);
719 }
720 else
721 {
723 }
724
725 word HP_str = "HP_" + name(NPrghmodes) + "_" + name(liftfieldT.size()) + "_" +
726 name(NTmodes);
727
728 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + HP_str))
729 {
730 ITHACAstream::ReadDenseMatrix(HP_matrix, "./ITHACAoutput/Matrices/", HP_str);
731 }
732 else
733 {
735 }
736
737 word W_str = "W_" + name(liftfieldT.size()) + "_" + name(NTmodes);
738
739 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + W_str))
740 {
741 ITHACAstream::ReadDenseMatrix(W_matrix, "./ITHACAoutput/Matrices/", W_str);
742 }
743 else
744 {
746 }
747
749 word Y_str = "Y_" + name(liftfieldT.size()) + "_" + name(NTmodes);
750
751 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + Y_str))
752 {
753 ITHACAstream::ReadDenseMatrix(Y_matrix, "./ITHACAoutput/Matrices/", Y_str);
754 }
755 else
756 {
758 }
759
760 word P_str = "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
761 NSUPmodes) + "_" + name(NPrghmodes);
762
763 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + P_str))
764 {
765 ITHACAstream::ReadDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/", P_str);
766 }
767 else
768 {
770 }
771 }
772 else
773 {
774 L_U_SUPmodes.resize(0);
775
776 if (liftfield.size() != 0)
777 {
778 for (label k = 0; k < liftfield.size(); k++)
779 {
780 L_U_SUPmodes.append(liftfield[k].clone());
781 }
782 }
783
784 if (NUmodes != 0)
785 {
786 for (label k = 0; k < NUmodes; k++)
787 {
788 L_U_SUPmodes.append(Umodes[k].clone());
789 }
790 }
791
792 if (NSUPmodes != 0)
793 {
794 for (label k = 0; k < NSUPmodes; k++)
795 {
796 L_U_SUPmodes.append(supmodes[k].clone());
797 }
798 }
799
800 L_T_modes.resize(0);
801
802 if (liftfieldT.size() != 0)
803 {
804 for (label k = 0; k < liftfieldT.size(); k++)
805 {
806 L_T_modes.append(liftfieldT[k].clone());
807 }
808 }
809
810 if (NTmodes != 0)
811 {
812 for (label k = 0; k < NTmodes; k++)
813 {
814 L_T_modes.append(Tmodes[k].clone());
815 }
816 }
817
833 }
834}
835
836// * * * * * * * * * * * * * * Momentum Eq. Methods * * * * * * * * * * * * * //
838 label NPrghmodes,
839 label NSUPmodes)
840{
841 label K1size = NUmodes + NSUPmodes + liftfield.size();
842 label K2size = NPrghmodes;
843 Eigen::MatrixXd K_matrix(K1size, K2size);
844 dimensionedVector g = _g();
845
846 // Project everything
847 for (label i = 0; i < K1size; i++)
848 {
849 for (label j = 0; j < K2size; j++)
850 {
851 K_matrix(i, j) = fvc::domainIntegrate(L_U_SUPmodes[i] &
852 fvc::reconstruct(fvc::snGrad(Prghmodes[j]) *
853 Prghmodes[j].mesh().magSf())).value();
854 }
855 }
856
857 if (Pstream::parRun())
858 {
859 reduce(K_matrix, sumOp<Eigen::MatrixXd>());
860 }
861
862 if (Pstream::master())
863 {
864 // Export the matrix
865 ITHACAstream::SaveDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/",
866 "K_" + name(liftfield.size()) + "_" + name(NUmodes)
867 + "_" + name(NSUPmodes) + "_" + name(NPrghmodes));
868 }
869
870 return K_matrix;
871}
872
873// * * * * * * * * * * * * * * Continuity Eq. Methods * * * * * * * * * * * * * //
874
875Eigen::MatrixXd UnsteadyBB::divergence_term(label NUmodes, label NPrghmodes,
876 label NSUPmodes)
877{
878 label P1size = NPrghmodes;
879 label P2size = NUmodes + NSUPmodes + liftfield.size();
880 Eigen::MatrixXd P_matrix(P1size, P2size);
881
882 // Project everything
883 for (label i = 0; i < P1size; i++)
884 {
885 for (label j = 0; j < P2size; j++)
886 {
887 P_matrix(i, j) = fvc::domainIntegrate(Prghmodes[i] * fvc::div (
888 L_U_SUPmodes[j])).value();
889 }
890 }
891
892 if (Pstream::parRun())
893 {
894 reduce(P_matrix, sumOp<Eigen::MatrixXd>());
895 }
896
897 if (Pstream::master())
898 {
899 //Export the matrix
900 ITHACAstream::SaveDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/",
901 "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_"
902 + name(NSUPmodes) + "_" + name(NPrghmodes));
903 }
904
905 return P_matrix;
906}
907
908Eigen::MatrixXd UnsteadyBB::buoyant_term(label NUmodes, label NTmodes,
909 label NSUPmodes)
910{
911 label H1size = NUmodes + liftfield.size() + NSUPmodes;
912 label H2size = NTmodes + liftfieldT.size() ;
913 Eigen::MatrixXd H_matrix(H1size, H2size);
914 dimensionedScalar beta = _beta();
915 dimensionedScalar TRef = _TRef();
916 dimensionedVector g = _g();
917 // volScalarField& gh = _gh();
918 surfaceScalarField& ghf = _ghf();
919
920 // Project everything
921 for (label i = 0; i < H1size; i++)
922 {
923 for (label j = 0; j < H2size; j++)
924 {
925 H_matrix(i, j) = fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::reconstruct(
926 ghf * fvc::snGrad(1.0 - (beta * (L_T_modes[j] - TRef)))
927 * L_T_modes[j].mesh().magSf())).value();
928 }
929 }
930
931 if (Pstream::parRun())
932 {
933 reduce(H_matrix, sumOp<Eigen::MatrixXd>());
934 }
935
936 if (Pstream::master())
937 {
938 //Export the matrix
939 ITHACAstream::SaveDenseMatrix(H_matrix, "./ITHACAoutput/Matrices/",
940 "H_" + name(liftfield.size()) + "_" + name(NUmodes) + "_"
941 + name(NSUPmodes) + "_" + name(liftfieldT.size()) + "_"
942 + name(NTmodes));
943 }
944
945 return H_matrix;
946}
947
948// * * * * * * * * * * * * * * Additional term for the PPE Method * * * * * * * * * * * * * //
950 label NTmodes)
951{
952 label H1size = NPrghmodes;
953 label H2size = NTmodes + liftfieldT.size() ;
954 Eigen::MatrixXd HP_matrix(H1size, H2size);
955 // Create PTRLIST with lift, velocities and temperatures
956 dimensionedScalar beta = _beta();
957 dimensionedScalar TRef = _TRef();
958 dimensionedVector g = _g();
959 volScalarField& gh = _gh();
960 surfaceScalarField& ghf = _ghf();
961
962 // Project everything
963 for (label i = 0; i < H1size; i++)
964 {
965 for (label j = 0; j < H2size; j++)
966 {
967 HP_matrix(i, j) = fvc::domainIntegrate(fvc::reconstruct(fvc::snGrad(
968 Prghmodes[i]) *
969 Prghmodes[i].mesh().magSf()) & fvc::reconstruct(
970 ghf * fvc::snGrad(-(beta * (L_T_modes[j])))
971 * L_T_modes[j].mesh().magSf())).value();
972 }
973 }
974
975 if (Pstream::parRun())
976 {
977 reduce(HP_matrix, sumOp<Eigen::MatrixXd>());
978 }
979
980 if (Pstream::master())
981 {
982 ITHACAstream::SaveDenseMatrix(HP_matrix, "./ITHACAoutput/Matrices/",
983 "HP_"
984 + name(NPrghmodes) + "_" + name(liftfieldT.size()) + "_"
985 + name(NTmodes));
986 }
987
988 return HP_matrix;
989}
990
991
992// * * * * * * * * * * * * * * Energy Eq. Methods * * * * * * * * * * * * * //
994 label NTmodes, label NSUPmodes)
995{
996 label Qsize = NUmodes + liftfield.size() + NSUPmodes;
997 label Qsizet = NTmodes + liftfieldT.size() ;
998 List <Eigen::MatrixXd> Q_matrix;
999 Q_matrix.setSize(Qsizet);
1000
1001 for (label j = 0; j < Qsizet; j++)
1002 {
1003 Q_matrix[j].resize(Qsize, Qsizet);
1004 }
1005
1006 for (label i = 0; i < Qsizet; i++)
1007 {
1008 for (label j = 0; j < Qsize; j++)
1009 {
1010 for (label k = 0; k < Qsizet; k++)
1011 {
1012 Q_matrix[i](j, k) = fvc::domainIntegrate(L_T_modes[i] * fvc::div(
1013 fvc::interpolate(L_U_SUPmodes[j]) & L_U_SUPmodes[j].mesh().Sf(),
1014 L_T_modes[k])).value();
1015 }
1016 }
1017
1018 if (Pstream::parRun())
1019 {
1020 reduce(Q_matrix[i], sumOp<Eigen::MatrixXd>());
1021 }
1022 }
1023
1024 if (Pstream::master())
1025 {
1026 // Export the matrix
1027 ITHACAstream::exportMatrix(Q_matrix, "Q", "python", "./ITHACAoutput/Matrices/");
1028 ITHACAstream::exportMatrix(Q_matrix, "Q", "matlab", "./ITHACAoutput/Matrices/");
1029 ITHACAstream::exportMatrix(Q_matrix, "Q", "eigen", "./ITHACAoutput/Matrices/Q");
1030 }
1031
1032 return Q_matrix;
1033}
1034
1035
1036
1038 label NTmodes, label NSUPmodes)
1039{
1040 label Ysize = NTmodes + liftfieldT.size();
1041 Eigen::MatrixXd Y_matrix(Ysize, Ysize);
1042
1043 for (label i = 0; i < Ysize; i++)
1044 {
1045 for (label j = 0; j < Ysize; j++)
1046 {
1047 Y_matrix(i, j) = fvc::domainIntegrate(L_T_modes[i] * fvc::laplacian(
1048 dimensionedScalar("1", dimless, 1), L_T_modes[j])).value();
1049 }
1050 }
1051
1052 if (Pstream::parRun())
1053 {
1054 reduce(Y_matrix, sumOp<Eigen::MatrixXd>());
1055 }
1056
1057 if (Pstream::master())
1058 {
1059 // Export the matrix
1060 ITHACAstream::SaveDenseMatrix(Y_matrix, "./ITHACAoutput/Matrices/",
1061 "Y_" + name(liftfieldT.size()) + "_" + name(NTmodes));
1062 }
1063
1064 return Y_matrix;
1065}
1066
1067
1069 label NSUPmodes)
1070{
1071 label Wsize = NTmodes + liftfieldT.size();
1072 Eigen::MatrixXd W_matrix(Wsize, Wsize);
1073
1074 for (label i = 0; i < Wsize; i++)
1075 {
1076 for (label j = 0; j < Wsize; j++)
1077 {
1078 W_matrix(i, j) = fvc::domainIntegrate(L_T_modes[i] * L_T_modes[j]).value();
1079 }
1080 }
1081
1082 if (Pstream::parRun())
1083 {
1084 reduce(W_matrix, sumOp<Eigen::MatrixXd>());
1085 }
1086
1087 if (Pstream::master())
1088 {
1089 // Export the matrix
1090 ITHACAstream::SaveDenseMatrix(W_matrix, "./ITHACAoutput/Matrices/",
1091 "W_" + name(liftfieldT.size()) + "_" + name(NTmodes));
1092 }
1093
1094 return W_matrix;
1095}
1096
1097
1098// Method to compute the lifting function for velocity
1100{
1101 for (label k = 0; k < inletIndex.rows(); k++)
1102 {
1103 Time& runTime = _runTime();
1104 surfaceScalarField& phi = _phi();
1105 fvMesh& mesh = _mesh();
1106 // volScalarField& p = _p();
1107 volVectorField& U = _U();
1108 // volScalarField& p_rgh = _p_rgh();
1109 volScalarField& UliftBC = _UliftBC();
1110 // volScalarField& T = _T();
1111 // dimensionedScalar& nu = _nu();
1112 IOMRFZoneList& MRF = _MRF();
1113 label BCind = inletIndex(k, 0);
1114 volVectorField Ulift("Ulift" + name(k), U);
1115 instantList Times = runTime.times();
1116 runTime.setTime(Times[1], 1);
1117 pisoControl potentialFlow(mesh, "potentialFlow");
1118 Info << "Solving a lifting Problem" << endl;
1119 Vector<double> v1(0, 0, 0);
1120 v1[inletIndex(k, 1)] = 1;
1121 Vector<double> v0(0, 0, 0);
1122
1123 for (label j = 0; j < U.boundaryField().size(); j++)
1124 {
1125 if (j == BCind)
1126 {
1127 assignBC(Ulift, j, v1);
1128 }
1129 else if (U.boundaryField()[BCind].type() == "fixedValue")
1130 {
1131 assignBC(Ulift, j, v0);
1132 }
1133 else
1134 {
1135 }
1136
1137 assignIF(Ulift, v0);
1138 phi = linearInterpolate(Ulift) & mesh.Sf();
1139 }
1140
1141 Info << "Constructing velocity potential field Phi\n" << endl;
1142 volScalarField Phi
1143 (
1144 IOobject
1145 (
1146 "Phi",
1147 runTime.timeName(),
1148 mesh,
1149 IOobject::READ_IF_PRESENT,
1150 IOobject::NO_WRITE
1151 ),
1152 mesh,
1153 dimensionedScalar("Phi", dimLength * dimVelocity, 0),
1154 UliftBC.boundaryField().types()
1155 );
1156 label PhiRefCell = 0;
1157 scalar PhiRefValue = 0.0;
1159 (
1160 Phi,
1161 potentialFlow.dict(),
1162 PhiRefCell,
1163 PhiRefValue
1164 );
1165 mesh.setFluxRequired(Phi.name());
1166 runTime.functionObjects().start();
1167 MRF.makeRelative(phi);
1168 adjustPhi(phi, Ulift, UliftBC);
1169
1170 while (potentialFlow.correctNonOrthogonal())
1171 {
1172 fvScalarMatrix PhiEqn
1173 (
1174 fvm::laplacian(dimensionedScalar("1", dimless, 1), Phi)
1175 ==
1176 fvc::div(phi)
1177 );
1178 PhiEqn.setReference(PhiRefCell, PhiRefValue);
1179 PhiEqn.solve();
1180
1181 if (potentialFlow.finalNonOrthogonalIter())
1182 {
1183 phi -= PhiEqn.flux();
1184 }
1185 }
1186
1187 MRF.makeAbsolute(phi);
1188 Info << "Continuity error = "
1189 << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
1190 << endl;
1191 Ulift = fvc::reconstruct(phi);
1192 Ulift.correctBoundaryConditions();
1193 Info << "Interpolated velocity error = "
1194 << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
1195 / sum(mesh.magSf())).value()
1196 << endl;
1197 Ulift.write();
1198 liftfield.append(Ulift.clone());
1199 }
1200}
1201
1203{
1204 const volScalarField& nu = _laminarTransport().nu();
1205 volScalarField& ciao = const_cast<volScalarField&>(nu);
1206 this->assignIF(ciao, mu);
1207
1208 for (label i = 0; i < ciao.boundaryFieldRef().size(); i++)
1209 {
1210 this->assignBC(ciao, i, mu);
1211 }
1212}
1213
1214
1215
1216
1217
1218
1219
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
fv::options & fvOptions
Definition NLsolve.H:25
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
TEqn solve()
Header file of the UnsteadyBB class.
autoPtr< volScalarField > _p_rgh
Shifted Pressure field.
Definition UnsteadyBB.H:131
Eigen::MatrixXd mass_term_temperature(label NUmodes, label NTmodes, label NSUPmodes)
Mass Term Energy Equation.
void liftSolveT()
Perform a lift solve for temperature.
Definition UnsteadyBB.C:423
void solvesupremizer(word type="snapshots")
solve the supremizer either with the use of the pressure snaphots or the pressure modes
Definition UnsteadyBB.C:282
Eigen::MatrixXd HP_matrix
Buoyancy term - PPE equation.
Definition UnsteadyBB.H:115
Eigen::MatrixXd Y_matrix
Diffusive term - energy equation.
Definition UnsteadyBB.H:118
autoPtr< volScalarField > _UliftBC
Definition UnsteadyBB.H:103
Eigen::MatrixXd diffusive_term_temperature(label NUmodes, label NTmodes, label NSUPmodes)
Diffusive Term Energy Equation.
PtrList< volScalarField > Prghmodes
List of pointers used to form the shifted pressure modes.
Definition UnsteadyBB.H:79
autoPtr< volScalarField > _T
Temperature field.
Definition UnsteadyBB.H:134
autoPtr< dimensionedScalar > _Prt
dimensionedScalar Prt;
Definition UnsteadyBB.H:146
Eigen::MatrixXd W_matrix
Mass Matrix - energy equation.
Definition UnsteadyBB.H:125
autoPtr< dimensionedScalar > _Pr
dimensionedScalar Pr;
Definition UnsteadyBB.H:143
List< Eigen::MatrixXd > convective_term_temperature(label NUmodes, label NTmodes, label NSUPmodes)
Convective Term Energy Equation.
Definition UnsteadyBB.C:993
Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPrghmodes, label NSUPmodes)
Gradient of pressure.
Definition UnsteadyBB.C:837
PtrList< volScalarField > Prghfield
List of pointers used to form the shifted pressure snapshots matrix.
Definition UnsteadyBB.H:73
Eigen::MatrixXd divergence_term(label NUmodes, label NPrghmodes, label NSUPmodes)
Divergence Term (supremizer approach)
Definition UnsteadyBB.C:875
autoPtr< dimensionedVector > _g
Definition UnsteadyBB.H:166
void change_viscosity(double mu)
Function to change the viscosity.
UnsteadyBB()
Null constructor.
Definition UnsteadyBB.C:41
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NTmodes, label NSUPmodes)
Project using a supremizer approach.
Definition UnsteadyBB.C:489
Eigen::MatrixXd H_matrix
Buoyancy term - momentum equation.
Definition UnsteadyBB.H:112
void liftSolve()
Perform a lift solve for velocity field.
autoPtr< volScalarField > _gh
List of pointers used to form the gravitational acceleration.
Definition UnsteadyBB.H:161
autoPtr< fvMesh > _mesh
Mesh.
Definition UnsteadyBB.H:128
label NPrghmodes
Number of pressure modes used for the projection.
Definition UnsteadyBB.H:109
Eigen::MatrixXd buoyant_term_poisson(label NPrghmodes, label NTmodes)
Buoyant Term PPE Equation.
Definition UnsteadyBB.C:949
PtrList< volScalarField > L_T_modes
List of pointers containing the lift for temperature and the temperature field.
Definition UnsteadyBB.H:100
List< Eigen::MatrixXd > Q_matrix
Non linear convective term - energy equation.
Definition UnsteadyBB.H:121
autoPtr< dimensionedScalar > _TRef
dimensionedScalar Tref;
Definition UnsteadyBB.H:140
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
Definition UnsteadyBB.H:76
bool checkWrite(Time &timeObject)
Function to check if the solution must be exported.
Definition UnsteadyBB.C:265
autoPtr< dimensionedScalar > _beta
dimensionedScalar beta;
Definition UnsteadyBB.H:137
autoPtr< volScalarField > _rhok
dimensionedScalar rhok;
Definition UnsteadyBB.H:155
autoPtr< volScalarField > _alphat
dimensionedScalar alphat;
Definition UnsteadyBB.H:149
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of lifting functions.
Definition UnsteadyBB.H:94
Eigen::MatrixXd buoyant_term(label NUmodes, label NTmodes, label NSUPmodes)
Buoyant Term Momentum Equation.
Definition UnsteadyBB.C:908
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
Definition UnsteadyBB.H:82
void projectPPE(fileName folder, label NUmodes, label NPrghmodes, label NTmodes, label NSUPmodes)
Project using a PPE approach.
Definition UnsteadyBB.C:653
autoPtr< surfaceScalarField > _ghf
List of pointers used to form the gravitational acceleration.
Definition UnsteadyBB.H:164
label NTmodes
Number of temperature modes used for the projection.
Definition UnsteadyBB.H:106
scalar startTime
Start Time (initial time to start storing the snapshots)
scalar writeEvery
Time step of the writing procedure.
scalar timeStep
Time step of the simulation.
scalar nextWrite
Auxiliary variable to store the next writing instant.
scalar finalTime
Final time (final time of the simulation and consequently of the acquisition of the snapshots)
Eigen::MatrixXi inletIndexT
void writeMu(List< scalar > mu_now)
Write out a list of scalar corresponding to the parameters used in the truthSolve.
label counter
Counter used for the output of the full order solutions.
void assignBC(volVectorField &s, label BC_ind, Vector< double > &value)
Assign Boundary Condition to a volVectorField.
void assignIF(T &s, G &value)
Assign internal field condition.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
Eigen::MatrixXd mu
Row matrix of parameters.
autoPtr< argList > _args
argList
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
label NPmodes
Number of pressure modes used for the projection.
Definition steadyNS.H:146
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:265
Eigen::MatrixXd BC1_matrix
PPE BC1.
Definition steadyNS.H:191
Eigen::MatrixXd diffusive_term(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Term.
Definition steadyNS.C:883
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:302
Eigen::MatrixXd BC3_matrix
PPE BC3.
Definition steadyNS.H:200
Eigen::Tensor< double, 3 > C_tensor
Diffusion term.
Definition steadyNS.H:170
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
autoPtr< fv::options > _fvOptions
fvOptions
Definition steadyNS.H:296
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:299
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
Eigen::MatrixXd pressure_BC1(label NPmodes, label NUmodes)
Term N° 1 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1309
List< Eigen::MatrixXd > G_matrix
Divergence of momentum PPE.
Definition steadyNS.H:185
Eigen::MatrixXd mass_term(label NUmodes, label NPmodes, label NSUPmodes)
Mass Term.
Definition steadyNS.C:1084
autoPtr< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model)
Definition steadyNS.H:311
label NUmodes
Number of velocity modes used for the projection.
Definition steadyNS.H:143
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
volVectorModes L_U_SUPmodes
List of pointers containing the total number of lift, supremizer and velocity modes.
Definition steadyNS.H:119
Eigen::MatrixXd pressure_BC3(label NPmodes, label NUmodes)
Term N° 3 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1420
Eigen::MatrixXd B_matrix
Diffusion term.
Definition steadyNS.H:160
Eigen::MatrixXd D_matrix
Laplacian term PPE.
Definition steadyNS.H:182
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:314
label NSUPmodes
Number of supremizer modes used for the projection.
Definition steadyNS.H:149
Eigen::MatrixXd K_matrix
Gradient of pressure matrix.
Definition steadyNS.H:166
Eigen::Tensor< double, 3 > convective_term_tens(label NUmodes, label NPmodes, label NSUPmodes)
Export convective term as a tensor.
Definition steadyNS.C:997
List< Eigen::MatrixXd > pressure_BC2(label NPmodes, label NUmodes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1343
PtrList< volVectorField > supfield
List of pointers used to form the supremizer snapshots matrix.
Definition steadyNS.H:92
Eigen::MatrixXd P_matrix
Div of velocity.
Definition steadyNS.H:173
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition steadyNS.H:163
List< Eigen::MatrixXd > BC2_matrix
PPE BC2.
Definition steadyNS.H:194
Eigen::MatrixXd laplacian_pressure(label NPmodes)
Laplacian of pressure term (PPE approach)
Definition steadyNS.C:1286
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:272
List< Eigen::MatrixXd > div_momentum(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach)
Definition steadyNS.C:1138
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:269
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition unsteadyNS.H:76
autoPtr< pimpleControl > _pimple
pimpleControl
Definition unsteadyNS.H:73
volScalarField & T
Definition createT.H:46
fvScalarMatrix & TEqn
Definition TEqn.H:15
volScalarField v0(v)
volScalarField & v
volScalarField & alphat
dimensionedScalar & nu
dimensionedScalar & Prt
dimensionedScalar & Pr
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.
void ReadDenseTensor(TensorType &Tensor, word folder, word MatrixName)
Read a dense tensor from file.
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
bool check_pod()
Check if the POD data folder "./ITHACAoutput/POD" exists.
bool check_off()
Check if the offline data folder "./ITHACAoutput/Offline" exists.
bool check_folder(word folder)
Checks if a folder exists.
bool check_file(std::string fileName)
Function that returns true if a file exists.
bool check_sup()
Check if the supremizer folder exists.
void changeBCtype(GeometricField< Type, fvPatchField, volMesh > &field, word BCtype, label BC_ind)
Change the boundary condition type for a GeometricField.
simpleControl simple(mesh)
surfaceScalarField & phi
volVectorField & U
volScalarField & p
volScalarField & gh
volScalarField & p_rgh
pimpleControl & pimple
volScalarField & rhok
volScalarField & UliftBC
dimensionedVector & g
dimensionedScalar & TRef
surfaceScalarField & ghf
dimensionedScalar & beta
singlePhaseTransportModel & laminarTransport
setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue)
adjustPhi(phiHbyA, U, p)
label i
Definition pEqn.H:46