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