Loading...
Searching...
No Matches
unsteadyNST.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
36#include "unsteadyNST.H"
37
38
40
41// Construct from zero
42unsteadyNST::unsteadyNST(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#include "fvOptions.H"
58 _piso = autoPtr<pisoControl>
59 (
60 new pisoControl
61 (
62 mesh
63 )
64 );
65#include "createFields.H"
66#include "createFvOptions.H"
70}
71
72// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73#include "fvCFD.H"
74
75void unsteadyNST::truthSolve(List<scalar> mu_now)
76
77{
78 Time& runTime = _runTime();
79 surfaceScalarField& phi = _phi();
80 fvMesh& mesh = _mesh();
81#include "initContinuityErrs.H"
82 pisoControl& piso = _piso();
83 volScalarField& p = _p();
84 volVectorField& U = _U();
85 volScalarField& T = _T();
86 dimensionedScalar DT = _DT();
87 dimensionedScalar nu = _nu();
88 instantList Times = runTime.times();
89 runTime.setEndTime(finalTime);
90 // Perform a TruthSolve
91 bool WRITE;
92 runTime.setTime(Times[1], 1);
93 runTime.setDeltaT(timeStep);
95
96 // Start the time loop
97 while (runTime.run())
98 {
99 runTime.setEndTime(finalTime);
100 runTime++;
101 Info << "Time = " << runTime.timeName() << nl << endl;
102#include "CourantNo.H"
103 // Momentum predictor
104 fvVectorMatrix UEqn
105 (
106 fvm::ddt(U)
107 + fvm::div(phi, U)
108 - fvm::laplacian(nu, U)
109 );
110
111 if (piso.momentumPredictor())
112 {
113 solve(UEqn == -fvc::grad(p));
114 }
115
116 // --- Pressure-velocity PISO corrector loop
117 while (piso.correct())
118 {
119 {
120#include "pEqn.H"
121 }
122 }
123
124 {
125#include "TEqn.H"
126 TEqn.solve();
127 }
128
129 //runTime.write();
130 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
131 << " ClockTime = " << runTime.elapsedClockTime() << " s"
132 << nl << endl;
133 //laminarTransport.correct();
134 //turbulence->correct();
135 WRITE = checkWrite(runTime);
136
137 if (WRITE)
138 {
139 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
140 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
141 ITHACAstream::exportSolution(T, name(counter), "./ITHACAoutput/Offline/");
142 std::ofstream of("./ITHACAoutput/Offline/" + name(counter) + "/" +
143 runTime.timeName());
144 Ufield.append(U.clone());
145 Pfield.append(p.clone());
146 Tfield.append(T.clone());
147 counter++;
149 writeMu(mu_now);
150 }
151 }
152}
153
154bool unsteadyNST::checkWrite(Time& timeObject)
155{
156 scalar diffnow = mag(nextWrite - atof(timeObject.timeName().c_str()));
157 scalar diffnext = mag(nextWrite - atof(timeObject.timeName().c_str()) -
158 timeObject.deltaTValue());
159
160 if ( diffnow < diffnext)
161 {
162 return true;
163 }
164 else
165 {
166 return false;
167 }
168}
169
170// Method to compute the lifting function for temperature
172{
173 for (label k = 0; k < inletIndexT.rows(); k++)
174 {
175 Time& runTime = _runTime();
176 fvMesh& mesh = _mesh();
177 volScalarField T = _T();
178 volScalarField p = _p();
179 volVectorField U = _U();
180 dimensionedScalar DT = _DT();
181 label BCind = inletIndexT(k, 0);
182 volScalarField Tlift("Tlift" + name(k), T);
183 volVectorField Ulift("Ulift" + name(k), U);
184 pisoControl potentialFlow(mesh, "potentialFlow");
185 instantList Times = runTime.times();
186 runTime.setTime(Times[1], 1);
187 simpleControl simple(mesh);
188 dimensionedScalar nu = _nu();
189 runTime.setEndTime(finalTime);
190 scalar t1 = 1;
191 scalar t0 = 0;
192
193 for (label j = 0; j < T.boundaryField().size(); j++)
194 {
195 if (j == BCind)
196 {
197 assignBC(Tlift, j, t1);
198 }
199 else if (T.boundaryField()[BCind].type() == "zeroGradient")
200 {
201 assignBC(Tlift, j, t0);
202 assignIF(Tlift, t1);
203 }
204 else if (T.boundaryField()[BCind].type() == "fixedValue")
205 {
206 assignBC(Tlift, j, t0);
207 assignIF(Tlift, t0);
208 }
209 else
210 {
211 }
212 }
213
214 while (simple.correctNonOrthogonal())
215 {
216 fvScalarMatrix TEqn
217 (
218 fvm::ddt(Tlift)
219 - fvm::laplacian(DT, Tlift)
220 );
221 TEqn.solve();
222 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
223 << " ClockTime = " << runTime.elapsedClockTime() << " s"
224 << nl << endl;
225 }
226
227 scalar totalTime = mesh.time().value();
228 scalar dt = mesh.time().deltaTValue();
229 forAll(Tlift, i)
230 {
231 Tlift[i] = (totalTime * Tlift[i] + dt * Tlift[i] ) / (totalTime + dt);
232 }
233 Tlift.write();
234 liftfieldT.append(Tlift.clone());
235 }
236}
237
238
239// * * * * * * * * * * * * * * Projection Methods * * * * * * * * * * * * * * //
240
241void unsteadyNST::projectSUP(fileName folder, label NU, label NP, label NT,
242 label NSUP)
243{
244 NUmodes = NU;
245 NPmodes = NP;
246 NTmodes = NT;
247 NSUPmodes = NSUP;
248
249 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
250 {
251 word M_str = "M_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
252 NSUPmodes);
253
254 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + M_str))
255 {
256 ITHACAstream::ReadDenseMatrix(M_matrix, "./ITHACAoutput/Matrices/", M_str);
257 }
258 else
259 {
261 }
262
263 word B_str = "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
264 NSUPmodes);
265
266 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + B_str))
267 {
268 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", B_str);
269 }
270 else
271 {
273 }
274
275 word C_str = "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
276 NSUPmodes) + "_t";
277
278 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + C_str))
279 {
280 ITHACAstream::ReadDenseTensor(C_tensor, "./ITHACAoutput/Matrices/", C_str);
281 }
282 else
283 {
285 }
286
287 word K_str = "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
288 NSUPmodes) + "_" + name(NPmodes);
289
290 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + K_str))
291 {
292 ITHACAstream::ReadDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/", K_str);
293 }
294 else
295 {
297 }
298
299 word MT_str = "MT_" + name(liftfieldT.size()) + "_" + name(NTmodes);
300
301 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + MT_str))
302 {
303 ITHACAstream::ReadDenseMatrix(MT_matrix, "./ITHACAoutput/Matrices/", MT_str);
304 }
305 else
306 {
308 }
309
311 word Y_str = "Y_" + name(liftfieldT.size()) + "_" + name(NTmodes);
312
313 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + Y_str))
314 {
315 ITHACAstream::ReadDenseMatrix(Y_matrix, "./ITHACAoutput/Matrices/", Y_str);
316 }
317 else
318 {
320 }
321
322 word P_str = "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
323 NSUPmodes) + "_" + name(NPmodes);
324
325 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + P_str))
326 {
327 ITHACAstream::ReadDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/", P_str);
328 }
329 else
330 {
332 }
333 }
334 else
335 {
336 L_U_SUPmodes.resize(0);
337
338 if (liftfield.size() != 0)
339 {
340 for (label k = 0; k < liftfield.size(); k++)
341 {
342 L_U_SUPmodes.append(liftfield[k].clone());
343 }
344 }
345
346 if (NUmodes != 0)
347 {
348 for (label k = 0; k < NUmodes; k++)
349 {
350 L_U_SUPmodes.append(Umodes[k].clone());
351 }
352 }
353
354 if (NSUPmodes != 0)
355 {
356 for (label k = 0; k < NSUPmodes; k++)
357 {
358 L_U_SUPmodes.append(supmodes[k].clone());
359 }
360 }
361
362 L_T_modes.resize(0);
363
364 if (liftfieldT.size() != 0)
365 {
366 for (label k = 0; k < liftfieldT.size(); k++)
367 {
368 L_T_modes.append(liftfieldT[k].clone());
369 }
370 }
371
372 if (NTmodes != 0)
373 {
374 for (label k = 0; k < NTmodes; k++)
375 {
376 L_T_modes.append(Tmodes[k].clone());
377 }
378 }
379
388 }
389}
390
391
392List< Eigen::MatrixXd > unsteadyNST::convective_term_temperature(label NUmodes,
393 label NTmodes, label NSUPmodes)
394{
395 label Qsize = NUmodes + liftfield.size() + NSUPmodes;
396 label Qsizet = NTmodes + liftfieldT.size() ;
397 List < Eigen::MatrixXd > Q_matrix;
398 Q_matrix.setSize(Qsizet);
399
400 for (label j = 0; j < Qsizet; j++)
401 {
402 Q_matrix[j].resize(Qsize, Qsizet);
403 }
404
405 {
406 for (label i = 0; i < Qsizet; i++)
407 {
408 for (label j = 0; j < Qsize; j++)
409 {
410 for (label k = 0; k < Qsizet; k++)
411 {
412 Q_matrix[i](j, k) = fvc::domainIntegrate(L_T_modes[i] * fvc::div(
413 fvc::interpolate(L_U_SUPmodes[j]) & L_U_SUPmodes[j].mesh().Sf(),
414 L_T_modes[k])).value();
415 }
416 }
417 }
418 }
419
420 // Export the matrix
421 ITHACAstream::exportMatrix(Q_matrix, "Q", "python", "./ITHACAoutput/Matrices/");
422 ITHACAstream::exportMatrix(Q_matrix, "Q", "matlab", "./ITHACAoutput/Matrices/");
423 ITHACAstream::exportMatrix(Q_matrix, "Q", "eigen", "./ITHACAoutput/Matrices/Q");
424 return Q_matrix;
425}
426
427Eigen::MatrixXd unsteadyNST::diffusive_term_temperature(label NUmodes,
428 label NTmodes, label NSUPmodes)
429{
430 label Ysize = NTmodes + liftfieldT.size();
431 Eigen::MatrixXd Y_matrix(Ysize, Ysize);
432 dimensionedScalar DT = _DT();
433 {
434 for (label i = 0; i < Ysize; i++)
435 {
436 for (label j = 0; j < Ysize; j++)
437 {
438 Y_matrix(i, j) = fvc::domainIntegrate(L_T_modes[i] * fvc::laplacian(
439 dimensionedScalar("1", dimless, 1), L_T_modes[j])).value();
440 }
441 }
442 }
443 // Export the matrix
444 ITHACAstream::SaveDenseMatrix(Y_matrix, "./ITHACAoutput/Matrices/",
445 "Y_" + name(liftfieldT.size()) + "_" + name(NTmodes));
446 return Y_matrix;
447}
448
449Eigen::MatrixXd unsteadyNST::mass_term_temperature(label NUmodes, label NTmodes,
450 label NSUPmodes)
451{
452 label Ysize = NTmodes + liftfieldT.size();
453 Eigen::MatrixXd MT_matrix(Ysize, Ysize);
454
455 for (label i = 0; i < Ysize; i++)
456 {
457 for (label j = 0; j < Ysize; j++)
458 {
459 MT_matrix(i, j) = fvc::domainIntegrate(L_T_modes[i] * L_T_modes[j]).value();
460 }
461 }
462
463 // Export the matrix
464 ITHACAstream::SaveDenseMatrix(MT_matrix, "./ITHACAoutput/Matrices/",
465 "MT_" + name(liftfieldT.size()) + "_" + name(NTmodes));
466 return MT_matrix;
467}
468// Calculate lifting function for velocity
470{
471 for (label k = 0; k < inletIndex.rows(); k++)
472 {
473 Time& runTime = _runTime();
474 surfaceScalarField& phi = _phi();
475 fvMesh& mesh = _mesh();
476 volScalarField p = _p();
477 volVectorField U = _U();
478 dimensionedScalar nu = _nu();
479 IOMRFZoneList& MRF = _MRF();
480 label BCind = inletIndex(k, 0);
481 volVectorField Ulift("Ulift" + name(k), U);
482 instantList Times = runTime.times();
483 runTime.setTime(Times[1], 1);
484 pisoControl potentialFlow(mesh, "potentialFlow");
485 Info << "Solving a lifting Problem" << endl;
486 Vector<double> v1(0, 0, 0);
487 v1[inletIndex(k, 1)] = 1;
488 Vector<double> v0(0, 0, 0);
489
490 for (label j = 0; j < U.boundaryField().size(); j++)
491 {
492 if (j == BCind)
493 {
494 assignBC(Ulift, j, v1);
495 }
496 else if (U.boundaryField()[BCind].type() == "fixedValue")
497 {
498 assignBC(Ulift, j, v0);
499 }
500 else
501 {
502 }
503
504 assignIF(Ulift, v0);
505 phi = linearInterpolate(Ulift) & mesh.Sf();
506 }
507
508 Info << "Constructing velocity potential field Phi\n" << endl;
509 volScalarField Phi
510 (
511 IOobject
512 (
513 "Phi",
514 runTime.timeName(),
515 mesh,
516 IOobject::READ_IF_PRESENT,
517 IOobject::NO_WRITE
518 ),
519 mesh,
520 dimensionedScalar("Phi", dimLength * dimVelocity, 0),
521 p.boundaryField().types()
522 );
523 label PhiRefCell = 0;
524 scalar PhiRefValue = 0;
526 (
527 Phi,
528 potentialFlow.dict(),
529 PhiRefCell,
530 PhiRefValue
531 );
532 mesh.setFluxRequired(Phi.name());
533 runTime.functionObjects().start();
534 MRF.makeRelative(phi);
535 adjustPhi(phi, Ulift, p);
536
537 while (potentialFlow.correctNonOrthogonal())
538 {
539 fvScalarMatrix PhiEqn
540 (
541 fvm::laplacian(dimensionedScalar("1", dimless, 1), Phi)
542 //fvm::laplacian(nu, Phi)
543 ==
544 fvc::div(phi)
545 );
546 PhiEqn.setReference(PhiRefCell, PhiRefValue);
547 PhiEqn.solve();
548
549 if (potentialFlow.finalNonOrthogonalIter())
550 {
551 phi -= PhiEqn.flux();
552 }
553 }
554
555 MRF.makeAbsolute(phi);
556 Info << "Continuity error = "
557 << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
558 << endl;
559 Ulift = fvc::reconstruct(phi);
560 Ulift.correctBoundaryConditions();
561 Info << "Interpolated velocity error = "
562 << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
563 / sum(mesh.magSf())).value()
564 << endl;
565 Ulift.write();
566 liftfield.append(Ulift.clone());
567 }
568}
569
570
571
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
TEqn solve()
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.
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 diffusive_term(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Term.
Definition steadyNS.C:861
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
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
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
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 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
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 > 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.
unsteadyNST()
Constructors.
Definition unsteadyNST.C:39
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
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NTmodes, label NSUPmodes)
Specific variable for the unstationary case.
bool checkWrite(Time &timeObject)
Function to check if the solution must be exported.
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
void liftSolve()
Perform a lift solve for velocity field.
Eigen::MatrixXd MT_matrix
Mass Matrix T.
autoPtr< dimensionedScalar > _DT
dimensionedScalar DT;
autoPtr< fvMesh > _mesh
Mesh.
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
PtrList< volScalarField > L_T_modes
List of pointers containing the lift for temperature and the temperature field.
void liftSolveT()
Perform a lift solve for temperature field.
volScalarField & T
Definition createT.H:46
fvScalarMatrix & TEqn
Definition TEqn.H:15
volScalarField v0(v)
dimensionedScalar & nu
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.
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.
simpleControl simple(mesh)
surfaceScalarField & phi
volVectorField & U
fvVectorMatrix & UEqn
Definition UEqn.H:37
volScalarField & p
setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue)
adjustPhi(phiHbyA, U, p)
label i
Definition pEqn.H:46
Header file of the unsteadyNST class.