Loading...
Searching...
No Matches
steadyNS.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13 License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
31
34
35#include "steadyNS.H"
36#include "viscosityModel.H"
37
38// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
39// Constructor
41
42steadyNS::steadyNS(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 _simple = autoPtr<simpleControl>
58 (
59 new simpleControl
60 (
61 mesh
62 )
63 );
64 simpleControl& simple = _simple();
65#include "createFields.H"
66#include "createFvOptions.H"
67 turbulence->validate();
68 ITHACAdict = new IOdictionary
69 (
70 IOobject
71 (
72 "ITHACAdict",
73 runTime.system(),
74 mesh,
75 IOobject::MUST_READ,
76 IOobject::NO_WRITE
77 )
78 );
79 tolerance = ITHACAdict->lookupOrDefault<scalar>("tolerance", 1e-5);
80 maxIter = ITHACAdict->lookupOrDefault<scalar>("maxIter", 1000);
81 bcMethod = ITHACAdict->lookupOrDefault<word>("bcMethod", "lift");
82 M_Assert(bcMethod == "lift" || bcMethod == "penalty" || bcMethod == "none",
83 "The BC method must be set to lift or penalty or none in ITHACAdict");
84 fluxMethod = ITHACAdict->lookupOrDefault<word>("fluxMethod", "inconsistent");
85 M_Assert(fluxMethod == "inconsistent" || bcMethod == "consistent",
86 "The flux method must be set to inconsistent or consistent in ITHACAdict");
87 para = ITHACAparameters::getInstance(mesh, runTime);
91}
92
93
94// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
95
96// Method to perform a truthSolve
97void steadyNS::truthSolve(List<scalar> mu_now)
98{
99 Time& runTime = _runTime();
100 fvMesh& mesh = _mesh();
101 volScalarField& p = _p();
102 volVectorField& U = _U();
103 surfaceScalarField& phi = _phi();
104 fv::options& fvOptions = _fvOptions();
105 simpleControl& simple = _simple();
106 IOMRFZoneList& MRF = _MRF();
107 singlePhaseTransportModel& laminarTransport = _laminarTransport();
108#include "NLsolvesteadyNS.H"
109 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
110 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
111 Ufield.append(U.clone());
112 Pfield.append(p.clone());
113 counter++;
114 writeMu(mu_now);
115 // --- Fill in the mu_samples with parameters (mu) to be used for the PODI sample points
116 mu_samples.conservativeResize(mu_samples.rows() + 1, mu_now.size());
117
118 for (label i = 0; i < mu_now.size(); i++)
119 {
120 mu_samples(mu_samples.rows() - 1, i) = mu_now[i];
121 }
122
123 // Resize to Unitary if not initialized by user (i.e. non-parametric problem)
124 if (mu.cols() == 0)
125 {
126 mu.resize(1, 1);
127 }
128
129 if (mu_samples.rows() == mu.cols())
130 {
131 ITHACAstream::exportMatrix(mu_samples, "mu_samples", "eigen",
132 "./ITHACAoutput/Offline");
133 }
134}
135
136// Method to solve the supremizer problem
138{
139 M_Assert(type == "modes"
140 || type == "snapshots",
141 "You must specify the variable type with either snapshots or modes");
142 PtrList<volScalarField> P_sup;
143
144 if (type == "snapshots")
145 {
146 P_sup = Pfield;
147 }
148 else
149 {
150 P_sup = Pmodes;
151 }
152
153 if (supex == 1)
154 {
155 volVectorField U = _U();
156 volVectorField Usup
157 (
158 IOobject
159 (
160 "Usup",
161 U.time().timeName(),
162 U.mesh(),
163 IOobject::NO_READ,
164 IOobject::AUTO_WRITE
165 ),
166 U.mesh(),
167 dimensionedVector("zero", U.dimensions(), vector::zero)
168 );
169
170 if (type == "snapshots")
171 {
172 ITHACAstream::read_fields(supfield, Usup, "./ITHACAoutput/supfield/");
173 }
174 else
175 {
176 ITHACAstream::read_fields(supmodes, Usup, "./ITHACAoutput/supremizer/");
177 }
178 }
179 else
180 {
181 volVectorField U = _U();
182 volVectorField Usup
183 (
184 IOobject
185 (
186 "Usup",
187 U.time().timeName(),
188 U.mesh(),
189 IOobject::NO_READ,
190 IOobject::AUTO_WRITE
191 ),
192 U.mesh(),
193 dimensionedVector("zero", U.dimensions(), vector::zero)
194 );
195 dimensionedScalar nu_fake
196 (
197 "nu_fake",
198 dimensionSet(0, 2, -1, 0, 0, 0, 0),
199 scalar(1)
200 );
201 Vector<double> v(0, 0, 0);
202
203 for (label i = 0; i < Usup.boundaryField().size(); i++)
204 {
205 if (Usup.boundaryField()[i].type() != "processor")
206 {
207 ITHACAutilities::changeBCtype(Usup, "fixedValue", i);
208 assignBC(Usup, i, v);
209 assignIF(Usup, v);
210 }
211 }
212
213 if (type == "snapshots")
214 {
215 for (label i = 0; i < P_sup.size(); i++)
216 {
217 fvVectorMatrix u_sup_eqn
218 (
219 - fvm::laplacian(nu_fake, Usup)
220 );
221 solve
222 (
223 u_sup_eqn == fvc::grad(P_sup[i])
224 );
225 supfield.append(Usup.clone());
226 ITHACAstream::exportSolution(Usup, name(i + 1), "./ITHACAoutput/supfield/");
227 }
228
229 ITHACAutilities::createSymLink("./ITHACAoutput/supfield");
230 }
231 else
232 {
233 for (label i = 0; i < P_sup.size(); i++)
234 {
235 fvVectorMatrix u_sup_eqn
236 (
237 - fvm::laplacian(nu_fake, Usup)
238 );
239 solve
240 (
241 u_sup_eqn == fvc::grad(P_sup[i])
242 );
243 supmodes.append(Usup.clone());
244 ITHACAstream::exportSolution(Usup, name(i + 1), "./ITHACAoutput/supremizer/");
245 }
246
247 ITHACAutilities::createSymLink("./ITHACAoutput/supremizer");
248 }
249 }
250}
251
252// Method to compute the lifting function
254{
255 for (label k = 0; k < inletIndex.rows(); k++)
256 {
257 Time& runTime = _runTime();
258 surfaceScalarField& phi = _phi();
259 fvMesh& mesh = _mesh();
260 volScalarField p = _p();
261 volVectorField U = _U();
262 IOMRFZoneList& MRF = _MRF();
263 label BCind = inletIndex(k, 0);
264 volVectorField Ulift("Ulift" + name(k), U);
265 instantList Times = runTime.times();
266 runTime.setTime(Times[1], 1);
267 pisoControl potentialFlow(mesh, "potentialFlow");
268 Info << "Solving a lifting Problem" << endl;
269 Vector<double> v1(0, 0, 0);
270 v1[inletIndex(k, 1)] = 1;
271 Vector<double> v0(0, 0, 0);
272
273 for (label j = 0; j < U.boundaryField().size(); j++)
274 {
275 if (j == BCind)
276 {
277 assignBC(Ulift, j, v1);
278 }
279 else if (U.boundaryField()[BCind].type() == "fixedValue")
280 {
281 assignBC(Ulift, j, v0);
282 }
283 else
284 {
285 }
286
287 assignIF(Ulift, v0);
288 phi = linearInterpolate(Ulift) & mesh.Sf();
289 }
290
291 Info << "Constructing velocity potential field Phi\n" << endl;
292 volScalarField Phi
293 (
294 IOobject
295 (
296 "Phi",
297 runTime.timeName(),
298 mesh,
299 IOobject::READ_IF_PRESENT,
300 IOobject::NO_WRITE
301 ),
302 mesh,
303 dimensionedScalar("Phi", dimLength * dimVelocity, 0),
304 p.boundaryField().types()
305 );
306 label PhiRefCell = 0;
307 scalar PhiRefValue = 0;
308 setRefCell
309 (
310 Phi,
311 potentialFlow.dict(),
312 PhiRefCell,
313 PhiRefValue
314 );
315 mesh.setFluxRequired(Phi.name());
316 runTime.functionObjects().start();
317 MRF.makeRelative(phi);
318 adjustPhi(phi, Ulift, p);
319
320 while (potentialFlow.correctNonOrthogonal())
321 {
322 fvScalarMatrix PhiEqn
323 (
324 fvm::laplacian(dimensionedScalar("1", dimless, 1), Phi)
325 ==
326 fvc::div(phi)
327 );
328 PhiEqn.setReference(PhiRefCell, PhiRefValue);
329 PhiEqn.solve();
330
331 if (potentialFlow.finalNonOrthogonalIter())
332 {
333 phi -= PhiEqn.flux();
334 }
335 }
336
337 MRF.makeAbsolute(phi);
338 Info << "Continuity error = "
339 << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
340 << endl;
341 Ulift = fvc::reconstruct(phi);
342 Ulift.correctBoundaryConditions();
343 Info << "Interpolated velocity error = "
344 << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
345 / sum(mesh.magSf())).value()
346 << endl;
347 Ulift.write();
348 liftfield.append(Ulift.clone());
349 }
350}
351
352// * * * * * * * * * * * * * * Projection Methods * * * * * * * * * * * * * * //
353
354void steadyNS::projectPPE(fileName folder, label NU, label NP, label NSUP)
355{
356 NUmodes = NU;
357 NPmodes = NP;
358 NSUPmodes = 0;
359 L_U_SUPmodes.resize(0);
360
361 if (liftfield.size() != 0)
362 {
363 for (label k = 0; k < liftfield.size(); k++)
364 {
365 L_U_SUPmodes.append(liftfield[k].clone());
366 }
367 }
368
369 if (NUmodes != 0)
370 {
371 for (label k = 0; k < NUmodes; k++)
372 {
373 L_U_SUPmodes.append(Umodes[k].clone());
374 }
375 }
376
377 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
378 {
379 word B_str = "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
380 NSUPmodes);
381
382 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + B_str))
383 {
384 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", B_str);
385 }
386 else
387 {
389 }
390
391 word K_str = "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
392 NSUPmodes) + "_" + name(NPmodes);
393
394 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + K_str))
395 {
396 ITHACAstream::ReadDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/", K_str);
397 }
398 else
399 {
401 }
402
403 word M_str = "M_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
404 NSUPmodes);
405
406 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + M_str))
407 {
408 ITHACAstream::ReadDenseMatrix(M_matrix, "./ITHACAoutput/Matrices/", M_str);
409 }
410 else
411 {
413 }
414
415 word D_str = "D_" + name(NPmodes);
416
417 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + D_str))
418 {
419 ITHACAstream::ReadDenseMatrix(D_matrix, "./ITHACAoutput/Matrices/", D_str);
420 }
421 else
422 {
424 }
425
426 word bc3_str = "BC3_" + name(liftfield.size()) + "_" + name(
427 NUmodes) + "_" + name(NPmodes);
428
429 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bc3_str))
430 {
431 ITHACAstream::ReadDenseMatrix(BC3_matrix, "./ITHACAoutput/Matrices/", bc3_str);
432 }
433 else
434 {
436 }
437
438 word bc4_str = "BC4_" + name(liftfield.size()) + "_" + name(
439 NUmodes) + "_" + name(NPmodes);
440
441 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bc4_str))
442 {
443 ITHACAstream::ReadDenseMatrix(BC4_matrix, "./ITHACAoutput/Matrices/", bc4_str);
444 }
445 else
446 {
448 }
449
450 word C_str = "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
451 NSUPmodes) + "_t";
452
453 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + C_str))
454 {
455 ITHACAstream::ReadDenseTensor(C_tensor, "./ITHACAoutput/Matrices/", C_str);
456 }
457 else
458 {
460 }
461
462 word G_str = "G_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
463 NSUPmodes) + "_" + name(NPmodes) + "_t";
464
465 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + G_str))
466 {
467 ITHACAstream::ReadDenseTensor(gTensor, "./ITHACAoutput/Matrices/", G_str);
468 }
469 else
470 {
472 }
473
474 if (bcMethod == "penalty")
475 {
478 }
479 }
480 else
481 {
490
491 if (bcMethod == "penalty")
492 {
495 }
496 }
497
498 // Export the matrices
499 if (para->exportPython)
500 {
501 ITHACAstream::exportMatrix(B_matrix, "B", "python", "./ITHACAoutput/Matrices/");
502 ITHACAstream::exportMatrix(K_matrix, "K", "python", "./ITHACAoutput/Matrices/");
503 ITHACAstream::exportMatrix(D_matrix, "D", "python", "./ITHACAoutput/Matrices/");
504 ITHACAstream::exportMatrix(M_matrix, "M", "python", "./ITHACAoutput/Matrices/");
505 ITHACAstream::exportMatrix(BC3_matrix, "BC3", "python",
506 "./ITHACAoutput/Matrices/");
507 ITHACAstream::exportMatrix(BC4_matrix, "BC4", "python",
508 "./ITHACAoutput/Matrices/");
509 ITHACAstream::exportTensor(C_tensor, "C", "python", "./ITHACAoutput/Matrices/");
510 ITHACAstream::exportTensor(gTensor, "G", "python", "./ITHACAoutput/Matrices/");
511 }
512
513 if (para->exportMatlab)
514 {
515 ITHACAstream::exportMatrix(B_matrix, "B", "matlab", "./ITHACAoutput/Matrices/");
516 ITHACAstream::exportMatrix(K_matrix, "K", "matlab", "./ITHACAoutput/Matrices/");
517 ITHACAstream::exportMatrix(D_matrix, "D", "matlab", "./ITHACAoutput/Matrices/");
518 ITHACAstream::exportMatrix(M_matrix, "M", "matlab", "./ITHACAoutput/Matrices/");
519 ITHACAstream::exportMatrix(BC3_matrix, "BC3", "matlab",
520 "./ITHACAoutput/Matrices/");
521 ITHACAstream::exportMatrix(BC4_matrix, "BC4", "matlab",
522 "./ITHACAoutput/Matrices/");
523 ITHACAstream::exportTensor(C_tensor, "C", "matlab", "./ITHACAoutput/Matrices/");
524 ITHACAstream::exportTensor(gTensor, "G", "matlab", "./ITHACAoutput/Matrices/");
525 }
526
527 if (para->exportTxt)
528 {
529 ITHACAstream::exportMatrix(B_matrix, "B", "eigen", "./ITHACAoutput/Matrices/");
530 ITHACAstream::exportMatrix(K_matrix, "K", "eigen", "./ITHACAoutput/Matrices/");
531 ITHACAstream::exportMatrix(D_matrix, "D", "eigen", "./ITHACAoutput/Matrices/");
532 ITHACAstream::exportMatrix(M_matrix, "M", "eigen", "./ITHACAoutput/Matrices/");
534 "./ITHACAoutput/Matrices/");
536 "./ITHACAoutput/Matrices/");
538 "./ITHACAoutput/Matrices/C");
540 "./ITHACAoutput/Matrices/G");
541 }
542
543 if (para->exportNpy)
544 {
545 cnpy::save(B_matrix, "./ITHACAoutput/Matrices/B.npy");
546 cnpy::save(K_matrix, "./ITHACAoutput/Matrices/K.npy");
547 cnpy::save(D_matrix, "./ITHACAoutput/Matrices/D.npy");
548 cnpy::save(M_matrix, "./ITHACAoutput/Matrices/M.npy");
549 cnpy::save(BC3_matrix, "./ITHACAoutput/Matrices/BC3.npy");
550 cnpy::save(BC4_matrix, "./ITHACAoutput/Matrices/BC4.npy");
551 cnpy::save(C_tensor, "./ITHACAoutput/Matrices/C.npy");
552 cnpy::save(gTensor, "./ITHACAoutput/Matrices/G.npy");
553 }
554}
555
556void steadyNS::projectSUP(fileName folder, label NU, label NP, label NSUP)
557{
558 NUmodes = NU;
559 NPmodes = NP;
560 NSUPmodes = NSUP;
561 L_U_SUPmodes.resize(0);
562
563 if (liftfield.size() != 0)
564 {
565 for (label k = 0; k < liftfield.size(); k++)
566 {
567 L_U_SUPmodes.append(liftfield[k].clone());
568 }
569 }
570
571 if (NUmodes != 0)
572 {
573 for (label k = 0; k < NUmodes; k++)
574 {
575 L_U_SUPmodes.append(Umodes[k].clone());
576 }
577 }
578
579 if (NSUPmodes != 0)
580 {
581 for (label k = 0; k < NSUPmodes; k++)
582 {
583 L_U_SUPmodes.append(supmodes[k].clone());
584 }
585 }
586
587 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
588 {
589 word B_str = "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
590 NSUPmodes);
591
592 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + B_str))
593 {
594 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", B_str);
595 }
596 else
597 {
599 }
600
601 word K_str = "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
602 NSUPmodes) + "_" + name(NPmodes);
603
604 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + K_str))
605 {
606 ITHACAstream::ReadDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/", K_str);
607 }
608 else
609 {
611 }
612
613 word P_str = "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
614 NSUPmodes) + "_" + name(NPmodes);
615
616 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + P_str))
617 {
618 ITHACAstream::ReadDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/", P_str);
619 }
620 else
621 {
623 }
624
625 word M_str = "M_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
626 NSUPmodes);
627
628 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + M_str))
629 {
630 ITHACAstream::ReadDenseMatrix(M_matrix, "./ITHACAoutput/Matrices/", M_str);
631 }
632 else
633 {
635 }
636
637 word C_str = "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
638 NSUPmodes) + "_t";
639
640 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + C_str))
641 {
642 ITHACAstream::ReadDenseTensor(C_tensor, "./ITHACAoutput/Matrices/", C_str);
643 }
644 else
645 {
647 }
648
649 if (bcMethod == "penalty")
650 {
653 }
654 }
655 else
656 {
662
663 if (bcMethod == "penalty")
664 {
667 }
668 }
669
670 // Export the matrices
671 if (para->exportPython)
672 {
673 ITHACAstream::exportMatrix(B_matrix, "B", "python", "./ITHACAoutput/Matrices/");
674 ITHACAstream::exportMatrix(K_matrix, "K", "python", "./ITHACAoutput/Matrices/");
675 ITHACAstream::exportMatrix(P_matrix, "P", "python", "./ITHACAoutput/Matrices/");
676 ITHACAstream::exportMatrix(M_matrix, "M", "python", "./ITHACAoutput/Matrices/");
677 ITHACAstream::exportTensor(C_tensor, "C", "python", "./ITHACAoutput/Matrices/");
678 }
679
680 if (para->exportMatlab)
681 {
682 ITHACAstream::exportMatrix(B_matrix, "B", "matlab", "./ITHACAoutput/Matrices/");
683 ITHACAstream::exportMatrix(K_matrix, "K", "matlab", "./ITHACAoutput/Matrices/");
684 ITHACAstream::exportMatrix(P_matrix, "P", "matlab", "./ITHACAoutput/Matrices/");
685 ITHACAstream::exportMatrix(M_matrix, "M", "matlab", "./ITHACAoutput/Matrices/");
686 ITHACAstream::exportTensor(C_tensor, "C", "python", "./ITHACAoutput/Matrices/");
687 }
688
689 if (para->exportTxt)
690 {
691 ITHACAstream::exportMatrix(B_matrix, "B", "eigen", "./ITHACAoutput/Matrices/");
692 ITHACAstream::exportMatrix(K_matrix, "K", "eigen", "./ITHACAoutput/Matrices/");
693 ITHACAstream::exportMatrix(P_matrix, "P", "eigen", "./ITHACAoutput/Matrices/");
694 ITHACAstream::exportMatrix(M_matrix, "M", "eigen", "./ITHACAoutput/Matrices/");
696 "./ITHACAoutput/Matrices/C");
697 }
698
699 if (para->exportNpy)
700 {
701 cnpy::save(B_matrix, "./ITHACAoutput/Matrices/B.npy");
702 cnpy::save(K_matrix, "./ITHACAoutput/Matrices/K.npy");
703 cnpy::save(P_matrix, "./ITHACAoutput/Matrices/P.npy");
704 cnpy::save(M_matrix, "./ITHACAoutput/Matrices/M.npy");
705 cnpy::save(C_tensor, "./ITHACAoutput/Matrices/C.npy");
706 }
707}
708
709void steadyNS::discretizeThenProject(fileName folder, label NU, label NP,
710 label NSUP)
711{
712 NUmodes = NU;
713 NPmodes = NP;
714 NSUPmodes = 0;
715 L_U_SUPmodes.resize(0);
716 Vector<double> inl(0, 0, 0);
717
718 if (NUmodes != 0)
719 {
720 for (label k = 0; k < NUmodes; k++)
721 {
722 L_U_SUPmodes.append(Umodes[k].clone());
723 // set homogenous boundary conditions
724 forAll(L_U_SUPmodes[k].mesh().boundary(), l)
725 {
726 assignBC(L_U_SUPmodes[k], l, inl);
727 }
728 }
729 }
730
731 if (fluxMethod == "consistent")
732 {
733 L_PHImodes.resize(0);
734
735 for (label k = 0; k < NUmodes; k++)
736 {
737 L_PHImodes.append(Phimodes[k].clone());
738 forAll(L_PHImodes[k].boundaryFieldRef()[0], l)
739 {
740 L_PHImodes[k].boundaryFieldRef()[0][l] = 0;
741 }
742 }
743 }
744
745 // Dummy field of which all boundary conditions converted to homogenous
746 // Dirichlet boundary conditions.
747 Uinl = autoPtr<volVectorField> (new volVectorField(_U()));
749 // dummy variables
750 dt_dummy = autoPtr<dimensionedScalar> (new dimensionedScalar
751 (
752 "dt_dummy",
753 dimensionSet(0, 0, 1, 0, 0, 0, 0),
754 scalar(1.0)
755 ));
756 nu_dummy = autoPtr<dimensionedScalar> (new dimensionedScalar
757 (
758 "nu_dummy",
759 dimensionSet(0, 2, -1, 0, 0, 0, 0),
760 scalar(1.0)
761 ));
762
763 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
764 {
765 word B_str = "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
766 NSUPmodes);
767
768 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + B_str))
769 {
770 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", B_str);
771 }
772 else
773 {
775 }
776
777 word K_str = "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
778 NSUPmodes) + "_" + name(NPmodes);
779
780 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + K_str))
781 {
782 ITHACAstream::ReadDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/", K_str);
783 }
784 else
785 {
787 }
788
789 word P_str = "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
790 NSUPmodes) + "_" + name(NPmodes);
791
792 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + P_str))
793 {
794 ITHACAstream::ReadDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/", P_str);
795 }
796 else
797 {
799 }
800
801 word C_str = "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
802 NSUPmodes) + "_t";
803
804 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + C_str))
805 {
806 ITHACAstream::ReadDenseTensor(C_tensor, "./ITHACAoutput/Matrices/", C_str);
807 }
808 else
809 {
811 }
812
813 word Cf_str = "Cf_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
814 NSUPmodes) + "_t";
815
816 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + Cf_str))
817 {
818 ITHACAstream::ReadDenseTensor(Cf_tensor, "./ITHACAoutput/Matrices/", Cf_str);
819 }
820 else
821 {
823 }
824
825 word BP_str = "BP_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
826 NSUPmodes);
827
828 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + BP_str))
829 {
830 ITHACAstream::ReadDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/", BP_str);
831 }
832 else
833 {
835 }
836
842
843 if (fluxMethod == "consistent")
844 {
852 }
853 }
854 else
855 {
867
868 if (fluxMethod == "consistent")
869 {
877 }
878 }
879}
880
881// * * * * * * * * * * * * * * Momentum Eq. Methods * * * * * * * * * * * * * //
882
883Eigen::MatrixXd steadyNS::diffusive_term(label NUmodes, label NPmodes,
884 label NSUPmodes)
885{
886 label Bsize = NUmodes + NSUPmodes + liftfield.size();
887 Eigen::MatrixXd B_matrix;
888 B_matrix.resize(Bsize, Bsize);
889
890 // Project everything
891 for (label i = 0; i < Bsize; i++)
892 {
893 for (label j = 0; j < Bsize; j++)
894 {
895 B_matrix(i, j) = fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::laplacian(
896 dimensionedScalar("1", dimless, 1), L_U_SUPmodes[j])).value();
897 }
898 }
899
900 if (Pstream::master())
901 {
902 ITHACAstream::SaveDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/",
903 "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(NSUPmodes));
904 }
905
906 return B_matrix;
907}
908
909Eigen::MatrixXd steadyNS::diffusive_term_sym(label NUmodes, label NPmodes,
910 label NSUPmodes)
911{
912 label Bsize = NUmodes + NSUPmodes + liftfield.size();
913 Eigen::MatrixXd B_matrix;
914 B_matrix.resize(Bsize, Bsize);
915
916 // Project everything
917 for (label i = 0; i < Bsize; i++)
918 {
919 for (label j = 0; j < Bsize; j++)
920 {
921 B_matrix(i, j) = - fvc::domainIntegrate(fvc::grad(L_U_SUPmodes[i])
922 && fvc::grad(L_U_SUPmodes[j])).value();
923 }
924 }
925
926 if (Pstream::master())
927 {
928 ITHACAstream::SaveDenseMatrix(B_matrix, "./ITHACAoutput/Matrices/",
929 "B_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(NSUPmodes));
930 }
931
932 return B_matrix;
933}
934
936 label NSUPmodes)
937{
938 label K1size = NUmodes + NSUPmodes + liftfield.size();
939 label K2size = NPmodes + liftfieldP.size();
940 Eigen::MatrixXd K_matrix(K1size, K2size);
941
942 // Project everything
943 for (label i = 0; i < K1size; i++)
944 {
945 for (label j = 0; j < K2size; j++)
946 {
947 K_matrix(i, j) = fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::grad(
948 Pmodes[j])).value();
949 }
950 }
951
952 if (Pstream::master())
953 {
954 ITHACAstream::SaveDenseMatrix(K_matrix, "./ITHACAoutput/Matrices/",
955 "K_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
956 NSUPmodes) + "_" + name(NPmodes));
957 }
958
959 return K_matrix;
960}
961
962List <Eigen::MatrixXd> steadyNS::convective_term(label NUmodes, label NPmodes,
963 label NSUPmodes)
964{
965 label Csize = NUmodes + NSUPmodes + liftfield.size();
966 List <Eigen::MatrixXd> C_matrix;
967 C_matrix.setSize(Csize);
968
969 for (label j = 0; j < Csize; j++)
970 {
971 C_matrix[j].resize(Csize, Csize);
972 }
973
974 for (label i = 0; i < Csize; i++)
975 {
976 for (label j = 0; j < Csize; j++)
977 {
978 for (label k = 0; k < Csize; k++)
979 {
980 C_matrix[i](j, k) = fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::div(
981 linearInterpolate(L_U_SUPmodes[j]) & L_U_SUPmodes[j].mesh().Sf(),
982 L_U_SUPmodes[k])).value();
983 }
984 }
985 }
986
987 if (Pstream::master())
988 {
989 // Export the matrix
990 ITHACAstream::exportMatrix(C_matrix, "C", "python", "./ITHACAoutput/Matrices/");
991 ITHACAstream::exportMatrix(C_matrix, "C", "matlab", "./ITHACAoutput/Matrices/");
992 ITHACAstream::exportMatrix(C_matrix, "C", "eigen", "./ITHACAoutput/Matrices/C");
993 return C_matrix;
994 }
995}
996
997Eigen::Tensor<double, 3> steadyNS::convective_term_tens(label NUmodes,
998 label NPmodes,
999 label NSUPmodes)
1000{
1001 label Csize = NUmodes + NSUPmodes + liftfield.size();
1002 Eigen::Tensor<double, 3> C_tensor;
1003 C_tensor.resize(Csize, Csize, Csize);
1004
1005 for (label i = 0; i < Csize; i++)
1006 {
1007 for (label j = 0; j < Csize; j++)
1008 {
1009 for (label k = 0; k < Csize; k++)
1010 {
1011 if (fluxMethod == "consistent")
1012 {
1013 C_tensor(i, j, k) = fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::div(
1014 L_PHImodes[j],
1015 L_U_SUPmodes[k])).value();
1016 }
1017 else
1018 {
1019 C_tensor(i, j, k) = fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::div(
1020 linearInterpolate(L_U_SUPmodes[j]) & L_U_SUPmodes[j].mesh().Sf(),
1021 L_U_SUPmodes[k])).value();
1022 }
1023 }
1024 }
1025 }
1026
1027 if (Pstream::master())
1028 {
1029 // Export the tensor
1030 ITHACAstream::SaveDenseTensor(C_tensor, "./ITHACAoutput/Matrices/",
1031 "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
1032 NSUPmodes) + "_t");
1033 }
1034
1035 return C_tensor;
1036}
1037
1038Eigen::Tensor<double, 3> steadyNS::convective_term_tens_cache(label NUmodes,
1039 label NPmodes,
1040 label NSUPmodes)
1041{
1042 label Csize = NUmodes + NSUPmodes + liftfield.size();
1043 Eigen::Tensor<double, 3> C_tensor(Csize, Csize, Csize);
1044
1045 for (label j = 0; j < Csize; ++j)
1046 {
1047 surfaceScalarField SfUj = linearInterpolate(L_U_SUPmodes[j]) &
1048 L_U_SUPmodes[j].mesh().Sf();
1049 // Cache all div fields for fixed j, all k
1050 List<tmp<volVectorField>> divRow(Csize);
1051
1052 for (label k = 0; k < Csize; ++k)
1053 {
1054 if (fluxMethod == "consistent")
1055 {
1056 divRow[k] = fvc::div(L_PHImodes[j], L_U_SUPmodes[k]);
1057 }
1058 else
1059 {
1060 divRow[k] = fvc::div(SfUj, L_U_SUPmodes[k]);
1061 }
1062 }
1063
1064 for (label i = 0; i < Csize; ++i)
1065 {
1066 for (label k = 0; k < Csize; ++k)
1067 {
1068 const volVectorField& divField = divRow[k]();
1069 C_tensor(i, j, k) = fvc::domainIntegrate(L_U_SUPmodes[i] & divField).value();
1070 }
1071 }
1072 }
1073
1074 if (Pstream::master())
1075 {
1076 // Export the tensor
1077 ITHACAstream::SaveDenseTensor(C_tensor, "./ITHACAoutput/Matrices/",
1078 "C_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
1079 NSUPmodes) + "_t");
1080 }
1081
1082 return C_tensor;
1083}
1084
1085Eigen::MatrixXd steadyNS::mass_term(label NUmodes, label NPmodes,
1086 label NSUPmodes)
1087{
1088 label Msize = NUmodes + NSUPmodes + liftfield.size();
1089 Eigen::MatrixXd M_matrix(Msize, Msize);
1090
1091 // Project everything
1092 for (label i = 0; i < Msize; i++)
1093 {
1094 for (label j = 0; j < Msize; j++)
1095 {
1096 M_matrix(i, j) = fvc::domainIntegrate(L_U_SUPmodes[i] &
1097 L_U_SUPmodes[j]).value();
1098 }
1099 }
1100
1101 if (Pstream::master())
1102 {
1103 ITHACAstream::SaveDenseMatrix(M_matrix, "./ITHACAoutput/Matrices/",
1104 "M_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(NSUPmodes));
1105 }
1106
1107 return M_matrix;
1108}
1109
1110// * * * * * * * * * * * * * * Continuity Eq. Methods * * * * * * * * * * * * * //
1111
1112Eigen::MatrixXd steadyNS::divergence_term(label NUmodes, label NPmodes,
1113 label NSUPmodes)
1114{
1115 label P1size = NPmodes;
1116 label P2size = NUmodes + NSUPmodes + liftfield.size();
1117 Eigen::MatrixXd P_matrix(P1size, P2size);
1118
1119 // Project everything
1120 for (label i = 0; i < P1size; i++)
1121 {
1122 for (label j = 0; j < P2size; j++)
1123 {
1124 P_matrix(i, j) = fvc::domainIntegrate(Pmodes[i] * fvc::div (
1125 L_U_SUPmodes[j])).value();
1126 }
1127 }
1128
1129 if (Pstream::master())
1130 {
1131 ITHACAstream::SaveDenseMatrix(P_matrix, "./ITHACAoutput/Matrices/",
1132 "P_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
1133 NSUPmodes) + "_" + name(NPmodes));
1134 }
1135
1136 return P_matrix;
1137}
1138
1139
1140List <Eigen::MatrixXd> steadyNS::div_momentum(label NUmodes, label NPmodes)
1141{
1142 label G1size = NPmodes;
1143 label G2size = NUmodes + NSUPmodes + liftfield.size();
1144 List <Eigen::MatrixXd> G_matrix;
1145 G_matrix.setSize(G1size);
1146
1147 for (label j = 0; j < G1size; j++)
1148 {
1149 G_matrix[j].resize(G2size, G2size);
1150 }
1151
1152 for (label i = 0; i < G1size; i++)
1153 {
1154 for (label j = 0; j < G2size; j++)
1155 {
1156 for (label k = 0; k < G2size; k++)
1157 {
1158 G_matrix[i](j, k) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & (fvc::div(
1159 fvc::interpolate(L_U_SUPmodes[j]) & L_U_SUPmodes[j].mesh().Sf(),
1160 L_U_SUPmodes[k]))).value();
1161 }
1162 }
1163 }
1164
1165 if (Pstream::master())
1166 {
1167 // Export the matrix
1168 ITHACAstream::exportMatrix(G_matrix, "G", "python", "./ITHACAoutput/Matrices/");
1169 ITHACAstream::exportMatrix(G_matrix, "G", "matlab", "./ITHACAoutput/Matrices/");
1170 ITHACAstream::exportMatrix(G_matrix, "G", "eigen", "./ITHACAoutput/Matrices/G");
1171 }
1172
1173 return G_matrix;
1174}
1175
1176Eigen::Tensor<double, 3> steadyNS::divMomentum(label NUmodes, label NPmodes)
1177{
1178 label g1Size = NPmodes + liftfieldP.size();
1179 label g2Size = NUmodes + NSUPmodes + liftfield.size();
1180 Eigen::Tensor<double, 3> gTensor;
1181 gTensor.resize(g1Size, g2Size, g2Size);
1182
1183 for (label i = 0; i < g1Size; i++)
1184 {
1185 for (label j = 0; j < g2Size; j++)
1186 {
1187 for (label k = 0; k < g2Size; k++)
1188 {
1189 gTensor(i, j, k) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & (fvc::div(
1190 fvc::interpolate(L_U_SUPmodes[j]) & L_U_SUPmodes[j].mesh().Sf(),
1191 L_U_SUPmodes[k]))).value();
1192 }
1193 }
1194 }
1195
1196 if (Pstream::master())
1197 {
1198 // Export the tensor
1199 ITHACAstream::SaveDenseTensor(gTensor, "./ITHACAoutput/Matrices/",
1200 "G_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
1201 NSUPmodes) + "_" + name(NPmodes) + "_t");
1202 }
1203
1204 return gTensor;
1205}
1206
1207Eigen::Tensor<double, 3> steadyNS::divMomentum_cache(label NUmodes,
1208 label NPmodes)
1209{
1210 label g1Size = NPmodes + liftfieldP.size();
1211 label g2Size = NUmodes + NSUPmodes + liftfield.size();
1212 Eigen::Tensor<double, 3> gTensor (g1Size, g2Size, g2Size);
1213 List<tmp<volVectorField>> PmodesGrad(g1Size);
1214
1215 for (label i = 0; i < g1Size; i++)
1216 {
1217 PmodesGrad[i] = fvc::grad(Pmodes[i]);
1218 }
1219
1220 for (label j = 0; j < g2Size; j++)
1221 {
1222 surfaceScalarField interpUj = fvc::interpolate(L_U_SUPmodes[j]) &
1223 L_U_SUPmodes[j].mesh().Sf();
1224 // Cache only one row (j fixed)
1225 List<tmp<volVectorField>> divRow(g2Size);
1226
1227 for (label k = 0; k < g2Size; ++k)
1228 {
1229 divRow[k] = fvc::div(interpUj, L_U_SUPmodes[k]);
1230 }
1231
1232 for (label i = 0; i < g1Size; i++)
1233 {
1234 const volVectorField& gradPi = PmodesGrad[i]();
1235
1236 for (label k = 0; k < g2Size; k++)
1237 {
1238 const volVectorField& divField = divRow[k]();
1239 gTensor(i, j, k) = fvc::domainIntegrate(gradPi & divField).value();
1240 }
1241 }
1242 }
1243
1244 if (Pstream::master())
1245 {
1246 // Export the tensor
1247 ITHACAstream::SaveDenseTensor(gTensor, "./ITHACAoutput/Matrices/",
1248 "G_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
1249 NSUPmodes) + "_" + name(NPmodes) + "_t");
1250 }
1251
1252 return gTensor;
1253}
1254
1255// large scale convection (or background convection)
1257 volVectorField vls)
1258{
1259 label Lsize = NUmodes + liftfield.size();
1260 Eigen::MatrixXd L_matrix(Lsize, Lsize);
1261
1262 for (label i = 0; i < Lsize; i++)
1263 {
1264 for (label j = 0; j < Lsize; j++)
1265 {
1266 L_matrix(i, j) = - fvc::domainIntegrate(L_U_SUPmodes[i] & fvc::div(
1267 fvc::interpolate(vls) & vls.mesh().Sf(),
1268 L_U_SUPmodes[j])).value();
1269 }
1270 }
1271
1272 return L_matrix;
1273}
1274
1276 label NUmodes, volVectorField vls)
1277{
1278 label LDsize1 = NPmodes + liftfieldP.size();
1279 label LDsize2 = NUmodes + liftfield.size();
1280 Eigen::MatrixXd L_D_matrix(LDsize1, LDsize2);
1281
1282 for (label i = 0; i < LDsize1; i++)
1283 {
1284 for (label j = 0; j < LDsize2; j++)
1285 {
1286 L_D_matrix(i, j) = - fvc::domainIntegrate(fvc::grad(Pmodes[i]) & fvc::div(
1287 fvc::interpolate(vls) & vls.mesh().Sf(),
1288 L_U_SUPmodes[j])).value();
1289 }
1290 }
1291
1292 return L_D_matrix;
1293}
1294
1296{
1297 label Dsize = NPmodes + liftfieldP.size();
1298 Eigen::MatrixXd D_matrix(Dsize, Dsize);
1299
1300 // Project everything
1301 for (label i = 0; i < Dsize; i++)
1302 {
1303 for (label j = 0; j < Dsize; j++)
1304 {
1305 D_matrix(i, j) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & fvc::grad(
1306 Pmodes[j])).value();
1307 }
1308 }
1309
1310 if (Pstream::master())
1311 {
1312 ITHACAstream::SaveDenseMatrix(D_matrix, "./ITHACAoutput/Matrices/",
1313 "D_" + name(NPmodes));
1314 }
1315
1316 return D_matrix;
1317}
1318
1319Eigen::MatrixXd steadyNS::pressure_BC1(label NUmodes, label NPmodes)
1320{
1321 label P_BC1size = NPmodes;
1322 label P_BC2size = NUmodes + liftfield.size();
1323 Eigen::MatrixXd BC1_matrix(P_BC1size, P_BC2size);
1324 const fvMesh& mesh = L_U_SUPmodes[0].mesh();
1325
1326 for (label i = 0; i < P_BC1size; i++)
1327 {
1328 for (label j = 0; j < P_BC2size; j++)
1329 {
1330 surfaceScalarField lpl((fvc::interpolate(fvc::laplacian(
1331 L_U_SUPmodes[j])) & mesh.Sf()) * fvc::interpolate(Pmodes[i]));
1332 double s = 0;
1333
1334 for (label k = 0; k < lpl.boundaryField().size(); k++)
1335 {
1336 if (lpl.boundaryField()[k].type() != "processor")
1337 {
1338 s += gSum(lpl.boundaryField()[k]);
1339 }
1340 }
1341
1342 BC1_matrix(i, j) = s;
1343 }
1344 }
1345
1346 if (Pstream::master())
1347 {
1348 ITHACAstream::SaveDenseMatrix(BC1_matrix, "./ITHACAoutput/Matrices/",
1349 "BC1_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(NPmodes));
1350 }
1351
1352 return BC1_matrix;
1353}
1354
1355
1356List <Eigen::MatrixXd> steadyNS::pressure_BC2(label NUmodes, label NPmodes)
1357{
1358 label P2_BC1size = NPmodes;
1359 label P2_BC2size = NUmodes + NSUPmodes + liftfield.size();
1360 List <Eigen::MatrixXd> BC2_matrix;
1361 const fvMesh& mesh = L_U_SUPmodes[0].mesh();
1362 BC2_matrix.setSize(P2_BC1size);
1363
1364 for (label j = 0; j < P2_BC1size; j++)
1365 {
1366 BC2_matrix[j].resize(P2_BC2size, P2_BC2size);
1367 }
1368
1369 for (label i = 0; i < P2_BC1size; i++)
1370 {
1371 for (label j = 0; j < P2_BC2size; j++)
1372 {
1373 for (label k = 0; k < P2_BC2size; k++)
1374 {
1375 surfaceScalarField div_m(fvc::interpolate(fvc::div(fvc::interpolate(
1376 L_U_SUPmodes[j]) & mesh.Sf(),
1377 L_U_SUPmodes[k])) & mesh.Sf() * fvc::interpolate(Pmodes[i]));
1378 double s = 0;
1379
1380 for (label k = 0; k < div_m.boundaryField().size(); k++)
1381 {
1382 if (div_m.boundaryField()[k].type() != "processor")
1383 {
1384 s += gSum(div_m.boundaryField()[k]);
1385 }
1386 }
1387
1388 BC2_matrix[i](j, k) = s;
1389 }
1390 }
1391 }
1392
1393 return BC2_matrix;
1394}
1395
1396Eigen::Tensor<double, 3> steadyNS::pressureBC2(label NUmodes, label NPmodes)
1397{
1398 label pressureBC1Size = NPmodes;
1399 label pressureBC2Size = NUmodes + NSUPmodes + liftfield.size();
1400 Eigen::Tensor<double, 3> bc2Tensor;
1401 const fvMesh& mesh = L_U_SUPmodes[0].mesh();
1402 bc2Tensor.resize(pressureBC1Size, pressureBC2Size, pressureBC2Size);
1403
1404 for (label i = 0; i < pressureBC1Size; i++)
1405 {
1406 for (label j = 0; j < pressureBC2Size; j++)
1407 {
1408 for (label k = 0; k < pressureBC2Size; k++)
1409 {
1410 surfaceScalarField div_m(fvc::interpolate(fvc::div(fvc::interpolate(
1411 L_U_SUPmodes[j]) & mesh.Sf(),
1412 L_U_SUPmodes[k])) & mesh.Sf() * fvc::interpolate(Pmodes[i]));
1413 double s = 0;
1414
1415 for (label k = 0; k < div_m.boundaryField().size(); k++)
1416 {
1417 if (div_m.boundaryField()[k].type() != "processor")
1418 {
1419 s += gSum(div_m.boundaryField()[k]);
1420 }
1421 }
1422
1423 bc2Tensor(i, j, k) = s;
1424 }
1425 }
1426 }
1427
1428 if (Pstream::master())
1429 {
1430 // Export the tensor
1431 ITHACAstream::SaveDenseTensor(bc2Tensor, "./ITHACAoutput/Matrices/",
1432 "BC2_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
1433 NSUPmodes) + "_" + name(NPmodes) + "_t");
1434 }
1435
1436 return bc2Tensor;
1437}
1438
1439Eigen::MatrixXd steadyNS::pressure_BC3(label NUmodes, label NPmodes)
1440{
1441 label P3_BC1size = NPmodes;
1442 label P3_BC2size = NUmodes + liftfield.size();
1443 Eigen::MatrixXd BC3_matrix(P3_BC1size, P3_BC2size);
1444 const fvMesh& mesh = L_U_SUPmodes[0].mesh();
1445 surfaceVectorField n(mesh.Sf() / mesh.magSf());
1446
1447 for (label i = 0; i < P3_BC1size; i++)
1448 {
1449 for (label j = 0; j < P3_BC2size; j++)
1450 {
1451 surfaceVectorField BC3 = fvc::interpolate(fvc::curl(L_U_SUPmodes[j])).ref();
1452 surfaceVectorField BC4 = (n ^ fvc::interpolate(fvc::grad(Pmodes[i]))).ref();
1453 surfaceScalarField BC5 = ((BC3 & BC4) * mesh.magSf()).ref();
1454 double s = 0;
1455
1456 for (label k = 0; k < BC5.boundaryField().size(); k++)
1457 {
1458 if (BC5.boundaryField()[k].type() != "processor")
1459 {
1460 s += gSum(BC5.boundaryField()[k]);
1461 }
1462 }
1463
1464 BC3_matrix(i, j) = s;
1465 }
1466 }
1467
1468 if (Pstream::master())
1469 {
1470 ITHACAstream::SaveDenseMatrix(BC3_matrix, "./ITHACAoutput/Matrices/",
1471 "BC3_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(NPmodes));
1472 }
1473
1474 return BC3_matrix;
1475}
1476
1477Eigen::MatrixXd steadyNS::pressure_BC4(label NUmodes, label NPmodes)
1478{
1479 label P4_BC1size = NPmodes;
1480 label P4_BC2size = NUmodes + liftfield.size();
1481 Eigen::MatrixXd BC4_matrix(P4_BC1size, P4_BC2size);
1482 const fvMesh& mesh = L_U_SUPmodes[0].mesh();
1483 surfaceVectorField n(mesh.Sf() / mesh.magSf());
1484
1485 for (label i = 0; i < P4_BC1size; i++)
1486 {
1487 for (label j = 0; j < P4_BC2size; j++)
1488 {
1489 surfaceScalarField BC3 = fvc::interpolate(Pmodes[i]).ref();
1490 surfaceScalarField BC4 = (n & fvc::interpolate(L_U_SUPmodes[j])).ref();
1491 surfaceScalarField BC5 = ((BC3 * BC4) * mesh.magSf()).ref();
1492 double s = 0;
1493
1494 for (label k = 0; k < BC5.boundaryField().size(); k++)
1495 {
1496 if (BC5.boundaryField()[k].type() != "processor")
1497 {
1498 s += gSum(BC5.boundaryField()[k]);
1499 }
1500 }
1501
1502 BC4_matrix(i, j) = s;
1503 }
1504 }
1505
1506 if (Pstream::master())
1507 {
1508 ITHACAstream::SaveDenseMatrix(BC4_matrix, "./ITHACAoutput/Matrices/",
1509 "BC4_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(NPmodes));
1510 }
1511
1512 return BC4_matrix;
1513}
1514
1515List<Eigen::MatrixXd> steadyNS::bcVelocityVec(label NUmodes,
1516 label NSUPmodes)
1517{
1518 label BCsize = NUmodes + NSUPmodes;
1519 List <Eigen::MatrixXd> bcVelVec(inletIndex.rows());
1520
1521 for (label j = 0; j < inletIndex.rows(); j++)
1522 {
1523 bcVelVec[j].resize(BCsize, 1);
1524 }
1525
1526 for (label k = 0; k < inletIndex.rows(); k++)
1527 {
1528 label BCind = inletIndex(k, 0);
1529 label BCcomp = inletIndex(k, 1);
1530
1531 for (label i = 0; i < BCsize; i++)
1532 {
1533 bcVelVec[k](i, 0) = gSum(L_U_SUPmodes[i].boundaryField()[BCind].component(
1534 BCcomp));
1535 }
1536 }
1537
1538 if (Pstream::master())
1539 {
1540 ITHACAstream::exportMatrix(bcVelVec, "bcVelVec", "eigen",
1541 "./ITHACAoutput/Matrices/bcVelVec");
1542 }
1543
1544 return bcVelVec;
1545}
1546
1547List<Eigen::MatrixXd> steadyNS::bcVelocityMat(label NUmodes,
1548 label NSUPmodes)
1549{
1550 label BCsize = NUmodes + NSUPmodes;
1551 label BCUsize = inletIndex.rows();
1552 List <Eigen::MatrixXd> bcVelMat(BCUsize);
1553
1554 for (label j = 0; j < inletIndex.rows(); j++)
1555 {
1556 bcVelMat[j].resize(BCsize, BCsize);
1557 }
1558
1559 for (label k = 0; k < inletIndex.rows(); k++)
1560 {
1561 label BCind = inletIndex(k, 0);
1562 label BCcomp = inletIndex(k, 1);
1563
1564 for (label i = 0; i < BCsize; i++)
1565 {
1566 for (label j = 0; j < BCsize; j++)
1567 {
1568 bcVelMat[k](i, j) = gSum(L_U_SUPmodes[i].boundaryField()[BCind].component(
1569 BCcomp) *
1570 L_U_SUPmodes[j].boundaryField()[BCind].component(BCcomp));
1571 }
1572 }
1573 }
1574
1575 if (Pstream::master())
1576 {
1577 ITHACAstream::exportMatrix(bcVelMat, "bcVelMat", "eigen",
1578 "./ITHACAoutput/Matrices/bcVelMat");
1579 }
1580
1581 return bcVelMat;
1582}
1583
1585 label NPmodes, label NSUPmodes)
1586{
1587 label BPsize1 = NPmodes;
1588 label BPsize2 = NUmodes + NSUPmodes + liftfield.size();
1589 Eigen::MatrixXd BP_matrix(BPsize1, BPsize2);
1590 volVectorField L_U_SUPmodesaux(L_U_SUPmodes[0]);
1591
1592 for (label i = 0; i < BPsize1; i++)
1593 {
1594 for (label j = 0; j < BPsize2; j++)
1595 {
1596 L_U_SUPmodesaux = dt_dummy * fvc::laplacian(
1597 nu_dummy(), L_U_SUPmodes[j]);
1598 BP_matrix(i, j) = fvc::domainIntegrate(Pmodes[i] *
1599 fvc::div(L_U_SUPmodesaux)).value();
1600 }
1601 }
1602
1603 if (Pstream::master())
1604 {
1605 ITHACAstream::SaveDenseMatrix(BP_matrix, "./ITHACAoutput/Matrices/",
1606 "BP_" + name(NPmodes));
1607 }
1608
1609 return BP_matrix;
1610}
1611
1613 label NPmodes, label NSUPmodes)
1614{
1615 Eigen::VectorXd ModeVector;
1616 label BCsize = inletIndex.rows();
1617 label RDsize = NUmodes + NSUPmodes;
1618 List <Eigen::MatrixXd> RD_matrix(BCsize);
1619 Eigen::SparseMatrix<double> A;
1620 Eigen::VectorXd b;
1621
1622 for (label i = 0; i < BCsize; i++)
1623 {
1624 label BCind = inletIndex(i, 0);
1625 Vector<double> v(0, 0, 0);
1626 v[inletIndex(i, 1)] = 1;
1627 volVectorField Upara(Uinl());
1628 assignBC(Upara, BCind, v);
1629 // Determine boundary vector
1630 fvVectorMatrix UEqn
1631 (
1632 -fvm::laplacian(nu_dummy(), Upara)
1633 );
1634 Foam2Eigen::fvMatrix2Eigen(UEqn, A, b);
1635 RD_matrix[i].resize(RDsize, 1);
1636
1637 for (label j = 0; j < RDsize; j++)
1638 {
1639 ModeVector = Foam2Eigen::field2Eigen(L_U_SUPmodes[j]);
1640 RD_matrix[i](j, 0) = ModeVector.dot(b.col(0));
1641 }
1642
1643 if (Pstream::parRun())
1644 {
1645 reduce(RD_matrix[i], sumOp<Eigen::MatrixXd>());
1646 }
1647
1648 if (Pstream::master())
1649 {
1650 ITHACAstream::SaveDenseMatrix(RD_matrix[i], "./ITHACAoutput/Matrices/RD/",
1651 "RD" + name(i) + "_" + name(NUmodes) + "_" + name(NSUPmodes));
1652 }
1653 }
1654
1655 return RD_matrix;
1656}
1657
1659 label NPmodes, label NSUPmodes)
1660{
1661 Eigen::VectorXd ModeVector;
1662 label BCsize = inletIndex.rows();
1663 label RCsize = NUmodes + NSUPmodes;
1664 List <Eigen::MatrixXd> RC_matrix(BCsize);
1665 Eigen::SparseMatrix<double> A;
1666 Eigen::VectorXd b;
1667
1668 for (label i = 0; i < BCsize; i++)
1669 {
1670 label BCind = inletIndex(i, 0);
1671 Vector<double> v(0, 0, 0);
1672 v[inletIndex(i, 1)] = 1;
1673 volVectorField Upara(Uinl());
1674 assignBC(Upara, BCind, v);
1675 // Determine boundary vector
1676 fvVectorMatrix UEqn
1677 (
1678 fvm::div(fvc::flux(Upara), Upara)
1679 );
1680 Foam2Eigen::fvMatrix2Eigen(UEqn, A, b);
1681 RC_matrix[i].resize(RCsize, 1);
1682
1683 for (label j = 0; j < RCsize; j++)
1684 {
1685 ModeVector = Foam2Eigen::field2Eigen(L_U_SUPmodes[j]);
1686 RC_matrix[i](j, 0) = ModeVector.dot(b.col(0));
1687 }
1688
1689 if (Pstream::parRun())
1690 {
1691 reduce(RC_matrix[i], sumOp<Eigen::MatrixXd>());
1692 }
1693
1694 if (Pstream::master())
1695 {
1696 ITHACAstream::SaveDenseMatrix(RC_matrix[i], "./ITHACAoutput/Matrices/RC/",
1697 "RC" + name(i) + "_" + name(NUmodes) + "_" + name(NSUPmodes));
1698 }
1699 }
1700
1701 return RC_matrix;
1702}
1703
1704Eigen::Tensor<double, 3> steadyNS::convective_term_flux_tens(label NUmodes,
1705 label NPmodes, label NSUPmodes)
1706{
1707 label Csize1 = NUmodes + NSUPmodes + liftfield.size();
1708 label Csize2 = NPmodes;
1709 Eigen::Tensor<double, 3> Cf_tensor;
1710 volVectorField L_U_SUPmodesaux(L_U_SUPmodes[0]);
1711 Cf_tensor.resize(Csize2, Csize1, Csize1);
1712
1713 for (label i = 0; i < Csize2; i++)
1714 {
1715 for (label j = 0; j < Csize1; j++)
1716 {
1717 for (label k = 0; k < Csize1; k++)
1718 {
1719 if (fluxMethod == "consistent")
1720 {
1721 L_U_SUPmodesaux = dt_dummy * (fvc::div(
1722 L_PHImodes[j],
1723 L_U_SUPmodes[k]));
1724 }
1725 else
1726 {
1727 L_U_SUPmodesaux = dt_dummy * (fvc::div(
1728 fvc::flux(L_U_SUPmodes[j]),
1729 L_U_SUPmodes[k]));
1730 }
1731
1732 Cf_tensor(i, j, k) = fvc::domainIntegrate(Pmodes[i]
1733 * fvc::div(L_U_SUPmodesaux)).value();
1734 }
1735 }
1736 }
1737
1738 if (Pstream::master())
1739 {
1740 ITHACAstream::SaveDenseTensor(Cf_tensor, "./ITHACAoutput/Matrices/",
1741 "Cf_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
1742 NSUPmodes) + "_t");
1743 }
1744
1745 return Cf_tensor;
1746}
1747
1749{
1750 label BCsize = inletIndex.rows();
1751 List<Eigen::MatrixXd> LinSysDivDummy;
1752 LinSysDiv.resize(BCsize + 1);
1753 volScalarField p(_p());
1754
1755 for (label i = 0; i < BCsize; i++)
1756 {
1757 label BCind = inletIndex(i, 0);
1758 Vector<double> v(0, 0, 0);
1759 v[inletIndex(i, 1)] = 1;
1760 volVectorField Upara(Uinl());
1761 assignBC(Upara, BCind, v);
1762 fvScalarMatrix pEqn
1763 (
1764 fvm::laplacian(p) == (1.0 / dt_dummy) * fvc::div(Upara)
1765 );
1766 pEqn.setReference(0, 0);
1767 LinSysDivDummy = Pmodes.project(pEqn, NPmodes);
1768
1769 if (i == 0)
1770 {
1771 LinSysDiv[0] = LinSysDivDummy[0];
1772 }
1773
1774 LinSysDiv[i + 1] = LinSysDivDummy[1];
1775 }
1776
1777 return LinSysDiv;
1778}
1779
1781 label NPmodes)
1782{
1783 label BCsize = inletIndex.rows();
1784 List<Eigen::MatrixXd> LinSysConvDummy;
1785 LinSysConv.resize(BCsize + 1);
1786 volScalarField p(_p());
1787
1788 for (label i = 0; i < BCsize; i++)
1789 {
1790 label BCind = inletIndex(i, 0);
1791 Vector<double> v(0, 0, 0);
1792 v[inletIndex(i, 1)] = 1;
1793 volVectorField Upara(Uinl());
1794 assignBC(Upara, BCind, v);
1795 volVectorField Caux(L_U_SUPmodes[0]);
1796 Caux = dt_dummy * (-fvc::div(fvc::flux(Upara), Upara));
1797 fvScalarMatrix pEqn
1798 (
1799 fvm::laplacian(p) == (1.0 / dt_dummy) * fvc::div(Caux)
1800 );
1801 pEqn.setReference(0, 0);
1802 LinSysConvDummy = Pmodes.project(pEqn, NPmodes);
1803
1804 if (i == 0)
1805 {
1806 LinSysConv[0] = LinSysConvDummy[0];
1807 }
1808
1809 LinSysConv[i + 1] = LinSysConvDummy[1];
1810 }
1811
1812 return LinSysConv;
1813}
1814
1816 label NPmodes)
1817{
1818 label BCsize = inletIndex.rows();
1819 List<Eigen::MatrixXd> LinSysDiffDummy;
1820 LinSysDiff.resize(BCsize + 1);
1821 volScalarField p(_p());
1822
1823 for (label i = 0; i < BCsize; i++)
1824 {
1825 label BCind = inletIndex(i, 0);
1826 Vector<double> v(0, 0, 0);
1827 v[inletIndex(i, 1)] = 1;
1828 volVectorField Upara(Uinl());
1829 assignBC(Upara, BCind, v);
1830 volVectorField Daux(L_U_SUPmodes[0]);
1831 Daux = dt_dummy * fvc::laplacian(nu_dummy(), Upara);
1832 fvScalarMatrix pEqn
1833 (
1834 fvm::laplacian(p) == (1.0 / dt_dummy) * fvc::div(Daux)
1835 );
1836 pEqn.setReference(0, 0);
1837 LinSysDiffDummy = Pmodes.project(pEqn, NPmodes);
1838
1839 if (i == 0)
1840 {
1841 LinSysDiff[0] = LinSysDiffDummy[0];
1842 }
1843
1844 LinSysDiff[i + 1] = LinSysDiffDummy[1];
1845 }
1846
1847 return LinSysDiff;
1848}
1849
1851 label NPmodes,
1852 label NSUPmodes)
1853{
1854 label Isize = NUmodes + NSUPmodes + liftfield.size();
1855 Eigen::MatrixXd I_matrix;
1856 I_matrix.resize(NUmodes, Isize);
1857
1858 for (label i = 0; i < NUmodes; i++)
1859 {
1860 volVectorField CoeffA = (fvc::reconstruct(L_PHImodes[i])).ref();
1861
1862 for (label j = 0; j < Isize; j++)
1863 {
1864 surfaceScalarField B = fvc::flux(L_U_SUPmodes[j]).ref();
1865 volVectorField CoeffB = fvc::reconstruct(B).ref();
1866 I_matrix(i, j) = fvc::domainIntegrate(CoeffA & CoeffB).value();
1867 }
1868 }
1869
1870 if (Pstream::master())
1871 {
1872 ITHACAstream::SaveDenseMatrix(I_matrix, "./ITHACAoutput/Matrices/",
1873 "I_" + name(NUmodes));
1874 }
1875
1876 return I_matrix;
1877}
1878
1880 label NPmodes,
1881 label NSUPmodes)
1882{
1883 label DFsize = NUmodes + NSUPmodes + liftfield.size();
1884 Eigen::MatrixXd DF_matrix;
1885 DF_matrix.resize(DFsize, NUmodes);
1886 surfaceScalarField phi_tmp("Phi_tmp", _phi());
1887
1888 for (label i = 0; i < NUmodes; i++)
1889 {
1890 for (label j = 0; j < DFsize; j++)
1891 {
1892 phi_tmp = dt_dummy * nu_dummy() * fvc::flux(fvc::laplacian(
1893 dimensionedScalar("1", dimless, 1),
1894 L_U_SUPmodes[j]));
1895 volVectorField CoeffB = fvc::reconstruct(phi_tmp).ref();
1896 volVectorField CoeffA = fvc::reconstruct(L_PHImodes[i]).ref();
1897 DF_matrix(i, j) = fvc::domainIntegrate(CoeffA & CoeffB).value();
1898 }
1899 }
1900
1901 if (Pstream::master())
1902 {
1903 ITHACAstream::SaveDenseMatrix(DF_matrix, "./ITHACAoutput/Matrices/",
1904 "DF_" + name(NUmodes) + "_" + name(
1905 NSUPmodes));
1906 }
1907
1908 return DF_matrix;
1909}
1910
1912 label NPmodes,
1913 label NSUPmodes)
1914{
1915 label KF1size = NUmodes ;
1916 label KF2size = NPmodes;
1917 Eigen::MatrixXd KF_matrix(KF1size, KF2size);
1918
1919 for (label i = 0; i < KF1size; i++)
1920 {
1921 for (label j = 0; j < KF2size; j++)
1922 {
1923 volVectorField CoeffA = (fvc::reconstruct(dt_dummy * fvc::snGrad(
1924 Pmodes[j]) * mag(Pmodes[j].mesh().magSf()))).ref();
1925 volVectorField CoeffB = fvc::reconstruct(L_PHImodes[i]).ref();
1926 KF_matrix(i, j) = fvc::domainIntegrate(CoeffA & CoeffB).value();
1927 }
1928 }
1929
1930 if (Pstream::master())
1931 {
1932 ITHACAstream::SaveDenseMatrix(KF_matrix, "./ITHACAoutput/Matrices/",
1933 "KF_" + name(NUmodes) + "_" + name(
1934 NSUPmodes));
1935 }
1936
1937 return KF_matrix;
1938}
1939
1941 label NUmodes,
1942 label NPmodes, label NSUPmodes)
1943{
1944 label Csize1 = NUmodes + NSUPmodes + liftfield.size();
1945 Eigen::Tensor<double, 3> Ci_tensor;
1946 volVectorField L_U_SUPmodesaux(L_U_SUPmodes[0]);
1947 Ci_tensor.resize(NUmodes, Csize1, Csize1);
1948 surfaceScalarField phi_tmp("Phi_tmp", _phi());
1949
1950 for (label i = 0; i < NUmodes; i++)
1951 {
1952 for (label j = 0; j < Csize1; j++)
1953 {
1954 for (label k = 0; k < Csize1; k++)
1955 {
1956 phi_tmp = dt_dummy * fvc::flux(fvc::div(L_PHImodes[i], L_U_SUPmodes[k]));
1957 volVectorField CoeffA = fvc::reconstruct(phi_tmp).ref();
1958 volVectorField CoeffB = fvc::reconstruct(L_PHImodes[j]).ref();
1959 Ci_tensor(i, j, k) = fvc::domainIntegrate(CoeffB & CoeffA).value();
1960 }
1961 }
1962 }
1963
1964 if (Pstream::master())
1965 {
1966 ITHACAstream::SaveDenseTensor(Ci_tensor, "./ITHACAoutput/Matrices/",
1967 "Ci_" + name(liftfield.size()) + "_" + name(NUmodes) + "_" + name(
1968 NSUPmodes) + "_t");
1969 }
1970
1971 return Ci_tensor;
1972}
1973
1975 label NUmodes,
1976 label NSUPmodes)
1977{
1978 label BCsize = inletIndex.rows();
1979 label SDsize = NUmodes + NSUPmodes;
1980 List <Eigen::MatrixXd> SD_matrix(BCsize);
1981 surfaceScalarField phi_tmp("Phi_tmp", _phi());
1982
1983 for (label i = 0; i < BCsize; i++)
1984 {
1985 label BCind = inletIndex(i, 0);
1986 Vector<double> v(0, 0, 0);
1987 v[inletIndex(i, 1)] = 1;
1988 volVectorField Upara(Uinl());
1989 assignBC(Upara, BCind, v);
1990 SD_matrix[i].resize(SDsize, 1);
1991
1992 // Project everything
1993 for (label j = 0; j < SDsize; j++)
1994 {
1995 volVectorField CoeffB = fvc::reconstruct(L_PHImodes[j]).ref();
1996 phi_tmp = dt_dummy * fvc::flux(fvc::laplacian(nu_dummy(), Upara));
1997 volVectorField CoeffA = fvc::reconstruct(phi_tmp).ref();
1998 SD_matrix[i](j, 0) = fvc::domainIntegrate(CoeffA & CoeffB).value();
1999 }
2000
2001 ITHACAstream::SaveDenseMatrix(SD_matrix[i], "./ITHACAoutput/Matrices/SD/",
2002 "SD" + name(i) + "_" + name(NUmodes) + "_" + name(NSUPmodes));
2003 }
2004
2005 return SD_matrix;
2006}
2007
2009 label NUmodes,
2010 label NSUPmodes)
2011{
2012 label BCsize = inletIndex.rows();
2013 label SCsize = NUmodes + NSUPmodes;
2014 List <Eigen::MatrixXd> SC_matrix(BCsize);
2015 surfaceScalarField phi_tmp("Phi_tmp", _phi());
2016
2017 for (label i = 0; i < BCsize; i++)
2018 {
2019 label BCind = inletIndex(i, 0);
2020 Vector<double> v(0, 0, 0);
2021 v[inletIndex(i, 1)] = 1;
2022 volVectorField Upara(Uinl());
2023 assignBC(Upara, BCind, v);
2024 SC_matrix[i].resize(SCsize, 1);
2025
2026 // Project everything
2027 for (label j = 0; j < SCsize; j++)
2028 {
2029 volVectorField CoeffB = fvc::reconstruct(L_PHImodes[j]).ref();
2030 phi_tmp = dt_dummy * fvc::flux(fvc::div(fvc::flux(Upara), Upara));
2031 volVectorField CoeffA = fvc::reconstruct(phi_tmp).ref();
2032 SC_matrix[i](j, 0) = fvc::domainIntegrate(CoeffA & CoeffB).value();
2033 }
2034
2035 ITHACAstream::SaveDenseMatrix(SC_matrix[i], "./ITHACAoutput/Matrices/SC/",
2036 "SD" + name(i) + "_" + name(NUmodes) + "_" + name(NSUPmodes));
2037 }
2038
2039 return SC_matrix;
2040}
2041
2043 label NPmodes,
2044 label NSUPmodes)
2045{
2046 Eigen::MatrixXd W_matrix;
2047 W_matrix.resize(NUmodes, NUmodes);
2048
2049 for (label i = 0; i < NUmodes; i++)
2050 {
2051 volVectorField A = fvc::reconstruct(L_PHImodes[i]).ref();
2052
2053 for (label j = 0; j < NUmodes; j++)
2054 {
2055 volVectorField B = fvc::reconstruct(L_PHImodes[j]).ref();
2056 W_matrix(i, j) = fvc::domainIntegrate(A & B).value();
2057 }
2058 }
2059
2060 ITHACAstream::SaveDenseMatrix(W_matrix, "./ITHACAoutput/Matrices/",
2061 "W_" + name(NUmodes));
2062 return W_matrix;
2063}
2064
2065
2067{
2068 const volScalarField& nu = _laminarTransport().nu();
2069 volScalarField& ciao = const_cast<volScalarField&>(nu);
2070 this->assignIF(ciao, mu);
2071
2072 for (label i = 0; i < ciao.boundaryFieldRef().size(); i++)
2073 {
2074 this->assignBC(ciao, i, mu);
2075 }
2076}
2077
2078
2080{
2081 tauMatrix.resize(L_U_SUPmodes.size(), 3);
2082 nMatrix.resize(NPmodes, 3);
2083 tauMatrix = tauMatrix * 0;
2084 nMatrix = nMatrix * 0;
2085 Time& runTime = _runTime();
2086 instantList Times = runTime.times();
2087 fvMesh& mesh = _mesh();
2088 volScalarField& p = _p();
2089 volVectorField& U = _U();
2090 //Read FORCESdict
2091 IOdictionary FORCESdict
2092 (
2093 IOobject
2094 (
2095 "FORCESdict",
2096 runTime.system(),
2097 mesh,
2098 IOobject::MUST_READ,
2099 IOobject::NO_WRITE
2100 )
2101 );
2102 IOdictionary transportProperties
2103 (
2104 IOobject
2105 (
2106 "transportProperties",
2107 runTime.constant(),
2108 mesh,
2109 IOobject::MUST_READ,
2110 IOobject::NO_WRITE
2111 )
2112 );
2113 word pName(FORCESdict.lookup("pName"));
2114 word UName(FORCESdict.lookup("UName"));
2115 functionObjects::ITHACAforces f("Forces", mesh, FORCESdict);
2116
2117 for (label i = 0; i < L_U_SUPmodes.size(); i++)
2118 {
2119 U = L_U_SUPmodes[i];
2120 p = Pmodes[0];
2121 mesh.readUpdate();
2122 f.write();
2123 f.calcForcesMoment();
2124
2125 for (label j = 0; j < 3; j++)
2126 {
2127 tauMatrix(i, j) = f.forceTau()[j];
2128 }
2129 }
2130
2131 for (label i = 0; i < NPmodes; i++)
2132 {
2133 U = L_U_SUPmodes[0];
2134 p = Pmodes[i];
2135 mesh.readUpdate();
2136 f.write();
2137 f.calcForcesMoment();
2138
2139 for (label j = 0; j < 3; j++)
2140 {
2141 nMatrix(i, j) = f.forcePressure()[j];
2142 }
2143 }
2144
2145 if (para->exportPython)
2146 {
2147 ITHACAstream::exportMatrix(tauMatrix, "tau", "python",
2148 "./ITHACAoutput/Matrices/");
2149 ITHACAstream::exportMatrix(nMatrix, "n", "python", "./ITHACAoutput/Matrices/");
2150 }
2151
2152 if (para->exportMatlab)
2153 {
2154 ITHACAstream::exportMatrix(tauMatrix, "tau", "matlab",
2155 "./ITHACAoutput/Matrices/");
2156 ITHACAstream::exportMatrix(nMatrix, "n", "matlab", "./ITHACAoutput/Matrices/");
2157 }
2158
2159 if (para->exportTxt)
2160 {
2161 ITHACAstream::exportMatrix(tauMatrix, "tau", "eigen",
2162 "./ITHACAoutput/Matrices/");
2163 ITHACAstream::exportMatrix(nMatrix, "n", "eigen", "./ITHACAoutput/Matrices/");
2164 }
2165}
2166
2168{
2169 tauMatrix.resize(nModes, 3);
2170 nMatrix.resize(nModes, 3);
2171 tauMatrix = tauMatrix * 0;
2172 nMatrix = nMatrix * 0;
2173 Time& runTime = _runTime();
2174 instantList Times = runTime.times();
2175 fvMesh& mesh = _mesh();
2176 volScalarField& p = _p();
2177 volVectorField& U = _U();
2178 //Read FORCESdict
2179 IOdictionary FORCESdict
2180 (
2181 IOobject
2182 (
2183 "FORCESdict",
2184 runTime.system(),
2185 mesh,
2186 IOobject::MUST_READ,
2187 IOobject::NO_WRITE
2188 )
2189 );
2190 IOdictionary transportProperties
2191 (
2192 IOobject
2193 (
2194 "transportProperties",
2195 runTime.constant(),
2196 mesh,
2197 IOobject::MUST_READ,
2198 IOobject::NO_WRITE
2199 )
2200 );
2201 word pName(FORCESdict.lookup("pName"));
2202 word UName(FORCESdict.lookup("UName"));
2203 functionObjects::ITHACAforces f("Forces", mesh, FORCESdict);
2204
2205 for (label i = 0; i < nModes; i++)
2206 {
2207 U = Umodes[i];
2208 p = Pmodes[0];
2209 mesh.readUpdate();
2210 f.write();
2211 f.calcForcesMoment();
2212
2213 for (label j = 0; j < 3; j++)
2214 {
2215 tauMatrix(i, j) = f.forceTau()[j];
2216 }
2217 }
2218
2219 for (label i = 0; i < nModes; i++)
2220 {
2221 U = Umodes[0];
2222 p = Pmodes[i];
2223 mesh.readUpdate();
2224 f.write();
2225 f.calcForcesMoment();
2226
2227 for (label j = 0; j < 3; j++)
2228 {
2229 nMatrix(i, j) = f.forcePressure()[j];
2230 }
2231 }
2232
2233 if (para->exportPython)
2234 {
2235 ITHACAstream::exportMatrix(tauMatrix, "tau", "python",
2236 "./ITHACAoutput/Matrices/");
2237 ITHACAstream::exportMatrix(nMatrix, "n", "python", "./ITHACAoutput/Matrices/");
2238 }
2239
2240 if (para->exportMatlab)
2241 {
2242 ITHACAstream::exportMatrix(tauMatrix, "tau", "matlab",
2243 "./ITHACAoutput/Matrices/");
2244 ITHACAstream::exportMatrix(nMatrix, "n", "matlab", "./ITHACAoutput/Matrices/");
2245 }
2246
2247 if (para->exportTxt)
2248 {
2249 ITHACAstream::exportMatrix(tauMatrix, "tau", "eigen",
2250 "./ITHACAoutput/Matrices/");
2251 ITHACAstream::exportMatrix(nMatrix, "n", "eigen", "./ITHACAoutput/Matrices/");
2252 }
2253}
2254
2255void steadyNS::reconstructLiftAndDrag(const Eigen::MatrixXd& velCoeffs,
2256 const Eigen::MatrixXd& pressureCoeffs, fileName folder)
2257{
2258 M_Assert(velCoeffs.cols() == tauMatrix.rows(),
2259 "The number of velocity modes in the coefficients matrix is not equal to the number of modes in the viscous forces matrix.");
2260 M_Assert(pressureCoeffs.cols() == nMatrix.rows(),
2261 "The number of pressure modes in the coefficients matrix is not equal to the number of modes in the pressure forces matrix.");
2262 mkDir(folder);
2263 system("ln -s ../../constant " + folder + "/constant");
2264 system("ln -s ../../0 " + folder + "/0");
2265 system("ln -s ../../system " + folder + "/system");
2266 //Read FORCESdict
2267 IOdictionary FORCESdict
2268 (
2269 IOobject
2270 (
2271 "FORCESdict",
2272 "./system",
2273 Umodes[0].mesh(),
2274 IOobject::MUST_READ,
2275 IOobject::NO_WRITE
2276 )
2277 );
2278 Eigen::MatrixXd fTau;
2279 Eigen::MatrixXd fN;
2280 fTau.setZero(velCoeffs.rows(), 3);
2281 fN.setZero(pressureCoeffs.rows(), 3);
2282 fTau = velCoeffs * tauMatrix;
2283 fN = pressureCoeffs * nMatrix;
2284
2285 // Export the matrices
2286 if (para->exportPython)
2287 {
2288 ITHACAstream::exportMatrix(fTau, "fTau", "python", folder);
2289 ITHACAstream::exportMatrix(fN, "fN", "python", folder);
2290 }
2291
2292 if (para->exportMatlab)
2293 {
2294 ITHACAstream::exportMatrix(fTau, "fTau", "matlab", folder);
2295 ITHACAstream::exportMatrix(fN, "fN", "matlab", folder);
2296 }
2297
2298 if (para->exportTxt)
2299 {
2300 ITHACAstream::exportMatrix(fTau, "fTau", "eigen", folder);
2301 ITHACAstream::exportMatrix(fN, "fN", "eigen", folder);
2302 }
2303}
2304
2306{
2307 //_runTime().objectRegistry::clear();
2308 //_mesh().objectRegistry::clear();
2309 // _mesh.clear();
2310 // _runTime.clear();
2311 //_simple.clear();
2312 _p.clear();
2313 _U.clear();
2314 _phi.clear();
2315 turbulence.clear();
2316 _fvOptions.clear();
2317 argList& args = _args();
2318 Time& runTime = _runTime();
2319 runTime.setTime(0, 1);
2320 Foam::fvMesh& mesh = _mesh();
2321 // _simple = autoPtr<simpleControl>
2322 // (
2323 // new simpleControl
2324 // (
2325 // mesh
2326 // )
2327 // );
2328 simpleControl& simple = _simple();
2329 Info << "ReReading field p\n" << endl;
2330 _p = autoPtr<volScalarField>
2331 (
2332 new volScalarField
2333 (
2334 IOobject
2335 (
2336 "p",
2337 runTime.timeName(),
2338 mesh,
2339 IOobject::MUST_READ,
2340 IOobject::AUTO_WRITE
2341 ),
2342 mesh
2343 )
2344 );
2345 volScalarField& p = _p();
2346 Info << "ReReading field U\n" << endl;
2347 _U = autoPtr<volVectorField>
2348 (
2349 new volVectorField
2350 (
2351 IOobject
2352 (
2353 "U",
2354 runTime.timeName(),
2355 mesh,
2356 IOobject::MUST_READ,
2357 IOobject::AUTO_WRITE
2358 ),
2359 mesh
2360 )
2361 );
2362 volVectorField& U = _U();
2363 Info << "ReReading/calculating face flux field phi\n" << endl;
2364 _phi = autoPtr<surfaceScalarField>
2365 (
2366 new surfaceScalarField
2367 (
2368 IOobject
2369 (
2370 "phi",
2371 runTime.timeName(),
2372 mesh,
2373 IOobject::READ_IF_PRESENT,
2374 IOobject::AUTO_WRITE
2375 ),
2376 linearInterpolate(U) & mesh.Sf()
2377 )
2378 );
2379 surfaceScalarField& phi = _phi();
2380 pRefCell = 0;
2381 pRefValue = 0.0;
2382 setRefCell(p, simple.dict(), pRefCell, pRefValue);
2383 _laminarTransport = autoPtr<singlePhaseTransportModel>
2384 (
2385 new singlePhaseTransportModel( U, phi )
2386 );
2387 singlePhaseTransportModel& laminarTransport = _laminarTransport();
2388 turbulence = autoPtr<incompressible::turbulenceModel>
2389 (
2390 incompressible::turbulenceModel::New(U, phi, laminarTransport)
2391 );
2392 _MRF = autoPtr<IOMRFZoneList>
2393 (
2394 new IOMRFZoneList(mesh)
2395 );
2396 _fvOptions = autoPtr<fv::options>(new fv::options(mesh));
2397 turbulence->validate();
2398}
static Eigen::MatrixXd field2Eigen(GeometricField< Type, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
static void fvMatrix2Eigen(fvMatrix< type_foam_matrix > foam_matrix, type_A &A, type_B &b)
Convert a FvMatrix OpenFOAM matrix (Linear System) into a Eigen Matrix A and a source vector b.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
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.
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.
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.
void change_viscosity(double mu)
Function to change the viscosity.
Definition steadyNS.C:2066
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
void restart()
set U and P back to the values into the 0 folder
Definition steadyNS.C:2305
List< Eigen::MatrixXd > pressure_gradient_term_linsys_div(label NPmodes)
Laplacian of pressure Linear System - Divergence term.
Definition steadyNS.C:1748
label NPmodes
Number of pressure modes used for the projection.
Definition steadyNS.H:146
List< Eigen::MatrixXd > pressure_gradient_term_linsys_conv(label NPmodes)
Laplacian of pressure Linear System - Convection term.
Definition steadyNS.C:1780
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:265
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Poisson Equation for pressure.
Definition steadyNS.C:354
Eigen::MatrixXd nMatrix
Pressure forces.
Definition steadyNS.H:224
Eigen::MatrixXd BC1_matrix
PPE BC1.
Definition steadyNS.H:191
Eigen::MatrixXd diffusive_term(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Term.
Definition steadyNS.C:883
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:302
Eigen::MatrixXd BC4_matrix
PPE BC4.
Definition steadyNS.H:203
Eigen::MatrixXd tauMatrix
Viscous forces.
Definition steadyNS.H:221
Eigen::MatrixXd BC3_matrix
PPE BC3.
Definition steadyNS.H:200
Eigen::Tensor< double, 3 > C_tensor
Diffusion term.
Definition steadyNS.H:170
Eigen::Tensor< double, 3 > convective_term_tens_cache(label NUmodes, label NPmodes, label NSUPmodes)
Export convective term as a tensor using the cached procedure.
Definition steadyNS.C:1038
void forcesMatrices(label NUmodes, label NPmodes, label NSUPmodes)
Compute lift and drag matrices.
Definition steadyNS.C:2079
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:293
surfaceScalarModes L_PHImodes
List of pointers containing the total number of flux modes.
Definition steadyNS.H:122
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:1515
List< Eigen::MatrixXd > boundary_vector_convection_consistent(label NUmodes, label NSUPmodes)
Boundary vector convection term - Consistent Flux Method.
Definition steadyNS.C:2008
Eigen::MatrixXd W_matrix
Mass Matrix New Time Step - Consistent Flux Method.
Definition steadyNS.H:206
Eigen::MatrixXd mass_matrix_oldtime_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Mass Matrix old time step (consistent flux method).
Definition steadyNS.C:1850
List< Eigen::MatrixXd > pressure_gradient_term_linsys_diff(label NPmodes)
Laplacian of pressure Linear System - Diffusion term.
Definition steadyNS.C:1815
autoPtr< fv::options > _fvOptions
fvOptions
Definition steadyNS.H:296
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:1477
void solvesupremizer(word type="snapshots")
solve the supremizer either with the use of the pressure snaphots or the pressure modes
Definition steadyNS.C:137
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:299
Eigen::MatrixXd I_matrix
Mass Matrix Old Time Step - Consistent Flux Method.
Definition steadyNS.H:209
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition steadyNS.H:101
Eigen::MatrixXd diffusive_term_sym(label NUmodes, label NPmodes, label NSUPmodes)
Symetric diffusive Term.
Definition steadyNS.C:909
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition steadyNS.H:89
List< Eigen::MatrixXd > boundary_vector_diffusion_consistent(label NUmodes, label NSUPmodes)
Boundary vector diffusion term (consistent flux method).
Definition steadyNS.C:1974
List< Eigen::MatrixXd > LinSysDiff
Projection Peqn onto Pressure modes - Diffusion term.
Definition steadyNS.H:257
autoPtr< dimensionedScalar > dt_dummy
Dummy time step including unit.
Definition steadyNS.H:284
Eigen::MatrixXd diffusive_term_flux_method(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Flux Method.
Definition steadyNS.C:1584
Eigen::MatrixXd diffusive_term_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Diffusion Term (consistent flux method).
Definition steadyNS.C:1879
Eigen::MatrixXd pressure_BC1(label NPmodes, label NUmodes)
Term N° 1 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1319
steadyNS()
Null constructor.
Definition steadyNS.C:40
Eigen::MatrixXd divergence_term(label NUmodes, label NPmodes, label NSUPmodes)
Divergence Term (supremizer approach).
Definition steadyNS.C:1112
List< Eigen::MatrixXd > G_matrix
Divergence of momentum PPE.
Definition steadyNS.H:185
scalar tolerance
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and p...
Definition steadyNS.H:125
Eigen::MatrixXd mass_term(label NUmodes, label NPmodes, label NSUPmodes)
Mass Term.
Definition steadyNS.C:1085
Eigen::Tensor< double, 3 > divMomentum_cache(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach) using the cached procedure.
Definition steadyNS.C:1207
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
autoPtr< dimensionedScalar > nu_dummy
Dummy viscocity including unit.
Definition steadyNS.H:287
List< Eigen::MatrixXd > boundary_vector_diffusion(label NUmodes, label NPmodes, label NSUPmodes)
Boundary vector diffusion term.
Definition steadyNS.C:1612
List< Eigen::MatrixXd > LinSysDiv
Projection Peqn onto Pressure modes - Divergence term.
Definition steadyNS.H:254
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
Eigen::Tensor< double, 3 > convective_term_flux_tens(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term.
Definition steadyNS.C:1704
Eigen::Tensor< double, 3 > convective_term_consistent_tens(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term (consistent flux method).
Definition steadyNS.C:1940
scalar pRefValue
Reference pressure value.
Definition steadyNS.H:320
label pRefCell
Reference pressure cell.
Definition steadyNS.H:317
word fluxMethod
Flux Method.
Definition steadyNS.H:329
volVectorModes L_U_SUPmodes
List of pointers containing the total number of lift, supremizer and velocity modes.
Definition steadyNS.H:119
Eigen::MatrixXd pressure_BC3(label NPmodes, label NUmodes)
Term N° 3 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1439
Eigen::MatrixXd L_matrix
Convective background / Large scale advection term.
Definition steadyNS.H:176
List< Eigen::MatrixXd > boundary_vector_convection(label NUmodes, label NPmodes, label NSUPmodes)
Boundary vector convection term.
Definition steadyNS.C:1658
Eigen::MatrixXd B_matrix
Diffusion term.
Definition steadyNS.H:160
List< Eigen::MatrixXd > RD_matrix
Boundary term for diffusion term.
Definition steadyNS.H:236
Eigen::MatrixXd D_matrix
Laplacian term PPE.
Definition steadyNS.H:182
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:314
List< Eigen::MatrixXd > SD_matrix
Boundary term for diffusion term - Consistent Flux Method.
Definition steadyNS.H:242
label NSUPmodes
Number of supremizer modes used for the projection.
Definition steadyNS.H:149
Eigen::Tensor< double, 3 > gTensor
Divergence of momentum PPE.
Definition steadyNS.H:188
Eigen::MatrixXd KF_matrix
Pressure Gradient Term - Consistent Flux Method.
Definition steadyNS.H:215
Eigen::MatrixXd K_matrix
Gradient of pressure matrix.
Definition steadyNS.H:166
List< Eigen::MatrixXd > C_matrix
Non linear term.
Definition steadyNS.H:169
Eigen::Tensor< double, 3 > convective_term_tens(label NUmodes, label NPmodes, label NSUPmodes)
Export convective term as a tensor.
Definition steadyNS.C:997
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition steadyNS.H:308
List< Eigen::MatrixXd > convective_term(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term.
Definition steadyNS.C:962
Eigen::MatrixXd pressure_gradient_term_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Pressure Gradient Term (consistent flux method).
Definition steadyNS.C:1911
List< Eigen::MatrixXd > pressure_BC2(label NPmodes, label NUmodes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1356
Eigen::MatrixXd divergent_convective_background(label NPmodes, label NUmodes, volVectorField vls)
Divergent of Large Scale / Background Advection.
Definition steadyNS.C:1275
List< Eigen::MatrixXd > LinSysConv
Projection Peqn onto Pressure modes - Convection term.
Definition steadyNS.H:260
PtrList< volVectorField > supfield
List of pointers used to form the supremizer snapshots matrix.
Definition steadyNS.H:92
Eigen::Tensor< double, 3 > divMomentum(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach).
Definition steadyNS.C:1176
List< Eigen::MatrixXd > SC_matrix
Boundary term for convection term - Consistent Flux Method.
Definition steadyNS.H:245
Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPmodes, label NSUPmodes)
Gradient of pressure.
Definition steadyNS.C:935
Eigen::MatrixXd DF_matrix
Diffusion Term - Consistent Flux Method.
Definition steadyNS.H:212
List< Eigen::MatrixXd > bcVelMat
Boundary term for penalty method - matrix.
Definition steadyNS.H:230
Eigen::MatrixXd P_matrix
Div of velocity.
Definition steadyNS.H:173
List< Eigen::MatrixXd > RC_matrix
Boundary vector for convection term.
Definition steadyNS.H:239
Eigen::MatrixXd BP_matrix
Diffusion term for flux method PPE.
Definition steadyNS.H:233
Eigen::Tensor< double, 3 > Cf_tensor
Convection term for flux method.
Definition steadyNS.H:248
Eigen::MatrixXd L_D_matrix
Divergent convective background / Large scale advection term.
Definition steadyNS.H:179
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition steadyNS.H:163
List< Eigen::MatrixXd > BC2_matrix
PPE BC2.
Definition steadyNS.H:194
autoPtr< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
Definition steadyNS.H:281
Eigen::MatrixXd convective_background(label NUmodes, volVectorField vls)
Large Scale / Background Advection.
Definition steadyNS.C:1256
void discretizeThenProject(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Discretize-then-project approach.
Definition steadyNS.C:709
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes)
Project using a supremizer approach.
Definition steadyNS.C:556
Eigen::MatrixXd laplacian_pressure(label NPmodes)
Laplacian of pressure term (PPE approach).
Definition steadyNS.C:1295
surfaceScalarModes Phimodes
List of pointers used to form the flux modes.
Definition steadyNS.H:107
Eigen::MatrixXd mass_matrix_newtime_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Mass Matrix new time step (consistent flux method).
Definition steadyNS.C:2042
List< Eigen::MatrixXd > bcVelocityMat(label NUmodes, label NSUPmodes)
Boundary integral modes on boundary used by the penaly method.
Definition steadyNS.C:1547
Eigen::Tensor< double, 3 > bc2Tensor
PPE BC2.
Definition steadyNS.H:197
PtrList< volScalarField > liftfieldP
List of pointer used to form the list of lifting functions for the pressure.
Definition steadyNS.H:113
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:272
Eigen::Tensor< double, 3 > pressureBC2(label NPmodes, label NUmodes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1396
Eigen::Tensor< double, 3 > Ci_tensor
Convection term - Consistent Flux Method.
Definition steadyNS.H:251
void liftSolve()
Perform a lift solve.
Definition steadyNS.C:253
List< Eigen::MatrixXd > div_momentum(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach).
Definition steadyNS.C:1140
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
word bcMethod
Boundary Method.
Definition steadyNS.H:326
void reconstructLiftAndDrag(const Eigen::MatrixXd &velCoeffs, const Eigen::MatrixXd &pressureCoeffs, fileName folder)
Method to reconstruct the forces using velocity and pressure coefficients.
Definition steadyNS.C:2255
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:269
void ReadDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Read a dense matrix from a binary format file.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
void SaveDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Save a dense matrix to a binary format file.
void ReadDenseTensor(TensorType &Tensor, word folder, word MatrixName)
Read a dense tensor from file.
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
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.
void assignZeroDirichlet(GeometricField< Type, fvPatchField, volMesh > &field)
Assign zero internal field.
bool check_folder(word folder)
Checks if a folder exists.
bool check_file(std::string fileName)
Function that returns true if a file exists.
bool check_sup()
Check if the supremizer folder exists.
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
void changeBCtype(GeometricField< Type, fvPatchField, volMesh > &field, word BCtype, label BC_ind)
Change the boundary condition type for a GeometricField.
Header file of the steadyNS class.