Loading...
Searching...
No Matches
msrProblem.C
Go to the documentation of this file.
1#include "msrProblem.H"
2
3
4// Constructor
6
7msrProblem::msrProblem(int argc, char* argv[])
8{
9#include "setRootCase.H"
10#include "createTime.H"
11#include "createMesh.H"
12 _simple = autoPtr<simpleControl>
13 (
14 new simpleControl
15 (
16 mesh
17 )
18 );
19 simpleControl& simple = _simple();
20#include "createFields.H"
23#include "createConstants.H"
24#include "createFvOptions.H"
25 turbulence->validate();
26 ITHACAdict = new IOdictionary
27 (
28 IOobject
29 (
30 "ITHACAdict",
31 runTime.system(),
32 mesh,
33 IOobject::MUST_READ,
34 IOobject::NO_WRITE
35 )
36 );
38 tolerance = ITHACAdict->lookupOrDefault<scalar>("tolerance", 1e-5);
39 maxIter = ITHACAdict->lookupOrDefault<scalar>("maxIter", 1000);
42}
43
44void msrProblem::truthSolve(List<scalar> mu_now)
45{
46 Time& runTime = _runTime();
47 fvMesh& mesh = _mesh();
48 volScalarField& p = _p();
49 volVectorField& U = _U();
50 surfaceScalarField& phi = _phi();
51 fv::options& fvOptions = _fvOptions();
52 simpleControl& simple = _simple();
53 IOMRFZoneList& MRF = _MRF();
54 singlePhaseTransportModel& laminarTransport = _laminarTransport();
55 dimensionedScalar& IV1 = _IV1();
56 dimensionedScalar& Keff = _Keff();
57 volScalarField& flux = _flux();
58 volScalarField flux_old = _flux();
59 dimensionedScalar& betaTot = _betaTot();
60 volScalarField& SP = _SP();
61 dimensionedScalar& SP1_0 = _SP1_0();
62 dimensionedScalar& alfa_SP1 = _alfa_SP1();
63 volScalarField& D = _D();
64 dimensionedScalar& D1_0 = _D1_0();
65 dimensionedScalar& alfa_D1 = _alfa_D1();
66 volScalarField& NSF = _NSF();
67 dimensionedScalar& NSF1_0 = _NSF1_0();
68 dimensionedScalar& alfa_NSF1 = _alfa_NSF1();
69 volScalarField& A = _A();
70 dimensionedScalar& A1_0 = _A1_0();
71 dimensionedScalar& alfa_A1 = _alfa_A1();
72 volScalarField& prec1 = _prec1();
73 dimensionedScalar& Sc = _Sc();
74 dimensionedScalar& Sct = _Sct();
75 dimensionedScalar& lam1 = _lam1();
76 dimensionedScalar& beta1 = _beta1();
77 volScalarField& prec2 = _prec2();
78 dimensionedScalar& lam2 = _lam2();
79 dimensionedScalar& beta2 = _beta2();
80 volScalarField& prec3 = _prec3();
81 dimensionedScalar& lam3 = _lam3();
82 dimensionedScalar& beta3 = _beta3();
83 volScalarField& prec4 = _prec4();
84 dimensionedScalar& lam4 = _lam4();
85 dimensionedScalar& beta4 = _beta4();
86 volScalarField& prec5 = _prec5();
87 dimensionedScalar& lam5 = _lam5();
88 dimensionedScalar& beta5 = _beta5();
89 volScalarField& prec6 = _prec6();
90 dimensionedScalar& lam6 = _lam6();
91 dimensionedScalar& beta6 = _beta6();
92 volScalarField& prec7 = _prec7();
93 dimensionedScalar& lam7 = _lam7();
94 dimensionedScalar& beta7 = _beta7();
95 volScalarField& prec8 = _prec8();
96 dimensionedScalar& lam8 = _lam8();
97 dimensionedScalar& beta8 = _beta8();
98 volScalarField& T = _T();
99 dimensionedScalar& Pr = _Pr();
100 dimensionedScalar& Prt = _Prt();
101 volScalarField& dec1 = _dec1();
102 dimensionedScalar& decLam1 = _decLam1();
103 dimensionedScalar& decBeta1 = _decBeta1();
104 volScalarField& dec2 = _dec2();
105 dimensionedScalar& decLam2 = _decLam2();
106 dimensionedScalar& decBeta2 = _decBeta2();
107 volScalarField& dec3 = _dec3();
108 dimensionedScalar& decLam3 = _decLam3();
109 dimensionedScalar& decBeta3 = _decBeta3();
110 dimensionedScalar& decbetaTot = _decbetaTot();
111 dimensionedScalar& rhoRef = _rhoRef();
112 dimensionedScalar& CpRef = _CpRef();
113 volScalarField v = _v();
114 volScalarField TXS = _TXS();
115 dimensionedScalar& nu = _nu();
116 dimensionedScalar& betaTE = _betaTE();
117 dimensionedScalar& Tref = _Tref();
118 dimensionedScalar& TrefXS = _TrefXS();
119 volScalarField& logT = _logT();
120 volScalarField& alphat = _alphat();
121 volScalarField& difft = _difft();
122 volScalarField powerDens = ((1 - decbetaTot) * flux * SP +
123 (decLam1 * dec1 + decLam2 * dec2 + decLam3 * dec3)).ref();
124 powerDens.rename("powerDens");
125#include "NLmsrProblem.H"
126 counter++;
127 writeMu(mu_now);
128 // --- Fill in the mu_samples with parameters (mu) to be used for the PODI sample points
129 mu_samples.conservativeResize(mu_samples.rows() + 1, mu_now.size());
130
131 for (label i = 0; i < mu_now.size(); i++)
132 {
133 mu_samples(mu_samples.rows() - 1, i) = mu_now[i];
134 }
135
136 // Resize to Unitary if not initialized by user (i.e. non-parametric problem)
137 if (mu.cols() == 0)
138 {
139 mu.resize(1, 1);
140 }
141
142 if (mu_samples.rows() == mu.cols())
143 {
144 ITHACAstream::exportMatrix(mu_samples, "mu_samples", "eigen",
145 "./ITHACAoutput/Offline");
146 }
147}
148
149
150// Get the modes function SVD
152{
153 label NU = para->ITHACAdict->lookupOrDefault<label>("NUout", 10);
154 label NP = para->ITHACAdict->lookupOrDefault<label>("NPout", 10);
155 label NF = para->ITHACAdict->lookupOrDefault<label>("NFluxout", 10);
156 label NPrec1 = para->ITHACAdict->lookupOrDefault<label>("NPrecout1", 10);
157 label NPrec2 = para->ITHACAdict->lookupOrDefault<label>("NPrecout2", 10);
158 label NPrec3 = para->ITHACAdict->lookupOrDefault<label>("NPrecout3", 10);
159 label NPrec4 = para->ITHACAdict->lookupOrDefault<label>("NPrecout4", 10);
160 label NPrec5 = para->ITHACAdict->lookupOrDefault<label>("NPrecout5", 10);
161 label NPrec6 = para->ITHACAdict->lookupOrDefault<label>("NPrecout6", 10);
162 label NPrec7 = para->ITHACAdict->lookupOrDefault<label>("NPrecout7", 10);
163 label NPrec8 = para->ITHACAdict->lookupOrDefault<label>("NPrecout8", 10);
164 label NT = para->ITHACAdict->lookupOrDefault<label>("NTout", 10);
165 label NDec1 = para->ITHACAdict->lookupOrDefault<label>("NDecout1", 10);
166 label NDec2 = para->ITHACAdict->lookupOrDefault<label>("NDecout2", 10);
167 label NDec3 = para->ITHACAdict->lookupOrDefault<label>("NDecout3", 10);
168 label NC = para->ITHACAdict->lookupOrDefault<label>("NCout", 10);
169
170 if (homboolU == true)
171 {
172 ITHACAPOD::getModesSVD(Uomfield, Umodes, _U().name(), podex, 0, 0, NU);
173 }
174 else
175 {
176 ITHACAPOD::getModesSVD(Ufield, Umodes, _U().name(), podex, 0, 0, NU);
177 }
178
179 ITHACAPOD::getModesSVD(Pfield, Pmodes, _p().name(), podex, 0, 0, NP);
182 NPrec1);
184 NPrec2);
186 NPrec3);
188 NPrec4);
190 NPrec5);
192 NPrec6);
194 NPrec7);
196 NPrec8);
197
198 if (homboolT == true)
199 {
200 ITHACAPOD::getModesSVD(Tomfield, Tmodes, _T().name(), podex, 0, 0, NT);
201 }
202 else
203 {
204 ITHACAPOD::getModesSVD(Tfield, Tmodes, _T().name(), podex, 0, 0, NT);
205 }
206
208 NDec1);
210 NDec2);
212 NDec3);
213 ITHACAPOD::getModesSVD(vFields, vmodes, _v().name(), podex, 0, 0, NC);
214 ITHACAPOD::getModesSVD(DFields, Dmodes, _D().name(), podex, 0, 0, NC);
215 ITHACAPOD::getModesSVD(NSFFields, NSFmodes, _NSF().name(), podex, 0, 0, NC);
216 ITHACAPOD::getModesSVD(AFields, Amodes, _A().name(), podex, 0, 0, NC);
217 ITHACAPOD::getModesSVD(SPFields, SPmodes, _SP().name(), podex, 0, 0, NC);
218 ITHACAPOD::getModesSVD(TXSFields, TXSmodes, _TXS().name(), podex, 0, 0, NC);
219 Info << "End\n" << endl;
220}
221
222// Get the modes function Compliance matrix method
224{
225 label NU = para->ITHACAdict->lookupOrDefault<label>("NUout", 10);
226 label NP = para->ITHACAdict->lookupOrDefault<label>("NPout", 10);
227 label NF = para->ITHACAdict->lookupOrDefault<label>("NFluxout", 10);
228 label NPrec1 = para->ITHACAdict->lookupOrDefault<label>("NPrecout1", 10);
229 label NPrec2 = para->ITHACAdict->lookupOrDefault<label>("NPrecout2", 10);
230 label NPrec3 = para->ITHACAdict->lookupOrDefault<label>("NPrecout3", 10);
231 label NPrec4 = para->ITHACAdict->lookupOrDefault<label>("NPrecout4", 10);
232 label NPrec5 = para->ITHACAdict->lookupOrDefault<label>("NPrecout5", 10);
233 label NPrec6 = para->ITHACAdict->lookupOrDefault<label>("NPrecout6", 10);
234 label NPrec7 = para->ITHACAdict->lookupOrDefault<label>("NPrecout7", 10);
235 label NPrec8 = para->ITHACAdict->lookupOrDefault<label>("NPrecout8", 10);
236 label NT = para->ITHACAdict->lookupOrDefault<label>("NTout", 10);
237 label NDec1 = para->ITHACAdict->lookupOrDefault<label>("NDecout1", 10);
238 label NDec2 = para->ITHACAdict->lookupOrDefault<label>("NDecout2", 10);
239 label NDec3 = para->ITHACAdict->lookupOrDefault<label>("NDecout3", 10);
240 label NC = para->ITHACAdict->lookupOrDefault<label>("NCout", 10);
241
242 if (homboolU == true)
243 {
244 ITHACAPOD::getModes(Uomfield, Umodes, _U().name(), podex, 0, 0, NU);
245 }
246 else
247 {
248 ITHACAPOD::getModes(Ufield, Umodes, _U().name(), podex, 0, 0, NU);
249 }
250
251 ITHACAPOD::getModes(Pfield, Pmodes, _U().name(), podex, 0, 0, NP);
252 ITHACAPOD::getModes(Fluxfield, Fluxmodes, _U().name(), podex, 0, 0, NF);
254 NPrec1);
256 NPrec2);
258 NPrec3);
260 NPrec4);
262 NPrec5);
264 NPrec6);
266 NPrec7);
268 NPrec8);
269
270 if (homboolT == true)
271 {
272 ITHACAPOD::getModes(Tomfield, Tmodes, _T().name(), podex, 0, 0, NT);
273 }
274 else
275 {
276 ITHACAPOD::getModes(Tfield, Tmodes, _T().name(), podex, 0, 0, NT);
277 }
278
279 ITHACAPOD::getModes(Dec1field, Dec1modes, _dec1().name(), podex, 0, 0, NDec1);
280 ITHACAPOD::getModes(Dec2field, Dec2modes, _dec2().name(), podex, 0, 0, NDec2);
281 ITHACAPOD::getModes(Dec3field, Dec3modes, _dec3().name(), podex, 0, 0, NDec3);
282 ITHACAPOD::getModes(vFields, vmodes, _v().name(), podex, 0, 0, NC);
283 ITHACAPOD::getModes(DFields, Dmodes, _D().name(), podex, 0, 0, NC);
284 ITHACAPOD::getModes(NSFFields, NSFmodes, _NSF().name(), podex, 0, 0, NC);
285 ITHACAPOD::getModes(AFields, Amodes, _A().name(), podex, 0, 0, NC);
286 ITHACAPOD::getModes(SPFields, SPmodes, _SP().name(), podex, 0, 0, NC);
287 ITHACAPOD::getModes(TXSFields, TXSmodes, _TXS().name(), podex, 0, 0, NC);
288 Info << "End\n" << endl;
289}
290
292{
293 for (label k = 0; k < inletIndex.rows(); k++)
294 {
295 Time& runTime = _runTime();
296 surfaceScalarField& phi = _phi();
297 fvMesh& mesh = _mesh();
298 volScalarField p = _p();
299 volVectorField U = _U();
300 IOMRFZoneList& MRF = _MRF();
301 label BCind = inletIndex(k, 0);
302 volVectorField Ulift("Ulift" + name(k), U);
303 instantList Times = runTime.times();
304 runTime.setTime(Times[1], 1);
305 pisoControl potentialFlow(mesh, "potentialFlow");
306 Info << "Solving a lifting Problem" << endl;
307 Vector<double> v1(0, 0, 0);
308 v1[inletIndex(k, 1)] = 1;
309 Vector<double> v0(0, 0, 0);
310
311 for (label j = 0; j < U.boundaryField().size(); j++)
312 {
313 if (j == BCind)
314 {
315 assignBC(Ulift, j, v1);
316 }
317 else if (U.boundaryField()[BCind].type() == "fixedValue")
318 {
319 assignBC(Ulift, j, v0);
320 }
321 else
322 {
323 }
324
325 assignIF(Ulift, v0);
326 phi = linearInterpolate(Ulift) & mesh.Sf();
327 }
328
329 Info << "Constructing velocity potential field Phi\n" << endl;
330 volScalarField Phi
331 (
332 IOobject
333 (
334 "Phi",
335 runTime.timeName(),
336 mesh,
337 IOobject::READ_IF_PRESENT,
338 IOobject::NO_WRITE
339 ),
340 mesh,
341 dimensionedScalar("Phi", dimLength * dimVelocity, 0),
342 p.boundaryField().types()
343 );
344 label PhiRefCell = 0;
345 scalar PhiRefValue = 0;
347 (
348 Phi,
349 potentialFlow.dict(),
350 PhiRefCell,
351 PhiRefValue
352 );
353 mesh.setFluxRequired(Phi.name());
354 runTime.functionObjects().start();
355 MRF.makeRelative(phi);
356 adjustPhi(phi, Ulift, p);
357
358 while (potentialFlow.correctNonOrthogonal())
359 {
360 fvScalarMatrix PhiEqn
361 (
362 fvm::laplacian(dimensionedScalar("1", dimless, 1), Phi)
363 ==
364 fvc::div(phi)
365 );
366 PhiEqn.setReference(PhiRefCell, PhiRefValue);
367 PhiEqn.solve();
368
369 if (potentialFlow.finalNonOrthogonalIter())
370 {
371 phi -= PhiEqn.flux();
372 }
373 }
374
375 MRF.makeAbsolute(phi);
376 Info << "Continuity error = "
377 << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
378 << endl;
379 Ulift = fvc::reconstruct(phi);
380 Ulift.correctBoundaryConditions();
381 Info << "Interpolated velocity error = "
382 << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
383 / sum(mesh.magSf())).value()
384 << endl;
385 Ulift.write();
386 liftfield.append(Ulift.clone());
387 }
388}
389
390
391
392
393// * * * * * * * * * * * * * * Projection Methods * * * * * * * * * * * * * * //
394
395void msrProblem::projectPPE(fileName folder, label NU, label NP, label NF,
396 Eigen::VectorXi NPrec, label NT, Eigen::VectorXi NDec, label NC)
397{
398 if (NPrec.size() != 8 || NDec.size() != 3)
399 {
400 std::cout <<
401 "The model assumes 8 groups of precursors and 3 of decay heat, check NDrec and NDec dimensions..."
402 << std::endl;
403 exit(0);
404 }
405
406 NUmodes = NU;
407 NPmodes = NP;
408 NFluxmodes = NF;
409 NPrecmodes = NPrec;
410 label NPrec1 = NPrecmodes(0);
411 label NPrec2 = NPrecmodes(1);
412 label NPrec3 = NPrecmodes(2);
413 label NPrec4 = NPrecmodes(3);
414 label NPrec5 = NPrecmodes(4);
415 label NPrec6 = NPrecmodes(5);
416 label NPrec7 = NPrecmodes(6);
417 label NPrec8 = NPrecmodes(7);
418 NTmodes = NT;
419 NDecmodes = NDec;
420 label NDec1 = NDecmodes(0);
421 label NDec2 = NDecmodes(1);
422 label NDec3 = NDecmodes(2);
423 NCmodes = NC;
424 Info << "\n Computing fluid-dynamics matrices\n" << endl;
434 Info << "\n End \n" << endl;
435 Info << "\n Computing neutronics matrices\n" << endl;
440 PS1_matrix = prec_source(NFluxmodes, NPrec1, 1);
441 PS2_matrix = prec_source(NFluxmodes, NPrec2, 2);
442 PS3_matrix = prec_source(NFluxmodes, NPrec3, 3);
443 PS4_matrix = prec_source(NFluxmodes, NPrec4, 4);
444 PS5_matrix = prec_source(NFluxmodes, NPrec5, 5);
445 PS6_matrix = prec_source(NFluxmodes, NPrec6, 6);
446 PS7_matrix = prec_source(NFluxmodes, NPrec7, 7);
447 PS8_matrix = prec_source(NFluxmodes, NPrec8, 8);
448 ST1_matrix = stream_term(NUmodes, NPrec1, 1);
449 ST2_matrix = stream_term(NUmodes, NPrec2, 2);
450 ST3_matrix = stream_term(NUmodes, NPrec3, 3);
451 ST4_matrix = stream_term(NUmodes, NPrec4, 4);
452 ST5_matrix = stream_term(NUmodes, NPrec5, 5);
453 ST6_matrix = stream_term(NUmodes, NPrec6, 6);
454 ST7_matrix = stream_term(NUmodes, NPrec7, 7);
455 ST8_matrix = stream_term(NUmodes, NPrec8, 8);
456 MP1_matrix = prec_mass(NPrec1, 1);
457 MP2_matrix = prec_mass(NPrec2, 2);
458 MP3_matrix = prec_mass(NPrec3, 3);
459 MP4_matrix = prec_mass(NPrec4, 4);
460 MP5_matrix = prec_mass(NPrec5, 5);
461 MP6_matrix = prec_mass(NPrec6, 6);
462 MP7_matrix = prec_mass(NPrec7, 7);
463 MP8_matrix = prec_mass(NPrec8, 8);
464 LP1_matrix = laplacian_prec(NPrec1, 1);
465 LP2_matrix = laplacian_prec(NPrec2, 2);
466 LP3_matrix = laplacian_prec(NPrec3, 3);
467 LP4_matrix = laplacian_prec(NPrec4, 4);
468 LP5_matrix = laplacian_prec(NPrec5, 5);
469 LP6_matrix = laplacian_prec(NPrec6, 6);
470 LP7_matrix = laplacian_prec(NPrec7, 7);
471 LP8_matrix = laplacian_prec(NPrec8, 8);
480 Info << "\n End \n" << endl;
481 Info << "\n Computing thermal matrices\n" << endl;
482 SD1_matrix = stream_dec(NUmodes, NDec1, 1);
483 SD2_matrix = stream_dec(NUmodes, NDec2, 2);
484 SD3_matrix = stream_dec(NUmodes, NDec3, 3);
485 MD1_matrix = dec_mass(NDec1, 1);
486 MD2_matrix = dec_mass(NDec2, 2);
487 MD3_matrix = dec_mass(NDec3, 3);
488 LD1_matrix = laplacian_dec(NDec1, 1);
489 LD2_matrix = laplacian_dec(NDec2, 2);
490 LD3_matrix = laplacian_dec(NDec3, 3);
501 Info << "\n End \n" << endl;
503}
504
505// * * * * * * * * * * * * * * Momentum Eq. Methods * * * * * * * * * * * * * //
506
507Eigen::MatrixXd msrProblem::diffusive_term(label NUmodes, label NPmodes)
508{
509 label Bsize = NUmodes + liftfield.size();
510 Eigen::MatrixXd B_matrix;
511 B_matrix.resize(Bsize, Bsize);
512 PtrList<volVectorField> Together(0);
513
514 if (liftfield.size() != 0)
515 {
516 for (label k = 0; k < liftfield.size(); k++)
517 {
518 Together.append(liftfield[k].clone());
519 }
520 }
521
522 if (NUmodes != 0)
523 {
524 for (label k = 0; k < NUmodes; k++)
525 {
526 Together.append(Umodes[k].clone());
527 }
528 }
529
530 // Project everything
531 for (label i = 0; i < Bsize; i++)
532 {
533 for (label j = 0; j < Bsize; j++)
534 {
535 B_matrix(i, j) = fvc::domainIntegrate(Together[i] & fvc::laplacian(
536 dimensionedScalar("1", dimless, 1), Together[j])).value();
537 }
538 }
539
540 // Export the matrix
542 "./ITHACAoutput/Matrices/fluid_dynamics/");
543 return B_matrix;
544}
545
546Eigen::MatrixXd msrProblem::pressure_gradient_term(label NUmodes, label NPmodes)
547{
548 label K1size = NUmodes + liftfield.size();
549 label K2size = NPmodes;
550 Eigen::MatrixXd K_matrix(K1size, K2size);
551 PtrList<volVectorField> Together(0);
552
553 if (liftfield.size() != 0)
554 {
555 for (label k = 0; k < liftfield.size(); k++)
556 {
557 Together.append(liftfield[k].clone());
558 }
559 }
560
561 if (NUmodes != 0)
562 {
563 for (label k = 0; k < NUmodes; k++)
564 {
565 Together.append(Umodes[k].clone());
566 }
567 }
568
569 // Project everything
570 for (label i = 0; i < K1size; i++)
571 {
572 for (label j = 0; j < K2size; j++)
573 {
574 K_matrix(i, j) = fvc::domainIntegrate(Together[i] & fvc::grad(
575 Pmodes[j])).value();
576 }
577 }
578
579 // Export the matrix
581 "./ITHACAoutput/Matrices/fluid_dynamics/");
582 return K_matrix;
583}
584
585List <Eigen::MatrixXd> msrProblem::convective_term(label NUmodes,
586 label NPmodes)
587{
588 label Csize = NUmodes + liftfield.size();
589 List <Eigen::MatrixXd> C_matrix;
590 C_matrix.setSize(Csize);
591
592 for (label j = 0; j < Csize; j++)
593 {
594 C_matrix[j].resize(Csize, Csize);
595 }
596
597 PtrList<volVectorField> Together(0);
598
599 if (liftfield.size() != 0)
600 {
601 for (label k = 0; k < liftfield.size(); k++)
602 {
603 Together.append(liftfield[k].clone());
604 }
605 }
606
607 if (NUmodes != 0)
608 {
609 for (label k = 0; k < NUmodes; k++)
610 {
611 Together.append(Umodes[k].clone());
612 }
613 }
614
615 for (label i = 0; i < Csize; i++)
616 {
617 for (label j = 0; j < Csize; j++)
618 {
619 for (label k = 0; k < Csize; k++)
620 {
621 C_matrix[i](j, k) = fvc::domainIntegrate(Together[i] & fvc::div(
622 linearInterpolate(Together[j]) & Together[j].mesh().Sf(), Together[k])).value();
623 }
624 }
625 }
626
627 // Export the matrix
629 "./ITHACAoutput/Matrices/fluid_dynamics/");
630 return C_matrix;
631}
632
633Eigen::MatrixXd msrProblem::mass_term(label NUmodes, label NPmodes)
634{
635 label Msize = NUmodes + liftfield.size();
636 Eigen::MatrixXd M_matrix(Msize, Msize);
637 PtrList<volVectorField> Together(0);
638
639 if (liftfield.size() != 0)
640 {
641 for (label k = 0; k < liftfield.size(); k++)
642 {
643 Together.append(liftfield[k].clone());
644 }
645 }
646
647 if (NUmodes != 0)
648 {
649 for (label k = 0; k < NUmodes; k++)
650 {
651 Together.append(Umodes[k].clone());
652 }
653 }
654
655 // Project everything
656 for (label i = 0; i < Msize; i++)
657 {
658 for (label j = 0; j < Msize; j++)
659 {
660 M_matrix(i, j) = fvc::domainIntegrate(Together[i] & Together[j]).value();
661 }
662 }
663
664 // Export the matrix
666 "./ITHACAoutput/Matrices/fluid_dynamics/");
667 return M_matrix;
668}
669
670// * * * * * * * * * * * * * * Continuity Eq. Methods * * * * * * * * * * * * * //
671
672Eigen::MatrixXd msrProblem::divergence_term(label NUmodes, label NPmodes)
673{
674 label P1size = NPmodes;
675 label P2size = NUmodes + liftfield.size();
676 Eigen::MatrixXd P_matrix(P1size, P2size);
677 PtrList<volVectorField> Together(0);
678
679 if (liftfield.size() != 0)
680 {
681 for (label k = 0; k < liftfield.size(); k++)
682 {
683 Together.append(liftfield[k].clone());
684 }
685 }
686
687 if (NUmodes != 0)
688 {
689 for (label k = 0; k < NUmodes; k++)
690 {
691 Together.append(Umodes[k].clone());
692 }
693 }
694
695 // Project everything
696 for (label i = 0; i < P1size; i++)
697 {
698 for (label j = 0; j < P2size; j++)
699 {
700 P_matrix(i, j) = fvc::domainIntegrate(Pmodes[i] * fvc::div (
701 Together[j])).value();
702 }
703 }
704
705 //Export the matrix
707 "./ITHACAoutput/Matrices/fluid_dynamics/");
708 return P_matrix;
709}
710
711
712List <Eigen::MatrixXd> msrProblem::div_momentum(label NUmodes, label NPmodes)
713{
714 label G1size = NPmodes;
715 label G2size = NUmodes + liftfield.size();
716 List <Eigen::MatrixXd> G_matrix;
717 G_matrix.setSize(G1size);
718
719 for (label j = 0; j < G1size; j++)
720 {
721 G_matrix[j].resize(G2size, G2size);
722 }
723
724 PtrList<volVectorField> Together(0);
725
726 if (liftfield.size() != 0)
727 {
728 for (label k = 0; k < liftfield.size(); k++)
729 {
730 Together.append(liftfield[k].clone());
731 }
732 }
733
734 if (NUmodes != 0)
735 {
736 for (label k = 0; k < NUmodes; k++)
737 {
738 Together.append(Umodes[k].clone());
739 }
740 }
741
742 for (label i = 0; i < G1size; i++)
743 {
744 for (label j = 0; j < G2size; j++)
745 {
746 for (label k = 0; k < G2size; k++)
747 {
748 G_matrix[i](j, k) = fvc::domainIntegrate(fvc::grad(Pmodes[i]) & (fvc::div(
749 fvc::interpolate(Together[j]) & Together[j].mesh().Sf(), Together[k]))).value();
750 }
751 }
752 }
753
754 // Export the matrix
756 "./ITHACAoutput/Matrices/fluid_dynamics/");
757 return G_matrix;
758}
759
760Eigen::MatrixXd msrProblem::laplacian_pressure(label NPmodes)
761{
762 label Dsize = NPmodes;
763 Eigen::MatrixXd D_matrix(Dsize, Dsize);
764
765 // Project everything
766 for (label i = 0; i < Dsize; i++)
767 {
768 for (label j = 0; j < Dsize; j++)
769 {
770 D_matrix(i, j) = fvc::domainIntegrate(fvc::grad(Pmodes[i])&fvc::grad(
771 Pmodes[j])).value();
772 }
773 }
774
775 //Export the matrix
777 "./ITHACAoutput/Matrices/fluid_dynamics/");
778 return D_matrix;
779}
780
781Eigen::MatrixXd msrProblem::pressure_BC1(label NUmodes, label NPmodes)
782{
783 label P_BC1size = NPmodes;
784 label P_BC2size = NUmodes + liftfield.size();
785 Eigen::MatrixXd BC1_matrix(P_BC1size, P_BC2size);
786 fvMesh& mesh = _mesh();
787 PtrList<volVectorField> Together(0);
788
789 if (liftfield.size() != 0)
790 {
791 for (label k = 0; k < liftfield.size(); k++)
792 {
793 Together.append(liftfield[k].clone());
794 }
795 }
796
797 if (NUmodes != 0)
798 {
799 for (label k = 0; k < NUmodes; k++)
800 {
801 Together.append(Umodes[k].clone());
802 }
803 }
804
805 for (label i = 0; i < P_BC1size; i++)
806 {
807 for (label j = 0; j < P_BC2size; j++)
808 {
809 surfaceScalarField lpl((fvc::interpolate(fvc::laplacian(
810 Together[j]))&mesh.Sf())*fvc::interpolate(Pmodes[i]));
811 double s = 0;
812
813 for (label k = 0; k < lpl.boundaryField().size(); k++)
814 {
815 s += gSum(lpl.boundaryField()[k]);
816 }
817
818 BC1_matrix(i, j) = s;
819 }
820 }
821
822 ITHACAstream::exportMatrix(BC1_matrix, "BC1", "matlab",
823 "./ITHACAoutput/Matrices/fluid_dynamics/");
824 return BC1_matrix;
825}
826
827
828List <Eigen::MatrixXd> msrProblem::pressure_BC2(label NUmodes, label NPmodes)
829{
830 label P2_BC1size = NPmodes;
831 label P2_BC2size = NUmodes + liftfield.size();
832 List <Eigen::MatrixXd> BC2_matrix;
833 fvMesh& mesh = _mesh();
834 BC2_matrix.setSize(P2_BC1size);
835
836 for (label j = 0; j < P2_BC1size; j++)
837 {
838 BC2_matrix[j].resize(P2_BC2size, P2_BC2size);
839 }
840
841 PtrList<volVectorField> Together(0);
842
843 if (liftfield.size() != 0)
844 {
845 for (label k = 0; k < liftfield.size(); k++)
846 {
847 Together.append(liftfield[k].clone());
848 }
849 }
850
851 if (NUmodes != 0)
852 {
853 for (label k = 0; k < NUmodes; k++)
854 {
855 Together.append(Umodes[k].clone());
856 }
857 }
858
859 for (label i = 0; i < P2_BC1size; i++)
860 {
861 for (label j = 0; j < P2_BC2size; j++)
862 {
863 for (label k = 0; k < P2_BC2size; k++)
864 {
865 surfaceScalarField div_m(fvc::interpolate(fvc::div(fvc::interpolate(
866 Together[j]) & mesh.Sf(), Together[k]))&mesh.Sf()*fvc::interpolate(Pmodes[i]));
867 double s = 0;
868
869 for (label k = 0; k < div_m.boundaryField().size(); k++)
870 {
871 s += gSum(div_m.boundaryField()[k]);
872 }
873
874 BC2_matrix[i](j, k) = s;
875 }
876 }
877 }
878
879 // Export the matrix
880 ITHACAstream::exportMatrix(BC2_matrix, "BC2", "matlab",
881 "./ITHACAoutput/Matrices/fluid_dynamics/");
882 return BC2_matrix;
883}
884
885Eigen::MatrixXd msrProblem::pressure_BC3(label NUmodes, label NPmodes)
886{
887 label P3_BC1size = NPmodes;
888 label P3_BC2size = NUmodes + liftfield.size();
889 Eigen::MatrixXd BC3_matrix(P3_BC1size, P3_BC2size);
890 fvMesh& mesh = _mesh();
891 PtrList<volVectorField> Together(0);
892
893 if (liftfield.size() != 0)
894 {
895 for (label k = 0; k < liftfield.size(); k++)
896 {
897 Together.append(liftfield[k].clone());
898 }
899 }
900
901 if (NUmodes != 0)
902 {
903 for (label k = 0; k < NUmodes; k++)
904 {
905 Together.append(Umodes[k].clone());
906 }
907 }
908
909 surfaceVectorField n(mesh.Sf() / mesh.magSf());
910
911 for (label i = 0; i < P3_BC1size; i++)
912 {
913 for (label j = 0; j < P3_BC2size; j++)
914 {
915 surfaceVectorField BC3 = fvc::interpolate(fvc::curl(Together[j])).ref();
916 surfaceVectorField BC4 = (n ^ fvc::interpolate(fvc::grad(Pmodes[i]))).ref();
917 surfaceScalarField BC5 = ((BC3 & BC4) * mesh.magSf()).ref();
918 double s = 0;
919
920 for (label k = 0; k < BC5.boundaryField().size(); k++)
921 {
922 s += gSum(BC5.boundaryField()[k]);
923 }
924
925 BC3_matrix(i, j) = s;
926 }
927 }
928
929 ITHACAstream::exportMatrix(BC3_matrix, "BC3", "matlab",
930 "./ITHACAoutput/Matrices/fluid_dynamics/");
931 return BC3_matrix;
932}
933
934
935// * * * * * * * * * * * * * * Diffusion Eq. Methods * * * * * * * * * * * * * //
936
937List<Eigen::MatrixXd> msrProblem::laplacian_flux(label NFluxmodes,
938 label NCmodes)
939{
940 label LFsize = NFluxmodes;
941 List<Eigen::MatrixXd> LF_matrix;
942 LF_matrix.setSize(LFsize);
943
944 for (label j = 0; j < LFsize; j++)
945 {
946 LF_matrix[j].resize(NCmodes, LFsize);
947 LF_matrix[j].setZero();
948 }
949
950 // Project everything
951 for (label i = 0; i < LFsize; i++)
952 {
953 for (label j = 0; j < NCmodes; j++)
954 {
955 for (label k = 0; k < LFsize; k++)
956 {
957 LF_matrix[i](j, k) = fvc::domainIntegrate(Fluxmodes[i] * fvc::laplacian(
958 Dmodes[j], Fluxmodes[k])).value();
959 }
960 }
961 }
962
963 // Export the matrix
965 "./ITHACAoutput/Matrices/neutronics/");
966 return LF_matrix;
967}
968
969Eigen::MatrixXd msrProblem::mass_flux(label NFluxmodes)
970{
971 label MFsize = NFluxmodes;
972 Eigen::MatrixXd MF_matrix;
973 MF_matrix.resize(MFsize, MFsize);
974
975 // Project everything
976 for (label i = 0; i < MFsize; i++)
977 {
978 for (label j = 0; j < MFsize; j++)
979 {
980 MF_matrix(i, j) = fvc::domainIntegrate(Fluxmodes[i] * Fluxmodes[j]).value();
981 }
982 }
983
984 // Export the matrix
986 "./ITHACAoutput/Matrices/neutronics/");
987 return MF_matrix;
988}
989
990List<Eigen::MatrixXd> msrProblem::prod_flux(label NFluxmodes, label NCmodes)
991{
992 label PFsize = NFluxmodes;
993 List<Eigen::MatrixXd> PF_matrix;
994 PF_matrix.setSize(PFsize);
995
996 for (label j = 0; j < PFsize; j++)
997 {
998 PF_matrix[j].resize(NCmodes, PFsize);
999 PF_matrix[j].setZero();
1000 }
1001
1002 //Project everything
1003 for (label i = 0; i < PFsize; i++)
1004 {
1005 for (label j = 0; j < NCmodes; j++)
1006 {
1007 for (label k = 0; k < PFsize; k++)
1008 {
1009 PF_matrix[i](j, k) = fvc::domainIntegrate(Fluxmodes[i] * NSFmodes[j] *
1010 Fluxmodes[k]).value();
1011 }
1012 }
1013 }
1014
1015 // Export the matrix
1016 ITHACAstream::exportMatrix(PF_matrix, "PF", "matlab",
1017 "./ITHACAoutput/Matrices/neutronics/");
1018 return PF_matrix;
1019}
1020
1021List<Eigen::MatrixXd> msrProblem::abs_flux(label NFluxmodes, label NCmodes)
1022{
1023 label AFsize = NFluxmodes;
1024 List<Eigen::MatrixXd> AF_matrix;
1025 AF_matrix.setSize(AFsize);
1026
1027 for (label j = 0; j < AFsize; j++)
1028 {
1029 AF_matrix[j].resize(NCmodes, AFsize);
1030 AF_matrix[j].setZero();
1031 }
1032
1033 //Project everything
1034 for (label i = 0; i < AFsize; i++)
1035 {
1036 for (label j = 0; j < NCmodes; j++)
1037 {
1038 for (label k = 0; k < AFsize; k++)
1039 {
1040 AF_matrix[i](j, k) = fvc::domainIntegrate(Fluxmodes[i] * Amodes[j] *
1041 Fluxmodes[k]).value();
1042 }
1043 }
1044 }
1045
1046 // Export the matrix
1047 ITHACAstream::exportMatrix(AF_matrix, "AF", "matlab",
1048 "./ITHACAoutput/Matrices/neutronics/");
1049 return AF_matrix;
1050}
1051
1052Eigen::MatrixXd msrProblem::prec_source(label NFluxmodes, label NPrecmodes,
1053 label family)
1054{
1055 label p = family;
1056 PtrList<volScalarField> Precmodes = choose_group("prec", p);
1057 label PS1size = NFluxmodes;
1058 label PS2size = NPrecmodes;
1059 Eigen::MatrixXd PS_matrix;
1060 PS_matrix.resize(PS1size, PS2size);
1061
1062 // Project everything
1063 for (label i = 0; i < PS1size; i++)
1064 {
1065 for (label j = 0; j < PS2size; j++)
1066 {
1067 PS_matrix(i, j) = fvc::domainIntegrate(Fluxmodes[i] * Precmodes[j]).value();
1068 }
1069 }
1070
1071 savegroupMatrix("PS", p, "./ITHACAoutput/Matrices/neutronics/", PS_matrix);
1072 return PS_matrix;
1073}
1074
1075// * * * * * * * * * * * * * * Precursor Eq. Methods * * * * * * * * * * * * * //
1076
1077List<Eigen::MatrixXd> msrProblem::stream_term(label NUmodes, label NPrecmodes,
1078 label family)
1079{
1080 label p = family;
1081 PtrList<volScalarField> Precmodes = choose_group("prec", p);
1082 label ST1size = NPrecmodes; //Qsizet
1083 label ST2size = NUmodes + liftfield.size();
1084 List<Eigen::MatrixXd> ST_matrix;
1085 ST_matrix.setSize(ST1size);
1086
1087 for (label j = 0; j < ST1size; j++)
1088 {
1089 ST_matrix[j].resize(ST2size, ST1size);
1090 ST_matrix[j].setZero();
1091 }
1092
1093 PtrList<volVectorField> Together(0);
1094
1095 if (liftfield.size() != 0)
1096 {
1097 for (label k = 0; k < liftfield.size(); k++)
1098 {
1099 Together.append(liftfield[k].clone());
1100 }
1101 }
1102
1103 if (NUmodes != 0)
1104 {
1105 for (label k = 0; k < NUmodes; k++)
1106 {
1107 Together.append(Umodes[k].clone());
1108 }
1109 }
1110
1111 // Project everything
1112 for (label i = 0; i < ST1size; i++)
1113 {
1114 for (label j = 0; j < ST2size; j++)
1115 {
1116 for (label k = 0; k < ST1size; k++)
1117 {
1118 ST_matrix[i](j, k) = fvc::domainIntegrate(Precmodes[i] * fvc::div(
1119 fvc::interpolate(Together[j]) & Together[j].mesh().Sf(), Precmodes[k])).value();
1120 }
1121 }
1122 }
1123
1124 savegroupMatrix("ST", p, "./ITHACAoutput/Matrices/neutronics/", ST_matrix);
1125 return ST_matrix;
1126}
1127
1128Eigen::MatrixXd msrProblem::prec_mass(label NPrecmodes, label family)
1129{
1130 label p = family;
1131 PtrList<volScalarField> Precmodes = choose_group("prec", p);
1132 label MPsize = NPrecmodes;
1133 Eigen::MatrixXd MP_matrix;
1134 MP_matrix.resize(MPsize, MPsize);
1135
1136 // Project everything
1137 for (label i = 0; i < MPsize; i++)
1138 {
1139 for (label j = 0; j < MPsize; j++)
1140 {
1141 MP_matrix(i, j) = fvc::domainIntegrate(Precmodes[i] * Precmodes[j]).value();
1142 }
1143 }
1144
1145 savegroupMatrix("MP", p, "./ITHACAoutput/Matrices/neutronics/", MP_matrix);
1146 return MP_matrix;
1147}
1148
1149Eigen::MatrixXd msrProblem::laplacian_prec(label NPrecmodes, label family)
1150{
1151 label p = family;
1152 PtrList<volScalarField> Precmodes = choose_group("prec", p);
1153 label LPsize = NPrecmodes;
1154 Eigen::MatrixXd LP_matrix;
1155 LP_matrix.resize(LPsize, LPsize);
1156
1157 for (label i = 0; i < LPsize; i++)
1158 {
1159 for (label j = 0; j < LPsize; j++)
1160 {
1161 LP_matrix(i, j) = fvc::domainIntegrate(Precmodes[i] * fvc::laplacian(
1162 dimensionedScalar("1", dimless, 1), Precmodes[j])).value();
1163 }
1164 }
1165
1166 savegroupMatrix("LP", p, "./ITHACAoutput/Matrices/neutronics/", LP_matrix);
1167 return LP_matrix;
1168}
1169
1170List<Eigen::MatrixXd> msrProblem::flux_source(label NFluxmodes,
1171 label NPrecmodes, label NCmodes, label family)
1172{
1173 label p = family;
1174 PtrList<volScalarField> Precmodes = choose_group("prec", p);
1175 label FSsize = NPrecmodes;
1176 List<Eigen::MatrixXd> FS_matrix;
1177 FS_matrix.setSize(FSsize);
1178
1179 for (label j = 0; j < FSsize; j++)
1180 {
1181 FS_matrix[j].resize(NCmodes, NFluxmodes);
1182 FS_matrix[j].setZero();
1183 }
1184
1185 // Project everything
1186 for (label i = 0; i < FSsize; i++)
1187 {
1188 for (label j = 0; j < NCmodes; j++)
1189 {
1190 for (label k = 0; k < NFluxmodes; k++)
1191 {
1192 FS_matrix[i](j, k) = fvc::domainIntegrate(Precmodes[i] * NSFmodes[j] *
1193 Fluxmodes[k]).value();
1194 }
1195 }
1196 }
1197
1198 savegroupMatrix("FS", p, "./ITHACAoutput/Matrices/neutronics/", FS_matrix);
1199 return FS_matrix;
1200}
1201
1202
1203
1204// * * * * * * * * * * Decay heat Eq. Methods * * * * * * * * * * * * * //
1205
1206
1207
1208List<Eigen::MatrixXd> msrProblem::stream_dec(label NUmodes, label NDecmodes,
1209 label decgroup)
1210{
1211 label g = decgroup;
1212 PtrList<volScalarField> Decmodes = choose_group("dec", g);
1213 label SD1size = NDecmodes; //Qsizet
1214 label SD2size = NUmodes + liftfield.size();
1215 List<Eigen::MatrixXd> SD_matrix;
1216 SD_matrix.setSize(SD1size);
1217
1218 for (label j = 0; j < SD1size; j++)
1219 {
1220 SD_matrix[j].resize(SD2size, SD1size);
1221 SD_matrix[j].setZero();
1222 }
1223
1224 PtrList<volVectorField> Together(0);
1225
1226 if (liftfield.size() != 0)
1227 {
1228 for (label k = 0; k < liftfield.size(); k++)
1229 {
1230 Together.append(liftfield[k].clone());
1231 }
1232 }
1233
1234 if (NUmodes != 0)
1235 {
1236 for (label k = 0; k < NUmodes; k++)
1237 {
1238 Together.append(Umodes[k].clone());
1239 }
1240 }
1241
1242 // Project everything
1243 for (label i = 0; i < SD1size; i++)
1244 {
1245 for (label j = 0; j < SD2size; j++)
1246 {
1247 for (label k = 0; k < SD1size; k++)
1248 {
1249 SD_matrix[i](j, k) = fvc::domainIntegrate(Decmodes[i] * fvc::div(
1250 fvc::interpolate(Together[j]) & Together[j].mesh().Sf(), Decmodes[k])).value();
1251 }
1252 }
1253 }
1254
1255 savegroupMatrix("SD", g, "./ITHACAoutput/Matrices/thermal/", SD_matrix);
1256 return SD_matrix;
1257}
1258
1259Eigen::MatrixXd msrProblem::dec_mass(label NDecmodes, label decgroup)
1260{
1261 label g = decgroup;
1262 PtrList<volScalarField> Decmodes = choose_group("dec", g);
1263 label MDsize = NDecmodes;
1264 Eigen::MatrixXd MD_matrix;
1265 MD_matrix.resize(MDsize, MDsize);
1266
1267 // Project everything
1268 for (label i = 0; i < MDsize; i++)
1269 {
1270 for (label j = 0; j < MDsize; j++)
1271 {
1272 MD_matrix(i, j) = fvc::domainIntegrate(Decmodes[i] * Decmodes[j]).value();
1273 }
1274 }
1275
1276 savegroupMatrix("MD", g, "./ITHACAoutput/Matrices/thermal/", MD_matrix);
1277 return MD_matrix;
1278}
1279
1280
1281Eigen::MatrixXd msrProblem::laplacian_dec(label NDecmodes, label decgroup)
1282{
1283 label g = decgroup;
1284 PtrList<volScalarField> Decmodes = choose_group("dec", g);
1285 label LDsize = NDecmodes;
1286 Eigen::MatrixXd LD_matrix;
1287 LD_matrix.resize(LDsize, LDsize);
1288
1289 // Project everything
1290 for (label i = 0; i < LDsize; i++)
1291 {
1292 for (label j = 0; j < LDsize; j++)
1293 {
1294 LD_matrix(i, j) = fvc::domainIntegrate(Decmodes[i] * fvc::laplacian(
1295 dimensionedScalar("1", dimless, 1), Decmodes[j])).value();
1296 }
1297 }
1298
1299 savegroupMatrix("LD", g, "./ITHACAoutput/Matrices/thermal/", LD_matrix);
1300 return LD_matrix;
1301}
1302
1303List<Eigen::MatrixXd> msrProblem::dec_fluxsource(label NFluxmodes,
1304 label NDecmodes, label NCmodes, label decgroup)
1305{
1306 label g = decgroup;
1307 PtrList<volScalarField> Decmodes = choose_group("dec", g);
1308 label DFSsize = NDecmodes;
1309 List<Eigen::MatrixXd> DFS_matrix;
1310 DFS_matrix.setSize(DFSsize);
1311
1312 for (label j = 0; j < DFSsize; j++)
1313 {
1314 DFS_matrix[j].resize(NCmodes, NFluxmodes);
1315 DFS_matrix[j].setZero();
1316 }
1317
1318 // Project everything
1319 for (label i = 0; i < DFSsize; i++)
1320 {
1321 for (label j = 0; j < NCmodes; j++)
1322 {
1323 for (label k = 0; k < NFluxmodes; k++)
1324 {
1325 DFS_matrix[i](j, k) = fvc::domainIntegrate(Decmodes[i] * SPmodes[j] *
1326 Fluxmodes[k]).value();
1327 }
1328 }
1329 }
1330
1331 savegroupMatrix("DFS", g, "./ITHACAoutput/Matrices/thermal/", DFS_matrix);
1332 return DFS_matrix;
1333}
1334
1335
1336
1337// * * * * * * * * * * * Temperature Eq. Methods * * * * * * * * * * * * * //
1338
1339Eigen::MatrixXd msrProblem::mass_temp(label NTmodes)
1340{
1341 label TMsize = NTmodes + liftfieldT.size();
1342 Eigen::MatrixXd TM_matrix;
1343 TM_matrix.resize(TMsize, TMsize);
1344 PtrList<volScalarField> TogetherT(0);
1345
1346 if (liftfieldT.size() != 0)
1347 {
1348 for (label k = 0; k < liftfieldT.size(); k++)
1349 {
1350 TogetherT.append(liftfieldT[k].clone());
1351 }
1352 }
1353
1354 if (NTmodes != 0)
1355 {
1356 for (label k = 0; k < NTmodes; k++)
1357 {
1358 TogetherT.append(Tmodes[k].clone());
1359 }
1360 }
1361
1362 // Project everything
1363 for (label i = 0; i < TMsize; i++)
1364 {
1365 for (label j = 0; j < TMsize; j++)
1366 {
1367 TM_matrix(i, j) = fvc::domainIntegrate(TogetherT[i] * TogetherT[j]).value();
1368 }
1369 }
1370
1371 // Export the matrix
1372 ITHACAstream::exportMatrix(TM_matrix, "TM", "matlab",
1373 "./ITHACAoutput/Matrices/thermal/");
1374 return TM_matrix;
1375}
1376
1377List<Eigen::MatrixXd> msrProblem::temp_stream(label NUmodes, label NTmodes)
1378{
1379 label TS1size = NTmodes + liftfieldT.size(); //Qsizet
1380 label TS2size = NUmodes + liftfield.size();
1381 List<Eigen::MatrixXd> TS_matrix;
1382 TS_matrix.setSize(TS1size);
1383
1384 for (label j = 0; j < TS1size; j++)
1385 {
1386 TS_matrix[j].resize(TS2size, TS1size);
1387 TS_matrix[j].setZero();
1388 }
1389
1390 PtrList<volVectorField> Together(0);
1391
1392 if (liftfield.size() != 0)
1393 {
1394 for (label k = 0; k < liftfield.size(); k++)
1395 {
1396 Together.append(liftfield[k].clone());
1397 }
1398 }
1399
1400 if (NUmodes != 0)
1401 {
1402 for (label k = 0; k < NUmodes; k++)
1403 {
1404 Together.append(Umodes[k].clone());
1405 }
1406 }
1407
1408 PtrList<volScalarField> TogetherT(0);
1409
1410 if (liftfieldT.size() != 0)
1411 {
1412 for (label k = 0; k < liftfieldT.size(); k++)
1413 {
1414 TogetherT.append(liftfieldT[k].clone());
1415 }
1416 }
1417
1418 if (NTmodes != 0)
1419 {
1420 for (label k = 0; k < NTmodes; k++)
1421 {
1422 TogetherT.append(Tmodes[k].clone());
1423 }
1424 }
1425
1426 // Project everything
1427 for (label i = 0; i < TS1size; i++)
1428 {
1429 for (label j = 0; j < TS2size; j++)
1430 {
1431 for (label k = 0; k < TS1size; k++)
1432 {
1433 TS_matrix[i](j, k) = fvc::domainIntegrate(TogetherT[i] * fvc::div(
1434 fvc::interpolate(Together[j]) & Together[j].mesh().Sf(), TogetherT[k])).value();
1435 }
1436 }
1437 }
1438
1439 // Export the matrix
1440 ITHACAstream::exportMatrix(TS_matrix, "TS", "matlab",
1441 "./ITHACAoutput/Matrices/thermal/");
1442 return TS_matrix;
1443}
1444
1445
1446Eigen::MatrixXd msrProblem::laplacian_temp(label NTmodes)
1447{
1448 label LTsize = NTmodes + liftfieldT.size();
1449 Eigen::MatrixXd LT_matrix;
1450 LT_matrix.resize(LTsize, LTsize);
1451 PtrList<volScalarField> TogetherT(0);
1452
1453 if (liftfieldT.size() != 0)
1454 {
1455 for (label k = 0; k < liftfieldT.size(); k++)
1456 {
1457 TogetherT.append(liftfieldT[k].clone());
1458 }
1459 }
1460
1461 if (NTmodes != 0)
1462 {
1463 for (label k = 0; k < NTmodes; k++)
1464 {
1465 TogetherT.append(Tmodes[k].clone());
1466 }
1467 }
1468
1469 // Project everything
1470 for (label i = 0; i < LTsize; i++)
1471 {
1472 for (label j = 0; j < LTsize; j++)
1473 {
1474 LT_matrix(i, j) = fvc::domainIntegrate(TogetherT[i] * fvc::laplacian(
1475 dimensionedScalar("1", dimless, 1), TogetherT[j])).value();
1476 }
1477 }
1478
1479 // Export the matrix
1480 ITHACAstream::exportMatrix(LT_matrix, "LT", "matlab",
1481 "./ITHACAoutput/Matrices/thermal/");
1482 return LT_matrix;
1483}
1484
1485List<Eigen::MatrixXd> msrProblem::temp_XSfluxsource(label NTmodes,
1486 label NFluxmodes, label NCmodes)
1487{
1488 label TXSsize = NTmodes + liftfieldT.size();
1489 List<Eigen::MatrixXd> TXS_matrix;
1490 TXS_matrix.setSize(TXSsize);
1491
1492 for (label j = 0; j < TXSsize; j++)
1493 {
1494 TXS_matrix[j].resize(NCmodes, NFluxmodes);
1495 TXS_matrix[j].setZero();
1496 }
1497
1498 PtrList<volScalarField> TogetherT(0);
1499
1500 if (liftfieldT.size() != 0)
1501 {
1502 for (label k = 0; k < liftfieldT.size(); k++)
1503 {
1504 TogetherT.append(liftfieldT[k].clone());
1505 }
1506 }
1507
1508 if (NTmodes != 0)
1509 {
1510 for (label k = 0; k < NTmodes; k++)
1511 {
1512 TogetherT.append(Tmodes[k].clone());
1513 }
1514 }
1515
1516 // Project everything
1517 for (label i = 0; i < TXSsize; i++)
1518 {
1519 for (label j = 0; j < NCmodes; j++)
1520 {
1521 for (label k = 0; k < NFluxmodes; k++)
1522 {
1523 TXS_matrix[i](j, k) = fvc::domainIntegrate(TogetherT[i] * TXSmodes[j] *
1524 Fluxmodes[k]).value();
1525 }
1526 }
1527 }
1528
1529 ITHACAstream::exportMatrix(TXS_matrix, "TXS", "matlab",
1530 "./ITHACAoutput/Matrices/thermal/");
1531 return TXS_matrix;
1532}
1533
1534List<Eigen::MatrixXd> msrProblem::temp_heatsource(label NTmodes,
1535 label NDecmodes, label NCmodes, label decgroup)
1536{
1537 label g = decgroup;
1538 PtrList<volScalarField> Decmodes = choose_group("dec", g);
1539 label THSsize = NTmodes + liftfieldT.size();
1540 List<Eigen::MatrixXd> THS_matrix;
1541 THS_matrix.setSize(THSsize);
1542
1543 for (label j = 0; j < THSsize; j++)
1544 {
1545 THS_matrix[j].resize(NCmodes, NDecmodes);
1546 THS_matrix[j].setZero();
1547 }
1548
1549 PtrList<volScalarField> TogetherT(0);
1550
1551 if (liftfieldT.size() != 0)
1552 {
1553 for (label k = 0; k < liftfieldT.size(); k++)
1554 {
1555 TogetherT.append(liftfieldT[k].clone());
1556 }
1557 }
1558
1559 if (NTmodes != 0)
1560 {
1561 for (label k = 0; k < NTmodes; k++)
1562 {
1563 TogetherT.append(Tmodes[k].clone());
1564 }
1565 }
1566
1567 // Project everything
1568 for (label i = 0; i < THSsize; i++)
1569 {
1570 for (label j = 0; j < NCmodes; j++)
1571 {
1572 for (label k = 0; k < NDecmodes; k++)
1573 {
1574 THS_matrix[i](j, k) = fvc::domainIntegrate(TogetherT[i] * vmodes[j] *
1575 Decmodes[k]).value();
1576 }
1577 }
1578 }
1579
1580 savegroupMatrix("THS", g, "./ITHACAoutput/Matrices/thermal/", THS_matrix);
1581 return THS_matrix;
1582}
1583
1584PtrList<volScalarField> msrProblem::choose_group(string field, label ith)
1585{
1586 if (field == "prec")
1587 {
1588 switch (ith)
1589 {
1590 case 1:
1591 return Prec1modes;
1592 break;
1593
1594 case 2:
1595 return Prec2modes;
1596 break;
1597
1598 case 3:
1599 return Prec3modes;
1600 break;
1601
1602 case 4:
1603 return Prec4modes;
1604 break;
1605
1606 case 5:
1607 return Prec5modes;
1608 break;
1609
1610 case 6:
1611 return Prec6modes;
1612 break;
1613
1614 case 7:
1615 return Prec7modes;
1616 break;
1617
1618 case 8:
1619 return Prec8modes;
1620 break;
1621 }
1622 }
1623
1624 if (field == "dec")
1625 {
1626 switch (ith)
1627 {
1628 case 1:
1629 return Dec1modes;
1630 break;
1631
1632 case 2:
1633 return Dec2modes;
1634 break;
1635
1636 case 3:
1637 return Dec3modes;
1638 break;
1639 }
1640 }
1641}
1642
1643
1644template<typename M>
1645void msrProblem::savegroupMatrix(string nome, label n, word folder, M matrice)
1646{
1647 nome.append(std::to_string(n));
1648 word name = nome;
1649 ITHACAstream::exportMatrix(matrice, name, "matlab", folder);
1650}
1651
1653{
1655 homboolU = true;
1656}
1657
1663
1665{
1666 volVectorField& U = _U();
1667 volScalarField& p = _p();
1668 volScalarField& flux = _flux();
1669 volScalarField& prec1 = _prec1();
1670 volScalarField& prec2 = _prec2();
1671 volScalarField& prec3 = _prec3();
1672 volScalarField& prec4 = _prec4();
1673 volScalarField& prec5 = _prec5();
1674 volScalarField& prec6 = _prec6();
1675 volScalarField& prec7 = _prec7();
1676 volScalarField& prec8 = _prec8();
1677 volScalarField& T = _T();
1678 volScalarField& dec1 = _dec1();
1679 volScalarField& dec2 = _dec2();
1680 volScalarField& dec3 = _dec3();
1681 volScalarField& v = _v();
1682 volScalarField& D = _D();
1683 volScalarField& NSF = _NSF();
1684 volScalarField& A = _A();
1685 volScalarField& SP = _SP();
1686 volScalarField& TXS = _TXS();
1687 volScalarField powerDens = (flux * (1 - _decbetaTot()) * _SP() + _decLam1() *
1688 dec1 + _decLam2() * dec2 + _decLam3() * dec3).ref();
1689 powerDens.rename("powerDens");
1690 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
1691 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
1692 ITHACAstream::read_fields(Fluxfield, flux, "./ITHACAoutput/Offline/");
1693 ITHACAstream::read_fields(Tfield, T, "./ITHACAoutput/Offline/");
1694 ITHACAstream::read_fields(Prec1field, prec1, "./ITHACAoutput/Offline/");
1695 ITHACAstream::read_fields(Prec2field, prec2, "./ITHACAoutput/Offline/");
1696 ITHACAstream::read_fields(Prec3field, prec3, "./ITHACAoutput/Offline/");
1697 ITHACAstream::read_fields(Prec4field, prec4, "./ITHACAoutput/Offline/");
1698 ITHACAstream::read_fields(Prec5field, prec5, "./ITHACAoutput/Offline/");
1699 ITHACAstream::read_fields(Prec6field, prec6, "./ITHACAoutput/Offline/");
1700 ITHACAstream::read_fields(Prec7field, prec7, "./ITHACAoutput/Offline/");
1701 ITHACAstream::read_fields(Prec8field, prec8, "./ITHACAoutput/Offline/");
1702 ITHACAstream::read_fields(Dec1field, dec1, "./ITHACAoutput/Offline/");
1703 ITHACAstream::read_fields(Dec2field, dec2, "./ITHACAoutput/Offline/");
1704 ITHACAstream::read_fields(Dec3field, dec3, "./ITHACAoutput/Offline/");
1705 ITHACAstream::read_fields(PowerDensfield, powerDens, "./ITHACAoutput/Offline/");
1706 ITHACAstream::read_fields(vFields, v, "./ITHACAoutput/Offline/");
1707 ITHACAstream::read_fields(DFields, D, "./ITHACAoutput/Offline/");
1708 ITHACAstream::read_fields(NSFFields, NSF, "./ITHACAoutput/Offline/");
1709 ITHACAstream::read_fields(AFields, A, "./ITHACAoutput/Offline/");
1710 ITHACAstream::read_fields(SPFields, SP, "./ITHACAoutput/Offline/");
1711 ITHACAstream::read_fields(TXSFields, TXS, "./ITHACAoutput/Offline/");
1712 mu_samples =
1713 ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt");
1714}
1715
1716void msrProblem::readMSRfields(std::string& dir)
1717{
1718 volVectorField& U = _U();
1719 volScalarField& p = _p();
1720 volScalarField& flux = _flux();
1721 volScalarField& prec1 = _prec1();
1722 volScalarField& prec2 = _prec2();
1723 volScalarField& prec3 = _prec3();
1724 volScalarField& prec4 = _prec4();
1725 volScalarField& prec5 = _prec5();
1726 volScalarField& prec6 = _prec6();
1727 volScalarField& prec7 = _prec7();
1728 volScalarField& prec8 = _prec8();
1729 volScalarField& T = _T();
1730 volScalarField& dec1 = _dec1();
1731 volScalarField& dec2 = _dec2();
1732 volScalarField& dec3 = _dec3();
1733 volScalarField& v = _v();
1734 volScalarField& D = _D();
1735 volScalarField& NSF = _NSF();
1736 volScalarField& A = _A();
1737 volScalarField& SP = _SP();
1738 volScalarField& TXS = _TXS();
1739 volScalarField powerDens = (flux * (1 - _decbetaTot()) * _SP() + _decLam1() *
1740 dec1 + _decLam2() * dec2 + _decLam3() * dec3).ref();
1741 powerDens.rename("powerDens");
1742 std::string folder = dir;
1743
1744 if (ITHACAutilities::check_folder(folder) == true)
1745 {
1746 for (label i = 0; i < mu.cols(); i++)
1747 {
1748 folder.append(std::to_string(i));
1749 folder.append("/");
1772 folder = dir;
1773 }
1774 }
1775 else
1776 {
1777 std::cout << "Error" << std::endl;
1778 }
1779}
1780
1781
1783{
1784 const volScalarField& nu = _laminarTransport().nu();
1785 volScalarField& ciao = const_cast<volScalarField&>(nu);
1786 this->assignIF(ciao, mu);
1787
1788 for (label i = 0; i < ciao.boundaryFieldRef().size(); i++)
1789 {
1790 this->assignBC(ciao, i, mu);
1791 }
1792}
1793
1795{
1796 NCmodes = NC;
1797 Eigen::MatrixXd Ncoeff_v = ITHACAutilities::getCoeffs(vFields, vmodes);
1798 ITHACAstream::exportMatrix(Ncoeff_v, "Ncoeff_v", "matlab",
1799 "./ITHACAoutput/Matrices/");
1800 label Ncol = Ncoeff_v.cols();
1801 SAMPLES_v.resize(NCmodes);
1802 rbfsplines_v.resize(NCmodes);
1803
1804 for (label i = 0; i < NCmodes; i++)
1805 {
1806 std::cout << "Constructing v RadialBasisFunction for mode " << i + 1 <<
1807 std::endl;
1808 SAMPLES_v[i] = new SPLINTER::DataTable(1, 1);
1809
1810 for (label j = 0; j < Ncol; j++)
1811 {
1812 SAMPLES_v[i]->addSample(mu_samples.row(j), Ncoeff_v(i, j));
1813 }
1814
1815 rbfsplines_v[i] = new SPLINTER::RBFSpline(*SAMPLES_v[i],
1816 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
1817 }
1818
1819 Eigen::MatrixXd Ncoeff_D = ITHACAutilities::getCoeffs(DFields, Dmodes);
1820 ITHACAstream::exportMatrix(Ncoeff_D, "Ncoeff_D", "matlab",
1821 "./ITHACAoutput/Matrices/");
1822 SAMPLES_D.resize(NCmodes);
1823 rbfsplines_D.resize(NCmodes);
1824
1825 for (label i = 0; i < NCmodes; i++)
1826 {
1827 std::cout << "Constructing D RadialBasisFunction for mode " << i + 1 <<
1828 std::endl;
1829 SAMPLES_D[i] = new SPLINTER::DataTable(1, 1);
1830
1831 for (label j = 0; j < Ncol; j++)
1832 {
1833 SAMPLES_D[i]->addSample(mu_samples.row(j), Ncoeff_D(i, j));
1834 }
1835
1836 rbfsplines_D[i] = new SPLINTER::RBFSpline(*SAMPLES_D[i],
1837 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
1838 }
1839
1840 Eigen::MatrixXd Ncoeff_NSF = ITHACAutilities::getCoeffs(NSFFields,
1841 NSFmodes);
1842 ITHACAstream::exportMatrix(Ncoeff_NSF, "Ncoeff_NSF", "matlab",
1843 "./ITHACAoutput/Matrices/");
1844 SAMPLES_NSF.resize(NCmodes);
1845 rbfsplines_NSF.resize(NCmodes);
1846
1847 for (label i = 0; i < NCmodes; i++)
1848 {
1849 std::cout << "Constructing NSF RadialBasisFunction for mode " << i + 1 <<
1850 std::endl;
1851 SAMPLES_NSF[i] = new SPLINTER::DataTable(1, 1);
1852
1853 for (label j = 0; j < Ncol; j++)
1854 {
1855 SAMPLES_NSF[i]->addSample(mu_samples.row(j), Ncoeff_NSF(i, j));
1856 }
1857
1858 rbfsplines_NSF[i] = new SPLINTER::RBFSpline(*SAMPLES_NSF[i],
1859 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
1860 }
1861
1862 Eigen::MatrixXd Ncoeff_A = ITHACAutilities::getCoeffs(AFields, Amodes);
1863 ITHACAstream::exportMatrix(Ncoeff_A, "Ncoeff_A", "matlab",
1864 "./ITHACAoutput/Matrices/");
1865 SAMPLES_A.resize(NCmodes);
1866 rbfsplines_A.resize(NCmodes);
1867
1868 for (label i = 0; i < NCmodes; i++)
1869 {
1870 std::cout << "Constructing A RadialBasisFunction for mode " << i + 1 <<
1871 std::endl;
1872 SAMPLES_A[i] = new SPLINTER::DataTable(1, 1);
1873
1874 for (label j = 0; j < Ncol; j++)
1875 {
1876 SAMPLES_A[i]->addSample(mu_samples.row(j), Ncoeff_A(i, j));
1877 }
1878
1879 rbfsplines_A[i] = new SPLINTER::RBFSpline(*SAMPLES_A[i],
1880 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
1881 }
1882
1883 Eigen::MatrixXd Ncoeff_SP = ITHACAutilities::getCoeffs(SPFields,
1884 SPmodes);
1885 ITHACAstream::exportMatrix(Ncoeff_SP, "Ncoeff_SP", "matlab",
1886 "./ITHACAoutput/Matrices/");
1887 SAMPLES_SP.resize(NCmodes);
1888 rbfsplines_SP.resize(NCmodes);
1889
1890 for (label i = 0; i < NCmodes; i++)
1891 {
1892 std::cout << "Constructing SP RadialBasisFunction for mode " << i + 1 <<
1893 std::endl;
1894 SAMPLES_SP[i] = new SPLINTER::DataTable(1, 1);
1895
1896 for (label j = 0; j < Ncol; j++)
1897 {
1898 SAMPLES_SP[i]->addSample(mu_samples.row(j), Ncoeff_SP(i, j));
1899 }
1900
1901 rbfsplines_SP[i] = new SPLINTER::RBFSpline(*SAMPLES_SP[i],
1902 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
1903 }
1904
1905 Eigen::MatrixXd Ncoeff_TXS = ITHACAutilities::getCoeffs(TXSFields,
1906 TXSmodes);
1907 ITHACAstream::exportMatrix(Ncoeff_TXS, "Ncoeff_TXS", "matlab",
1908 "./ITHACAoutput/Matrices/");
1909 SAMPLES_TXS.resize(NCmodes);
1910 rbfsplines_TXS.resize(NCmodes);
1911
1912 for (label i = 0; i < NCmodes; i++)
1913 {
1914 std::cout << "Constructing TXS RadialBasisFunction for mode " << i + 1 <<
1915 std::endl;
1916 SAMPLES_TXS[i] = new SPLINTER::DataTable(1, 1);
1917
1918 for (label j = 0; j < Ncol; j++)
1919 {
1920 SAMPLES_TXS[i]->addSample(mu_samples.row(j), Ncoeff_TXS(i, j));
1921 }
1922
1923 rbfsplines_TXS[i] = new SPLINTER::RBFSpline(*SAMPLES_TXS[i],
1924 SPLINTER::RadialBasisFunctionType::GAUSSIAN);
1925 }
1926}
1927
1928
1930{
1931 for (label k = 0; k < inletIndexT.rows(); k++)
1932 {
1933 Time& runTime = _runTime();
1934 fvMesh& mesh = _mesh();
1935 volScalarField& T = _T();
1936 volVectorField& U = _U();
1937 surfaceScalarField& phi = _phi();
1938 phi = linearInterpolate(U) & mesh.Sf();
1939 //pisoControl potentialFlow(mesh,"potentialFLow");
1940 simpleControl simple(mesh);
1941 IOMRFZoneList& MRF = _MRF();
1942 singlePhaseTransportModel& laminarTransport = _laminarTransport();
1943 turbulence = autoPtr<incompressible::turbulenceModel>
1944 (
1945 incompressible::turbulenceModel::New(U, phi, laminarTransport)
1946 );
1947 dimensionedScalar& nu = _nu();
1948 dimensionedScalar& Pr = _Pr();
1949 dimensionedScalar& Prt = _Prt();
1950 volScalarField& alphat = _alphat();
1951 volScalarField& v = _v();
1952 dimensionedScalar& cp = _CpRef();
1953 dimensionedScalar& rho = _rhoRef();
1954 dimensionedScalar& sp = _SP1_0();
1955 dimensionedScalar& dbtot = _decbetaTot();
1956 label BCind = inletIndexT(k, 0);
1957 volScalarField Tlift("Tlift" + name(k), T);
1958 instantList Times = runTime.times();
1959 runTime.setTime(Times[1], 1);
1960 Info << "Solving a lifting Problem" << endl;
1961 scalar t1 = 1;
1962 scalar t0 = 0;
1963 alphat = turbulence->nut() / Prt;
1964 alphat.correctBoundaryConditions();
1965 volScalarField source
1966 (
1967 IOobject
1968 (
1969 "source",
1970 runTime.timeName(),
1971 mesh,
1972 IOobject::NO_READ,
1973 IOobject::AUTO_WRITE
1974 ),
1975 mesh,
1976 dimensionedScalar("source", dimensionSet(0, -2, -1, 0, 0, 0, 0), 1)
1977 );
1978
1979 for (label j = 0; j < T.boundaryField().size(); j++)
1980 {
1981 if (j == BCind)
1982 {
1983 assignBC(Tlift, j, t1);
1984 assignIF(Tlift, t0);
1985 }
1986 else if (T.boundaryField()[BCind].type() == "fixedValue")
1987 {
1988 assignBC(Tlift, j, t0);
1989 assignIF(Tlift, t0);
1990 }
1991 else
1992 {
1993 }
1994 }
1995
1996 while (simple.correctNonOrthogonal())
1997 {
1998 fvScalarMatrix TEqn
1999 (
2000 fvm::div(phi, Tlift) == fvm::laplacian(turbulence->nu() / Pr + alphat, Tlift)
2001 + sp * (1 - dbtot) / (cp * rho)*source
2002 );
2003 TEqn.solve();
2004 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
2005 << " ClockTime = " << runTime.elapsedClockTime() << " s"
2006 << nl << endl;
2007 }
2008
2009 Tlift.write();
2010 liftfieldT.append(Tlift.clone());
2011 }
2012}
2013
2014
2016{
2017 volScalarField& p = _p();
2018 volScalarField& p0 = _p0();
2019 volVectorField& U = _U();
2020 volVectorField& U0 = _U0();
2021 surfaceScalarField& phi = _phi();
2022 surfaceScalarField& phi0 = _phi0();
2023 volScalarField& flux = _flux();
2024 volScalarField& flux0 = _flux0();
2025 volScalarField& prec1 = _prec1();
2026 volScalarField& prec10 = _prec10();
2027 volScalarField& prec2 = _prec2();
2028 volScalarField& prec20 = _prec20();
2029 volScalarField& prec3 = _prec3();
2030 volScalarField& prec30 = _prec30();
2031 volScalarField& prec4 = _prec4();
2032 volScalarField& prec40 = _prec40();
2033 volScalarField& prec5 = _prec5();
2034 volScalarField& prec50 = _prec50();
2035 volScalarField& prec6 = _prec6();
2036 volScalarField& prec60 = _prec60();
2037 volScalarField& prec7 = _prec7();
2038 volScalarField& prec70 = _prec70();
2039 volScalarField& prec8 = _prec8();
2040 volScalarField& prec80 = _prec80();
2041 volScalarField& T = _T();
2042 volScalarField& T0 = _T0();
2043 volScalarField& dec1 = _dec1();
2044 volScalarField& dec10 = _dec10();
2045 volScalarField& dec2 = _dec2();
2046 volScalarField& dec20 = _dec20();
2047 volScalarField& dec3 = _dec3();
2048 volScalarField& dec30 = _dec30();
2049 dimensionedScalar& Keff = _Keff();
2050 dimensionedScalar& K0 = _K0();
2051 volScalarField& v = _v();
2052 volScalarField& v0 = _v0();
2053 volScalarField& NSF = _NSF();
2054 volScalarField& NSF0 = _NSF0();
2055 volScalarField& A = _A();
2056 volScalarField& A0 = _A0();
2057 volScalarField& D = _D();
2058 volScalarField& D0 = _D0();
2059 volScalarField& SP = _SP();
2060 volScalarField& SP0 = _SP0();
2061 volScalarField& TXS = _TXS();
2062 volScalarField& TXS0 = _TXS0();
2063 fvMesh& mesh = _mesh();
2064 p = p0;
2065 U = U0;
2066 phi = phi0;
2067 turbulence.reset(
2068 (incompressible::turbulenceModel::New(U, phi, _laminarTransport())).ptr()
2069 );
2070 flux = flux0;
2071 prec1 = prec10;
2072 prec2 = prec20;
2073 prec3 = prec30;
2074 prec4 = prec40;
2075 prec5 = prec50;
2076 prec6 = prec60;
2077 prec7 = prec70;
2078 prec8 = prec80;
2079
2080 if (precInBool == true)
2081 {
2082 prec1.boundaryFieldRef().set(precinIndex,
2083 fvPatchField<scalar>::New(prec10.boundaryField()[precinIndex].type(),
2084 mesh.boundary()[precinIndex], prec1));
2085 prec2.boundaryFieldRef().set(precinIndex,
2086 fvPatchField<scalar>::New(prec20.boundaryField()[precinIndex].type(),
2087 mesh.boundary()[precinIndex], prec2));
2088 prec3.boundaryFieldRef().set(precinIndex,
2089 fvPatchField<scalar>::New(prec30.boundaryField()[precinIndex].type(),
2090 mesh.boundary()[precinIndex], prec3));
2091 prec4.boundaryFieldRef().set(precinIndex,
2092 fvPatchField<scalar>::New(prec40.boundaryField()[precinIndex].type(),
2093 mesh.boundary()[precinIndex], prec4));
2094 prec5.boundaryFieldRef().set(precinIndex,
2095 fvPatchField<scalar>::New(prec50.boundaryField()[precinIndex].type(),
2096 mesh.boundary()[precinIndex], prec5));
2097 prec6.boundaryFieldRef().set(precinIndex,
2098 fvPatchField<scalar>::New(prec60.boundaryField()[precinIndex].type(),
2099 mesh.boundary()[precinIndex], prec6));
2100 prec7.boundaryFieldRef().set(precinIndex,
2101 fvPatchField<scalar>::New(prec70.boundaryField()[precinIndex].type(),
2102 mesh.boundary()[precinIndex], prec7));
2103 prec8.boundaryFieldRef().set(precinIndex,
2104 fvPatchField<scalar>::New(prec80.boundaryField()[precinIndex].type(),
2105 mesh.boundary()[precinIndex], prec8));
2106 }
2107
2108 T = T0;
2109 dec1 = dec10;
2110 dec2 = dec20;
2111 dec3 = dec30;
2112 Keff = K0;
2113 v = v0;
2114 NSF = NSF0;
2115 A = A0;
2116 D = D0;
2117 SP = SP0;
2118 TXS = TXS0;
2119}
2120
2121
2122
2123
2124
2125
2126
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
IOdictionary * ITHACAdict
Dictionary for input objects from file.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
List< Eigen::MatrixXd > stream_term(label NUmodes, label NPrecmodes, label family)
precursor eq. methods:
Eigen::MatrixXd MP7_matrix
precursor mass term-7
Definition msrProblem.H:418
autoPtr< dimensionedScalar > _betaTot
Definition msrProblem.H:130
autoPtr< volScalarField > _NSF
Definition msrProblem.H:109
Eigen::MatrixXd MD1_matrix
decay heat mass term-1
Definition msrProblem.H:464
void homogenizeT()
Method to compute the homogenized temperature field, it also sets homboolT=true.
std::vector< SPLINTER::DataTable * > SAMPLES_NSF
Definition msrProblem.H:534
List< Eigen::MatrixXd > ST5_matrix
precursor stream term-5
Definition msrProblem.H:398
PtrList< volScalarField > Prec6field
List of pointers used to form the prec6 snapshots matrix.
Definition msrProblem.H:219
autoPtr< volScalarField > _prec8
Definition msrProblem.H:147
List< Eigen::MatrixXd > THS2_matrix
temperature decay heat source term-2
Definition msrProblem.H:493
Eigen::MatrixXd LP3_matrix
precursor laplacian term-3
Definition msrProblem.H:426
std::vector< SPLINTER::DataTable * > SAMPLES_v
Definition msrProblem.H:530
List< Eigen::MatrixXd > SD1_matrix
decay heat stream term-1
Definition msrProblem.H:458
autoPtr< incompressible::turbulenceModel > turbulence
Definition msrProblem.H:58
Eigen::MatrixXd BC1_matrix
PPE BC1.
Definition msrProblem.H:350
std::vector< SPLINTER::DataTable * > SAMPLES_A
Definition msrProblem.H:536
PtrList< volScalarField > NSFmodes
List of pointers used to form the NSF snapshosts matrix.
Definition msrProblem.H:312
autoPtr< volScalarField > _dec3
Definition msrProblem.H:163
PtrList< volScalarField > Prec5modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:282
autoPtr< volScalarField > _SP
Definition msrProblem.H:112
autoPtr< dimensionedScalar > _rhoRef
Definition msrProblem.H:66
List< Eigen::MatrixXd > PF_matrix
production flux
Definition msrProblem.H:366
autoPtr< dimensionedScalar > _betaTE
Definition msrProblem.H:139
autoPtr< dimensionedScalar > _Pr
Definition msrProblem.H:64
autoPtr< volScalarField > _T
Definition msrProblem.H:159
autoPtr< volScalarField > _v0
Definition msrProblem.H:90
PtrList< volScalarField > choose_group(string field, label ith)
method to choose one field among precs & decs field can be "prec" or "dec" only if field==prec then i...
void savegroupMatrix(string nome, label n, word folder, M matrice)
method to save matrices for precs and decs M can be an Eigen::MatrixXd or List<Eigen::MatrixXd> nome ...
autoPtr< volScalarField > _dec30
Definition msrProblem.H:88
PtrList< volScalarField > Prec4modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:279
List< Eigen::MatrixXd > FS3_matrix
precursor flux source term-3
Definition msrProblem.H:442
List< Eigen::MatrixXd > BC2_matrix
PPE BC2.
Definition msrProblem.H:353
List< Eigen::MatrixXd > ST7_matrix
precursor stream term-7
Definition msrProblem.H:402
PtrList< volScalarField > Prec8field
List of pointers used to form the prec8 snapshots matrix.
Definition msrProblem.H:225
autoPtr< dimensionedScalar > _beta8
Definition msrProblem.H:129
autoPtr< volScalarField > _dec2
Definition msrProblem.H:162
List< Eigen::MatrixXd > FS1_matrix
precursor flux source term-1
Definition msrProblem.H:438
PtrList< volScalarField > Prec1modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:270
List< Eigen::MatrixXd > TXS_matrix
temperature flux source term TXS
Definition msrProblem.H:489
List< Eigen::MatrixXd > LF_matrix
laplacian_flux
Definition msrProblem.H:362
PtrList< volScalarField > Pmodes
List of pointers used to form the pressure modes.
Definition msrProblem.H:261
label NUmodes
Number of modes adopted during Galerkin projection.
Definition msrProblem.H:500
autoPtr< dimensionedScalar > _Tref
Definition msrProblem.H:69
Eigen::MatrixXd MP1_matrix
precursor mass term-1
Definition msrProblem.H:406
autoPtr< dimensionedScalar > _CpRef
Definition msrProblem.H:68
autoPtr< dimensionedScalar > _decBeta2
Definition msrProblem.H:136
autoPtr< dimensionedScalar > _beta5
Definition msrProblem.H:126
List< Eigen::MatrixXd > ST8_matrix
precursor stream term-8
Definition msrProblem.H:404
PtrList< volVectorField > Uomfield
List of pointers used to form the homogeneous velocity snapshots.
Definition msrProblem.H:517
List< Eigen::MatrixXd > ST4_matrix
precursor stream term-4
Definition msrProblem.H:396
autoPtr< dimensionedScalar > _lam7
Definition msrProblem.H:120
PtrList< volScalarField > Prec7modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:288
PtrList< volScalarField > vFields
List of pointers used to form the v snapshosts matrix.
Definition msrProblem.H:243
ITHACAparameters * para
Definition msrProblem.H:370
autoPtr< dimensionedScalar > _SP1_0
Definition msrProblem.H:111
Eigen::MatrixXd PS3_matrix
prec_source 3
Definition msrProblem.H:377
Eigen::MatrixXd prec_mass(label NPrecmodes, label family)
List< Eigen::MatrixXd > abs_flux(label NFluxmodes, label NCmodes)
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
Definition msrProblem.H:294
label precinIndex
indexes of inlet and outlet to adopt for precursors boundary conditions
Definition msrProblem.H:184
Eigen::MatrixXd K_matrix
Gradient of pressure matrix.
Definition msrProblem.H:335
autoPtr< volScalarField > _prec2
Definition msrProblem.H:141
List< Eigen::MatrixXd > G_matrix
Divergence of momentum PPE.
Definition msrProblem.H:347
List< Eigen::MatrixXd > AF_matrix
absorption flux
Definition msrProblem.H:368
std::vector< SPLINTER::RBFSpline * > rbfsplines_A
Definition msrProblem.H:548
Eigen::MatrixXd MP2_matrix
precursor mass term-2
Definition msrProblem.H:408
List< Eigen::MatrixXd > FS2_matrix
precursor flux source term-2
Definition msrProblem.H:440
Eigen::MatrixXd LP8_matrix
precursor laplacian term-8
Definition msrProblem.H:436
Eigen::MatrixXd PS6_matrix
prec_source 6
Definition msrProblem.H:383
autoPtr< dimensionedScalar > _lam3
Definition msrProblem.H:116
PtrList< volScalarField > Prec8modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:291
autoPtr< dimensionedScalar > _Keff
Definition msrProblem.H:100
List< Eigen::MatrixXd > FS4_matrix
precursor flux source term-4
Definition msrProblem.H:444
PtrList< volScalarField > Fluxfield
List of pointers used to form the flux snapshots matrix.
Definition msrProblem.H:201
Eigen::MatrixXd BC3_matrix
PPE BC3.
Definition msrProblem.H:356
Eigen::MatrixXd LP4_matrix
precursor laplacian term-4
Definition msrProblem.H:428
Eigen::MatrixXd LP7_matrix
precursor laplacian term-7
Definition msrProblem.H:434
Eigen::MatrixXd TM_matrix
temperature mass term
Definition msrProblem.H:483
autoPtr< volScalarField > _prec4
Definition msrProblem.H:143
std::vector< SPLINTER::RBFSpline * > rbfsplines_NSF
Definition msrProblem.H:546
label NCmodes
Definition msrProblem.H:506
autoPtr< dimensionedScalar > _K0
Definition msrProblem.H:89
autoPtr< volScalarField > _alphat
Definition msrProblem.H:150
Eigen::MatrixXd laplacian_prec(label NPrecmodes, label family)
List< Eigen::MatrixXd > SD3_matrix
decay heat stream term-3
Definition msrProblem.H:462
PtrList< volScalarField > Prec2modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:273
Eigen::MatrixXd LP6_matrix
precursor laplacian term-6
Definition msrProblem.H:432
autoPtr< volScalarField > _NSF0
Definition msrProblem.H:91
List< Eigen::MatrixXd > TS_matrix
temperature stream term
Definition msrProblem.H:485
Eigen::MatrixXd MF_matrix
mass flux
Definition msrProblem.H:364
label NPmodes
Definition msrProblem.H:501
autoPtr< dimensionedScalar > _beta2
Definition msrProblem.H:123
PtrList< volScalarField > Prec4field
List of pointers used to form the prec4 snapshots matrix.
Definition msrProblem.H:213
List< Eigen::MatrixXd > temp_heatsource(label NTmodes, label NDecmodes, label NCmodes, label decgroup)
autoPtr< dimensionedScalar > _Sc
Definition msrProblem.H:148
Eigen::MatrixXd MP3_matrix
precursor mass term-3
Definition msrProblem.H:410
PtrList< volScalarField > DFields
List of pointers used to form the D snapshosts matrix.
Definition msrProblem.H:246
Eigen::MatrixXd D_matrix
Laplacian term PPE.
Definition msrProblem.H:344
List< Eigen::MatrixXd > DFS2_matrix
decay heat flux source term-2
Definition msrProblem.H:478
void change_viscosity(double mu)
method to change the viscosity in UEqn
autoPtr< dimensionedScalar > _lam6
Definition msrProblem.H:119
List< Eigen::MatrixXd > laplacian_flux(label NFluxmodes, label NCmodes)
diffusion eq. methods:
Definition msrProblem.C:937
autoPtr< volScalarField > _flux0
Definition msrProblem.H:76
autoPtr< IOMRFZoneList > _MRF
Definition msrProblem.H:59
PtrList< volScalarField > AFields
List of pointers used to form the A snapshosts matrix.
Definition msrProblem.H:252
Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPmodes)
Definition msrProblem.C:546
autoPtr< dimensionedScalar > _beta3
Definition msrProblem.H:124
Eigen::MatrixXd PS7_matrix
prec_source 7
Definition msrProblem.H:385
autoPtr< volScalarField > _T0
Definition msrProblem.H:85
List< Eigen::MatrixXd > FS7_matrix
precursor flux source term-7
Definition msrProblem.H:450
List< Eigen::MatrixXd > ST6_matrix
precursor stream term-6
Definition msrProblem.H:400
void readMSRfields()
Method to read all the fieds of the MSR problem in the offline folder, it also reads mu_samples matri...
Eigen::MatrixXd LP1_matrix
precursor laplacian term-1
Definition msrProblem.H:422
autoPtr< dimensionedScalar > _TrefXS
Definition msrProblem.H:70
PtrList< volScalarField > Dec1modes
List of pointers used to form the dec1 modes.
Definition msrProblem.H:297
List< Eigen::MatrixXd > THS1_matrix
temperature decay heat source term-1
Definition msrProblem.H:491
PtrList< volScalarField > Dec3field
List of pointers used to form the dec3 snapshots matrix.
Definition msrProblem.H:237
autoPtr< dimensionedScalar > _Prt
Definition msrProblem.H:65
List< Eigen::MatrixXd > ST1_matrix
precursor stream term-1
Definition msrProblem.H:390
void liftSolve()
Perform a lift solve for the velocity field.
Definition msrProblem.C:291
autoPtr< dimensionedScalar > _IV1
Definition msrProblem.H:101
autoPtr< volVectorField > _U
Definition msrProblem.H:55
Eigen::MatrixXd LT_matrix
temperature laplacian term
Definition msrProblem.H:487
autoPtr< dimensionedScalar > _D1_0
Definition msrProblem.H:102
PtrList< volScalarField > Prec6modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:285
autoPtr< dimensionedScalar > _beta7
Definition msrProblem.H:128
autoPtr< volScalarField > _prec70
Definition msrProblem.H:83
void liftSolveT()
Perform a lift solve for the temperature.
autoPtr< dimensionedScalar > _lam4
Definition msrProblem.H:117
List< Eigen::MatrixXd > FS6_matrix
precursor flux source term-6
Definition msrProblem.H:448
Eigen::VectorXi NPrecmodes
Definition msrProblem.H:503
autoPtr< dimensionedScalar > _beta1
Definition msrProblem.H:122
autoPtr< volScalarField > _difft
Definition msrProblem.H:151
Eigen::MatrixXd mass_temp(label NTmodes)
temperature eq. methods
autoPtr< dimensionedScalar > _decLam3
Definition msrProblem.H:134
autoPtr< volScalarField > _prec5
Definition msrProblem.H:144
std::vector< SPLINTER::RBFSpline * > rbfsplines_v
Definition msrProblem.H:542
autoPtr< volScalarField > _TXS0
Definition msrProblem.H:95
PtrList< volScalarField > vmodes
List of pointers used to form the v modes.
Definition msrProblem.H:306
Eigen::MatrixXd P_matrix
Div of velocity.
Definition msrProblem.H:341
autoPtr< volScalarField > _v
Definition msrProblem.H:67
autoPtr< dimensionedScalar > _beta4
Definition msrProblem.H:125
autoPtr< fv::options > _fvOptions
fvOptions
Definition msrProblem.H:168
scalar tolerance
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and p...
Definition msrProblem.H:175
autoPtr< fvMesh > _mesh
Definition msrProblem.H:51
PtrList< volScalarField > Dec1field
List of pointers used to form the dec1 snapshots matrix.
Definition msrProblem.H:231
scalar maxIter
Number of maximum iterations to be done for the computation of the truth solution.
Definition msrProblem.H:177
autoPtr< dimensionedScalar > _alfa_NSF1
Definition msrProblem.H:110
PtrList< volScalarField > TXSmodes
List of pointers used to form the SP snapshosts matrix.
Definition msrProblem.H:321
autoPtr< dimensionedScalar > _lam5
Definition msrProblem.H:118
Eigen::MatrixXd divergence_term(label NUmodes, label NPmodes)
continuity eq. methods:
Definition msrProblem.C:672
Eigen::VectorXi NDecmodes
Definition msrProblem.H:504
autoPtr< dimensionedScalar > _lam1
Definition msrProblem.H:114
Eigen::MatrixXd dec_mass(label NDecmodes, label decgroup)
Eigen::MatrixXd LD1_matrix
decay heat laplacian term-1
Definition msrProblem.H:470
List< Eigen::MatrixXd > SD2_matrix
decay heat stream term-2
Definition msrProblem.H:460
autoPtr< simpleControl > _simple
simpleControl
Definition msrProblem.H:173
std::vector< SPLINTER::DataTable * > SAMPLES_SP
Definition msrProblem.H:538
autoPtr< dimensionedScalar > _Sct
Definition msrProblem.H:149
Eigen::MatrixXd MP5_matrix
precursor mass term-5
Definition msrProblem.H:414
autoPtr< volScalarField > _prec3
Definition msrProblem.H:142
PtrList< volScalarField > PowerDensfield
List of pointers used to form the powerDens snapshots matrix.
Definition msrProblem.H:240
bool homboolT
Definition msrProblem.H:523
autoPtr< volScalarField > _prec40
Definition msrProblem.H:80
autoPtr< volScalarField > _A0
Definition msrProblem.H:92
autoPtr< dimensionedScalar > _NSF1_0
Definition msrProblem.H:108
PtrList< volScalarField > Prec3modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:276
Eigen::MatrixXd LP2_matrix
precursor laplacian term-2
Definition msrProblem.H:424
autoPtr< volScalarField > _p0
Initial fields (for restart purposes)
Definition msrProblem.H:73
autoPtr< volVectorField > _U0
Definition msrProblem.H:74
Eigen::MatrixXd diffusive_term(label NUmodes, label NPmodes)
sub-functions needed by projectPPE
Definition msrProblem.C:507
autoPtr< dimensionedScalar > _decBeta3
Definition msrProblem.H:137
autoPtr< volScalarField > _prec6
Definition msrProblem.H:145
Eigen::MatrixXd prec_source(label NFluxmodes, label NPrecmodes, label family)
autoPtr< volScalarField > _D
Definition msrProblem.H:103
bool homboolU
boolean variables to check if the homogenization of U and T is performed (true) or not (false)
Definition msrProblem.H:522
autoPtr< volScalarField > _SP0
Definition msrProblem.H:94
Eigen::MatrixXd LD2_matrix
decay heat laplacian term-2
Definition msrProblem.H:472
autoPtr< volScalarField > _dec20
Definition msrProblem.H:87
Eigen::MatrixXd pressure_BC1(label NUmodes, label NPmodes)
Definition msrProblem.C:781
autoPtr< volScalarField > _prec60
Definition msrProblem.H:82
PtrList< volScalarField > Prec7field
List of pointers used to form the prec7 snapshots matrix.
Definition msrProblem.H:222
Eigen::MatrixXd PS8_matrix
prec_source 8
Definition msrProblem.H:387
autoPtr< dimensionedScalar > _alfa_A1
Definition msrProblem.H:107
Eigen::MatrixXd MP8_matrix
precursor mass term-8
Definition msrProblem.H:420
autoPtr< volScalarField > _prec7
Definition msrProblem.H:146
List< Eigen::MatrixXd > div_momentum(label NUmodes, label NPmodes)
Definition msrProblem.C:712
autoPtr< dimensionedScalar > _lam8
Definition msrProblem.H:121
label NFluxmodes
Definition msrProblem.H:502
autoPtr< dimensionedScalar > _nu
Definition msrProblem.H:63
autoPtr< volScalarField > _TXS
Definition msrProblem.H:71
autoPtr< dimensionedScalar > _decbetaTot
Definition msrProblem.H:138
autoPtr< volScalarField > _prec30
Definition msrProblem.H:79
autoPtr< volScalarField > _p
Definition msrProblem.H:60
Eigen::MatrixXd MD2_matrix
decay heat mass term-2
Definition msrProblem.H:466
Eigen::MatrixXd laplacian_pressure(label NPmodes)
Definition msrProblem.C:760
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of lifting functions.
Definition msrProblem.H:515
List< Eigen::MatrixXd > dec_fluxsource(label NFluxmodes, label NDecmodes, label NCmodes, label decgroup)
PtrList< volScalarField > Dec2modes
List of pointers used to form the dec2 modes.
Definition msrProblem.H:300
void restart()
method to set all fields back to values in 0 folder
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition msrProblem.H:513
Eigen::MatrixXd PS4_matrix
prec_source 4
Definition msrProblem.H:379
std::vector< SPLINTER::DataTable * > SAMPLES_TXS
Definition msrProblem.H:540
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition msrProblem.H:195
autoPtr< dimensionedScalar > _decBeta1
Definition msrProblem.H:135
autoPtr< dimensionedScalar > _decLam2
Definition msrProblem.H:133
Eigen::MatrixXd MP6_matrix
precursor mass term-6
Definition msrProblem.H:416
autoPtr< dimensionedScalar > _lam2
Definition msrProblem.H:115
PtrList< volScalarField > NSFFields
List of pointers used to form the NSF snapshosts matrix.
Definition msrProblem.H:249
PtrList< volScalarField > Prec3field
List of pointers used to form the prec3 snapshots matrix.
Definition msrProblem.H:210
autoPtr< dimensionedScalar > _decLam1
Definition msrProblem.H:132
PtrList< volVectorField > Umodes
List of pointers used to form the velocity modes.
Definition msrProblem.H:264
autoPtr< Time > _runTime
Definition msrProblem.H:47
PtrList< volScalarField > Dec3modes
List of pointers used to form the dec3 modes.
Definition msrProblem.H:303
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
Definition msrProblem.H:228
Eigen::MatrixXd B_matrix
Diffusion term.
Definition msrProblem.H:329
PtrList< volScalarField > Tomfield
List of pointers used to form the homogeneous velocity snapshots.
Definition msrProblem.H:519
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition msrProblem.H:198
PtrList< volScalarField > Prec1field
List of pointers used to form the prec1 snapshots matrix.
Definition msrProblem.H:204
std::vector< SPLINTER::RBFSpline * > rbfsplines_SP
Definition msrProblem.H:550
Eigen::MatrixXd laplacian_temp(label NTmodes)
autoPtr< dimensionedScalar > _alfa_D1
Definition msrProblem.H:104
List< Eigen::MatrixXd > pressure_BC2(label NUmodes, label NPmodes)
Definition msrProblem.C:828
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition msrProblem.H:332
Eigen::MatrixXd mass_flux(label NFluxmodes)
Definition msrProblem.C:969
Eigen::MatrixXd MP4_matrix
precursor mass term-4
Definition msrProblem.H:412
autoPtr< surfaceScalarField > _phi0
Definition msrProblem.H:75
autoPtr< volScalarField > _prec1
Definition msrProblem.H:140
PtrList< volScalarField > SPmodes
List of pointers used to form the SP snapshosts matrix.
Definition msrProblem.H:318
autoPtr< volScalarField > _prec10
Definition msrProblem.H:77
List< Eigen::MatrixXd > ST2_matrix
precursor stream term-2
Definition msrProblem.H:392
std::vector< SPLINTER::DataTable * > SAMPLES_D
Definition msrProblem.H:532
List< Eigen::MatrixXd > FS8_matrix
precursor flux source term-8
Definition msrProblem.H:452
PtrList< volScalarField > SPFields
List of pointers used to form the SP snapshosts matrix.
Definition msrProblem.H:255
List< Eigen::MatrixXd > convective_term(label NUmodes, label NPmodes)
Definition msrProblem.C:585
void msrcoeff(label &NC)
method to apply RBF interpolation procedure NC is the number of modes to adopt for construncting the ...
Eigen::MatrixXd PS1_matrix
prec_source 1
Definition msrProblem.H:373
List< Eigen::MatrixXd > DFS1_matrix
decay heat flux source term-1
Definition msrProblem.H:476
List< Eigen::MatrixXd > DFS3_matrix
decay heat flux source term-3
Definition msrProblem.H:480
autoPtr< dimensionedScalar > _A1_0
Definition msrProblem.H:105
List< Eigen::MatrixXd > flux_source(label NFluxmodes, label NPrecmodes, label NCmodes, label family)
Eigen::MatrixXd pressure_BC3(label NUmodes, label NPmodes)
Definition msrProblem.C:885
void homogenizeU()
Method to compute the homogenized velocity field, it also sets homboolU=true.
List< Eigen::MatrixXd > stream_dec(label NUmodes, label NDecmodes, label decgroup)
decay heat eq. methods:
List< Eigen::MatrixXd > FS5_matrix
precursor flux source term-5
Definition msrProblem.H:446
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NFluxmodes, Eigen::VectorXi Nprecmodes, label NTmodes, Eigen::VectorXi Ndecmodes, label NCmodes)
Project using the Poisson Equation for pressure.
Definition msrProblem.C:395
List< Eigen::MatrixXd > THS3_matrix
temperature decay heat source term-3
Definition msrProblem.H:495
List< Eigen::MatrixXd > temp_stream(label NUmodes, label NTmodes)
PtrList< volScalarField > Dec2field
List of pointers used to form the dec2 snapshots matrix.
Definition msrProblem.H:234
PtrList< volScalarField > Prec5field
List of pointers used to form the prec5 snapshots matrix.
Definition msrProblem.H:216
Eigen::MatrixXd LD3_matrix
decay heat laplacian term-3
Definition msrProblem.H:474
bool precInBool
boolean variable to decide if apply prec inlet BC
Definition msrProblem.H:182
List< Eigen::MatrixXd > C_matrix
Non linear term.
Definition msrProblem.H:338
PtrList< volScalarField > TXSFields
List of pointers used to form the SP snapshosts matrix.
Definition msrProblem.H:258
Eigen::MatrixXd LP5_matrix
precursor laplacian term-5
Definition msrProblem.H:430
void msrgetModesEVD()
Method to compute the modes for all the fields in the MSR problem, if hombool==false the velocity mod...
Definition msrProblem.C:223
Eigen::MatrixXd MD3_matrix
decay heat mass term-3
Definition msrProblem.H:468
autoPtr< dimensionedScalar > _alfa_SP1
Definition msrProblem.H:113
PtrList< volScalarField > Amodes
List of pointers used to form the A snapshosts matrix.
Definition msrProblem.H:315
autoPtr< volScalarField > _prec20
Definition msrProblem.H:78
Eigen::MatrixXd PS2_matrix
prec_source 2
Definition msrProblem.H:375
autoPtr< volScalarField > _D0
Definition msrProblem.H:93
PtrList< volScalarField > Fluxmodes
List of pointers used to form the flux modes.
Definition msrProblem.H:267
autoPtr< surfaceScalarField > _phi
Definition msrProblem.H:56
Eigen::MatrixXd laplacian_dec(label NDecmodes, label decgroup)
PtrList< volScalarField > Prec2field
List of pointers used to form the prec2 snapshots matrix.
Definition msrProblem.H:207
autoPtr< volScalarField > _prec80
Definition msrProblem.H:84
autoPtr< volScalarField > _logT
Definition msrProblem.H:160
PtrList< volScalarField > Dmodes
List of pointers used to form the D modes.
Definition msrProblem.H:309
void msrgetModesSVD()
Method to compute the modes for all the fields in the MSR problem, if hombool==false the velocity mod...
Definition msrProblem.C:151
Eigen::MatrixXd mass_term(label NUmodes, label NPmodes)
Definition msrProblem.C:633
autoPtr< volScalarField > _dec10
Definition msrProblem.H:86
List< Eigen::MatrixXd > temp_XSfluxsource(label NTmodes, label NFluxmodes, label NCmodes)
autoPtr< volScalarField > _flux
Definition msrProblem.H:131
autoPtr< volScalarField > _A
Definition msrProblem.H:106
autoPtr< volScalarField > _prec50
Definition msrProblem.H:81
autoPtr< dimensionedScalar > _beta6
Definition msrProblem.H:127
List< Eigen::MatrixXd > prod_flux(label NFluxmodes, label NCmodes)
Definition msrProblem.C:990
label NTmodes
Definition msrProblem.H:505
autoPtr< singlePhaseTransportModel > _laminarTransport
Definition msrProblem.H:57
std::vector< SPLINTER::RBFSpline * > rbfsplines_D
Definition msrProblem.H:544
std::vector< SPLINTER::RBFSpline * > rbfsplines_TXS
Definition msrProblem.H:552
List< Eigen::MatrixXd > ST3_matrix
precursor stream term-3
Definition msrProblem.H:394
autoPtr< volScalarField > _dec1
Definition msrProblem.H:161
Eigen::MatrixXd PS5_matrix
prec_source 5
Definition msrProblem.H:381
Eigen::MatrixXi inletIndexT
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.
void computeLift(T &Lfield, T &liftfield, T &omfield)
Homogenize the snapshot matrix, it works with PtrList of volVectorField and volScalarField.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
void computeLiftT(T &Lfield, T &liftfield, T &omfield)
Virtual function to compute the lifting function.
IOdictionary * ITHACAdict
dictionary to store input output infos
Eigen::MatrixXd mu
Row matrix of parameters.
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.
volScalarField & T
Definition createT.H:46
volScalarField & powerDens
fvScalarMatrix & TEqn
Definition TEqn.H:15
volScalarField & logT
volScalarField & SP
volScalarField NSF0(NSF)
volScalarField D0(D)
dimensionedScalar & decbetaTot
volScalarField & NSF
volScalarField A0(A)
volScalarField & A
dimensionedScalar & betaTot
volScalarField & D
volScalarField v0(v)
volScalarField SP0(SP)
volScalarField & v
volScalarField & TXS
volScalarField & difft
volScalarField & alphat
volScalarField TXS0(TXS)
volScalarField prec20(prec2)
volScalarField & prec5
volScalarField & prec6
volScalarField prec40(prec4)
volScalarField prec10(prec1)
volScalarField prec70(prec7)
volScalarField prec50(prec5)
volScalarField & prec8
volScalarField & flux
volScalarField & prec3
volScalarField prec60(prec6)
volScalarField flux0(flux)
volScalarField prec30(prec3)
volScalarField & prec2
volScalarField & prec7
volScalarField & prec4
volScalarField prec80(prec8)
volScalarField & prec1
volScalarField dec10(dec1)
volScalarField T0(T)
volScalarField & dec3
volScalarField dec20(dec2)
volScalarField & dec1
volScalarField & dec2
volScalarField dec30(dec3)
dimensionedScalar & Sc
dimensionedScalar & Tref
dimensionedScalar & rhoRef
dimensionedScalar & CpRef
dimensionedScalar & betaTE
dimensionedScalar & nu
dimensionedScalar & Prt
dimensionedScalar & Pr
dimensionedScalar & Sct
dimensionedScalar & beta1
dimensionedScalar & beta3
dimensionedScalar & D1_0
dimensionedScalar & A1_0
dimensionedScalar & beta2
dimensionedScalar & Keff
dimensionedScalar & lam2
dimensionedScalar & beta5
dimensionedScalar & beta4
dimensionedScalar & beta8
dimensionedScalar & lam6
dimensionedScalar & IV1
dimensionedScalar & lam4
dimensionedScalar & lam8
dimensionedScalar & alfa_NSF1
dimensionedScalar & lam1
dimensionedScalar & alfa_A1
dimensionedScalar & beta6
dimensionedScalar & NSF1_0
dimensionedScalar & TrefXS
dimensionedScalar & lam5
dimensionedScalar & lam7
dimensionedScalar & lam3
dimensionedScalar & beta7
dimensionedScalar K0(Keff)
dimensionedScalar & alfa_D1
dimensionedScalar & alfa_SP1
dimensionedScalar & SP1_0
dimensionedScalar & decLam1
dimensionedScalar & decBeta3
dimensionedScalar & decBeta1
dimensionedScalar & decBeta2
dimensionedScalar & decLam3
dimensionedScalar & decLam2
void getModes(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Computes the bases or reads them for a field.
Definition ITHACAPOD.C:93
void getModesSVD(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Gets the bases for a scalar field using SVD instead of the method of snapshots.
Definition ITHACAPOD.C:471
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 ...
List< Eigen::MatrixXd > readMatrix(word folder, word mat_name)
Read a three dimensional matrix from a txt file in Eigen format.
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.
Eigen::VectorXd getCoeffs(GeometricField< Type, PatchField, GeoMesh > &snapshot, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Projects a snapshot on a basis function and gets the coefficients of the projection.
bool check_pod()
Check if the POD data folder "./ITHACAoutput/POD" exists.
bool check_off()
Check if the offline data folder "./ITHACAoutput/Offline" exists.
bool check_folder(word folder)
Checks if a folder exists.
simpleControl simple(mesh)
surfaceScalarField & phi
surfaceScalarField phi0(phi)
volVectorField & U
volVectorField U0(U)
dimensionedVector & g
volScalarField & p
singlePhaseTransportModel & laminarTransport
volScalarField p0(p)
setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue)
adjustPhi(phiHbyA, U, p)
label i
Definition pEqn.H:46