Loading...
Searching...
No Matches
UnsteadyNSTurbIntrusive.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
32
35
36// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
37
38// Construct Null
40
41// Construct from zero
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#include "createFields.H"
65#include "createFvOptions.H"
66 ITHACAdict = new IOdictionary
67 (
68 IOobject
69 (
70 "ITHACAdict",
71 runTime.system(),
72 mesh,
73 IOobject::MUST_READ,
74 IOobject::NO_WRITE
75 )
76 );
77 tolerance = ITHACAdict->lookupOrDefault<scalar>("tolerance", 1e-5);
78 maxIter = ITHACAdict->lookupOrDefault<scalar>("maxIter", 1000);
79 bcMethod = ITHACAdict->lookupOrDefault<word>("bcMethod", "lift");
81 ITHACAdict->lookupOrDefault<word>("timeDerivativeSchemeOrder", "second");
82 M_Assert(bcMethod == "lift" || bcMethod == "penalty",
83 "The BC method must be set to lift or penalty in ITHACAdict");
85 || timeDerivativeSchemeOrder == "second",
86 "The time derivative approximation must be set to either first or second order scheme in ITHACAdict");
91}
92
93// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94
95void UnsteadyNSTurbIntrusive::truthSolve(List<scalar> mu_now)
96{
97 Time& runTime = _runTime();
98 surfaceScalarField& phi = _phi();
99 fvMesh& mesh = _mesh();
100#include "initContinuityErrs.H"
101 fv::options& fvOptions = _fvOptions();
102 pimpleControl& pimple = _pimple();
103 volScalarField p = _p();
104 volVectorField U = _U();
105 volScalarField nut = _nut();
106 IOMRFZoneList& MRF = _MRF();
107 singlePhaseTransportModel& laminarTransport = _laminarTransport();
108 instantList Times = runTime.times();
109 runTime.setEndTime(finalTime);
110 // Perform a TruthSolve
111 runTime.setTime(Times[1], 1);
112 runTime.setDeltaT(timeStep);
114 // Initialize Nsnapshots
115 label nsnapshots = 0;
116
117 // Start the time loop
118 while (runTime.run())
119 {
120#include "readTimeControls.H"
121#include "CourantNo.H"
122#include "setDeltaT.H"
123 runTime.setEndTime(finalTime + timeStep);
124 Info << "Time = " << runTime.timeName() << nl << endl;
125
126 // --- Pressure-velocity PIMPLE corrector loop
127 while (pimple.loop())
128 {
129#include "UEqn.H"
130
131 // --- Pressure corrector loop
132 while (pimple.correct())
133 {
134#include "pEqn.H"
135 }
136
137 if (pimple.turbCorr())
138 {
139 laminarTransport.correct();
140 turbulence->correct();
141 }
142 }
143
144 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
145 << " ClockTime = " << runTime.elapsedClockTime() << " s"
146 << nl << endl;
147
148 if (checkWrite(runTime))
149 {
150 nsnapshots += 1;
151 volScalarField nut = turbulence->nut().ref();
152 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
153 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
154 ITHACAstream::exportSolution(nut, name(counter), "./ITHACAoutput/Offline/");
155 std::ofstream of("./ITHACAoutput/Offline/" + name(counter) + "/" +
156 runTime.timeName());
157 Ufield.append(U.clone());
158 Pfield.append(p.clone());
159 nutFields.append(nut.clone());
160 counter++;
162 writeMu(mu_now);
163 // --- Fill in the mu_samples with parameters (time, mu) to be used for the PODI sample points
164 mu_samples.conservativeResize(mu_samples.rows() + 1, mu_now.size() + 1);
165 mu_samples(mu_samples.rows() - 1, 0) = atof(runTime.timeName().c_str());
166
167 for (label i = 0; i < mu_now.size(); i++)
168 {
169 mu_samples(mu_samples.rows() - 1, i + 1) = mu_now[i];
170 }
171 }
172
173 runTime++;
174 }
175
176 // Resize to Unitary if not initialized by user (i.e. non-parametric problem)
177 if (mu.cols() == 0)
178 {
179 mu.resize(1, 1);
180 }
181
182 if (mu_samples.rows() == nsnapshots * mu.cols())
183 {
184 ITHACAstream::exportMatrix(mu_samples, "mu_samples", "eigen",
185 "./ITHACAoutput/Offline");
186 }
187}
188
190 label nModes)
191{
192 Eigen::Tensor<double, 3> ct1Tensor;
193 ct1Tensor.resize(nModes, nModes, nModes);
194
195 for (label i = 0; i < nModes; i++)
196 {
197 for (label j = 0; j < nModes; j++)
198 {
199 for (label k = 0; k < nModes; k++)
200 {
201 ct1Tensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & fvc::laplacian(
202 nutModes[j], Umodes[k])).value();
203 }
204 }
205 }
206
207 if (Pstream::parRun())
208 {
209 reduce(ct1Tensor, sumOp<Eigen::Tensor<double, 3 >> ());
210 }
211
212 // Export the tensor
213 ITHACAstream::SaveDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/",
214 "ct1_" + name(nModes) + "_t");
215 return ct1Tensor;
216}
217
219 label nModes)
220{
221 Eigen::Tensor<double, 3> ct2Tensor;
222 ct2Tensor.resize(nModes, nModes, nModes);
223
224 for (label i = 0; i < nModes; i++)
225 {
226 for (label j = 0; j < nModes; j++)
227 {
228 for (label k = 0; k < nModes; k++)
229 {
230 ct2Tensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & (fvc::div(
231 nutModes[j] * dev((fvc::grad(Umodes[k]))().T())))).value();
232 }
233 }
234 }
235
236 if (Pstream::parRun())
237 {
238 reduce(ct2Tensor, sumOp<Eigen::Tensor<double, 3 >> ());
239 }
240
241 // Export the tensor
242 ITHACAstream::SaveDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/",
243 "ct2_" + name(nModes) + "_t");
244 return ct2Tensor;
245}
246
247Eigen::MatrixXd UnsteadyNSTurbIntrusive::btTurbulence(label nModes)
248{
249 Eigen::MatrixXd btMatrix(nModes, nModes);
250 btMatrix = btMatrix * 0;
251
252 // Project everything
253 for (label i = 0; i < nModes; i++)
254 {
255 for (label j = 0; j < nModes; j++)
256 {
257 btMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & (fvc::div(dev((T(
258 fvc::grad(
259 Umodes[j]))))))).value();
260 }
261 }
262
263 // Export the matrix
264 ITHACAstream::SaveDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/",
265 "bt_" + name(nModes) + "_t");
266 return btMatrix;
267}
268
269void UnsteadyNSTurbIntrusive::project(fileName folder, label nModes)
270{
271 nModesOnline = nModes;
272 L_U_SUPmodes.resize(0);
273 if (nModes != 0)
274 {
275 for (label k = 0; k < nModes; k++)
276 {
277 L_U_SUPmodes.append(Umodes[k].clone());
278 }
279 }
280
281 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
282 {
283 word bStr = "b_" + name(nModes);
284 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bStr))
285 {
286 ITHACAstream::ReadDenseMatrix(bMatrix, "./ITHACAoutput/Matrices/", bStr);
287 }
288
289 else
290 {
291 bMatrix = diffusiveTerm(nModes);
292 }
293 word btStr = "bt_" + name(nModes);
294 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + btStr))
295 {
296 ITHACAstream::ReadDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/", btStr);
297 }
298
299 else
300 {
301 btMatrix = btTurbulence(nModes);
302 }
303 word kStr = "k_" + name(nModes);
304 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + kStr))
305 {
306 ITHACAstream::ReadDenseMatrix(kMatrix, "./ITHACAoutput/Matrices/", kStr);
307 }
308
309 else
310 {
312 }
313 word cStr = "c_" + name(nModes) + "_t";
314 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + cStr))
315 {
316 ITHACAstream::ReadDenseTensor(convTensor, "./ITHACAoutput/Matrices/", cStr);
317 }
318
319 else
320 {
321 convTensor = convectiveTerm(nModes);
322 }
323 word ct1Str = "ct1_" + name(nModes) + "_t";
324 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct1Str))
325 {
326 ITHACAstream::ReadDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/", ct1Str);
327 }
328
329 else
330 {
332 }
333 word ct2Str = "ct2_" + name(nModes) + "_t";
334 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct2Str))
335 {
336 ITHACAstream::ReadDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/", ct2Str);
337 }
338
339 else
340 {
342 }
343 if (bcMethod == "penalty")
344 {
345 bcVelVec = bcVelocityVec(nModes);
346 bcVelMat = bcVelocityMat(nModes);
347 }
348 }
349
350 else
351 {
352 bMatrix = diffusiveTerm(nModes);
353 convTensor = convectiveTerm(nModes);
355 btMatrix = btTurbulence(nModes);
358
359 if (bcMethod == "penalty")
360 {
361 bcVelVec = bcVelocityVec(nModes);
362 bcVelMat = bcVelocityMat(nModes);
363 }
364 }
365
366 // Export the matrices
367 if (para->exportPython)
368 {
369 ITHACAstream::exportMatrix(bMatrix, "b", "python", "./ITHACAoutput/Matrices/");
370 ITHACAstream::exportMatrix(kMatrix, "k", "python", "./ITHACAoutput/Matrices/");
371 ITHACAstream::exportMatrix(btMatrix, "bt", "python",
372 "./ITHACAoutput/Matrices/");
374 "./ITHACAoutput/Matrices/");
375 ITHACAstream::exportTensor(ct1Tensor, "ct1", "python",
376 "./ITHACAoutput/Matrices/");
377 ITHACAstream::exportTensor(ct2Tensor, "ct2", "python",
378 "./ITHACAoutput/Matrices/");
379 }
380
381 if (para->exportMatlab)
382 {
383 ITHACAstream::exportMatrix(bMatrix, "b", "matlab", "./ITHACAoutput/Matrices/");
384 ITHACAstream::exportMatrix(kMatrix, "k", "matlab", "./ITHACAoutput/Matrices/");
385 ITHACAstream::exportMatrix(btMatrix, "bt", "matlab",
386 "./ITHACAoutput/Matrices/");
388 "./ITHACAoutput/Matrices/");
389 ITHACAstream::exportTensor(ct1Tensor, "ct1", "matlab",
390 "./ITHACAoutput/Matrices/");
391 ITHACAstream::exportTensor(ct2Tensor, "ct2", "matlab",
392 "./ITHACAoutput/Matrices/");
393 }
394
395 if (para->exportTxt)
396 {
397 ITHACAstream::exportMatrix(bMatrix, "b", "eigen", "./ITHACAoutput/Matrices/");
398 ITHACAstream::exportMatrix(kMatrix, "k", "eigen", "./ITHACAoutput/Matrices/");
399 ITHACAstream::exportMatrix(btMatrix, "bt", "eigen", "./ITHACAoutput/Matrices/");
401 "./ITHACAoutput/Matrices/c");
402 ITHACAstream::exportTensor(ct1Tensor, "ct1_", "eigen",
403 "./ITHACAoutput/Matrices/ct1");
404 ITHACAstream::exportTensor(ct2Tensor, "ct2_", "eigen",
405 "./ITHACAoutput/Matrices/ct2");
406 }
407
409 cTotalTensor.resize(nModes, nModes, nModes);
411}
412
413void UnsteadyNSTurbIntrusive::projectPPE(fileName folder, label nUModes,
414 label nPModes)
415{
416 NUmodes = nUModes;
417 NPmodes = nPModes;
418 L_U_SUPmodes.resize(0);
419 if (nUModes != 0)
420 {
421 for (label k = 0; k < nUModes; k++)
422 {
423 L_U_SUPmodes.append(Umodes[k].clone());
424 }
425 }
426
427 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
428 {
429 word bStr = "b_" + name(nUModes);
430 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bStr))
431 {
432 ITHACAstream::ReadDenseMatrix(bMatrix, "./ITHACAoutput/Matrices/", bStr);
433 }
434
435 else
436 {
437 bMatrix = diffusiveTerm(nUModes);
438 }
439 word btStr = "bt_" + name(nUModes);
440 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + btStr))
441 {
442 ITHACAstream::ReadDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/", btStr);
443 }
444
445 else
446 {
447 btMatrix = btTurbulence(nUModes);
448 }
449 word kStr = "k_" + name(nUModes) + "_" + name(nPModes);
450 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + kStr))
451 {
452 ITHACAstream::ReadDenseMatrix(kMatrix, "./ITHACAoutput/Matrices/", kStr);
453 }
454
455 else
456 {
457 kMatrix = pressureGradientTerm(nUModes, nPModes);
458 }
459 word D_str = "D_" + name(nPModes);
460 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + D_str))
461 {
462 ITHACAstream::ReadDenseMatrix(D_matrix, "./ITHACAoutput/Matrices/", D_str);
463 }
464
465 else
466 {
467 D_matrix = laplacianPressure(nPModes);
468 }
469 word bc1_str = "BC1_" + name(nUModes) + "_" + name(nPModes);
470 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bc1_str))
471 {
472 ITHACAstream::ReadDenseMatrix(BC1_matrix, "./ITHACAoutput/Matrices/", bc1_str);
473 }
474
475 else
476 {
477 BC1_matrix = pressureBC1(nUModes, nPModes);
478 }
479 word bc2_str = "BC2_" + name(nUModes) + "_" + name(nPModes) + "_t";
480 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bc2_str))
481 {
482 ITHACAstream::ReadDenseTensor(bc2Tensor, "./ITHACAoutput/Matrices/", bc2_str);
483 }
484
485 else
486 {
487 bc2Tensor = pressureBC2(nUModes, nPModes);
488 }
489 word bc3_str = "BC3_" + name(nUModes) + "_" + name(nPModes);
490 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bc3_str))
491 {
492 ITHACAstream::ReadDenseMatrix(BC3_matrix, "./ITHACAoutput/Matrices/", bc3_str);
493 }
494
495 else
496 {
497 BC3_matrix = pressureBC3(nUModes, nPModes);
498 }
499 word cStr = "c_" + name(nUModes) + "_t";
500 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + cStr))
501 {
502 ITHACAstream::ReadDenseTensor(convTensor, "./ITHACAoutput/Matrices/", cStr);
503 }
504
505 else
506 {
507 convTensor = convectiveTerm(nUModes);
508 }
509 word ct1Str = "ct1_" + name(nUModes) + "_t";
510 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct1Str))
511 {
512 ITHACAstream::ReadDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/", ct1Str);
513 }
514
515 else
516 {
517 ct1Tensor = turbulenceTensor1(nUModes);
518 }
519 word ct2Str = "ct2_" + name(nUModes) + "_t";
520 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct2Str))
521 {
522 ITHACAstream::ReadDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/", ct2Str);
523 }
524
525 else
526 {
527 ct2Tensor = turbulenceTensor2(nUModes);
528 }
529 word ct1PPEStr = "ct1PPE_" + name(nUModes) + "_" + name(nPModes) + "_t";
530 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct1PPEStr))
531 {
532 ITHACAstream::ReadDenseTensor(ct1PPETensor, "./ITHACAoutput/Matrices/",
533 ct1PPEStr);
534 }
535
536 else
537 {
538 ct1PPETensor = turbulencePPETensor1(nUModes, nPModes);
539 }
540 word ct2PPEStr = "ct2PPE_" + name(nUModes) + "_" + name(nPModes) + "_t";
541 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct2PPEStr))
542 {
543 ITHACAstream::ReadDenseTensor(ct2PPETensor, "./ITHACAoutput/Matrices/",
544 ct2PPEStr);
545 }
546
547 else
548 {
549 ct2PPETensor = turbulencePPETensor2(nUModes, nPModes);
550 }
551 word G_str = "G_" + name(nUModes) + "_" + name(nPModes) + "_t";
552 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + G_str))
553 {
554 ITHACAstream::ReadDenseTensor(gTensor, "./ITHACAoutput/Matrices/", G_str);
555 }
556
557 else
558 {
559 gTensor = divMomentum(nUModes, nPModes);
560 }
561 if (bcMethod == "penalty")
562 {
563 bcVelVec = bcVelocityVec(nUModes);
564 bcVelMat = bcVelocityMat(nUModes);
565 }
566 }
567
568 else
569 {
570 bMatrix = diffusiveTerm(nUModes);
571 convTensor = convectiveTerm(nUModes);
572 kMatrix = pressureGradientTerm(nUModes, nPModes);
573 btMatrix = btTurbulence(nUModes);
574 ct1Tensor = turbulenceTensor1(nUModes);
575 ct2Tensor = turbulenceTensor2(nUModes);
576 D_matrix = laplacianPressure(nPModes);
577 ct1PPETensor = turbulencePPETensor1(nUModes, nPModes);
578 ct2PPETensor = turbulencePPETensor2(nUModes, nPModes);
579 gTensor = divMomentum(nUModes, nPModes);
580 BC1_matrix = pressureBC1(nUModes, nPModes);
581 bc2Tensor = pressureBC2(nUModes, nPModes);
582 BC3_matrix = pressureBC3(nUModes, nPModes);
583
584 if (bcMethod == "penalty")
585 {
586 bcVelVec = bcVelocityVec(nUModes);
587 bcVelMat = bcVelocityMat(nUModes);
588 }
589 }
590
591 // Export the matrices
592 if (para->exportPython)
593 {
594 ITHACAstream::exportMatrix(bMatrix, "b", "python", "./ITHACAoutput/Matrices/");
595 ITHACAstream::exportMatrix(kMatrix, "k", "python", "./ITHACAoutput/Matrices/");
596 ITHACAstream::exportMatrix(btMatrix, "bt", "python",
597 "./ITHACAoutput/Matrices/");
598 ITHACAstream::exportMatrix(D_matrix, "D", "python", "./ITHACAoutput/Matrices/");
599 ITHACAstream::exportMatrix(BC1_matrix, "BC1", "python",
600 "./ITHACAoutput/Matrices/");
601 ITHACAstream::exportMatrix(BC3_matrix, "BC3", "python",
602 "./ITHACAoutput/Matrices/");
603 ITHACAstream::exportTensor(ct1PPETensor, "ct1PPE", "python",
604 "./ITHACAoutput/Matrices/");
605 ITHACAstream::exportTensor(ct2PPETensor, "ct2PPE", "python",
606 "./ITHACAoutput/Matrices/");
607 ITHACAstream::exportTensor(gTensor, "G", "python", "./ITHACAoutput/Matrices/");
608 ITHACAstream::exportTensor(bc2Tensor, "BC2", "python",
609 "./ITHACAoutput/Matrices/");
611 "./ITHACAoutput/Matrices/");
612 ITHACAstream::exportTensor(ct1Tensor, "ct1", "python",
613 "./ITHACAoutput/Matrices/");
614 ITHACAstream::exportTensor(ct2Tensor, "ct2", "python",
615 "./ITHACAoutput/Matrices/");
616 }
617
618 if (para->exportMatlab)
619 {
620 ITHACAstream::exportMatrix(bMatrix, "b", "matlab", "./ITHACAoutput/Matrices/");
621 ITHACAstream::exportMatrix(kMatrix, "k", "matlab", "./ITHACAoutput/Matrices/");
622 ITHACAstream::exportMatrix(btMatrix, "bt", "matlab",
623 "./ITHACAoutput/Matrices/");
624 ITHACAstream::exportMatrix(D_matrix, "D", "matlab", "./ITHACAoutput/Matrices/");
625 ITHACAstream::exportMatrix(BC1_matrix, "BC1", "matlab",
626 "./ITHACAoutput/Matrices/");
627 ITHACAstream::exportMatrix(BC3_matrix, "BC3", "matlab",
628 "./ITHACAoutput/Matrices/");
629 ITHACAstream::exportTensor(ct1PPETensor, "ct1PPE", "matlab",
630 "./ITHACAoutput/Matrices/");
631 ITHACAstream::exportTensor(ct2PPETensor, "ct2PPE", "matlab",
632 "./ITHACAoutput/Matrices/");
633 ITHACAstream::exportTensor(gTensor, "G", "matlab", "./ITHACAoutput/Matrices/");
634 ITHACAstream::exportTensor(bc2Tensor, "BC2", "matlab",
635 "./ITHACAoutput/Matrices/");
637 "./ITHACAoutput/Matrices/");
638 ITHACAstream::exportTensor(ct1Tensor, "ct1", "matlab",
639 "./ITHACAoutput/Matrices/");
640 ITHACAstream::exportTensor(ct2Tensor, "ct2", "matlab",
641 "./ITHACAoutput/Matrices/");
642 }
643
644 if (para->exportTxt)
645 {
646 ITHACAstream::exportMatrix(bMatrix, "b", "eigen", "./ITHACAoutput/Matrices/");
647 ITHACAstream::exportMatrix(kMatrix, "k", "eigen", "./ITHACAoutput/Matrices/");
648 ITHACAstream::exportMatrix(btMatrix, "bt", "eigen", "./ITHACAoutput/Matrices/");
649 ITHACAstream::exportMatrix(D_matrix, "D", "eigen", "./ITHACAoutput/Matrices/");
651 "./ITHACAoutput/Matrices/");
653 "./ITHACAoutput/Matrices/");
654 ITHACAstream::exportTensor(ct1PPETensor, "ct1PPE_", "eigen",
655 "./ITHACAoutput/Matrices/ct1PPE");
656 ITHACAstream::exportTensor(ct2PPETensor, "ct2PPE_", "eigen",
657 "./ITHACAoutput/Matrices/ct2PPE");
659 "./ITHACAoutput/Matrices/G");
660 ITHACAstream::exportTensor(bc2Tensor, "BC2_", "eigen",
661 "./ITHACAoutput/Matrices/BC2");
663 "./ITHACAoutput/Matrices/c");
664 ITHACAstream::exportTensor(ct1Tensor, "ct1_", "eigen",
665 "./ITHACAoutput/Matrices/ct1");
666 ITHACAstream::exportTensor(ct2Tensor, "ct2_", "eigen",
667 "./ITHACAoutput/Matrices/ct2");
668 }
669
671 cTotalTensor.resize(nUModes, nUModes, nUModes);
673 cTotalPPETensor.resize(nPModes, nUModes, nUModes);
675}
676
677Eigen::MatrixXd UnsteadyNSTurbIntrusive::diffusiveTerm(label nModes)
678{
679 Eigen::MatrixXd bMatrix;
680 bMatrix.resize(nModes, nModes);
681
682 // Project everything
683 for (label i = 0; i < nModes; i++)
684 {
685 for (label j = 0; j < nModes; j++)
686 {
687 bMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & fvc::laplacian(
688 dimensionedScalar("1", dimless, 1), Umodes[j])).value();
689 }
690 }
691
692 if (Pstream::parRun())
693 {
694 reduce(bMatrix, sumOp<Eigen::MatrixXd>());
695 }
696
697 ITHACAstream::SaveDenseMatrix(bMatrix, "./ITHACAoutput/Matrices/",
698 "b_" + name(nModes));
699 return bMatrix;
700}
701
702Eigen::Tensor<double, 3> UnsteadyNSTurbIntrusive::convectiveTerm(label nModes)
703{
704 Eigen::Tensor<double, 3> convTensor;
705 convTensor.resize(nModes, nModes, nModes);
706
707 for (label i = 0; i < nModes; i++)
708 {
709 for (label j = 0; j < nModes; j++)
710 {
711 for (label k = 0; k < nModes; k++)
712 {
713 convTensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & fvc::div(
714 linearInterpolate(Umodes[j]) & Umodes[j].mesh().Sf(),
715 Umodes[k])).value();
716 }
717 }
718 }
719
720 if (Pstream::parRun())
721 {
722 reduce(convTensor, sumOp<Eigen::Tensor<double, 3 >> ());
723 }
724
725 // Export the tensor
726 ITHACAstream::SaveDenseTensor(convTensor, "./ITHACAoutput/Matrices/",
727 "c_" + name(nModes) + "_t");
728 return convTensor;
729}
730
732{
733 Eigen::MatrixXd kMatrix(nModes, nModes);
734
735 // Project everything
736 for (label i = 0; i < nModes; i++)
737 {
738 for (label j = 0; j < nModes; j++)
739 {
740 kMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & fvc::grad(
741 Pmodes[j])).value();
742 }
743 }
744
745 if (Pstream::parRun())
746 {
747 reduce(kMatrix, sumOp<Eigen::MatrixXd>());
748 }
749
750 ITHACAstream::SaveDenseMatrix(kMatrix, "./ITHACAoutput/Matrices/",
751 "k_" + name(nModes));
752 return kMatrix;
753}
754
756 label nPModes)
757{
758 Eigen::MatrixXd kMatrix(nUModes, nPModes);
759
760 // Project everything
761 for (label i = 0; i < nUModes; i++)
762 {
763 for (label j = 0; j < nPModes; j++)
764 {
765 kMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & fvc::grad(
766 Pmodes[j])).value();
767 }
768 }
769
770 if (Pstream::parRun())
771 {
772 reduce(kMatrix, sumOp<Eigen::MatrixXd>());
773 }
774
775 ITHACAstream::SaveDenseMatrix(kMatrix, "./ITHACAoutput/Matrices/",
776 "k_" + name(nUModes) + "_" + name(nPModes));
777 return kMatrix;
778}
779
780Eigen::MatrixXd UnsteadyNSTurbIntrusive::laplacianPressure(label nPModes)
781{
782 Eigen::MatrixXd dMatrix(nPModes, nPModes);
783
784 // Project everything
785 for (label i = 0; i < nPModes; i++)
786 {
787 for (label j = 0; j < nPModes; j++)
788 {
789 dMatrix(i, j) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & fvc::grad(
790 Pmodes[j])).value();
791 }
792 }
793
794 if (Pstream::parRun())
795 {
796 reduce(dMatrix, sumOp<Eigen::MatrixXd>());
797 }
798
799 ITHACAstream::SaveDenseMatrix(dMatrix, "./ITHACAoutput/Matrices/",
800 "D_" + name(nPModes));
801 return dMatrix;
802}
803
804Eigen::Tensor<double, 3> UnsteadyNSTurbIntrusive::divMomentum(label nUModes,
805 label nPModes)
806{
807 label g1Size = nPModes;
808 label g2Size = nUModes;
809 Eigen::Tensor<double, 3> gTensor;
810 gTensor.resize(g1Size, g2Size, g2Size);
811
812 for (label i = 0; i < g1Size; i++)
813 {
814 for (label j = 0; j < g2Size; j++)
815 {
816 for (label k = 0; k < g2Size; k++)
817 {
818 gTensor(i, j, k) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & (fvc::div(
819 fvc::interpolate(Umodes[j]) & Umodes[j].mesh().Sf(),
820 Umodes[k]))).value();
821 }
822 }
823 }
824
825 if (Pstream::parRun())
826 {
827 reduce(gTensor, sumOp<Eigen::Tensor<double, 3 >> ());
828 }
829
830 // Export the tensor
831 ITHACAstream::SaveDenseTensor(gTensor, "./ITHACAoutput/Matrices/",
832 "G_" + name(nUModes) + "_" + name(
833 nPModes) + "_t");
834 return gTensor;
835}
836
837Eigen::MatrixXd UnsteadyNSTurbIntrusive::pressureBC1(label nUModes,
838 label nPModes)
839{
840 label P_BC1size = nPModes;
841 label P_BC2size = nUModes;
842 Eigen::MatrixXd BC1_matrix(P_BC1size, P_BC2size);
843 fvMesh& mesh = _mesh();
844 for (label i = 0; i < P_BC1size; i++)
845 {
846 for (label j = 0; j < P_BC2size; j++)
847 {
848 surfaceScalarField lpl((fvc::interpolate(fvc::laplacian(
849 Umodes[j])) & mesh.Sf()) * fvc::interpolate(Pmodes[i]));
850 double s = 0;
851 for (label k = 0; k < lpl.boundaryField().size(); k++)
852 {
853 s += gSum(lpl.boundaryField()[k]);
854 }
855
856 BC1_matrix(i, j) = s;
857 }
858 }
859
860 if (Pstream::parRun())
861 {
862 reduce(BC1_matrix, sumOp<Eigen::MatrixXd>());
863 }
864
865 ITHACAstream::SaveDenseMatrix(BC1_matrix, "./ITHACAoutput/Matrices/",
866 "BC1_" + name(nUModes) + "_" + name(nPModes));
867 return BC1_matrix;
868}
869
870Eigen::Tensor<double, 3 > UnsteadyNSTurbIntrusive::pressureBC2(label nUModes,
871 label nPModes)
872{
873 label pressureBC1Size = nPModes;
874 label pressureBC2Size = nUModes;
875 Eigen::Tensor<double, 3 > bc2Tensor;
876 fvMesh& mesh = _mesh();
877 bc2Tensor.resize(pressureBC1Size, pressureBC2Size, pressureBC2Size);
878
879 for (label i = 0; i < pressureBC1Size; i++)
880 {
881 for (label j = 0; j < pressureBC2Size; j++)
882 {
883 for (label k = 0; k < pressureBC2Size; k++)
884 {
885 surfaceScalarField div_m(fvc::interpolate(fvc::div(fvc::interpolate(
886 Umodes[j]) & mesh.Sf(),
887 Umodes[k])) & mesh.Sf() * fvc::interpolate(Pmodes[i]));
888 double s = 0;
889
890 for (label k = 0; k < div_m.boundaryField().size(); k++)
891 {
892 s += gSum(div_m.boundaryField()[k]);
893 }
894
895 bc2Tensor(i, j, k) = s;
896 }
897 }
898 }
899
900 if (Pstream::parRun())
901 {
902 reduce(bc2Tensor, sumOp<Eigen::Tensor<double, 3 >> ());
903 }
904
905 // Export the tensor
906 ITHACAstream::SaveDenseTensor(bc2Tensor, "./ITHACAoutput/Matrices/",
907 "BC2_" + name(nUModes) + "_" + name(nPModes) + "_t");
908 return bc2Tensor;
909}
910
911Eigen::MatrixXd UnsteadyNSTurbIntrusive::pressureBC3(label nUModes,
912 label nPModes)
913{
914 label P3_BC1size = nPModes;
915 label P3_BC2size = nUModes;
916 Eigen::MatrixXd BC3_matrix(P3_BC1size, P3_BC2size);
917 fvMesh& mesh = _mesh();
918 surfaceVectorField n(mesh.Sf() / mesh.magSf());
919 for (label i = 0; i < P3_BC1size; i++)
920 {
921 for (label j = 0; j < P3_BC2size; j++)
922 {
923 surfaceVectorField BC3 = fvc::interpolate(fvc::curl(Umodes[j])).ref();
924 surfaceVectorField BC4 = (n^ fvc::interpolate(fvc::grad(Pmodes[i]))).ref();
925 surfaceScalarField BC5 = ((BC3& BC4) * mesh.magSf()).ref();
926 double s = 0;
927 for (label k = 0; k < BC5.boundaryField().size(); k++)
928 {
929 s += gSum(BC5.boundaryField()[k]);
930 }
931
932 BC3_matrix(i, j) = s;
933 }
934 }
935
936 if (Pstream::parRun())
937 {
938 reduce(BC3_matrix, sumOp<Eigen::MatrixXd>());
939 }
940
941 ITHACAstream::SaveDenseMatrix(BC3_matrix, "./ITHACAoutput/Matrices/",
942 "BC3_" + name(nUModes) + "_" + name(nPModes));
943 return BC3_matrix;
944}
945
946Eigen::MatrixXd UnsteadyNSTurbIntrusive::pressureBC4(label nUModes,
947 label nPModes)
948{
949 label P4_BC1size = nPModes;
950 label P4_BC2size = nUModes;
951 Eigen::MatrixXd BC4_matrix(P4_BC1size, P4_BC2size);
952 fvMesh& mesh = _mesh();
953 surfaceVectorField n(mesh.Sf() / mesh.magSf());
954
955 for (label i = 0; i < P4_BC1size; i++)
956 {
957 for (label j = 0; j < P4_BC2size; j++)
958 {
959 surfaceScalarField BC3 = fvc::interpolate(Pmodes[i]).ref();
960 surfaceScalarField BC4 = (n & fvc::interpolate(Umodes[j])).ref();
961 surfaceScalarField BC5 = ((BC3 * BC4) * mesh.magSf()).ref();
962 double s = 0;
963
964 for (label k = 0; k < BC5.boundaryField().size(); k++)
965 {
966 s += gSum(BC5.boundaryField()[k]);
967 }
968
969 BC4_matrix(i, j) = s;
970 }
971 }
972
973 if (Pstream::parRun())
974 {
975 reduce(BC4_matrix, sumOp<Eigen::MatrixXd>());
976 }
977
978 ITHACAstream::SaveDenseMatrix(BC4_matrix, "./ITHACAoutput/Matrices/",
979 "BC4_" + name(nUModes) + "_" + name(nPModes));
980 return BC4_matrix;
981}
982
983List< Eigen::MatrixXd > UnsteadyNSTurbIntrusive::bcVelocityVec(label nModes)
984{
985 List < Eigen::MatrixXd > bcVelVec(inletIndex.rows());
986
987 for (label j = 0; j < inletIndex.rows(); j++)
988 {
989 bcVelVec[j].resize(nModes, 1);
990 }
991
992 for (label k = 0; k < inletIndex.rows(); k++)
993 {
994 label BCind = inletIndex(k, 0);
995 label BCcomp = inletIndex(k, 1);
996
997 for (label i = 0; i < nModes; i++)
998 {
999 bcVelVec[k](i, 0) = gSum(Umodes[i].boundaryField()[BCind].component(
1000 BCcomp));
1001 }
1002 }
1003
1004 ITHACAstream::exportMatrix(bcVelVec, "bcVelVec", "eigen",
1005 "./ITHACAoutput/Matrices/bcVelVec");
1006 return bcVelVec;
1007}
1008
1009List< Eigen::MatrixXd > UnsteadyNSTurbIntrusive::bcVelocityMat(label nModes)
1010{
1011 label BCUsize = inletIndex.rows();
1012 List < Eigen::MatrixXd > bcVelMat(BCUsize);
1013
1014 for (label j = 0; j < inletIndex.rows(); j++)
1015 {
1016 bcVelMat[j].resize(nModes, nModes);
1017 }
1018
1019 for (label k = 0; k < inletIndex.rows(); k++)
1020 {
1021 label BCind = inletIndex(k, 0);
1022 label BCcomp = inletIndex(k, 1);
1023
1024 for (label i = 0; i < nModes; i++)
1025 {
1026 for (label j = 0; j < nModes; j++)
1027 {
1028 bcVelMat[k](i, j) = gSum(Umodes[i].boundaryField()[BCind].component(
1029 BCcomp) *
1030 Umodes[j].boundaryField()[BCind].component(BCcomp));
1031 }
1032 }
1033 }
1034
1035 ITHACAstream::exportMatrix(bcVelMat, "bcVelMat", "eigen",
1036 "./ITHACAoutput/Matrices/bcVelMat");
1037 return bcVelMat;
1038}
1039
1041 label nUModes, label nPModes)
1042{
1043 Eigen::Tensor<double, 3> ct1PPETensor;
1044 ct1PPETensor.resize(nPModes, nUModes, nUModes);
1045
1046 for (label i = 0; i < nPModes; i++)
1047 {
1048 for (label j = 0; j < nUModes; j++)
1049 {
1050 for (label k = 0; k < nUModes; k++)
1051 {
1052 ct1PPETensor(i, j, k) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & (
1053 fvc::laplacian(
1054 nutModes[j], Umodes[k]))).value();
1055 }
1056 }
1057 }
1058
1059 if (Pstream::parRun())
1060 {
1061 reduce(ct1PPETensor, sumOp<Eigen::Tensor<double, 3 >> ());
1062 }
1063
1064 // Export the tensor
1065 ITHACAstream::SaveDenseTensor(ct1PPETensor, "./ITHACAoutput/Matrices/",
1066 "ct1PPE_" + name(nUModes) + "_" + name(nPModes) + "_t");
1067 return ct1PPETensor;
1068}
1069
1071 label nUModes, label nPModes)
1072{
1073 Eigen::Tensor<double, 3> ct2PPETensor;
1074 ct2PPETensor.resize(nPModes, nUModes, nUModes);
1075
1076 for (label i = 0; i < nPModes; i++)
1077 {
1078 for (label j = 0; j < nUModes; j++)
1079 {
1080 for (label k = 0; k < nUModes; k++)
1081 {
1082 ct2PPETensor(i, j, k) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & ((fvc::div(
1083 nutModes[j] * dev2((fvc::grad(Umodes[k]))().T()))))).value();
1084 }
1085 }
1086 }
1087
1088 if (Pstream::parRun())
1089 {
1090 reduce(ct2PPETensor, sumOp<Eigen::Tensor<double, 3 >> ());
1091 }
1092
1093 // Export the tensor
1094 ITHACAstream::SaveDenseTensor(ct2PPETensor, "./ITHACAoutput/Matrices/",
1095 "ct2PPE_" + name(nUModes) + "_" + name(nPModes) + "_t");
1096 return ct2PPETensor;
1097}
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
#define M_Assert(Expr, Msg)
fv::options & fvOptions
Definition NLsolve.H:25
auto nut
Definition NLsolve.H:195
Header file of the UnsteadyNSTurbIntrusive class.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
Eigen::Tensor< double, 3 > convectiveTerm(label nModes)
The method for computing the convective term tensor.
Eigen::MatrixXd pressureBC4(label nUModes, label nPModes)
Term N° 4 given by the additional boundary condition using a PPE approach.
Eigen::MatrixXd kMatrix
Pressure Gradient matrix.
Eigen::Tensor< double, 3 > cTotalTensor
Total Turbulent tensor.
Eigen::Tensor< double, 3 > turbulencePPETensor1(label nUModes, label nPModes)
Method to compute one of the turbulence eddy viscosity tensors in the PPE approach.
List< Eigen::MatrixXd > bcVelocityMat(label nModes)
Boundary integral modes on boundary used by the penaly method.
Eigen::MatrixXd pressureGradientTerm(label nModes)
The method for computing the pressure gradient term with number of modes of pressure being equal to t...
Eigen::Tensor< double, 3 > ct1PPETensor
Turbulent viscosity tensor PPE approach.
Eigen::Tensor< double, 3 > turbulenceTensor1(label nModes)
Method to compute one of the turbulence tensors.
Eigen::MatrixXd btTurbulence(label nModes)
bt added matrix for the turbulence treatement
Eigen::Tensor< double, 3 > divMomentum(label nUModes, label nPModes)
Divergence of convective term (PPE approach)
Eigen::Tensor< double, 3 > turbulenceTensor2(label nModes)
Method to compute one of the turbulence tensors.
Eigen::Tensor< double, 3 > ct2PPETensor
Turbulent viscosity tensor PPE approach.
Eigen::Tensor< double, 3 > pressureBC2(label nUModes, label nPModes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Eigen::Tensor< double, 3 > turbulencePPETensor2(label nUModes, label nPModes)
Method to compute one of the turbulence eddy viscosity tensors in the PPE approach.
Eigen::MatrixXd pressureBC3(label nUModes, label nPModes)
Term N° 3 given by the additional boundary condition using a PPE approach.
Eigen::Tensor< double, 3 > convTensor
Convective tensor.
void projectPPE(fileName folder, label nUModes, label nPModes)
The projection method for computing the reduced order matrices with the usage of Poisson equation for...
Eigen::Tensor< double, 3 > ct2Tensor
Turbulent viscosity tensor.
Eigen::Tensor< double, 3 > cTotalPPETensor
Total Turbulent tensor PPE approach.
Eigen::MatrixXd btMatrix
Turbulent viscosity matrix.
Eigen::MatrixXd bTotalMatrix
Total B Matrix.
Eigen::MatrixXd laplacianPressure(label nPModes)
Laplacian of pressure term (PPE approach)
Eigen::Tensor< double, 3 > ct1Tensor
Turbulent viscosity tensor.
UnsteadyNSTurbIntrusive()
Construct Null.
autoPtr< volScalarField > _nut
Eddy viscosity field.
Eigen::MatrixXd bMatrix
Diffusive matrix.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
label nModesOnline
Number of modes used in the online stage.
volScalarModes nutModes
List of POD modes for eddy viscosity.
Eigen::MatrixXd diffusiveTerm(label nModes)
Diffusive Term.
Eigen::MatrixXd pressureBC1(label nUModes, label nPModes)
Term N° 1 given by the additional boundary condition using a PPE approach.
List< Eigen::MatrixXd > bcVelocityVec(label nModes)
Boundary integral modes on boundary used by the penaly method.
bool checkWrite(Time &timeObject)
Function to check if the solution must be exported.
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)
void project()
General projection operation.
void writeMu(List< scalar > mu_now)
Write out a list of scalar corresponding to the parameters used in the truthSolve.
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
label counter
Counter used for the output of the full order solutions.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
IOdictionary * ITHACAdict
dictionary to store input output infos
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.
List< Eigen::MatrixXd > bcVelVec
Boundary term for penalty method - vector.
Definition steadyNS.H:227
scalar maxIter
Number of maximum iterations to be done for the computation of the truth solution.
Definition steadyNS.H:128
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
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:302
Eigen::MatrixXd BC4_matrix
PPE BC4.
Definition steadyNS.H:203
Eigen::MatrixXd BC3_matrix
PPE BC3.
Definition steadyNS.H:200
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
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
scalar tolerance
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and p...
Definition steadyNS.H:125
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:290
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
volVectorModes L_U_SUPmodes
List of pointers containing the total number of lift, supremizer and velocity modes.
Definition steadyNS.H:119
ITHACAparameters * para
Definition steadyNS.H:82
Eigen::MatrixXd D_matrix
Laplacian term PPE.
Definition steadyNS.H:182
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:314
Eigen::Tensor< double, 3 > gTensor
Divergence of momentum PPE.
Definition steadyNS.H:188
List< Eigen::MatrixXd > bcVelMat
Boundary term for penalty method - matrix.
Definition steadyNS.H:230
Eigen::Tensor< double, 3 > bc2Tensor
PPE BC2.
Definition steadyNS.H:197
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:272
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
word bcMethod
Boundary Method.
Definition steadyNS.H:326
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:269
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition unsteadyNS.H:76
word timeDerivativeSchemeOrder
Definition unsteadyNS.H:99
autoPtr< pimpleControl > _pimple
pimpleControl
Definition unsteadyNS.H:73
volScalarField & T
Definition createT.H:46
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 exportTensor(Eigen::Tensor< T, 3 > tensor, word Name, word type, word folder)
Export the reduced tensor in numpy (tipo=python), matlab (tipo=matlab) and txt (tipo=eigen) format.
void SaveDenseTensor(TensorType &Tensor, word folder, word MatrixName)
Save a dense tensor to 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.
surfaceScalarField & phi
volVectorField & U
volScalarField & p
pimpleControl & pimple
singlePhaseTransportModel & laminarTransport
label i
Definition pEqn.H:46