Loading...
Searching...
No Matches
msrProblem.H
1#ifndef msrProblem_H
2#define msrProblem_H
3#include "fvCFD.H"
4#include "singlePhaseTransportModel.H"
5#include "turbulentTransportModel.H"
6#include "simpleControl.H"
7#include "pisoControl.H"
8#include "fvOptions.H"
9#include "reductionProblem.H"
10#include "ITHACAstream.H"
11#include "ITHACAforces.H"
12#include "volFields.H"
13#include <iostream>
14#include "IOmanip.H"
15#include "IFstream.H"
16#include "primitiveFields.H"
17#include "FieldFields.H"
18#include "scalarMatrices.H"
19#include "SortableList.H"
20#include "volFieldsFwd.H"
21#include "forces.H"
22#include "forceCoeffs.H"
23#include "volFields.H"
24#include <fstream>
25#include <sstream>
26#include <vector>
27#include <string>
28#include <stdio.h>
29#include "ITHACAPOD.H"
30#include <math.h>
31
32
33
35class msrProblem: public reductionProblem
36{
37
38 public:
39 //Constructors
40 msrProblem();
42 msrProblem(int argc, char* argv[]);
43 ~msrProblem() {};
44
45 // Dummy variables to transform msrFoam into a class, createTime.H
46
47 autoPtr<Time> _runTime;
48
49 // Dummy variables to transform msrFoam into a class, createMesh.H
50
51 mutable autoPtr<fvMesh> _mesh;
52
53 // Dummy variables to transform msrFoam into a class, createFields.H
54
55 autoPtr<volVectorField> _U;
56 autoPtr<surfaceScalarField> _phi;
57 autoPtr<singlePhaseTransportModel> _laminarTransport;
58 autoPtr<incompressible::turbulenceModel> turbulence;
59 autoPtr<IOMRFZoneList> _MRF;
60 autoPtr<volScalarField> _p;
61 label pRefCell;
62 scalar pRefValue;
63 autoPtr<dimensionedScalar> _nu;
64 autoPtr<dimensionedScalar> _Pr;
65 autoPtr<dimensionedScalar> _Prt;
66 autoPtr<dimensionedScalar> _rhoRef;
67 autoPtr<volScalarField> _v;
68 autoPtr<dimensionedScalar> _CpRef;
69 autoPtr<dimensionedScalar> _Tref;
70 autoPtr<dimensionedScalar> _TrefXS;
71 autoPtr<volScalarField> _TXS;
73 autoPtr<volScalarField> _p0;
74 autoPtr<volVectorField> _U0;
75 autoPtr<surfaceScalarField> _phi0;
76 autoPtr<volScalarField> _flux0;
77 autoPtr<volScalarField> _prec10;
78 autoPtr<volScalarField> _prec20;
79 autoPtr<volScalarField> _prec30;
80 autoPtr<volScalarField> _prec40;
81 autoPtr<volScalarField> _prec50;
82 autoPtr<volScalarField> _prec60;
83 autoPtr<volScalarField> _prec70;
84 autoPtr<volScalarField> _prec80;
85 autoPtr<volScalarField> _T0;
86 autoPtr<volScalarField> _dec10;
87 autoPtr<volScalarField> _dec20;
88 autoPtr<volScalarField> _dec30;
89 autoPtr<dimensionedScalar> _K0;
90 autoPtr<volScalarField> _v0;
91 autoPtr<volScalarField> _NSF0;
92 autoPtr<volScalarField> _A0;
93 autoPtr<volScalarField> _D0;
94 autoPtr<volScalarField> _SP0;
95 autoPtr<volScalarField> _TXS0;
96
97
98 // Dummy variables to transform msrFoam into a class, createFields_neutronics.H
99
100 autoPtr<dimensionedScalar> _Keff;
101 autoPtr<dimensionedScalar> _IV1;
102 autoPtr<dimensionedScalar> _D1_0;
103 autoPtr<volScalarField> _D;
104 autoPtr<dimensionedScalar> _alfa_D1;
105 autoPtr<dimensionedScalar> _A1_0;
106 autoPtr<volScalarField> _A;
107 autoPtr<dimensionedScalar> _alfa_A1;
108 autoPtr<dimensionedScalar> _NSF1_0;
109 autoPtr<volScalarField> _NSF;
110 autoPtr<dimensionedScalar> _alfa_NSF1;
111 autoPtr<dimensionedScalar> _SP1_0;
112 autoPtr<volScalarField> _SP;
113 autoPtr<dimensionedScalar> _alfa_SP1;
114 autoPtr<dimensionedScalar> _lam1;
115 autoPtr<dimensionedScalar> _lam2;
116 autoPtr<dimensionedScalar> _lam3;
117 autoPtr<dimensionedScalar> _lam4;
118 autoPtr<dimensionedScalar> _lam5;
119 autoPtr<dimensionedScalar> _lam6;
120 autoPtr<dimensionedScalar> _lam7;
121 autoPtr<dimensionedScalar> _lam8;
122 autoPtr<dimensionedScalar> _beta1;
123 autoPtr<dimensionedScalar> _beta2;
124 autoPtr<dimensionedScalar> _beta3;
125 autoPtr<dimensionedScalar> _beta4;
126 autoPtr<dimensionedScalar> _beta5;
127 autoPtr<dimensionedScalar> _beta6;
128 autoPtr<dimensionedScalar> _beta7;
129 autoPtr<dimensionedScalar> _beta8;
130 autoPtr<dimensionedScalar> _betaTot;
131 autoPtr<volScalarField> _flux;
132 autoPtr<dimensionedScalar> _decLam1;
133 autoPtr<dimensionedScalar> _decLam2;
134 autoPtr<dimensionedScalar> _decLam3;
135 autoPtr<dimensionedScalar> _decBeta1;
136 autoPtr<dimensionedScalar> _decBeta2;
137 autoPtr<dimensionedScalar> _decBeta3;
138 autoPtr<dimensionedScalar> _decbetaTot;
139 autoPtr<dimensionedScalar> _betaTE;
140 autoPtr<volScalarField> _prec1;
141 autoPtr<volScalarField> _prec2;
142 autoPtr<volScalarField> _prec3;
143 autoPtr<volScalarField> _prec4;
144 autoPtr<volScalarField> _prec5;
145 autoPtr<volScalarField> _prec6;
146 autoPtr<volScalarField> _prec7;
147 autoPtr<volScalarField> _prec8;
148 autoPtr<dimensionedScalar> _Sc;
149 autoPtr<dimensionedScalar> _Sct;
150 autoPtr<volScalarField> _alphat;
151 autoPtr<volScalarField> _difft;
152
153 //recirculation time
154 double tau;
155
156
157 // Dummy variables to transform msrFoam into a class, createFields_T.H
158
159 autoPtr<volScalarField> _T;
160 autoPtr<volScalarField> _logT;
161 autoPtr<volScalarField> _dec1;
162 autoPtr<volScalarField> _dec2;
163 autoPtr<volScalarField> _dec3;
164
165 // Dummy variables to transform msrFoam into a class, createFvOptions.H
166
168 autoPtr<fv::options> _fvOptions;
169
170 // Dummy variables to transform msrFoam into a class, remaining in .C
171
173 autoPtr<simpleControl> _simple;
175 scalar tolerance;
177 scalar maxIter;
180
182 bool precInBool = false;
185 label precoutIndex;
187 Eigen::MatrixXd bc_prec;
188
189
190
191
192 // Fields containing the offline solution i.e. the snapshot matrix
193
195 PtrList<volScalarField> Pfield;
196
198 PtrList<volVectorField> Ufield;
199
201 PtrList<volScalarField> Fluxfield;
202
204 PtrList<volScalarField> Prec1field;
205
207 PtrList<volScalarField> Prec2field;
208
210 PtrList<volScalarField> Prec3field;
211
213 PtrList<volScalarField> Prec4field;
214
216 PtrList<volScalarField> Prec5field;
217
219 PtrList<volScalarField> Prec6field;
220
222 PtrList<volScalarField> Prec7field;
223
225 PtrList<volScalarField> Prec8field;
226
228 PtrList<volScalarField> Tfield;
229
231 PtrList<volScalarField> Dec1field;
232
234 PtrList<volScalarField> Dec2field;
235
237 PtrList<volScalarField> Dec3field;
238
240 PtrList<volScalarField> PowerDensfield;
241
243 PtrList<volScalarField> vFields;
244
246 PtrList<volScalarField> DFields;
247
249 PtrList<volScalarField> NSFFields;
250
252 PtrList<volScalarField> AFields;
253
255 PtrList<volScalarField> SPFields;
256
258 PtrList<volScalarField> TXSFields;
259
261 PtrList<volScalarField> Pmodes;
262
264 PtrList<volVectorField> Umodes;
265
267 PtrList<volScalarField> Fluxmodes;
268
270 PtrList<volScalarField> Prec1modes;
271
273 PtrList<volScalarField> Prec2modes;
274
276 PtrList<volScalarField> Prec3modes;
277
279 PtrList<volScalarField> Prec4modes;
280
282 PtrList<volScalarField> Prec5modes;
283
285 PtrList<volScalarField> Prec6modes;
286
288 PtrList<volScalarField> Prec7modes;
289
291 PtrList<volScalarField> Prec8modes;
292
294 PtrList<volScalarField> Tmodes;
295
297 PtrList<volScalarField> Dec1modes;
298
300 PtrList<volScalarField> Dec2modes;
301
303 PtrList<volScalarField> Dec3modes;
304
306 PtrList<volScalarField> vmodes;
307
309 PtrList<volScalarField> Dmodes;
310
312 PtrList<volScalarField> NSFmodes;
313
315 PtrList<volScalarField> Amodes;
316
318 PtrList<volScalarField> SPmodes;
319
321 PtrList<volScalarField> TXSmodes;
322
323 //----------------------------------------------------------------------//
324 // Coefficient matrices:
325
326 //Fluid-dynamics matrices:
327
329 Eigen::MatrixXd B_matrix;
330
332 Eigen::MatrixXd M_matrix;
333
335 Eigen::MatrixXd K_matrix;
336
338 List <Eigen::MatrixXd> C_matrix;
339
341 Eigen::MatrixXd P_matrix;
342
344 Eigen::MatrixXd D_matrix;
345
347 List <Eigen::MatrixXd> G_matrix;
348
350 Eigen::MatrixXd BC1_matrix;
351
353 List <Eigen::MatrixXd> BC2_matrix;
354
356 Eigen::MatrixXd BC3_matrix;
357
358
359 //Neutronics matrices:
360
362 List<Eigen::MatrixXd> LF_matrix;
364 Eigen::MatrixXd MF_matrix;
366 List<Eigen::MatrixXd> PF_matrix;
368 List<Eigen::MatrixXd> AF_matrix;
369
370 ITHACAparameters* para;
371
373 Eigen::MatrixXd PS1_matrix;
375 Eigen::MatrixXd PS2_matrix;
377 Eigen::MatrixXd PS3_matrix;
379 Eigen::MatrixXd PS4_matrix;
381 Eigen::MatrixXd PS5_matrix;
383 Eigen::MatrixXd PS6_matrix;
385 Eigen::MatrixXd PS7_matrix;
387 Eigen::MatrixXd PS8_matrix;
388
390 List<Eigen::MatrixXd> ST1_matrix;
392 List<Eigen::MatrixXd> ST2_matrix;
394 List<Eigen::MatrixXd> ST3_matrix;
396 List<Eigen::MatrixXd> ST4_matrix;
398 List<Eigen::MatrixXd> ST5_matrix;
400 List<Eigen::MatrixXd> ST6_matrix;
402 List<Eigen::MatrixXd> ST7_matrix;
404 List<Eigen::MatrixXd> ST8_matrix;
406 Eigen::MatrixXd MP1_matrix;
408 Eigen::MatrixXd MP2_matrix;
410 Eigen::MatrixXd MP3_matrix;
412 Eigen::MatrixXd MP4_matrix;
414 Eigen::MatrixXd MP5_matrix;
416 Eigen::MatrixXd MP6_matrix;
418 Eigen::MatrixXd MP7_matrix;
420 Eigen::MatrixXd MP8_matrix;
422 Eigen::MatrixXd LP1_matrix;
424 Eigen::MatrixXd LP2_matrix;
426 Eigen::MatrixXd LP3_matrix;
428 Eigen::MatrixXd LP4_matrix;
430 Eigen::MatrixXd LP5_matrix;
432 Eigen::MatrixXd LP6_matrix;
434 Eigen::MatrixXd LP7_matrix;
436 Eigen::MatrixXd LP8_matrix;
438 List<Eigen::MatrixXd> FS1_matrix;
440 List<Eigen::MatrixXd> FS2_matrix;
442 List<Eigen::MatrixXd> FS3_matrix;
444 List<Eigen::MatrixXd> FS4_matrix;
446 List<Eigen::MatrixXd> FS5_matrix;
448 List<Eigen::MatrixXd> FS6_matrix;
450 List<Eigen::MatrixXd> FS7_matrix;
452 List<Eigen::MatrixXd> FS8_matrix;
453
454
455 //Energy matrices:
456
458 List<Eigen::MatrixXd> SD1_matrix;
460 List<Eigen::MatrixXd> SD2_matrix;
462 List<Eigen::MatrixXd> SD3_matrix;
464 Eigen::MatrixXd MD1_matrix;
466 Eigen::MatrixXd MD2_matrix;
468 Eigen::MatrixXd MD3_matrix;
470 Eigen::MatrixXd LD1_matrix;
472 Eigen::MatrixXd LD2_matrix;
474 Eigen::MatrixXd LD3_matrix;
476 List<Eigen::MatrixXd> DFS1_matrix;
478 List<Eigen::MatrixXd> DFS2_matrix;
480 List<Eigen::MatrixXd> DFS3_matrix;
481
483 Eigen::MatrixXd TM_matrix;
485 List<Eigen::MatrixXd> TS_matrix;
487 Eigen::MatrixXd LT_matrix;
489 List<Eigen::MatrixXd> TXS_matrix;
491 List<Eigen::MatrixXd> THS1_matrix;
493 List<Eigen::MatrixXd> THS2_matrix;
495 List<Eigen::MatrixXd> THS3_matrix;
496 //--------------------------------------------------------------------------------//
497
499
500 label NUmodes;
501 label NPmodes;
502 label NFluxmodes;
503 Eigen::VectorXi NPrecmodes;
504 Eigen::VectorXi NDecmodes;
505 label NTmodes;
506 label NCmodes;
507
508
509 //--------------------------------------------------------------------------------//
510 //Lift-functions and homogenization related terms:
511
513 PtrList<volVectorField> liftfield;
515 PtrList<volScalarField> liftfieldT;
517 PtrList<volVectorField> Uomfield;
519 PtrList<volScalarField> Tomfield;
522 bool homboolU = false;
523 bool homboolT = false;
524
525 //-----------------------------------------------------------------------------//
526 //RBF interpolation procedure quantities, i.e. constants changing with
527 //temperature
528
529 // Create a SAMPLES for interpolation
530 std::vector<SPLINTER::DataTable*> SAMPLES_v;
531 // Create a SAMPLES for interpolation
532 std::vector<SPLINTER::DataTable*> SAMPLES_D;
533 // Create a SAMPLES for interpolation
534 std::vector<SPLINTER::DataTable*> SAMPLES_NSF;
535 // Create a SAMPLES for interpolation
536 std::vector<SPLINTER::DataTable*> SAMPLES_A;
537 // Create a SAMPLES for interpolation
538 std::vector<SPLINTER::DataTable*> SAMPLES_SP;
539 // Create a SAMPLES for interpolation
540 std::vector<SPLINTER::DataTable*> SAMPLES_TXS;
541 // Create a RBF splines for interpolation
542 std::vector<SPLINTER::RBFSpline*> rbfsplines_v;
543 // Create a RBF splines for interpolation
544 std::vector<SPLINTER::RBFSpline*> rbfsplines_D;
545 // Create a RBF splines for interpolation
546 std::vector<SPLINTER::RBFSpline*> rbfsplines_NSF;
547 // Create a RBF splines for interpolation
548 std::vector<SPLINTER::RBFSpline*> rbfsplines_A;
549 // Create a RBF splines for interpolation
550 std::vector<SPLINTER::RBFSpline*> rbfsplines_SP;
551 // Create a RBF splines for interpolation
552 std::vector<SPLINTER::RBFSpline*> rbfsplines_TXS;
553
554 //-----------------------------------------------------------------------------------------//
555 // Methods:
556
558 void truthSolve(List<scalar> mu_now);
559 //----------------------------------------------------------------------------------
560
563 void msrcoeff(label& NC);
564 //----------------------------------------------------------------------------------
565
567 void liftSolve();
568 //----------------------------------------------------------------------------------
569
571 void liftSolveT();
572 //----------------------------------------------------------------------------------
573
578 void readMSRfields();
579 void readMSRfields(std::string& dir);
580 //-----------------------------------------------------------------------------------
581
583 void homogenizeU();
584 //-----------------------------------------------------------------------------------
585
587 void homogenizeT();
588 //--------------------------------------------------------------------------
589
595 void msrgetModesSVD();
596 //--------------------------------------------------------------------------
597
603 void msrgetModesEVD();
604 //--------------------------------------------------------------------------
605
616 void projectPPE(fileName folder, label NUmodes, label NPmodes, label NFluxmodes,
617 Eigen::VectorXi Nprecmodes, label NTmodes, Eigen::VectorXi Ndecmodes,
618 label NCmodes);
619 //--------------------------------------------------------------------------
620
622
624 Eigen::MatrixXd diffusive_term(label NUmodes, label NPmodes);
625 Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPmodes);
626 List < Eigen::MatrixXd > convective_term(label NUmodes, label NPmodes);
627 Eigen::MatrixXd mass_term(label NUmodes, label NPmodes);
628
630 Eigen::MatrixXd divergence_term(label NUmodes, label NPmodes);
631 List < Eigen::MatrixXd > div_momentum(label NUmodes, label NPmodes);
632 Eigen::MatrixXd laplacian_pressure(label NPmodes);
633 Eigen::MatrixXd pressure_BC1(label NUmodes, label NPmodes);
634 List < Eigen::MatrixXd > pressure_BC2(label NUmodes, label NPmodes);
635 Eigen::MatrixXd pressure_BC3(label NUmodes, label NPmodes);
636
638 List<Eigen::MatrixXd> laplacian_flux(label NFluxmodes, label NCmodes);
639 Eigen::MatrixXd mass_flux(label NFluxmodes);
640 List<Eigen::MatrixXd> prod_flux(label NFluxmodes, label NCmodes);
641 List<Eigen::MatrixXd> abs_flux(label NFluxmodes, label NCmodes);
642 Eigen::MatrixXd prec_source(label NFluxmodes, label NPrecmodes, label family);
643
645 List<Eigen::MatrixXd> stream_term(label NUmodes, label NPrecmodes,
646 label family);
647 Eigen::MatrixXd prec_mass(label NPrecmodes, label family);
648 Eigen::MatrixXd laplacian_prec(label NPrecmodes, label family);
649 List<Eigen::MatrixXd> flux_source(label NFluxmodes, label NPrecmodes,
650 label NCmodes, label family);
651
653 List<Eigen::MatrixXd> stream_dec(label NUmodes, label NDecmodes,
654 label decgroup);
655 Eigen::MatrixXd dec_mass(label NDecmodes, label decgroup);
656 Eigen::MatrixXd laplacian_dec(label NDecmodes, label decgroup);
657 List<Eigen::MatrixXd> dec_fluxsource(label NFluxmodes, label NDecmodes,
658 label NCmodes, label decgroup);
659
661 Eigen::MatrixXd mass_temp(label NTmodes);
662 List<Eigen::MatrixXd> temp_stream(label NUmodes, label NTmodes);
663 Eigen::MatrixXd laplacian_temp(label NTmodes);
664 List<Eigen::MatrixXd> temp_heatsource(label NTmodes, label NDecmodes,
665 label NCmodes, label decgroup);
666 List<Eigen::MatrixXd> temp_XSfluxsource(label NTmodes, label NFluxmodes,
667 label NCmodes);
668
669 //--------------------------------------------------------------------------
671 void change_viscosity(double mu);
672
673 //--------------------------------------------------------------------------
675 void restart();
676
677
678
679 protected:
680
681 //--------------------------------------------------------------------------
686
687 PtrList<volScalarField> choose_group(string field, label ith);
688
689 //--------------------------------------------------------------------------
695
696 template<typename M>
697 void savegroupMatrix(string nome, label n, word folder, M matrice);
698
702
703
704};
705
706#endif
707
708
709
710
711
712
713
714
715
716
717
Header file of the ITHACAPOD class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
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
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.
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
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
List< Eigen::MatrixXd > SD1_matrix
decay heat stream term-1
Definition msrProblem.H:458
Eigen::MatrixXd BC1_matrix
PPE BC1.
Definition msrProblem.H:350
PtrList< volScalarField > NSFmodes
List of pointers used to form the NSF snapshosts matrix.
Definition msrProblem.H:312
PtrList< volScalarField > Prec5modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:282
List< Eigen::MatrixXd > PF_matrix
production flux
Definition msrProblem.H:366
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 ...
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
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
Eigen::MatrixXd MP1_matrix
precursor mass term-1
Definition msrProblem.H:406
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
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
Eigen::MatrixXd PS3_matrix
prec_source 3
Definition msrProblem.H:377
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
List< Eigen::MatrixXd > G_matrix
Divergence of momentum PPE.
Definition msrProblem.H:347
List< Eigen::MatrixXd > AF_matrix
absorption flux
Definition msrProblem.H:368
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
PtrList< volScalarField > Prec8modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:291
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
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
List< Eigen::MatrixXd > TS_matrix
temperature stream term
Definition msrProblem.H:485
Eigen::MatrixXd MF_matrix
mass flux
Definition msrProblem.H:364
PtrList< volScalarField > Prec4field
List of pointers used to form the prec4 snapshots matrix.
Definition msrProblem.H:213
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
List< Eigen::MatrixXd > laplacian_flux(label NFluxmodes, label NCmodes)
diffusion eq. methods:
Definition msrProblem.C:938
PtrList< volScalarField > AFields
List of pointers used to form the A snapshosts matrix.
Definition msrProblem.H:252
Eigen::MatrixXd PS7_matrix
prec_source 7
Definition msrProblem.H:385
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
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
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
Eigen::MatrixXd LT_matrix
temperature laplacian term
Definition msrProblem.H:487
PtrList< volScalarField > Prec6modes
List of pointers used to form the prec1 modes.
Definition msrProblem.H:285
void liftSolveT()
Perform a lift solve for the temperature.
List< Eigen::MatrixXd > FS6_matrix
precursor flux source term-6
Definition msrProblem.H:448
Eigen::MatrixXd mass_temp(label NTmodes)
temperature eq. methods
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
scalar cumulativeContErr
continuity error
Definition msrProblem.H:179
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
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
PtrList< volScalarField > TXSmodes
List of pointers used to form the SP snapshosts matrix.
Definition msrProblem.H:321
Eigen::MatrixXd divergence_term(label NUmodes, label NPmodes)
continuity eq. methods:
Definition msrProblem.C:672
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
Eigen::MatrixXd MP5_matrix
precursor mass term-5
Definition msrProblem.H:414
PtrList< volScalarField > PowerDensfield
List of pointers used to form the powerDens snapshots matrix.
Definition msrProblem.H:240
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
Eigen::MatrixXd diffusive_term(label NUmodes, label NPmodes)
sub-functions needed by projectPPE
Definition msrProblem.C:507
bool homboolU
boolean variables to check if the homogenization of U and T is performed (true) or not (false)
Definition msrProblem.H:522
Eigen::MatrixXd LD2_matrix
decay heat laplacian term-2
Definition msrProblem.H:472
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
Eigen::MatrixXd MP8_matrix
precursor mass term-8
Definition msrProblem.H:420
Eigen::MatrixXd MD2_matrix
decay heat mass term-2
Definition msrProblem.H:466
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of lifting functions.
Definition msrProblem.H:515
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
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition msrProblem.H:195
Eigen::MatrixXd MP6_matrix
precursor mass term-6
Definition msrProblem.H:416
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
PtrList< volVectorField > Umodes
List of pointers used to form the velocity modes.
Definition msrProblem.H:264
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
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition msrProblem.H:332
Eigen::MatrixXd MP4_matrix
precursor mass term-4
Definition msrProblem.H:412
PtrList< volScalarField > SPmodes
List of pointers used to form the SP snapshosts matrix.
Definition msrProblem.H:318
List< Eigen::MatrixXd > ST2_matrix
precursor stream term-2
Definition msrProblem.H:392
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
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
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
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
PtrList< volScalarField > Amodes
List of pointers used to form the A snapshosts matrix.
Definition msrProblem.H:315
Eigen::MatrixXd PS2_matrix
prec_source 2
Definition msrProblem.H:375
PtrList< volScalarField > Fluxmodes
List of pointers used to form the flux modes.
Definition msrProblem.H:267
PtrList< volScalarField > Prec2field
List of pointers used to form the prec2 snapshots matrix.
Definition msrProblem.H:207
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 bc_prec
matrix to store the values of precursors BC inlet conditions
Definition msrProblem.H:187
List< Eigen::MatrixXd > ST3_matrix
precursor stream term-3
Definition msrProblem.H:394
Eigen::MatrixXd PS5_matrix
prec_source 5
Definition msrProblem.H:381
reductionProblem()
Construct Null.
Eigen::MatrixXd mu
Row matrix of parameters.
void truthSolve()
Perform a TruthSolve.
Header file of the reductionProblem class.