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