Loading...
Searching...
No Matches
SteadyNSTurb.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
13License
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 "steadyNS.H"
32#include "SteadyNSTurb.H"
33#include "viscosityModel.H"
34
35// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
36
37// Constructor
39SteadyNSTurb::SteadyNSTurb(int argc, char* argv[])
40{
41 _args = autoPtr<argList>
42 (
43 new argList(argc, argv)
44 );
45
46 if (!_args->checkRootCase())
47 {
48 Foam::FatalError.exit();
49 }
50
51 argList& args = _args();
52#include "createTime.H"
53#include "createMesh.H"
54 _simple = autoPtr<simpleControl>
55 (
56 new simpleControl
57 (
58 mesh
59 )
60 );
61 simpleControl& simple = _simple();
62#include "createFields.H"
63#include "createFvOptions.H"
64 //
65#pragma GCC diagnostic push
66#pragma GCC diagnostic ignored "-Wunused-variable"
67#include "initContinuityErrs.H"
68#pragma GCC diagnostic pop
69 //
70 turbulence->validate();
71 ITHACAdict = new IOdictionary
72 (
73 IOobject
74 (
75 "ITHACAdict",
76 runTime.system(),
77 mesh,
78 IOobject::MUST_READ,
79 IOobject::NO_WRITE
80 )
81 );
82 tolerance = ITHACAdict->lookupOrDefault<scalar>("tolerance", 1e-5);
83 maxIter = ITHACAdict->lookupOrDefault<scalar>("maxIter", 1000);
84 bcMethod = ITHACAdict->lookupOrDefault<word>("bcMethod", "lift");
85 M_Assert(bcMethod == "lift" || bcMethod == "penalty",
86 "The BC method must be set to lift or penalty in ITHACAdict");
87 viscCoeff = ITHACAdict->lookupOrDefault<word>("viscCoeff", "RBF");
92}
93
94
95// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
96
97// Method to performa a truthSolve
98void SteadyNSTurb::truthSolve(List<scalar> mu_now)
99{
100 Time& runTime = _runTime();
101 fvMesh& mesh = _mesh();
102 volScalarField& p = _p();
103 volVectorField& U = _U();
104 surfaceScalarField& phi = _phi();
105 fv::options& fvOptions = _fvOptions();
106 simpleControl& simple = _simple();
107 IOMRFZoneList& MRF = _MRF();
108 singlePhaseTransportModel& laminarTransport = _laminarTransport();
109#include "NLsolveSteadyNSTurb.H"
110 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
111 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
112 volScalarField _nut(turbulence->nut());
113 ITHACAstream::exportSolution(_nut, name(counter), "./ITHACAoutput/Offline/");
114 Ufield.append(U.clone());
115 Pfield.append(p.clone());
116 nutFields.append(_nut.clone());
117 counter++;
118 writeMu(mu_now);
119 // --- Fill in the mu_samples with parameters (mu) to be used for the PODI sample points
120 mu_samples.conservativeResize(mu_samples.rows() + 1, mu_now.size());
121
122 for (label i = 0; i < mu_now.size(); i++)
123 {
124 mu_samples(mu_samples.rows() - 1, i) = mu_now[i];
125 }
126
127 // Resize to Unitary if not initialized by user (i.e. non-parametric problem)
128 if (mu.cols() == 0)
129 {
130 mu.resize(1, 1);
131 }
132
133 if (mu_samples.rows() == mu.cols())
134 {
135 ITHACAstream::exportMatrix(mu_samples, "mu_samples", "eigen",
136 "./ITHACAoutput/Offline");
137 }
138}
139
140List < Eigen::MatrixXd > SteadyNSTurb::turbulenceTerm1(label NUmodes,
141 label NSUPmodes, label nNutModes)
142{
143 label cSize = NUmodes + NSUPmodes + liftfield.size();
144 List < Eigen::MatrixXd > ct1Matrix;
145 ct1Matrix.setSize(cSize);
146
147 for (label j = 0; j < cSize; j++)
148 {
149 ct1Matrix[j].resize(nNutModes, cSize);
150 ct1Matrix[j] = ct1Matrix[j] * 0;
151 }
152
153 for (label i = 0; i < cSize; i++)
154 {
155 Info << "Filling layer number " << i + 1 << " in the matrix ct1Matrix" << endl;
156
157 for (label j = 0; j < nNutModes; j++)
158 {
159 for (label k = 0; k < cSize; k++)
160 {
161 ct1Matrix[i](j, k) = fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::laplacian(
162 nutModes[j], L_U_SUPmodes[k])).value();
163 }
164 }
165 }
166
167 // Export the matrix
168 ITHACAstream::exportMatrix(ct1Matrix, "ct1Matrix", "python",
169 "./ITHACAoutput/Matrices/");
170 ITHACAstream::exportMatrix(ct1Matrix, "ct1Matrix", "matlab",
171 "./ITHACAoutput/Matrices/");
172 ITHACAstream::exportMatrix(ct1Matrix, "ct1Matrix", "eigen",
173 "./ITHACAoutput/Matrices/ct1");
174 return ct1Matrix;
175}
176
177Eigen::Tensor<double, 3> SteadyNSTurb::turbulenceTensor1(label NUmodes,
178 label NSUPmodes, label nNutModes)
179{
180 label cSize = NUmodes + NSUPmodes + liftfield.size();
181 Eigen::Tensor<double, 3> ct1Tensor;
182 ct1Tensor.resize(cSize, nNutModes, cSize);
183
184 for (label i = 0; i < cSize; i++)
185 {
186 for (label j = 0; j < nNutModes; j++)
187 {
188 for (label k = 0; k < cSize; k++)
189 {
190 ct1Tensor(i, j, k) = fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::laplacian(
191 nutModes[j], L_U_SUPmodes[k])).value();
192 }
193 }
194 }
195
196 // Export the tensor
197 ITHACAstream::SaveDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/",
198 "ct1_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
199 NSUPmodes) + "_" + name(nNutModes) + "_t");
200 return ct1Tensor;
201}
202
203
204
205
206
207List < Eigen::MatrixXd > SteadyNSTurb::turbulenceTerm2(label NUmodes,
208 label NSUPmodes, label nNutModes)
209{
210 label cSize = NUmodes + NSUPmodes + liftfield.size();
211 List < Eigen::MatrixXd > ct2Matrix;
212 ct2Matrix.setSize(cSize);
213
214 for (label j = 0; j < cSize; j++)
215 {
216 ct2Matrix[j].resize(nNutModes, cSize);
217 ct2Matrix[j] = ct2Matrix[j] * 0;
218 }
219
220 for (label i = 0; i < cSize; i++)
221 {
222 Info << "Filling layer number " << i + 1 << " in the matrix ct2Matrix" << endl;
223
224 for (label j = 0; j < nNutModes; j++)
225 {
226 for (label k = 0; k < cSize; k++)
227 {
228 ct2Matrix[i](j, k) = fvc::domainIntegrate(L_U_SUPmodes[i] & (fvc::div(
229 nutModes[j] * dev((fvc::grad(L_U_SUPmodes[k]))().T())))).value();
230 }
231 }
232 }
233
234 // Export the matrix
235 ITHACAstream::exportMatrix(ct2Matrix, "ct2Matrix", "python",
236 "./ITHACAoutput/Matrices/");
237 ITHACAstream::exportMatrix(ct2Matrix, "ct2Matrix", "matlab",
238 "./ITHACAoutput/Matrices/");
239 ITHACAstream::exportMatrix(ct2Matrix, "ct2Matrix", "eigen",
240 "./ITHACAoutput/Matrices/ct2");
241 return ct2Matrix;
242}
243
244Eigen::Tensor<double, 3> SteadyNSTurb::turbulenceTensor2(label NUmodes,
245 label NSUPmodes, label nNutModes)
246{
247 label cSize = NUmodes + NSUPmodes + liftfield.size();
248 Eigen::Tensor<double, 3> ct2Tensor;
249 ct2Tensor.resize(cSize, nNutModes, cSize);
250
251 for (label i = 0; i < cSize; i++)
252 {
253 for (label j = 0; j < nNutModes; j++)
254 {
255 for (label k = 0; k < cSize; k++)
256 {
257 ct2Tensor(i, j, k) = fvc::domainIntegrate(L_U_SUPmodes[i] & (fvc::div(
258 nutModes[j] * dev((fvc::grad(L_U_SUPmodes[k]))().T())))).value();
259 }
260 }
261 }
262
263 // Export the tensor
264 ITHACAstream::SaveDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/",
265 "ct2_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
266 NSUPmodes) + "_" + name(nNutModes) + "_t");
267 return ct2Tensor;
268}
269
270Eigen::Tensor<double, 3> SteadyNSTurb::turbulencePPETensor1(label NUmodes,
271 label NSUPmodes, label NPmodes, label nNutModes)
272{
273 label cSize = NUmodes + NSUPmodes + liftfield.size();
274 Eigen::Tensor<double, 3> ct1PPETensor;
275 ct1PPETensor.resize(NPmodes, nNutModes, cSize);
276
277 for (label i = 0; i < NPmodes; i++)
278 {
279 for (label j = 0; j < nNutModes; j++)
280 {
281 for (label k = 0; k < cSize; k++)
282 {
283 // ct1PPETensor(i, j, k) = fvc::domainIntegrate(2 * Pmodes[i] * (fvc::laplacian(
284 // L_U_SUPmodes[k]) & fvc::grad(nutModes[j]))).value();
285 // ct1PPETensor(i, j, k) = fvc::domainIntegrate(Pmodes[i] * (fvc::div(
286 // fvc::laplacian(
287 // nutModes[j], L_U_SUPmodes[k])))).value();
288 ct1PPETensor(i, j, k) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & (
289 fvc::laplacian(
290 nutModes[j], L_U_SUPmodes[k]))).value();
291 }
292 }
293 }
294
295 // Export the tensor
296 ITHACAstream::SaveDenseTensor(ct1PPETensor, "./ITHACAoutput/Matrices/",
297 "ct1PPE_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
298 NSUPmodes) + "_" + name(NPmodes) + "_" + name(nNutModes) + "_t");
299 return ct1PPETensor;
300}
301
302Eigen::Tensor<double, 3> SteadyNSTurb::turbulencePPETensor2(label NUmodes,
303 label NSUPmodes, label NPmodes, label nNutModes)
304{
305 label cSize = NUmodes + NSUPmodes + liftfield.size();
306 Eigen::Tensor<double, 3> ct2PPETensor;
307 ct2PPETensor.resize(NPmodes, nNutModes, cSize);
308
309 for (label i = 0; i < NPmodes; i++)
310 {
311 for (label j = 0; j < nNutModes; j++)
312 {
313 for (label k = 0; k < cSize; k++)
314 {
315 // ct2PPETensor(i, j, k) = fvc::domainIntegrate(Pmodes[i] * (fvc::grad(fvc::grad(
316 // nutModes[j])) && (dev2((fvc::grad(L_U_SUPmodes[k]))() + fvc::grad(
317 // L_U_SUPmodes[k]))().T()))).value();
318 // ct2PPETensor(i, j, k) = fvc::domainIntegrate(Pmodes[i] * ((fvc::div(fvc::div(
319 // nutModes[j] * dev2((fvc::grad(L_U_SUPmodes[k]))().T())))))).value();
320 ct2PPETensor(i, j, k) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & ((fvc::div(
321 nutModes[j] * dev2((fvc::grad(L_U_SUPmodes[k]))().T()))))).value();
322 }
323 }
324 }
325
326 // Export the tensor
327 ITHACAstream::SaveDenseTensor(ct2PPETensor, "./ITHACAoutput/Matrices/",
328 "ct2PPE_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
329 NSUPmodes) + "_" + name(NPmodes) + "_" + name(nNutModes) + "_t");
330 return ct2PPETensor;
331}
332
333Eigen::MatrixXd SteadyNSTurb::btTurbulence(label NUmodes, label NSUPmodes)
334{
335 label btSize = NUmodes + NSUPmodes + liftfield.size();
336 Eigen::MatrixXd btMatrix(btSize, btSize);
337 btMatrix = btMatrix * 0;
338
339 // Project everything
340 for (label i = 0; i < btSize; i++)
341 {
342 for (label j = 0; j < btSize; j++)
343 {
344 btMatrix(i, j) = fvc::domainIntegrate(L_U_SUPmodes[i] & (fvc::div(dev((T(
345 fvc::grad(
346 L_U_SUPmodes[j]))))))).value();
347 }
348 }
349
350 // Export the matrix
351 ITHACAstream::SaveDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/",
352 "bt_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(NSUPmodes));
353 return btMatrix;
354}
355
356void SteadyNSTurb::projectPPE(fileName folder, label NU, label NP, label NSUP,
357 label Nnut)
358{
359 NUmodes = NU;
360 NPmodes = NP;
361 NSUPmodes = 0;
362 nNutModes = Nnut;
363 L_U_SUPmodes.resize(0);
364
365 if (liftfield.size() != 0)
366 {
367 for (label k = 0; k < liftfield.size(); k++)
368 {
369 L_U_SUPmodes.append(liftfield[k].clone());
370 }
371 }
372
373 if (NUmodes != 0)
374 {
375 for (label k = 0; k < NUmodes; k++)
376 {
377 L_U_SUPmodes.append(Umodes[k].clone());
378 }
379 }
380
381 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
382 {
383 word bStr = "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
384 NSUPmodes);
385
386 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bStr))
387 {
388 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", bStr);
389 }
390 else
391 {
393 }
394
395 word btStr = "bt_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
396 NSUPmodes);
397
398 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + btStr))
399 {
400 ITHACAstream::ReadDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/", btStr);
401 }
402 else
403 {
405 }
406
407 word kStr = "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
408 NSUPmodes) + "_" + name(NPmodes);
409
410 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + kStr))
411 {
412 ITHACAstream::ReadDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/", kStr);
413 }
414 else
415 {
417 }
418
419 word pStr = "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
420 NSUPmodes) + "_" + name(NPmodes);
421
422 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + pStr))
423 {
424 ITHACAstream::ReadDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/", pStr);
425 }
426 else
427 {
429 }
430
431 word mStr = "M_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
432 NSUPmodes);
433
434 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + mStr))
435 {
436 ITHACAstream::ReadDenseMatrix(M_matrix, "./ITHACAoutput/Matrices/", mStr);
437 }
438 else
439 {
441 }
442
443 word D_str = "D_" + name(NPmodes);
444
445 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + D_str))
446 {
447 ITHACAstream::ReadDenseMatrix(D_matrix, "./ITHACAoutput/Matrices/", D_str);
448 }
449 else
450 {
452 }
453
454 word bc3_str = "BC3_" + name(liftfield.size()) + "_" + name(
455 NUmodes) + "_" + name(NPmodes);
456
457 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bc3_str))
458 {
459 ITHACAstream::ReadDenseMatrix(BC3_matrix, "./ITHACAoutput/Matrices/", bc3_str);
460 }
461 else
462 {
464 }
465
466 word bc4_str = "BC4_" + name(liftfield.size()) + "_" + name(
467 NUmodes) + "_" + name(NPmodes);
468
469 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bc4_str))
470 {
471 ITHACAstream::ReadDenseMatrix(BC4_matrix, "./ITHACAoutput/Matrices/", bc4_str);
472 }
473 else
474 {
476 }
477
478 word C_str = "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
479 NSUPmodes) + "_t";
480
481 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + C_str))
482 {
483 ITHACAstream::ReadDenseTensor(C_tensor, "./ITHACAoutput/Matrices/", C_str);
484 }
485 else
486 {
488 }
489
490 word ct1Str = "ct1_" + name(liftfield.size()) + "_" + name(
491 NUmodes) + "_" + name(
492 NSUPmodes) + "_" + name(nNutModes) + "_t";
493
494 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct1Str))
495 {
496 ITHACAstream::ReadDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/", ct1Str);
497 }
498 else
499 {
501 }
502
503 word ct2Str = "ct2_" + name(liftfield.size()) + "_" + name(
504 NUmodes) + "_" + name(
505 NSUPmodes) + "_" + name(nNutModes) + "_t";
506
507 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct2Str))
508 {
509 ITHACAstream::ReadDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/", ct2Str);
510 }
511 else
512 {
514 }
515
516 word G_str = "G_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
517 NSUPmodes) + "_" + name(NPmodes) + "_t";
518
519 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + G_str))
520 {
521 ITHACAstream::ReadDenseTensor(gTensor, "./ITHACAoutput/Matrices/", G_str);
522 }
523 else
524 {
526 }
527
528 if (bcMethod == "penalty")
529 {
532 }
533 }
534 else
535 {
536 L_U_SUPmodes.resize(0);
537
538 if (liftfield.size() != 0)
539 {
540 for (label k = 0; k < liftfield.size(); k++)
541 {
542 L_U_SUPmodes.append(liftfield[k].clone());
543 }
544 }
545
546 if (NUmodes != 0)
547 {
548 for (label k = 0; k < NUmodes; k++)
549 {
550 L_U_SUPmodes.append(Umodes[k].clone());
551 }
552 }
553
554 if (NSUPmodes != 0)
555 {
556 for (label k = 0; k < NSUPmodes; k++)
557 {
558 L_U_SUPmodes.append(supmodes[k].clone());
559 }
560 }
561
570 // BC4_matrix = pressure_BC4(NUmodes, NPmodes);
575
576 if (bcMethod == "penalty")
577 {
580 }
581 }
582
583 // Export the matrices
584 if (para->exportPython)
585 {
586 ITHACAstream::exportMatrix(B_matrix, "B", "python", "./ITHACAoutput/Matrices/");
587 ITHACAstream::exportMatrix(K_matrix, "K", "python", "./ITHACAoutput/Matrices/");
588 ITHACAstream::exportMatrix(D_matrix, "D", "python", "./ITHACAoutput/Matrices/");
589 ITHACAstream::exportMatrix(M_matrix, "M", "python", "./ITHACAoutput/Matrices/");
590 ITHACAstream::exportMatrix(BC3_matrix, "BC3", "python",
591 "./ITHACAoutput/Matrices/");
592 ITHACAstream::exportMatrix(BC4_matrix, "BC4", "python",
593 "./ITHACAoutput/Matrices/");
594 ITHACAstream::exportMatrix(btMatrix, "bt", "python",
595 "./ITHACAoutput/Matrices/");
596 ITHACAstream::exportTensor(C_tensor, "C", "python", "./ITHACAoutput/Matrices/");
597 ITHACAstream::exportTensor(gTensor, "G", "python", "./ITHACAoutput/Matrices/");
598 ITHACAstream::exportTensor(ct1Tensor, "ct1", "python",
599 "./ITHACAoutput/Matrices/");
600 ITHACAstream::exportTensor(ct2Tensor, "ct2", "python",
601 "./ITHACAoutput/Matrices/");
602 ITHACAstream::exportTensor(ct1PPETensor, "ct1PPE", "python",
603 "./ITHACAoutput/Matrices/");
604 ITHACAstream::exportTensor(ct2PPETensor, "ct2PPE", "python",
605 "./ITHACAoutput/Matrices/");
606 }
607
608 if (para->exportMatlab)
609 {
610 ITHACAstream::exportMatrix(B_matrix, "B", "matlab", "./ITHACAoutput/Matrices/");
611 ITHACAstream::exportMatrix(K_matrix, "K", "matlab", "./ITHACAoutput/Matrices/");
612 ITHACAstream::exportMatrix(D_matrix, "D", "matlab", "./ITHACAoutput/Matrices/");
613 ITHACAstream::exportMatrix(M_matrix, "M", "matlab", "./ITHACAoutput/Matrices/");
614 ITHACAstream::exportMatrix(BC3_matrix, "BC3", "matlab",
615 "./ITHACAoutput/Matrices/");
616 ITHACAstream::exportMatrix(BC4_matrix, "BC4", "matlab",
617 "./ITHACAoutput/Matrices/");
618 ITHACAstream::exportMatrix(btMatrix, "bt", "matlab",
619 "./ITHACAoutput/Matrices/");
620 ITHACAstream::exportTensor(gTensor, "G", "matlab", "./ITHACAoutput/Matrices/");
621 ITHACAstream::exportTensor(ct1Tensor, "ct1", "matlab",
622 "./ITHACAoutput/Matrices/");
623 ITHACAstream::exportTensor(ct2Tensor, "ct2", "matlab",
624 "./ITHACAoutput/Matrices/");
625 ITHACAstream::exportTensor(ct1PPETensor, "ct1PPE", "matlab",
626 "./ITHACAoutput/Matrices/");
627 ITHACAstream::exportTensor(ct2PPETensor, "ct2PPE", "matlab",
628 "./ITHACAoutput/Matrices/");
629 }
630
631 if (para->exportTxt)
632 {
633 ITHACAstream::exportMatrix(B_matrix, "B", "eigen", "./ITHACAoutput/Matrices/");
634 ITHACAstream::exportMatrix(K_matrix, "K", "eigen", "./ITHACAoutput/Matrices/");
635 ITHACAstream::exportMatrix(D_matrix, "D", "eigen", "./ITHACAoutput/Matrices/");
636 ITHACAstream::exportMatrix(M_matrix, "M", "eigen", "./ITHACAoutput/Matrices/");
638 "./ITHACAoutput/Matrices/");
640 "./ITHACAoutput/Matrices/");
641 ITHACAstream::exportMatrix(btMatrix, "bt", "eigen", "./ITHACAoutput/Matrices/");
642 ITHACAstream::exportTensor(gTensor, "G", "eigen", "./ITHACAoutput/Matrices/G");
643 ITHACAstream::exportTensor(C_tensor, "C", "eigen", "./ITHACAoutput/Matrices/C");
644 ITHACAstream::exportTensor(ct1Tensor, "ct1_", "eigen",
645 "./ITHACAoutput/Matrices/ct1");
646 ITHACAstream::exportTensor(ct2Tensor, "ct2_", "eigen",
647 "./ITHACAoutput/Matrices/ct2");
648 ITHACAstream::exportTensor(ct1PPETensor, "ct1PPE_", "eigen",
649 "./ITHACAoutput/Matrices/ct1PPE");
650 ITHACAstream::exportTensor(ct2PPETensor, "ct2PPE_", "eigen",
651 "./ITHACAoutput/Matrices/ct2PPE");
652 }
653
655 label cSize = NUmodes + NSUPmodes + liftfield.size();
656 cTotalTensor.resize(cSize, nNutModes, cSize);
658 cTotalPPETensor.resize(NPmodes, nNutModes, cSize);
660 // Get the coeffs for interpolation (the orthonormal one is used because basis are orthogonal)
663 ITHACAstream::exportMatrix(coeffL2, "coeffL2", "python",
664 "./ITHACAoutput/Matrices/");
665 ITHACAstream::exportMatrix(coeffL2, "coeffL2", "matlab",
666 "./ITHACAoutput/Matrices/");
667 ITHACAstream::exportMatrix(coeffL2, "coeffL2", "eigen",
668 "./ITHACAoutput/Matrices/");
669 // Export the matrix
670 ITHACAstream::SaveDenseMatrix(coeffL2, "./ITHACAoutput/Matrices/",
671 "coeffL2_nut_" + name(nNutModes));
672 samples.resize(nNutModes);
673 rbfSplines.resize(nNutModes);
674 Eigen::MatrixXd weights;
675
676 for (label i = 0; i < nNutModes; i++)
677 {
678 word weightName = "wRBF_N" + name(i + 1) + "_" + name(liftfield.size()) + "_"
679 + name(NUmodes) + "_" + name(NSUPmodes) ;
680
681 if (ITHACAutilities::check_file("./ITHACAoutput/weightsPPE/" + weightName))
682 {
683 samples[i] = new SPLINTER::DataTable(1, 1);
684
685 for (label j = 0; j < coeffL2.cols(); j++)
686 {
687 samples[i]->addSample(mu.row(j), coeffL2(i, j));
688 }
689
690 ITHACAstream::ReadDenseMatrix(weights, "./ITHACAoutput/weightsPPE/",
691 weightName);
692 rbfSplines[i] = new SPLINTER::RBFSpline(*samples[i],
693 SPLINTER::RadialBasisFunctionType::GAUSSIAN, weights);
694 std::cout << "Constructing RadialBasisFunction for mode " << i + 1 << std::endl;
695 }
696 else
697 {
698 samples[i] = new SPLINTER::DataTable(1, 1);
699
700 for (label j = 0; j < coeffL2.cols(); j++)
701 {
702 samples[i]->addSample(mu.row(j), coeffL2(i, j));
703 }
704
705 rbfSplines[i] = new SPLINTER::RBFSpline(*samples[i],
706 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
708 "./ITHACAoutput/weightsPPE/", weightName);
709 std::cout << "Constructing RadialBasisFunction for mode " << i + 1 << std::endl;
710 }
711 }
712}
713
714void SteadyNSTurb::projectSUP(fileName folder, label NU, label NP, label NSUP,
715 label Nnut)
716{
717 NUmodes = NU;
718 NPmodes = NP;
719 NSUPmodes = NSUP;
720 nNutModes = Nnut;
721 L_U_SUPmodes.resize(0);
722
723 if (liftfield.size() != 0)
724 {
725 for (label k = 0; k < liftfield.size(); k++)
726 {
727 L_U_SUPmodes.append(liftfield[k].clone());
728 }
729 }
730
731 if (NUmodes != 0)
732 {
733 for (label k = 0; k < NUmodes; k++)
734 {
735 L_U_SUPmodes.append(Umodes[k].clone());
736 }
737 }
738
739 if (NSUPmodes != 0)
740 {
741 for (label k = 0; k < NSUPmodes; k++)
742 {
743 L_U_SUPmodes.append(supmodes[k].clone());
744 }
745 }
746
747 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
748 {
749 word bStr = "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
750 NSUPmodes);
751
752 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bStr))
753 {
754 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", bStr);
755 }
756 else
757 {
759 }
760
761 word btStr = "bt_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
762 NSUPmodes);
763
764 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + btStr))
765 {
766 ITHACAstream::ReadDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/", btStr);
767 }
768 else
769 {
771 }
772
773 word kStr = "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
774 NSUPmodes) + "_" + name(NPmodes);
775
776 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + kStr))
777 {
778 ITHACAstream::ReadDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/", kStr);
779 }
780 else
781 {
783 }
784
785 word pStr = "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
786 NSUPmodes) + "_" + name(NPmodes);
787
788 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + pStr))
789 {
790 ITHACAstream::ReadDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/", pStr);
791 }
792 else
793 {
795 }
796
797 word mStr = "M_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
798 NSUPmodes);
799
800 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + mStr))
801 {
802 ITHACAstream::ReadDenseMatrix(M_matrix, "./ITHACAoutput/Matrices/", mStr);
803 }
804 else
805 {
807 }
808
809 word C_str = "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
810 NSUPmodes) + "_t";
811
812 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + C_str))
813 {
814 ITHACAstream::ReadDenseTensor(C_tensor, "./ITHACAoutput/Matrices/", C_str);
815 }
816 else
817 {
819 }
820
821 word ct1Str = "ct1_" + name(liftfield.size()) + "_" + name(
822 NUmodes) + "_" + name(
823 NSUPmodes) + "_" + name(nNutModes) + "_t";
824
825 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct1Str))
826 {
827 ITHACAstream::ReadDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/", ct1Str);
828 }
829 else
830 {
832 }
833
834 word ct2Str = "ct2_" + name(liftfield.size()) + "_" + name(
835 NUmodes) + "_" + name(
836 NSUPmodes) + "_" + name(nNutModes) + "_t";
837
838 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct2Str))
839 {
840 ITHACAstream::ReadDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/", ct2Str);
841 }
842 else
843 {
845 }
846
847 if (bcMethod == "penalty")
848 {
851 }
852 }
853 else
854 {
855 L_U_SUPmodes.resize(0);
856
857 if (liftfield.size() != 0)
858 {
859 for (label k = 0; k < liftfield.size(); k++)
860 {
861 L_U_SUPmodes.append(liftfield[k].clone());
862 }
863 }
864
865 if (NUmodes != 0)
866 {
867 for (label k = 0; k < NUmodes; k++)
868 {
869 L_U_SUPmodes.append(Umodes[k].clone());
870 }
871 }
872
873 if (NSUPmodes != 0)
874 {
875 for (label k = 0; k < NSUPmodes; k++)
876 {
877 L_U_SUPmodes.append(supmodes[k].clone());
878 }
879 }
880
889
890 if (bcMethod == "penalty")
891 {
894 }
895 }
896
897 // Export the matrices
898 if (para->exportPython)
899 {
900 ITHACAstream::exportMatrix(B_matrix, "B", "python", "./ITHACAoutput/Matrices/");
901 ITHACAstream::exportMatrix(K_matrix, "K", "python", "./ITHACAoutput/Matrices/");
902 ITHACAstream::exportMatrix(P_matrix, "P", "python", "./ITHACAoutput/Matrices/");
903 ITHACAstream::exportMatrix(M_matrix, "M", "python", "./ITHACAoutput/Matrices/");
904 ITHACAstream::exportMatrix(btMatrix, "bt", "python",
905 "./ITHACAoutput/Matrices/");
906 ITHACAstream::exportTensor(C_tensor, "C", "python", "./ITHACAoutput/Matrices/");
907 ITHACAstream::exportTensor(ct1Tensor, "ct1", "python",
908 "./ITHACAoutput/Matrices/");
909 ITHACAstream::exportTensor(ct2Tensor, "ct2", "python",
910 "./ITHACAoutput/Matrices/");
911 }
912
913 if (para->exportMatlab)
914 {
915 ITHACAstream::exportMatrix(B_matrix, "B", "matlab", "./ITHACAoutput/Matrices/");
916 ITHACAstream::exportMatrix(K_matrix, "K", "matlab", "./ITHACAoutput/Matrices/");
917 ITHACAstream::exportMatrix(P_matrix, "P", "matlab", "./ITHACAoutput/Matrices/");
918 ITHACAstream::exportMatrix(M_matrix, "M", "matlab", "./ITHACAoutput/Matrices/");
919 ITHACAstream::exportMatrix(btMatrix, "bt", "matlab",
920 "./ITHACAoutput/Matrices/");
921 ITHACAstream::exportTensor(C_tensor, "C", "matlab", "./ITHACAoutput/Matrices/");
922 ITHACAstream::exportTensor(ct1Tensor, "ct1", "matlab",
923 "./ITHACAoutput/Matrices/");
924 ITHACAstream::exportTensor(ct2Tensor, "ct2", "matlab",
925 "./ITHACAoutput/Matrices/");
926 }
927
928 if (para->exportTxt)
929 {
930 ITHACAstream::exportMatrix(B_matrix, "B", "eigen", "./ITHACAoutput/Matrices/");
931 ITHACAstream::exportMatrix(K_matrix, "K", "eigen", "./ITHACAoutput/Matrices/");
932 ITHACAstream::exportMatrix(P_matrix, "P", "eigen", "./ITHACAoutput/Matrices/");
933 ITHACAstream::exportMatrix(M_matrix, "M", "eigen", "./ITHACAoutput/Matrices/");
934 ITHACAstream::exportMatrix(btMatrix, "bt", "eigen", "./ITHACAoutput/Matrices/");
935 ITHACAstream::exportTensor(C_tensor, "C", "eigen", "./ITHACAoutput/Matrices/C");
936 ITHACAstream::exportTensor(ct1Tensor, "ct1_", "eigen",
937 "./ITHACAoutput/Matrices/ct1");
938 ITHACAstream::exportTensor(ct2Tensor, "ct2_", "eigen",
939 "./ITHACAoutput/Matrices/ct2");
940 }
941
943 label cSize = NUmodes + NSUPmodes + liftfield.size();
944 cTotalTensor.resize(cSize, nNutModes, cSize);
946 // Get the coeffs for interpolation (the orthonormal one is used because basis are orthogonal)
949 ITHACAstream::exportMatrix(coeffL2, "coeffL2", "python",
950 "./ITHACAoutput/Matrices/");
951 ITHACAstream::exportMatrix(coeffL2, "coeffL2", "matlab",
952 "./ITHACAoutput/Matrices/");
953 ITHACAstream::exportMatrix(coeffL2, "coeffL2", "eigen",
954 "./ITHACAoutput/Matrices/");
955 // Export the matrix
956 ITHACAstream::SaveDenseMatrix(coeffL2, "./ITHACAoutput/Matrices/",
957 "coeffL2_nut_" + name(nNutModes));
958 samples.resize(nNutModes);
959 rbfSplines.resize(nNutModes);
960 Eigen::MatrixXd weights;
961
962 for (label i = 0; i < nNutModes; i++)
963 {
964 word weightName = "wRBF_N" + name(i + 1) + "_" + name(liftfield.size()) + "_"
965 + name(NUmodes) + "_" + name(NSUPmodes) ;
966
967 if (ITHACAutilities::check_file("./ITHACAoutput/weightsSUP/" + weightName))
968 {
969 samples[i] = new SPLINTER::DataTable(1, 1);
970
971 for (label j = 0; j < coeffL2.cols(); j++)
972 {
973 samples[i]->addSample(mu.row(j), coeffL2(i, j));
974 }
975
976 ITHACAstream::ReadDenseMatrix(weights, "./ITHACAoutput/weightsSUP/",
977 weightName);
978 rbfSplines[i] = new SPLINTER::RBFSpline(*samples[i],
979 SPLINTER::RadialBasisFunctionType::GAUSSIAN, weights);
980 std::cout << "Constructing RadialBasisFunction for mode " << i + 1 << std::endl;
981 }
982 else
983 {
984 samples[i] = new SPLINTER::DataTable(1, 1);
985
986 for (label j = 0; j < coeffL2.cols(); j++)
987 {
988 samples[i]->addSample(mu.row(j), coeffL2(i, j));
989 }
990
991 rbfSplines[i] = new SPLINTER::RBFSpline(*samples[i],
992 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
994 "./ITHACAoutput/weightsSUP/", weightName);
996 "wRBF_" + name(i), "eigen",
997 "./ITHACAoutput/weightsSUP/"
998 );
999 std::cout << "Constructing RadialBasisFunction for mode " << i + 1 << std::endl;
1000 }
1001 }
1002}
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
#define M_Assert(Expr, Msg)
Header file of the SteadyNSTurb class.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
Eigen::MatrixXd coeffL2
The matrix of L2 projection coefficients for the eddy viscosity.
Eigen::MatrixXd btMatrix
Turbulent viscosity matrix.
Eigen::Tensor< double, 3 > cTotalPPETensor
Turbulent total viscosity tensor in the PPE equation.
label nNutModes
Number of viscoisty modes used for the projection.
Eigen::MatrixXd bTotalMatrix
Total B Matrix.
std::vector< SPLINTER::RBFSpline * > rbfSplines
Create a RBF splines for interpolation.
Eigen::MatrixXd btTurbulence(label NUmodes, label NSUPmodes)
bt added matrix for the turbulence treatement
autoPtr< volScalarField > _nut
Eddy viscosity field.
std::vector< SPLINTER::DataTable * > samples
Create a samples for interpolation.
List< Eigen::MatrixXd > turbulenceTerm1(label NUmodes, label NSUPmodes, label nNutModes)
ct1 added matrix for the turbulence treatement
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
List< Eigen::MatrixXd > ct2Matrix
Turbulent viscosity tensor.
Eigen::Tensor< double, 3 > ct1PPETensor
Turbulent viscosity tensor in the PPE equation.
Eigen::Tensor< double, 3 > ct1Tensor
List< Eigen::MatrixXd > ct1Matrix
Turbulent viscosity tensor.
List< Eigen::MatrixXd > turbulenceTerm2(label NUmodes, label NSUPmodes, label nNutModes)
Method to compute one of the turbulence eddy viscosity tensors.
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes)
Project using a supremizer approach.
Eigen::Tensor< double, 3 > ct2PPETensor
Turbulent viscosity tensor in the PPE equation.
Eigen::Tensor< double, 3 > turbulencePPETensor1(label NUmodes, label NSUPmodes, label NPmodes, label nNutModes)
ct1PPE added tensor for the turbulence treatement in the PPE method
Eigen::Tensor< double, 3 > cTotalTensor
Eigen::Tensor< double, 3 > turbulencePPETensor2(label NUmodes, label NSUPmodes, label NPmodes, label nNutModes)
ct2PPE added tensor for the turbulence treatement in the PPE method
Eigen::Tensor< double, 3 > turbulenceTensor1(label NUmodes, label NSUPmodes, label nNutModes)
ct1 tensor for the turbulence treatement
Eigen::Tensor< double, 3 > ct2Tensor
Eigen::Tensor< double, 3 > turbulenceTensor2(label NUmodes, label NSUPmodes, label nNutModes)
Method to compute one of the turbulence eddy viscosity tensors.
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes)
Project using a PPE approach.
word viscCoeff
The way to compute the eddy viscosity snapshots.
volScalarModes nutModes
List of POD modes for eddy viscosity.
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
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:218
scalar maxIter
Number of maximum iterations to be done for the computation of the truth solution.
Definition steadyNS.H:125
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::MatrixXd BC4_matrix
PPE BC4.
Definition steadyNS.H:194
Eigen::MatrixXd BC3_matrix
PPE BC3.
Definition steadyNS.H:191
Eigen::Tensor< double, 3 > C_tensor
Diffusion term.
Definition steadyNS.H:167
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:284
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
List< Eigen::MatrixXd > bcVelocityVec(label NUmodes, label NSUPmodes)
Boundary integral modes on boundary used by the penaly method.
Definition steadyNS.C:1418
autoPtr< fv::options > _fvOptions
fvOptions
Definition steadyNS.H:287
Eigen::MatrixXd pressure_BC4(label NPmodes, label NUmodes)
Term N° 4 given by the additional boundary condition using a PPE approach for time-dependent BCs.
Definition steadyNS.C:1378
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 divergence_term(label NUmodes, label NPmodes, label NSUPmodes)
Divergence Term (supremizer approach)
Definition steadyNS.C:1073
scalar tolerance
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and p...
Definition steadyNS.H:122
Eigen::MatrixXd mass_term(label NUmodes, label NPmodes, label NSUPmodes)
Mass Term.
Definition steadyNS.C:1041
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:281
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
ITHACAparameters * para
Definition steadyNS.H:82
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::Tensor< double, 3 > gTensor
Divergence of momentum PPE.
Definition steadyNS.H:179
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
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition steadyNS.H:299
Eigen::Tensor< double, 3 > divMomentum(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach)
Definition steadyNS.C:1147
Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPmodes, label NSUPmodes)
Gradient of pressure.
Definition steadyNS.C:923
List< Eigen::MatrixXd > bcVelMat
Boundary term for penalty method - matrix.
Definition steadyNS.H:221
Eigen::MatrixXd P_matrix
Div of velocity.
Definition steadyNS.H:170
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition steadyNS.H:160
Eigen::MatrixXd laplacian_pressure(label NPmodes)
Laplacian of pressure term (PPE approach)
Definition steadyNS.C:1183
List< Eigen::MatrixXd > bcVelocityMat(label NUmodes, label NSUPmodes)
Boundary integral modes on boundary used by the penaly method.
Definition steadyNS.C:1455
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
word bcMethod
Boundary Method.
Definition steadyNS.H:317
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
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.
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_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
volScalarField & p
singlePhaseTransportModel & laminarTransport
label i
Definition pEqn.H:46
Header file of the steadyNS class.