Loading...
Searching...
No Matches
15MSR_cavity.C
Go to the documentation of this file.
1#include "usmsrProblem.H"
3#include "ITHACAstream.H"
4#include "LRSensitivity.H"
5#include "Tm.H"
6#include "Ptot.H"
7#include <Eigen/Dense>
8#include <math.h>
9#include <iomanip>
10#include <chrono>
11#include <ctime>
12#include "Tm_time.H"
13#include "Ptot_time.H"
14
15class msr : public usmsrProblem
16{
17 public:
18 explicit msr(int argc, char* argv[])
19 :
20 usmsrProblem(argc, argv),
21 U(_U()),
22 p(_p()),
23 flux(_flux()),
24 prec1(_prec1()),
25 prec2(_prec2()),
26 prec3(_prec3()),
27 prec4(_prec4()),
28 prec5(_prec5()),
29 prec6(_prec6()),
30 prec7(_prec7()),
31 prec8(_prec8()),
32 T(_T()),
33 dec1(_dec1()),
34 dec2(_dec2()),
35 dec3(_dec3()),
36 v(_v()),
37 D(_D()),
38 NSF(_NSF()),
39 A(_A()),
40 SP(_SP()),
41 TXS(_TXS())
42 {}
43
44 volVectorField& U;
45 volScalarField& p;
46 volScalarField& flux;
47 volScalarField& prec1;
48 volScalarField& prec2;
49 volScalarField& prec3;
50 volScalarField& prec4;
51 volScalarField& prec5;
52 volScalarField& prec6;
53 volScalarField& prec7;
54 volScalarField& prec8;
55 volScalarField& T;
56 volScalarField& dec1;
57 volScalarField& dec2;
58 volScalarField& dec3;
59 volScalarField& v;
60 volScalarField& D;
61 volScalarField& NSF;
62 volScalarField& A;
63 volScalarField& SP;
64 volScalarField& TXS;
65
67 int NUproj = para->ITHACAdict->lookupOrDefault<int>("NUproj", 2);
68 int NPproj = para->ITHACAdict->lookupOrDefault<int>("NPproj", 2);
69 int NFluxproj = para->ITHACAdict->lookupOrDefault<int>("NFluxproj", 2);
70 int NPrecproj1 = para->ITHACAdict->lookupOrDefault<int>("NPrecproj1", 2);
71 int NPrecproj2 = para->ITHACAdict->lookupOrDefault<int>("NPrecproj2", 2);
72 int NPrecproj3 = para->ITHACAdict->lookupOrDefault<int>("NPrecproj3", 2);
73 int NPrecproj4 = para->ITHACAdict->lookupOrDefault<int>("NPrecproj4", 2);
74 int NPrecproj5 = para->ITHACAdict->lookupOrDefault<int>("NPrecproj5", 2);
75 int NPrecproj6 = para->ITHACAdict->lookupOrDefault<int>("NPrecproj6", 2);
76 int NPrecproj7 = para->ITHACAdict->lookupOrDefault<int>("NPrecproj7", 2);
77 int NPrecproj8 = para->ITHACAdict->lookupOrDefault<int>("NPrecproj8", 2);
78 int NTproj = para->ITHACAdict->lookupOrDefault<int>("NTproj", 2);
79 int NDecproj1 = para->ITHACAdict->lookupOrDefault<int>("NDecproj1", 2);
80 int NDecproj2 = para->ITHACAdict->lookupOrDefault<int>("NDecproj2", 2);
81 int NDecproj3 = para->ITHACAdict->lookupOrDefault<int>("NDecproj3", 2);
82 int NCproj = para->ITHACAdict->lookupOrDefault<int>("NCproj", 2);
83
84 double r1 = _beta1().value() / _betaTot().value();
85 double r2 = _beta2().value() / _betaTot().value();
86 double r3 = _beta3().value() / _betaTot().value();
87 double r4 = _beta4().value() / _betaTot().value();
88 double r5 = _beta5().value() / _betaTot().value();
89 double r6 = _beta6().value() / _betaTot().value();
90 double r7 = _beta7().value() / _betaTot().value();
91 double r8 = _beta8().value() / _betaTot().value();
92
94 {
95 if (offline == false)
96 {
97 List<scalar> mu_now(3);
98
99 for (int i = 0; i < mu.cols(); i++)
100 {
101 _nu().value() = mu(0, i);
102 change_viscosity(mu(0, i));
103 _betaTot().value() = mu(1, i);
104 _beta1().value() = r1 * mu(1, i);
105 _beta2().value() = r2 * mu(1, i);
106 _beta3().value() = r3 * mu(1, i);
107 _beta4().value() = r4 * mu(1, i);
108 _beta5().value() = r5 * mu(1, i);
109 _beta6().value() = r6 * mu(1, i);
110 _beta7().value() = r7 * mu(1, i);
111 _beta8().value() = r8 * mu(1, i);
112 _decLam3().value() = mu(2, i);
113 mu_now[0] = mu(0, i);
114 mu_now[1] = mu(1, i);
115 mu_now[2] = mu(2, i);
116 truthSolve(mu_now);
117 restart();
118 }
119 }
120 else
121 {
123 }
124 }
125
126};
127
128class Tmlocal : public Tm
129{
130 public:
131 explicit Tmlocal(int argc, char* argv[], int Nsampled)
132 :
133 Tm(argc, argv, Nsampled),
134 T(_T())
135 {}
136
137 volScalarField& T;
138};
139
140class Ploct: public Ptot_time
141{
142 public:
143 explicit Ploct(int argc, char* argv[], int Nsampled)
144 :
145 Ptot_time(argc, argv, Nsampled),
147 {}
148
149 volScalarField& powerDens;
150};
151class Plocal: public Ptot
152{
153 public:
154 explicit Plocal(int argc, char* argv[], int Nsampled)
155 :
156 Ptot(argc, argv, Nsampled),
158 {}
159
160 volScalarField& powerDens;
161};
162class Tmloct: public Tm_time
163{
164 public:
165 explicit Tmloct(int argc, char* argv[], int Nsampled)
166 :
167 Tm_time(argc, argv, Nsampled),
168 T(_T())
169 {}
170
171 volScalarField& T;
172};
173
174int main(int argc, char* argv[])
175{
176 //object initialization
177 msr prova(argc, argv);
178 //inletIndex set
179 prova.inletIndex.resize(1, 2);
180 prova.inletIndex(0, 0) = 0;
181 prova.inletIndex(0, 1) = 0;
182 prova.inletIndexT.resize(4, 1);
183 prova.inletIndexT(0, 0) = 0;
184 prova.inletIndexT(1, 0) = 1;
185 prova.inletIndexT(2, 0) = 2;
186 prova.inletIndexT(3, 0) = 3;
187 //set the number of parameters and values of each
188 prova.Pnumber = 3;
189 prova.Tnumber = 31;
190 prova.setParameters();
191 int central = static_cast<int>(prova.Tnumber / 2);
192 double nu0 = 2.46E-06;
193 double betatot0 = 321.8E-05;
194 double dlam30 = 3.58e-04;
195 //+-10% from default values
196 prova.mu_range(0, 0) = 1.1 * nu0; //nu range
197 prova.mu_range(0, 1) = 0.9 * nu0;
198 prova.mu_range(1, 0) = 1.1 * betatot0; //betaTot range
199 prova.mu_range(1, 1) = 0.9 * betatot0;
200 prova.mu_range(2, 0) = 1.1 * dlam30; //decLam3 range
201 prova.mu_range(2, 1) = 0.9 * dlam30;
202 prova.genRandPar();
203 //override central column of mu with reference values of the parameters
204 prova.mu(0, central) = nu0;
205 prova.mu(1, central) = betatot0;
206 prova.mu(2, central) = dlam30;
207 double tstart = std::time(0);
208 Eigen::MatrixXd arange(prova.Pnumber, 2);
209 Eigen::MatrixXd mu_f = prova.mu;
210
211 if (prova.offline == false)
212 {
213 for (int i = 0; i < prova.Pnumber; i++)
214 {
215 arange(i, 0) = prova.mu.row(i).minCoeff();
216 arange(i, 1) = prova.mu.row(i).maxCoeff();
217 }
218
219 ITHACAstream::exportMatrix(arange, "range", "eigen", "./ITHACAoutput");
220 ITHACAstream::exportMatrix(prova.mu, "mu_off", "eigen", "./ITHACAoutput");
221 }
222 else
223 {
224 mu_f = ITHACAstream::readMatrix("./ITHACAoutput/mu_off_mat.txt");
225 arange = ITHACAstream::readMatrix("./ITHACAoutput/range_mat.txt");
226 }
227
228 //perform the offline step
229 prova.offlineSolve();
230 prova._nu().value() = nu0;
231 prova.change_viscosity(nu0);
232 double tend = std::time(0);
233 double deltat_offline;
234 deltat_offline = difftime(tend, tstart);
235
236 //save the time needed to perform offline stage
237 if (prova.offline == false)
238 {
239 std::ofstream diff_t("./ITHACAoutput/t_offline=" + name(deltat_offline));
240 }
241
242 std::cout << "Offline stage time= " << deltat_offline << std::endl;
243 //compute the liftfunction for the velocity
244 prova.liftSolve();
245 //compute the liftfuncion for the temperature
246 prova.liftSolveT();
247 //homogenize the velocity field
248 prova.homogenizeU();
249 //homogenize the temperature field
250 prova.homogenizeT();
251 //compute the spatial modes
252 prova.msrgetModesEVD();
253 // perform the projection
254 Eigen::VectorXi NPrec(8);
255 NPrec(0) = prova.NPrecproj1;
256 NPrec(1) = prova.NPrecproj2;
257 NPrec(2) = prova.NPrecproj3;
258 NPrec(3) = prova.NPrecproj4;
259 NPrec(4) = prova.NPrecproj5;
260 NPrec(5) = prova.NPrecproj6;
261 NPrec(6) = prova.NPrecproj7;
262 NPrec(7) = prova.NPrecproj8;
263 Eigen::VectorXi NDec(3);
264 NDec(0) = prova.NDecproj1;
265 NDec(1) = prova.NDecproj2;
266 NDec(2) = prova.NDecproj3;
267 prova.projectPPE("./Matrices", prova.NUproj, prova.NPproj, prova.NFluxproj,
268 NPrec, prova.NTproj, NDec, prova.NCproj);
269 //define the moving wall bc condition for U
270 double umedio = 2.46E-03;
271 //define the temperature at the moving wall
272 Eigen::MatrixXd Twall(4, 2);
273 Twall.setZero();
274 Twall(0, 0) = 850;
275 Twall(1, 0) = 950;
276 Twall(2, 0) = 950;
277 Twall(3, 0) = 950;
278 //define the ROM object and its time settings
279 reducedusMSR ridotto(prova);
280 ridotto.tstart = 0;
281 ridotto.finalTime = 10;
282 ridotto.dt = 0.1;
283 //vel_now for the ROM object
284 Eigen::MatrixXd vel_now(1, 1);
285 vel_now(0, 0) = umedio;
286 //temp_now for the ROM object
287 Eigen::MatrixXd temp_now = Twall;
288 // set the online value of mu equal to mu(:,central)
289 Eigen::VectorXd mu_on(3);
290 ridotto.nu = mu_f(0, central);
291 ridotto.btot = mu_f(1, central);
292 ridotto.b1 = prova.r1 * ridotto.btot;
293 ridotto.b2 = prova.r2 * ridotto.btot;
294 ridotto.b3 = prova.r3 * ridotto.btot;
295 ridotto.b4 = prova.r4 * ridotto.btot;
296 ridotto.b5 = prova.r5 * ridotto.btot;
297 ridotto.b6 = prova.r6 * ridotto.btot;
298 ridotto.b7 = prova.r7 * ridotto.btot;
299 ridotto.b8 = prova.r8 * ridotto.btot;
300 ridotto.dl3 = mu_f(2, central);
301 // index corresponding at central solution
302 int ref_start = 465;
303 mu_on(0) = ridotto.nu;
304 mu_on(1) = ridotto.btot;
305 mu_on(2) = ridotto.dl3;
306 ridotto.solveOnline(vel_now, temp_now, mu_on, ref_start);
307 ridotto.reconstructAP("./ITHACAoutput/Reconstruction/", 10);
308 //ridotto.reconstructAP("./ITHACAoutput/Reconstruction/",50);
309 //final snapshots for ROM and FOM respectively
310 int tf = 21;
311 int tf_fom = 31;
312 //Some post process calculations
313 int c = 0;
314 Eigen::MatrixXd err(10, tf);
315 Eigen::MatrixXd TmFOM(prova.Tnumber, tf);
316 Eigen::MatrixXd PtotFOM(prova.Tnumber, tf);
317 Eigen::MatrixXd TmROM(prova.Tnumber, tf);
318 Eigen::MatrixXd PtotROM(prova.Tnumber, tf);
319 std::cout << "Computing errors on reference case..." << std::endl;
320
321 for (int i = ref_start; i < ref_start + tf; i++)
322 {
323 err(0, c) = ITHACAutilities::errorL2Rel(prova.Ufield[i], ridotto.UREC[c]);
324 err(1, c) = ITHACAutilities::errorL2Rel(prova.Fluxfield[i],
325 ridotto.FLUXREC[c]);
326 err(2, c) = ITHACAutilities::errorL2Rel(prova.Prec1field[i],
327 ridotto.PREC1REC[c]);
328 err(3, c) = ITHACAutilities::errorL2Rel(prova.Prec4field[i],
329 ridotto.PREC4REC[c]);
330 err(4, c) = ITHACAutilities::errorL2Rel(prova.Prec8field[i],
331 ridotto.PREC8REC[c]);
332 err(5, c) = ITHACAutilities::errorL2Rel(prova.Tfield[i], ridotto.TREC[c]);
333 err(6, c) = ITHACAutilities::errorL2Rel(prova.Dec1field[i],
334 ridotto.DEC1REC[c]);
335 err(7, c) = ITHACAutilities::errorL2Rel(prova.Dec3field[i],
336 ridotto.DEC3REC[c]);
338 ridotto.POWERDENSREC[c]);
339 err(9, c) = ITHACAutilities::errorL2Rel(prova.TXSFields[i],
340 ridotto.TXSREC[c]);
341 c++;
342 }
343
344 std::cout << "End" << std::endl;
345 ITHACAstream::exportMatrix(err, "err", "eigen",
346 "./ITHACAoutput");
347 std::cout << "computing FOM output..." << std::endl;
348 double tmpt;
349 double tmpp;
350 int rig = 0;
351 int col = 0;
352
353 for (int i = 0; i < prova.Tfield.size(); i++)
354 {
355 tmpt = prova.Tfield[i].weightedAverage(prova._mesh().V()).value();
356 tmpp = fvc::domainIntegrate(prova.PowerDensfield[i]).value();
357
358 if (col < tf)
359 {
360 TmFOM(rig, col) = tmpt;
361 PtotFOM(rig, col) = tmpp;
362 col++;
363 }
364 else if (col == tf_fom)
365 {
366 col = 0;
367 rig++;
368 }
369 else
370 {
371 col++;
372 }
373 }
374
375 std::cout << "End" << std::endl;
376 ITHACAstream::exportMatrix(TmFOM, "fom_resT", "eigen",
377 "./ITHACAoutput");
378 ITHACAstream::exportMatrix(PtotFOM, "fom_resP", "eigen",
379 "./ITHACAoutput");
380 // perform online solve for all the values of mu_offline
381 std::string folder = {"./ITHACAoutput/ROMOutput/"};
382
383 for (int i = 0; i < prova.Tnumber; i++)
384 {
385 ridotto.nu = mu_f(0, i);
386 ridotto.btot = mu_f(1, i);
387 ridotto.b1 = prova.r1 * ridotto.btot;
388 ridotto.b2 = prova.r2 * ridotto.btot;
389 ridotto.b3 = prova.r3 * ridotto.btot;
390 ridotto.b4 = prova.r4 * ridotto.btot;
391 ridotto.b5 = prova.r5 * ridotto.btot;
392 ridotto.b6 = prova.r6 * ridotto.btot;
393 ridotto.b7 = prova.r7 * ridotto.btot;
394 ridotto.b8 = prova.r8 * ridotto.btot;
395 ridotto.dl3 = mu_f(2, i);
396 mu_on(0) = ridotto.nu;
397 mu_on(1) = ridotto.btot;
398 mu_on(2) = ridotto.dl3;
399 folder.append(std::to_string(i));
400 ridotto.solveOnline(vel_now, temp_now, mu_on, ref_start);
401 ridotto.reconstructAP(folder, 50);
402 folder = {"./ITHACAoutput/ROMOutput/"};
403 ridotto.clearFields();
404 }
405
406 //compute the figures of merit for every "snapshot" time
407 Tmloct Tmedia(argc, argv, prova.Tnumber);
408 Ploct Ptotale(argc, argv, prova.Tnumber);
409 c = 0;
410
411 for (int i = 0; i < tf; i++)
412 {
413 Tmedia.buildMO(folder, i);
414 Ptotale.buildMO(folder, i);
415 TmROM.col(c) = Tmedia.modelOutput;
416 PtotROM.col(c) = Ptotale.modelOutput;
417 c++;
418 }
419
420 ITHACAstream::exportMatrix(TmROM, "rom_resT", "eigen",
421 "./ITHACAoutput");
422 ITHACAstream::exportMatrix(PtotROM, "rom_resP", "eigen",
423 "./ITHACAoutput");
424 //Sensitivity Analysis, two figures of merit considered:
425 // average temperature and total power at the last time instant, i.e rom.finalTime
426 //define the folder where save the output
427 folder = {"./ITHACAoutput/Sensitivity/ModelOutput/"};
428 int Nparameters = prova.Pnumber;
429 std::vector<std::string> pdflist = {"normal", "normal", "normal"};
430 int samplingPoints = 1000;
431 LRSensitivity analisi(Nparameters, samplingPoints);
432 //define the training range, values sampled outside are rejected
433 analisi.trainingRange = arange;
434 //set the parameters of each distributions
435 Eigen::MatrixXd Mp(Nparameters, 2);
436 Mp(0, 0) = nu0;
437 Mp(0, 1) = nu0 * 0.1 / 3;
438 Mp(1, 0) = betatot0;
439 Mp(1, 1) = betatot0 * 0.1 / 3;
440 Mp(2, 0) = dlam30;
441 Mp(2, 1) = dlam30 * 0.1 / 3;
442 //build and export the sampling sets
443 analisi.buildSamplingSet(pdflist, Mp);
444 ITHACAstream::exportMatrix(analisi.MatX, "sampled", "eigen", "./ITHACAoutput");
445 //compute parameters statistics
446 analisi.getXstats();
447 //perform the simulations with the reduced object
448 int counter = 0;
449 tstart = std::time(0);
450 ridotto.clearFields();
451
452 for (int j = 0; j < samplingPoints; j++)
453 {
454 ridotto.nu = analisi.MatX(j, 0);
455 ridotto.btot = analisi.MatX(j, 1);
456 ridotto.b1 = prova.r1 * ridotto.btot;
457 ridotto.b2 = prova.r2 * ridotto.btot;
458 ridotto.b3 = prova.r3 * ridotto.btot;
459 ridotto.b4 = prova.r4 * ridotto.btot;
460 ridotto.b5 = prova.r5 * ridotto.btot;
461 ridotto.b6 = prova.r6 * ridotto.btot;
462 ridotto.b7 = prova.r7 * ridotto.btot;
463 ridotto.b8 = prova.r8 * ridotto.btot;
464 ridotto.dl3 = analisi.MatX(j, 2);
465 mu_on(0) = ridotto.nu;
466 mu_on(1) = ridotto.btot;
467 mu_on(2) = ridotto.dl3;
468 folder.append(std::to_string(counter));
469 ridotto.solveOnline(vel_now, temp_now, mu_on, ref_start);
470 ridotto.reconstructAP(folder, 1000);
471 folder = {"./ITHACAoutput/Sensitivity/ModelOutput/"};
472 ridotto.clearFields();
473 counter++;
474 }
475
476 tend = std::time(0);
477 int Ntot = samplingPoints;
478 double deltat_online = difftime(tend, tstart);
479 // save SA online time
480 std::ofstream diff_ton("./ITHACAoutput/t_onlineSA_" + name(Ntot) + "=" + name(
481 deltat_online));
482 //initialize the figure of merit objects
483 Tmlocal TmediaSA(argc, argv, Ntot);
484 Plocal PtotSA(argc, argv, Ntot);
485 //build the Model Output
486 TmediaSA.buildMO(folder);
487 PtotSA.buildMO(folder);
488 //assign FofM to LRSensitivity object, first figure of merit considered is the average temperature
489 analisi.M = autoPtr<FofM>(&TmediaSA);
490 //load the model output in the LRSensitivity object
491 analisi.load_output();
492 //compute output statistics
493 analisi.getYstat();
494 std::cout << "Mean value and variance of the output:" << std::endl;
495 std::cout << analisi.Ey << "\t" << analisi.Vy << std::endl;
496 std::cout << "-------------" << std::endl;
497 // compute linear regression coefficients
498 analisi.getBetas();
499 std::cout << "analisi regression coefficients:" << std::endl;
500 std::cout << analisi.betas << std::endl;
501 std::cout << "-------------" << std::endl;
502 // assess the quality of the linear regression, i.e. R^2
503 analisi.assessQuality();
504 std::cout << "quality: " << analisi.QI << std::endl;
505 // save the ROM output, linear output, statistics
506 Eigen::MatrixXd saveyout(samplingPoints, 1);
507 Eigen::MatrixXd saveymodel(samplingPoints, 1);
508 saveyout.col(0) = analisi.y;
509 saveymodel.col(0) = analisi.ylin;
510 ITHACAstream::exportMatrix(saveyout, "T_out", "eigen", "./ITHACAoutput");
511 ITHACAstream::exportMatrix(saveymodel, "T_model", "eigen", "./ITHACAoutput");
512 Eigen::MatrixXd coeffsT(4, 3);
513 coeffsT.setZero();
514 coeffsT(0, 0) = analisi.Ey;
515 coeffsT(0, 1) = analisi.Vy;
516 coeffsT(0, 2) = analisi.QI;
517 coeffsT.row(1) = analisi.EX;
518 coeffsT.row(2) = analisi.VX;
519 coeffsT.row(3) = analisi.betas;
520 ITHACAstream::exportMatrix(coeffsT, "coeffsT", "eigen", "./ITHACAoutput");
521 //assign FofM object to LRSensitivity object
522 analisi.M = autoPtr<FofM>(&PtotSA);
523 //load the model output in the LRSensitivity object, repeat all the previous step
524 // for the total power
525 analisi.load_output();
526 analisi.getYstat();
527 std::cout << "Mean value and variance of the output:" << std::endl;
528 std::cout << analisi.Ey << "\t" << analisi.Vy << std::endl;
529 std::cout << "-------------" << std::endl;
530 analisi.getBetas();
531 std::cout << "analisi regression coefficients:" << std::endl;
532 std::cout << analisi.betas << std::endl;
533 std::cout << "-------------" << std::endl;
534 analisi.assessQuality();
535 std::cout << "quality: " << analisi.QI << std::endl;
536 Eigen::MatrixXd saveyoutP(samplingPoints, 1);
537 Eigen::MatrixXd saveymodelP(samplingPoints, 1);
538 saveyoutP.col(0) = analisi.y;
539 saveymodelP.col(0) = analisi.ylin;
540 ITHACAstream::exportMatrix(saveyoutP, "P_out", "eigen", "./ITHACAoutput");
541 ITHACAstream::exportMatrix(saveymodelP, "P_model", "eigen", "./ITHACAoutput");
542 Eigen::MatrixXd coeffsP(4, 3);
543 coeffsP.setZero();
544 coeffsP(0, 0) = analisi.Ey;
545 coeffsP(0, 1) = analisi.Vy;
546 coeffsP(0, 2) = analisi.QI;
547 coeffsP.row(1) = analisi.EX;
548 coeffsP.row(2) = analisi.VX;
549 coeffsP.row(3) = analisi.betas;
550 ITHACAstream::exportMatrix(coeffsP, "coeffsP", "eigen", "./ITHACAoutput");
551 return 0;
552}
553
int main(int argc, char *argv[])
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Eigen::VectorXd modelOutput
Vector to store the desired output of the model, i.e. the figure of merit.
Definition FofM.H:49
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
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.
void getYstat()
Method to compute Ey and Vy, it sets Ydone=true.
void getXstats()
Method to compute EX and VX, it sets Xdone=true.
double Vy
Variance of the output.
Eigen::VectorXd y
Vectors to store the output of the model.
double QI
double to quantify the goodness of linear regression from 0 (no linearity) to 1 (perfect linearity) a...
Eigen::MatrixXd MatX
Matrices storing the independent variables' values.
Eigen::VectorXd VX
Vector to store the variances of the parameters.
double Ey
Mean value of the output.
void getBetas()
Method to compute SRCs, it sets bdone=true.
void load_output()
Method to load the output of FOM/ROM inside member LRsensitivity member y.
Eigen::VectorXd ylin
Vector to store the result of linear approximation.
Eigen::MatrixXd trainingRange
Matrix to store the range used in training/offline stage.
Eigen::VectorXd betas
Vector to store the standardized regression coefficient (SRC)
void assessQuality()
Method to compute the quality of linear regression approximation.
void buildSamplingSet(std::vector< std::string > &pdflist, Eigen::MatrixXd plist)
Method to build MatX, it internally calls ITHACAsampling::samplingMC pdflist is the list with the nam...
Eigen::VectorXd EX
Vector to store the meanvalues of the parameters.
autoPtr< FofM > M
Figure of merit object.
volScalarField & powerDens
Plocal(int argc, char *argv[], int Nsampled)
Ploct(int argc, char *argv[], int Nsampled)
volScalarField & powerDens
autoPtr< volScalarField > _powerDens
List of pointers to power density field.
Definition Ptot_time.H:51
void buildMO(std::string dir, label t)
Method that computes the total power for the snapshot number t (starting from 0), output are sought i...
Definition Ptot_time.C:17
Definition Ptot.H:42
void buildMO(std::string dir)
Method that computes the total power at the last time instant of the simulation, output are sought in...
Definition Ptot.C:17
autoPtr< volScalarField > _powerDens
List of pointers to power density field.
Definition Ptot.H:52
autoPtr< volScalarField > _T
List of pointers to temperature field.
Definition Tm_time.H:50
void buildMO(std::string dir, label t)
Method that computes the average temperature (integral average) for the snapshot number t (starting f...
Definition Tm_time.C:17
Definition Tm.H:42
void buildMO(std::string dir)
Method that computes the average temperature (integral average) at the last time instant of the simul...
Definition Tm.C:17
autoPtr< volScalarField > _T
List of pointers to temperature field.
Definition Tm.H:51
Tmlocal(int argc, char *argv[], int Nsampled)
volScalarField & T
volScalarField & T
Tmloct(int argc, char *argv[], int Nsampled)
autoPtr< dimensionedScalar > _betaTot
Definition msrProblem.H:130
autoPtr< volScalarField > _NSF
Definition msrProblem.H:109
void homogenizeT()
Method to compute the homogenized temperature field, it also sets homboolT=true.
autoPtr< volScalarField > _prec8
Definition msrProblem.H:147
autoPtr< volScalarField > _dec3
Definition msrProblem.H:163
autoPtr< volScalarField > _SP
Definition msrProblem.H:112
autoPtr< volScalarField > _T
Definition msrProblem.H:159
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
autoPtr< dimensionedScalar > _beta5
Definition msrProblem.H:126
autoPtr< volScalarField > _prec2
Definition msrProblem.H:141
PtrList< volScalarField > Fluxfield
List of pointers used to form the flux snapshots matrix.
Definition msrProblem.H:201
autoPtr< volScalarField > _prec4
Definition msrProblem.H:143
autoPtr< dimensionedScalar > _beta2
Definition msrProblem.H:123
PtrList< volScalarField > Prec4field
List of pointers used to form the prec4 snapshots matrix.
Definition msrProblem.H:213
void change_viscosity(double mu)
method to change the viscosity in UEqn
autoPtr< dimensionedScalar > _beta3
Definition msrProblem.H:124
void readMSRfields()
Method to read all the fieds of the MSR problem in the offline folder, it also reads mu_samples matri...
PtrList< volScalarField > Dec3field
List of pointers used to form the dec3 snapshots matrix.
Definition msrProblem.H:237
void liftSolve()
Perform a lift solve for the velocity field.
Definition msrProblem.C:291
autoPtr< volVectorField > _U
Definition msrProblem.H:55
autoPtr< dimensionedScalar > _beta7
Definition msrProblem.H:128
void liftSolveT()
Perform a lift solve for the temperature.
autoPtr< dimensionedScalar > _beta1
Definition msrProblem.H:122
autoPtr< dimensionedScalar > _decLam3
Definition msrProblem.H:134
autoPtr< volScalarField > _prec5
Definition msrProblem.H:144
autoPtr< volScalarField > _v
Definition msrProblem.H:67
autoPtr< dimensionedScalar > _beta4
Definition msrProblem.H:125
autoPtr< fvMesh > _mesh
Definition msrProblem.H:51
PtrList< volScalarField > Dec1field
List of pointers used to form the dec1 snapshots matrix.
Definition msrProblem.H:231
autoPtr< volScalarField > _prec3
Definition msrProblem.H:142
PtrList< volScalarField > PowerDensfield
List of pointers used to form the powerDens snapshots matrix.
Definition msrProblem.H:240
autoPtr< volScalarField > _prec6
Definition msrProblem.H:145
autoPtr< volScalarField > _D
Definition msrProblem.H:103
autoPtr< volScalarField > _prec7
Definition msrProblem.H:146
autoPtr< dimensionedScalar > _nu
Definition msrProblem.H:63
autoPtr< volScalarField > _TXS
Definition msrProblem.H:71
autoPtr< volScalarField > _p
Definition msrProblem.H:60
void restart()
method to set all fields back to values in 0 folder
autoPtr< Time > _runTime
Definition msrProblem.H:47
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
Definition msrProblem.H:228
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
autoPtr< volScalarField > _prec1
Definition msrProblem.H:140
void homogenizeU()
Method to compute the homogenized velocity field, it also sets homboolU=true.
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
PtrList< volScalarField > TXSFields
List of pointers used to form the SP snapshosts matrix.
Definition msrProblem.H:258
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
autoPtr< volScalarField > _flux
Definition msrProblem.H:131
autoPtr< volScalarField > _A
Definition msrProblem.H:106
autoPtr< dimensionedScalar > _beta6
Definition msrProblem.H:127
autoPtr< volScalarField > _dec1
Definition msrProblem.H:161
int NDecproj1
int NTproj
volScalarField & p
double r7
int NPrecproj8
volScalarField & prec3
ITHACAparameters * para
double r6
volScalarField & SP
volScalarField & prec4
volScalarField & prec6
volScalarField & dec1
volScalarField & prec1
int NCproj
volScalarField & v
void offlineSolve()
volScalarField & dec2
int NPrecproj2
volScalarField & A
volScalarField & TXS
int NDecproj2
int NPproj
int NDecproj3
int NPrecproj6
int NPrecproj4
int NFluxproj
int NPrecproj7
int NPrecproj1
volScalarField & prec7
volScalarField & dec3
volScalarField & D
double r2
volScalarField & flux
volScalarField & prec2
int NUproj
int NPrecproj5
int NPrecproj3
double r1
double r5
volScalarField & prec5
double r3
volVectorField & U
volScalarField & prec8
volScalarField & NSF
double r4
msr(int argc, char *argv[])
volScalarField & T
double r8
PtrList< volScalarField > DEC3REC
Definition ReducedMSR.H:320
scalar b4
Definition ReducedMSR.H:233
PtrList< volScalarField > FLUXREC
Definition ReducedMSR.H:308
scalar nu
characteristic constants of the problem
Definition ReducedMSR.H:211
PtrList< volScalarField > TXSREC
Definition ReducedMSR.H:327
PtrList< volScalarField > PREC4REC
Definition ReducedMSR.H:312
PtrList< volVectorField > UREC
Recontructed fields.
Definition ReducedMSR.H:306
scalar b2
Definition ReducedMSR.H:231
scalar b1
Definition ReducedMSR.H:230
scalar b8
Definition ReducedMSR.H:237
scalar b3
Definition ReducedMSR.H:232
PtrList< volScalarField > DEC1REC
Definition ReducedMSR.H:318
PtrList< volScalarField > PREC1REC
Definition ReducedMSR.H:309
scalar btot
Definition ReducedMSR.H:238
PtrList< volScalarField > POWERDENSREC
Definition ReducedMSR.H:321
PtrList< volScalarField > TREC
Definition ReducedMSR.H:317
void clearFields()
Method to clear all the fields of MSR (sets the size to zero)
scalar dl3
Definition ReducedMSR.H:243
scalar b6
Definition ReducedMSR.H:235
scalar b5
Definition ReducedMSR.H:234
scalar b7
Definition ReducedMSR.H:236
PtrList< volScalarField > PREC8REC
Definition ReducedMSR.H:316
void solveOnline(Eigen::MatrixXd vel_now, Eigen::MatrixXd temp_now, Eigen::VectorXd mu_online, int startSnap=0)
void reconstructAP(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Eigen::MatrixXi inletIndexT
label Pnumber
Number of parameters.
label Tnumber
Dimension of the training set (used only when gerating parameters without input)
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
void genRandPar()
Generate Random Numbers.
Eigen::MatrixXd mu
Row matrix of parameters.
Eigen::MatrixXd mu_range
Range of the parameter spaces.
void setParameters()
Set Parameters Problems.
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
void truthSolve()
Perform a TruthSolve.
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.
double errorL2Rel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in L2 norm.
singVal col(i)
label i
Definition pEqn.H:46