Loading...
Searching...
No Matches
UnsteadyNSTTurb.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#include "UnsteadyNSTTurb.H"
32#include "viscosityModel.H"
33#include "alphatJayatillekeWallFunctionFvPatchScalarField.H"
34
35
38
39// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
40
41// Construct Null
43
44// Construct from zero
45UnsteadyNSTTurb::UnsteadyNSTTurb(int argc, char* argv[])
46{
47#include "setRootCase.H"
48#include "createTime.H"
49#include "createMesh.H"
50 _piso = autoPtr<pisoControl>
51 (
52 new pisoControl
53 (
54 mesh
55 )
56 );
57#pragma GCC diagnostic push
58#pragma GCC diagnostic ignored "-Wunused-variable"
59#include "createFields.H"
60#include "createFvOptions.H"
61#pragma GCC diagnostic pop
63}
64
65// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66
67void UnsteadyNSTTurb::truthSolve(List<scalar> mu_now)
68{
69 Time& runTime = _runTime();
70 surfaceScalarField& phi = _phi();
71 fvMesh& mesh = _mesh();
72#include "initContinuityErrs.H"
73 fv::options& fvOptions = _fvOptions();
74 pisoControl& piso = _piso();
75 volScalarField p = _p();
76 volVectorField U = _U();
77 volScalarField T = _T();
78 volScalarField _nut(turbulence->nut());
79 dimensionedScalar nu = _nu();
80 dimensionedScalar Pr = _Pr();
81 dimensionedScalar Prt = _Prt();
82 volScalarField alphat = _alphat();
83 singlePhaseTransportModel& laminarTransport = _laminarTransport();
84 //dimensionedScalar alphaEff = _alphaEff();
85 IOMRFZoneList& MRF = _MRF();
86 instantList Times = runTime.times();
87 runTime.setEndTime(finalTime);
88 // Perform a TruthSolve
89 bool WRITE;
90 runTime.setTime(Times[1], 1);
91 runTime.setDeltaT(timeStep);
93 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94 turbulence->validate();
95
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97
98 // Start the time loop
99 while (runTime.run())
100 {
101#include "CourantNo.H"
102 runTime.setEndTime(finalTime);
103 runTime++;
104 Info << "Time = " << runTime.timeName() << nl << endl;
105 // --- Pressure-velocity PIMPLE corrector loop
106 {
107#include "UEqn.H"
108
109 // --- Pressure corrector loop
110 while (piso.correct())
111 {
112#include "pEqn.H"
113 }
114 }
115 {
116#include "TEqn.H"
117 //TEqn.solve();
118 }
119 volScalarField _nut(turbulence->nut());
120 laminarTransport.correct();
121 turbulence->correct();
122 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
123 << " ClockTime = " << runTime.elapsedClockTime() << " s"
124 << nl << endl;
125 WRITE = checkWrite(runTime);
126
127 if (WRITE)
128 {
129 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
130 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
131 ITHACAstream::exportSolution(_nut, name(counter), "./ITHACAoutput/Offline/");
132 ITHACAstream::exportSolution(T, name(counter), "./ITHACAoutput/Offline/");
133 ITHACAstream::exportSolution(alphat, name(counter), "./ITHACAoutput/Offline/");
134 std::ofstream of("./ITHACAoutput/Offline/" + name(counter) + "/" +
135 runTime.timeName());
136 Ufield.append(U.clone());
137 Pfield.append(p.clone());
138 nutFields.append(_nut.clone());
139 Tfield.append(T.clone());
140 counter++;
142 writeMu(mu_now);
143 }
144 }
145}
146
147bool UnsteadyNSTTurb::checkWrite(Time& timeObject)
148{
149 scalar diffnow = mag(nextWrite - atof(timeObject.timeName().c_str()));
150 scalar diffnext = mag(nextWrite - atof(timeObject.timeName().c_str()) -
151 timeObject.deltaTValue());
152
153 if ( diffnow < diffnext)
154 {
155 return true;
156 }
157 else
158 {
159 return false;
160 }
161}
162
164{
165 for (label k = 0; k < inletIndexT.rows(); k++)
166 {
167 Time& runTime = _runTime();
168 fvMesh& mesh = _mesh();
169 volScalarField& T = _T();
170 simpleControl simple(mesh);
171 //IOMRFZoneList& MRF = _MRF();
172 //singlePhaseTransportModel& laminarTransport = _laminarTransport();
173 //volScalarField& nut = _nut();
174 volScalarField& alphat = _alphat();
175 //dimensionedScalar& nu = _nu();
176 dimensionedScalar& Pr = _Pr();
177 dimensionedScalar& Prt = _Prt();
178 label BCind = inletIndexT(k, 0);
179 volScalarField Tlift("Tlift" + name(k), T);
180 instantList Times = runTime.times();
181 runTime.setTime(Times[1], 1);
182 Info << "Solving a lifting Problem" << endl;
183 scalar t1 = 1;
184 scalar t0 = 0;
185
186 for (label j = 0; j < T.boundaryField().size(); j++)
187 {
188 if (j == BCind)
189 {
190 assignBC(Tlift, j, t1);
191 }
192 else if (T.boundaryField()[BCind].type() == "zeroGradient")
193 {
194 assignBC(Tlift, j, t0);
195 assignIF(Tlift, t1);
196 }
197 else if (T.boundaryField()[BCind].type() == "fixedValue")
198 {
199 assignBC(Tlift, j, t0);
200 assignIF(Tlift, t0);
201 }
202 else
203 {
204 }
205 }
206
207 while (simple.correctNonOrthogonal())
208 {
209 alphat = turbulence->nut() / Prt;
210 alphat.correctBoundaryConditions();
211 volScalarField alphaEff("alphaEff", turbulence->nu() / Pr + alphat);
212 fvScalarMatrix TEqn
213 (
214 fvm::ddt(Tlift)
215 - fvm::laplacian(alphaEff, Tlift)
216 );
217 TEqn.solve();
218 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
219 << " ClockTime = " << runTime.elapsedClockTime() << " s"
220 << nl << endl;
221 }
222
223 scalar totalTime = mesh.time().value();
224 scalar dt = mesh.time().deltaTValue();
225 forAll(Tlift, i)
226 {
227 Tlift[i] = (totalTime * Tlift[i] + dt * Tlift[i] ) / (totalTime + dt);
228 }
229 Tlift.write();
230 liftfieldT.append(Tlift.clone());
231 }
232}
233
234List <Eigen::MatrixXd> UnsteadyNSTTurb::turbulenceTerm1(label NUmodes,
235 label NSUPmodes, label Nnutmodes)
236{
237 label Csize = NUmodes + NSUPmodes + liftfield.size();
238 List <Eigen::MatrixXd> CT1_matrix;
239 CT1_matrix.setSize(Csize);
240
241 for (label j = 0; j < Csize; j++)
242 {
243 CT1_matrix[j].resize(Nnutmodes, Csize);
244 CT1_matrix[j] = CT1_matrix[j] * 0;
245 }
246
247 PtrList<volVectorField> Together(0);
248
249 // Create PTRLIST with lift, velocities and supremizers
250
251 if (liftfield.size() != 0)
252 {
253 for (label k = 0; k < liftfield.size(); k++)
254 {
255 Together.append(liftfield[k].clone());
256 }
257 }
258
259 if (NUmodes != 0)
260 {
261 for (label k = 0; k < NUmodes; k++)
262 {
263 Together.append(Umodes[k].clone());
264 }
265 }
266
267 if (NSUPmodes != 0)
268 {
269 for (label k = 0; k < NSUPmodes; k++)
270 {
271 Together.append(supmodes[k].clone());
272 }
273 }
274
275 for (label i = 0; i < Csize; i++)
276 {
277 Info << "Filling layer number " << i + 1 << " in the matrix CT1_matrix" << endl;
278
279 for (label j = 0; j < Nnutmodes; j++)
280 {
281 for (label k = 0; k < Csize; k++)
282 {
283 CT1_matrix[i](j, k) = fvc::domainIntegrate(Together[i] & fvc::laplacian(
284 nuTmodes[j], Together[k])).value();
285 }
286 }
287 }
288
289 // Export the matrix
290 ITHACAstream::exportMatrix(CT1_matrix, "CT1_matrix", "python",
291 "./ITHACAoutput/Matrices/");
292 ITHACAstream::exportMatrix(CT1_matrix, "CT1_matrix", "matlab",
293 "./ITHACAoutput/Matrices/");
294 ITHACAstream::exportMatrix(CT1_matrix, "CT1_matrix", "eigen",
295 "./ITHACAoutput/Matrices/CT1");
296 return CT1_matrix;
297}
298
299List <Eigen::MatrixXd> UnsteadyNSTTurb::turbulenceTerm2(label NUmodes,
300 label NSUPmodes, label Nnutmodes)
301{
302 label Csize = NUmodes + NSUPmodes + liftfield.size();
303 List <Eigen::MatrixXd> CT2_matrix;
304 CT2_matrix.setSize(Csize);
305
306 for (label j = 0; j < Csize; j++)
307 {
308 CT2_matrix[j].resize(Nnutmodes, Csize);
309 CT2_matrix[j] = CT2_matrix[j] * 0;
310 }
311
312 PtrList<volVectorField> Together(0);
313
314 // Create PTRLIST with lift, velocities and supremizers
315
316 if (liftfield.size() != 0)
317 {
318 for (label k = 0; k < liftfield.size(); k++)
319 {
320 Together.append(liftfield[k].clone());
321 }
322 }
323
324 if (NUmodes != 0)
325 {
326 for (label k = 0; k < NUmodes; k++)
327 {
328 Together.append(Umodes[k].clone());
329 }
330 }
331
332 if (NSUPmodes != 0)
333 {
334 for (label k = 0; k < NSUPmodes; k++)
335 {
336 Together.append(supmodes[k].clone());
337 }
338 }
339
340 for (label i = 0; i < Csize; i++)
341 {
342 Info << "Filling layer number " << i + 1 << " in the matrix CT2_matrix" << endl;
343
344 for (label j = 0; j < Nnutmodes; j++)
345 {
346 for (label k = 0; k < Csize; k++)
347 {
348 CT2_matrix[i](j, k) = fvc::domainIntegrate(Together[i] & (fvc::div(
349 nuTmodes[j] * dev((fvc::grad(Together[k]))().T())))).value();
350 }
351 }
352 }
353
354 // Export the matrix
355 ITHACAstream::exportMatrix(CT2_matrix, "CT2_matrix", "python",
356 "./ITHACAoutput/Matrices/");
357 ITHACAstream::exportMatrix(CT2_matrix, "CT2_matrix", "matlab",
358 "./ITHACAoutput/Matrices/");
359 ITHACAstream::exportMatrix(CT2_matrix, "CT2_matrix", "eigen",
360 "./ITHACAoutput/Matrices/CT2");
361 return CT2_matrix;
362}
363
364Eigen::MatrixXd UnsteadyNSTTurb::BTturbulence(label NUmodes, label NSUPmodes)
365{
366 label BTsize = NUmodes + NSUPmodes + liftfield.size();
367 Eigen::MatrixXd BT_matrix(BTsize, BTsize);
368 BT_matrix = BT_matrix * 0;
369 // Create PTRLIST with lift, velocities and supremizers
370 PtrList<volVectorField> Together(0);
371
372 if (liftfield.size() != 0)
373 {
374 for (label k = 0; k < liftfield.size(); k++)
375 {
376 Together.append(liftfield[k].clone());
377 }
378 }
379
380 if (NUmodes != 0)
381 {
382 for (label k = 0; k < NUmodes; k++)
383 {
384 Together.append(Umodes[k].clone());
385 }
386 }
387
388 if (NSUPmodes != 0)
389 {
390 for (label k = 0; k < NSUPmodes; k++)
391 {
392 Together.append(supmodes[k].clone());
393 }
394 }
395
396 // Project everything
397 for (label i = 0; i < BTsize; i++)
398 {
399 for (label j = 0; j < BTsize; j++)
400 {
401 BT_matrix(i, j) = fvc::domainIntegrate(Together[i] & (fvc::div(dev((T(fvc::grad(
402 Together[j]))))))).value();
403 }
404 }
405
406 // Export the matrix
407 ITHACAstream::exportMatrix(BT_matrix, "BT_matrix", "python",
408 "./ITHACAoutput/Matrices/");
409 ITHACAstream::exportMatrix(BT_matrix, "BT_matrix", "matlab",
410 "./ITHACAoutput/Matrices/");
411 ITHACAstream::exportMatrix(BT_matrix, "BT_matrix", "eigen",
412 "./ITHACAoutput/Matrices/");
413 return BT_matrix;
414}
415
417 label NTmodes, label Nnutmodes)
418{
419 label Stsize = NTmodes + liftfieldT.size();
420 List <Eigen::MatrixXd> S_matrix;
421 S_matrix.setSize(Stsize);
422
423 for (label j = 0; j < Stsize; j++)
424 {
425 S_matrix[j].resize(Nnutmodes, Stsize);
426 }
427
428 PtrList<volScalarField> Together(0);
429
430 // Create PTRLIST with lift, velocities and supremizers
431
432 if (liftfieldT.size() != 0)
433 {
434 for (label k = 0; k < liftfieldT.size(); k++)
435 {
436 Together.append(liftfieldT[k].clone());
437 }
438 }
439
440 if (NTmodes != 0)
441 {
442 for (label k = 0; k < NTmodes; k++)
443 {
444 Together.append(Tmodes[k].clone());
445 }
446 }
447
448 for (label i = 0; i < Stsize; i++)
449 {
450 Info << "Filling layer number " << i + 1 << " in the matrix S_matrix" << endl;
451
452 for (label j = 0; j < Nnutmodes; j++)
453 {
454 for (label k = 0; k < Stsize; k++)
455 {
456 S_matrix[i](j, k) = fvc::domainIntegrate(Together[i] * fvc::laplacian(
457 nuTmodes[j], Together[k])).value();
458 }
459 }
460 }
461
462 // Export the matrix
463 ITHACAstream::exportMatrix(S_matrix, "S_matrix", "python",
464 "./ITHACAoutput/Matrices/");
465 ITHACAstream::exportMatrix(S_matrix, "S_matrix", "matlab",
466 "./ITHACAoutput/Matrices/");
467 ITHACAstream::exportMatrix(S_matrix, "S_matrix", "eigen",
468 "./ITHACAoutput/Matrices/S");
469 return S_matrix;
470}
471
472void UnsteadyNSTTurb::projectSUP(fileName folder, label NU, label NP,
473 label NSUP, label Nnut, label NT)
474{
475 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
476 {
477 B_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/B_mat.txt");
478 C_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/C", "C");
479 K_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/K_mat.txt");
480 P_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/P_mat.txt");
481 M_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/M_mat.txt");
482 BT_matrix =
483 ITHACAstream::readMatrix("./ITHACAoutput/Matrices/BT_matrix_mat.txt");
484 CT1_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/CT1",
485 "CT1_matrix");
486 CT2_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/CT2",
487 "CT2_matrix");
488 Q_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/Q", "Q");
489 Y_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/Y_mat.txt");
490 S_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/S_mat.txt");
491 MT_matrix = ITHACAstream::readMatrix("./ITHACAoutput/Matrices/MT_mat.txt");
492 }
493 else
494 {
495 NUmodes = NU;
496 NPmodes = NP;
497 NSUPmodes = NSUP;
498 Nnutmodes = Nnut;
499 NTmodes = NT;
512 }
513
515 C_total_matrix.setSize(C_matrix.size());
516
517 for (label i = 0; i < C_matrix.size(); i++)
518 {
520 }
521
522 NUmodes = NU;
523 NPmodes = NP;
524 NSUPmodes = NSUP;
525 Nnutmodes = Nnut;
526 // Get the coeffs for interpolation (the orthonormal one is used because basis are orthogonal)
527 Eigen::MatrixXd Ncoeff = ITHACAutilities::getCoeffs(nutFields, nuTmodes);
528 ITHACAstream::exportMatrix(Ncoeff, "Ncoeff", "python",
529 "./ITHACAoutput/Matrices/");
530 SAMPLES.resize(Nnutmodes);
531 rbfsplines.resize(Nnutmodes);
532
533 for (label i = 0; i < Nnutmodes; i++)
534 {
535 SAMPLES[i] = new SPLINTER::DataTable(1, 1);
536
537 for (label j = 0; j < Ncoeff.cols(); j++)
538 {
539 SAMPLES[i]->addSample(mu.row(j), Ncoeff(i, j));
540 }
541
542 rbfsplines[i] = new SPLINTER::RBFSpline(*SAMPLES[i],
543 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
544 std::cout << "Constructing RadialBasisFunction for mode " << i + 1 << std::endl;
545 }
546}
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Header file of the unsteadyNSTTurb class.
label Nnutmodes
Number of viscoisty modes used for the projection.
scalar Prt
Scalar to store the turbulent Prandtl number.
Eigen::MatrixXd BT_matrix
Turbulent viscosity term.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
List< Eigen::MatrixXd > CT1_matrix
Turbulent viscosity term.
bool checkWrite(Time &timeObject)
Function to check if the solution must be exported.
List< Eigen::MatrixXd > S_matrix
Turbulent diffusivity term.
scalar Pr
Scalar to store the Prandtl number.
List< Eigen::MatrixXd > C_total_matrix
Total C Matrix.
void projectSUP(fileName folder, label NU, label NP, label NSUP, label Nnut, label NT)
Project using a supremizer approach.
List< Eigen::MatrixXd > CT2_matrix
Turbulent viscosity term.
Eigen::MatrixXd B_total_matrix
Total B Matrix.
autoPtr< volScalarField > _alphat
Turbulent thermal diffusivity.
List< Eigen::MatrixXd > turbulenceTerm1(label NUmodes, label NSUPmodes, label Nnutmodes)
CT1 added matrix for the turbulence treatement.
List< Eigen::MatrixXd > turbulenceTerm2(label NUmodes, label NSUPmodes, label Nnutmodes)
CT2 added matrix for the turbulence treatement.
std::vector< SPLINTER::RBFSpline * > rbfsplines
Create a SAMPLES for interpolation.
Eigen::MatrixXd BTturbulence(label NUmodes, label NSUPmodes)
BT added matrix for the turbulence treatement.
std::vector< SPLINTER::DataTable * > SAMPLES
Create a Rbf splines for interpolation.
PtrList< volScalarField > nuTmodes
List of POD modes for eddy viscosity.
autoPtr< dimensionedScalar > _Prt
Turbulent Prandtl number.
autoPtr< volScalarField > _nut
Eddy viscosity field.
List< Eigen::MatrixXd > temperatureTurbulenceTerm(label NTmodes, label Nnutmodes)
S added matrix for the thermal turbulence treatement.
autoPtr< dimensionedScalar > _Pr
Prandtl number.
void liftSolveT()
Perform a lift solve for temperature field.
UnsteadyNSTTurb()
Construct Null.
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.
Eigen::MatrixXd mu
Row matrix of parameters.
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 diffusive_term(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Term.
Definition steadyNS.C:861
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
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
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 divergence_term(label NUmodes, label NPmodes, label NSUPmodes)
Divergence Term (supremizer approach)
Definition steadyNS.C:1073
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
Eigen::MatrixXd B_matrix
Diffusion term.
Definition steadyNS.H:157
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
List< Eigen::MatrixXd > C_matrix
Non linear term.
Definition steadyNS.H:166
List< Eigen::MatrixXd > convective_term(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term.
Definition steadyNS.C:955
Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPmodes, label NSUPmodes)
Gradient of pressure.
Definition steadyNS.C:923
Eigen::MatrixXd P_matrix
Div of velocity.
Definition steadyNS.H:170
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition steadyNS.H:160
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
Eigen::MatrixXd Y_matrix
Gradient of pressure matrix.
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
Definition unsteadyNST.H:99
List< Eigen::MatrixXd > convective_term_temperature(label NUmodes, label NTmodes, label NSUPmodes)
Convective Term for Temperature.
List< Eigen::MatrixXd > Q_matrix
Non linear term.
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of the temperature lifting functions.
scalar timeStep
Time step of the simulation.
Definition unsteadyNST.H:81
autoPtr< Time > _runTime
Time.
label NTmodes
Number of temperature modes used for the projection.
scalar finalTime
Final time (final time of the simulation and consequently of the acquisition of the snapshots)
Definition unsteadyNST.H:78
Eigen::MatrixXd diffusive_term_temperature(label NUmodes, label NTmodes, label NSUPmodes)
Diffusive Term for Temperature.
autoPtr< dimensionedScalar > _nu
dimensionedScalar nu;
Eigen::MatrixXd mass_term_temperature(label NUmodes, label NTmodes, label NSUPmodes)
Mass Term for Temperature.
scalar nextWrite
Auxiliary variable to store the next writing instant.
Definition unsteadyNST.H:87
scalar startTime
Start Time (initial time to start storing the snapshots)
Definition unsteadyNST.H:75
autoPtr< volScalarField > _T
Temperature field.
autoPtr< pisoControl > _piso
pisoControl
scalar writeEvery
Time step of the writing procedure.
Definition unsteadyNST.H:84
Eigen::MatrixXd MT_matrix
Mass Matrix T.
autoPtr< fvMesh > _mesh
Mesh.
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition unsteadyNS.H:76
volScalarField & T
Definition createT.H:46
fvScalarMatrix & TEqn
Definition TEqn.H:15
volScalarField & alphat
dimensionedScalar & nu
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 ...
List< Eigen::MatrixXd > readMatrix(word folder, word mat_name)
Read a three dimensional matrix from a txt file in 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.
bool check_sup()
Check if the supremizer folder exists.
simpleControl simple(mesh)
surfaceScalarField & phi
volVectorField & U
volScalarField & p
singlePhaseTransportModel & laminarTransport
label i
Definition pEqn.H:46