55 for (
int k = 0; k <
problem->liftfield.size(); k++)
60 for (
int k = 0; k <
problem->NUmodes; k++)
65 for (
int k = 0; k <
problem->NPmodes; k++)
70 for (
int k = 0; k <
problem->NFluxmodes; k++)
115 for (
int k = 0; k <
problem->liftfieldT.size(); k++)
120 for (
int k = 0; k <
problem->NTmodes; k++)
151 for (
int k = 0; k <
problem->Ufield.size(); k++)
186 Eigen::VectorXd& fvec)
const
188 Eigen::VectorXd a_tmp(
Nphi_u);
189 Eigen::VectorXd b_tmp(
Nphi_p);
194 Eigen::MatrixXd cc(1, 1);
195 Eigen::MatrixXd gg(1, 1);
196 Eigen::MatrixXd bb(1, 1);
198 Eigen::VectorXd M1 =
problem->B_matrix * a_tmp *
nu;
200 Eigen::VectorXd M2 =
problem->K_matrix * b_tmp;
202 Eigen::VectorXd M3 =
problem->D_matrix * b_tmp;
204 Eigen::VectorXd M6 =
problem->BC1_matrix * a_tmp *
nu;
206 Eigen::VectorXd M7 =
problem->BC3_matrix * a_tmp *
nu;
210 cc = a_tmp.transpose() *
problem->C_matrix[
i] * a_tmp;
211 fvec(
i) = M1(
i) - cc(0, 0) - M2(
i);
219 gg = a_tmp.transpose() *
problem->G_matrix[
i] * a_tmp;
221 fvec(k) = M3(
i, 0) + gg(0, 0) - M7(
i, 0);
224 for (
int j = 0; j <
N_BC; j++)
226 fvec(j) = x(j) -
BC(j);
234 Eigen::NumericalDiff<newton_msr_fd> numDiff(*
this);
240 Eigen::VectorXd& fvecn)
const
270 Eigen::MatrixXd lf(1, 1);
272 Eigen::MatrixXd pf(1, 1);
274 Eigen::MatrixXd af(1, 1);
276 Eigen::VectorXd F3_1 =
problem->PS1_matrix * d1_tmp *
l1;
277 Eigen::VectorXd F3_2 =
problem->PS2_matrix * d2_tmp *
l2;
278 Eigen::VectorXd F3_3 =
problem->PS3_matrix * d3_tmp *
l3;
279 Eigen::VectorXd F3_4 =
problem->PS4_matrix * d4_tmp *
l4;
280 Eigen::VectorXd F3_5 =
problem->PS5_matrix * d5_tmp *
l5;
281 Eigen::VectorXd F3_6 =
problem->PS6_matrix * d6_tmp *
l6;
282 Eigen::VectorXd F3_7 =
problem->PS7_matrix * d7_tmp *
l7;
283 Eigen::VectorXd F3_8 =
problem->PS8_matrix * d8_tmp *
l8;
285 Eigen::MatrixXd pp1(1, 1);
286 Eigen::MatrixXd pp2(1, 1);
287 Eigen::MatrixXd pp3(1, 1);
288 Eigen::MatrixXd pp4(1, 1);
289 Eigen::MatrixXd pp5(1, 1);
290 Eigen::MatrixXd pp6(1, 1);
291 Eigen::MatrixXd pp7(1, 1);
292 Eigen::MatrixXd pp8(1, 1);
294 Eigen::VectorXd P1_1 =
problem->LP1_matrix * d1_tmp * (
nu /
Sc);
295 Eigen::VectorXd P1_2 =
problem->LP2_matrix * d2_tmp * (
nu /
Sc);
296 Eigen::VectorXd P1_3 =
problem->LP3_matrix * d3_tmp * (
nu /
Sc);
297 Eigen::VectorXd P1_4 =
problem->LP4_matrix * d4_tmp * (
nu /
Sc);
298 Eigen::VectorXd P1_5 =
problem->LP5_matrix * d5_tmp * (
nu /
Sc);
299 Eigen::VectorXd P1_6 =
problem->LP6_matrix * d6_tmp * (
nu /
Sc);
300 Eigen::VectorXd P1_7 =
problem->LP7_matrix * d7_tmp * (
nu /
Sc);
301 Eigen::VectorXd P1_8 =
problem->LP8_matrix * d8_tmp * (
nu /
Sc);
303 Eigen::VectorXd P2_1 =
problem->MP1_matrix * d1_tmp *
l1;
304 Eigen::VectorXd P2_2 =
problem->MP2_matrix * d2_tmp *
l2;
305 Eigen::VectorXd P2_3 =
problem->MP3_matrix * d3_tmp *
l3;
306 Eigen::VectorXd P2_4 =
problem->MP4_matrix * d4_tmp *
l4;
307 Eigen::VectorXd P2_5 =
problem->MP5_matrix * d5_tmp *
l5;
308 Eigen::VectorXd P2_6 =
problem->MP6_matrix * d6_tmp *
l6;
309 Eigen::VectorXd P2_7 =
problem->MP7_matrix * d7_tmp *
l7;
310 Eigen::VectorXd P2_8 =
problem->MP8_matrix * d8_tmp *
l8;
312 Eigen::MatrixXd fs1(1, 1);
313 Eigen::MatrixXd fs2(1, 1);
314 Eigen::MatrixXd fs3(1, 1);
315 Eigen::MatrixXd fs4(1, 1);
316 Eigen::MatrixXd fs5(1, 1);
317 Eigen::MatrixXd fs6(1, 1);
318 Eigen::MatrixXd fs7(1, 1);
319 Eigen::MatrixXd fs8(1, 1);
323 lf =
d_c.transpose() *
problem->LF_matrix[
i] * c_tmp;
325 af =
a_c.transpose() *
problem->AF_matrix[
i] * c_tmp;
326 fvecn(
i) = lf(0, 0) + pf(0, 0) - af(0,
327 0) + F3_1(
i) + F3_2(
i) + F3_3(
i) + F3_4(
i) + F3_5(
i) + F3_6(
i) + F3_7(
i) + F3_8(
338 fvecn(k) = -pp1(0, 0) + P1_1(
i) - P2_1(
i) + fs1(0, 0);
348 fvecn(k) = -pp2(0, 0) + P1_2(
i) - P2_2(
i) + fs2(0, 0);
358 fvecn(k) = -pp3(0, 0) + P1_3(
i) - P2_3(
i) + fs3(0, 0);
368 fvecn(k) = -pp4(0, 0) + P1_4(
i) - P2_4(
i) + fs4(0, 0);
378 fvecn(k) = -pp5(0, 0) + P1_5(
i) - P2_5(
i) + fs5(0, 0);
388 fvecn(k) = -pp6(0, 0) + P1_6(
i) - P2_6(
i) + fs6(0, 0);
398 fvecn(k) = -pp7(0, 0) + P1_7(
i) - P2_7(
i) + fs7(0, 0);
408 fvecn(k) = -pp8(0, 0) + P1_8(
i) - P2_8(
i) + fs8(0, 0);
415 Eigen::MatrixXd& fjacn)
const
417 Eigen::NumericalDiff<newton_msr_n> numDiff(*
this);
418 numDiff.df(n, fjacn);
423 Eigen::VectorXd& fvect)
const
425 Eigen::VectorXd e_tmp(
Nphi_T);
438 Eigen::MatrixXd tt(1, 1);
440 Eigen::VectorXd T1 =
problem->LT_matrix * e_tmp *
nu /
Pr;
442 Eigen::MatrixXd xsf(1, 1);
444 Eigen::MatrixXd dhs1(1, 1);
445 Eigen::MatrixXd dhs2(1, 1);
446 Eigen::MatrixXd dhs3(1, 1);
448 Eigen::MatrixXd dh1(1, 1);
449 Eigen::MatrixXd dh2(1, 1);
450 Eigen::MatrixXd dh3(1, 1);
452 Eigen::VectorXd DH1_1 =
problem->LD1_matrix * f1_tmp *
nu /
Sc;
453 Eigen::VectorXd DH1_2 =
problem->LD2_matrix * f2_tmp *
nu /
Sc;
454 Eigen::VectorXd DH1_3 =
problem->LD3_matrix * f3_tmp *
nu /
Sc;
456 Eigen::VectorXd DH2_1 =
problem->MD1_matrix * f1_tmp *
dl1;
457 Eigen::VectorXd DH2_2 =
problem->MD2_matrix * f2_tmp *
dl2;
458 Eigen::VectorXd DH2_3 =
problem->MD3_matrix * f3_tmp *
dl3;
460 Eigen::MatrixXd dfs1(1, 1);
461 Eigen::MatrixXd dfs2(1, 1);
462 Eigen::MatrixXd dfs3(1, 1);
471 fvect(
i) = -tt(0, 0) + T1(
i) + xsf(0, 0) + dhs1(0, 0) + dhs2(0, 0) + dhs3(0, 0);
481 fvect(k) = -dh1(0, 0) + DH1_1(
i) - DH2_1(
i) + dfs1(0, 0);
491 fvect(k) = -dh2(0, 0) + DH1_2(
i) - DH2_2(
i) + dfs2(0, 0);
501 fvect(k) = -dh3(0, 0) + DH1_3(
i) - DH2_3(
i) + dfs3(0, 0);
513 Eigen::MatrixXd& fjact)
const
515 Eigen::NumericalDiff<newton_msr_t> numDiff(*
this);
516 numDiff.df(t, fjact);
522 Eigen::VectorXd mu_online)
524 Info <<
"\n Starting online stage...\n" << endl;
533 for (
int j = 0; j <
N_BC; j++)
535 y(j) = vel_now(j, 0);
538 for (
int j = 0; j <
N_BCt; j++)
540 z(j) = temp_now(j, 0);
569 for (
int j = 0; j <
N_BC; j++)
605 for (
int j = 0; j <
N_BCt; j++)
611 Eigen::VectorXd res_fd(
y);
614 Eigen::VectorXd res_n(
w);
618 Eigen::VectorXd res_t(
z);
623 " ##################" << std::endl;
625 if (res_fd.norm() /
y.norm() < 1e-5)
627 std::cout << green <<
"|F_fd(x)| = " << res_fd.norm() /
y.norm() <<
628 " - Minimun reached in " << hnls_fd.iter <<
" iterations " << def << std::endl
633 std::cout << red <<
"|F_fd(x)| = " << res_fd.norm() /
y.norm() <<
634 " - Minimun reached in " << hnls_fd.iter <<
" iterations " << def << std::endl
638 if (res_n.norm() /
w.norm() < 1e-5)
640 std::cout << green <<
"|F_n(x)| = " << res_n.norm() /
w.norm() <<
641 " - Minimun reached in " << hnls_n.iter <<
" iterations " << def << std::endl <<
646 std::cout << red <<
"|F_n(x)| = " << res_n.norm() /
w.norm() <<
647 " - Minimun reached in " << hnls_n.iter <<
" iterations " << def << std::endl <<
651 if (res_t.norm() /
z.norm() < 1e-5)
653 std::cout << green <<
"|F_t(x)| = " << res_t.norm() /
z.norm() <<
654 " - Minimun reached in " << hnls_t.iter <<
" iterations " << def << std::endl <<
659 std::cout << red <<
"|F_t(x)| = " << res_t.norm() /
z.norm() <<
660 " - Minimun reached in " << hnls_t.iter <<
" iterations " << def << std::endl <<
684 "./ITHACAoutput/red_coeff_fd");
686 "./ITHACAoutput/red_coeff_n");
688 "./ITHACAoutput/red_coeff_t");
690 "./ITHACAoutput/red_coeff_C");
714 Info <<
"Reconstructing online solution | fluid-dynamics" << endl;
721 if (counter == nextwrite)
723 volVectorField U_rec(
"U",
Umodes[0] * 0);
725 for (
int j = 0; j <
Nphi_u; j++)
731 volScalarField P_rec(
"p",
Pmodes[0] * 0);
733 for (
int j = 0; j <
Nphi_p; j++)
739 nextwrite += printevery;
741 UREC.append(U_rec.clone());
742 PREC.append(P_rec.clone());
748 Info <<
"End" << endl;
759 Info <<
"Reconstructing online solution | neutronics" << endl;
766 if (counter == nextwrite)
768 volScalarField Flux_rec(
"flux",
Fluxmodes[0] * 0);
777 volScalarField Prec1_rec(
"prec1",
Prec1modes[0] * 0);
786 volScalarField Prec2_rec(
"prec2",
Prec2modes[0] * 0);
795 volScalarField Prec3_rec(
"prec3",
Prec3modes[0] * 0);
804 volScalarField Prec4_rec(
"prec4",
Prec4modes[0] * 0);
813 volScalarField Prec5_rec(
"prec5",
Prec5modes[0] * 0);
822 volScalarField Prec6_rec(
"prec6",
Prec6modes[0] * 0);
831 volScalarField Prec7_rec(
"prec7",
Prec7modes[0] * 0);
840 volScalarField Prec8_rec(
"prec8",
Prec8modes[0] * 0);
848 nextwrite += printevery;
850 FLUXREC.append(Flux_rec.clone());
864 Info <<
"End" << endl;
875 Info <<
"Reconstructing online solution | thermal" << endl;
879 dimensionedScalar
decLam1(
"decLam1", dimensionSet(0, 0, -1, 0, 0, 0, 0),
dl1);
880 dimensionedScalar
decLam2(
"decLam2", dimensionSet(0, 0, -1, 0, 0, 0, 0),
dl2);
881 dimensionedScalar
decLam3(
"decLam3", dimensionSet(0, 0, -1, 0, 0, 0, 0),
dl3);
885 if (counter == nextwrite)
887 volScalarField T_rec(
"T",
Tmodes[0] * 0);
889 for (
int j = 0; j <
Nphi_T; j++)
897 volScalarField Dec1_rec(
"dec1",
Dec1modes[0] * 0);
907 volScalarField Dec2_rec(
"dec2",
Dec2modes[0] * 0);
917 volScalarField Dec3_rec(
"dec3",
Dec3modes[0] * 0);
928 nextwrite += printevery;
930 TREC.append(T_rec.clone());
931 DEC1REC.append(Dec1_rec.clone());
932 DEC2REC.append(Dec2_rec.clone());
933 DEC3REC.append(Dec3_rec.clone());
940 Info <<
"End" << endl;
951 Info <<
"Reconstructing temperature changing constants" << endl;
958 if (counter == nextwrite)
960 volScalarField v_rec(
"v",
vmodes[0] * 0);
969 volScalarField D_rec(
"D",
Dmodes[0] * 0);
978 volScalarField NSF_rec(
"NSF",
NSFmodes[0] * 0);
987 volScalarField A_rec(
"A",
Amodes[0] * 0);
996 volScalarField SP_rec(
"SP",
SPmodes[0] * 0);
1005 volScalarField TXS_rec(
"TXS",
TXSmodes[0] * 0);
1013 std::ofstream of(folder +
"/" + name(counter2) +
"/" + name(
1015 nextwrite += printevery;
1017 vREC.append(v_rec.clone());
1018 DREC.append(D_rec.clone());
1019 NSFREC.append(NSF_rec.clone());
1020 AREC.append(A_rec.clone());
1021 SPREC.append(SP_rec.clone());
1022 TXSREC.append(TXS_rec.clone());
1028 Info <<
"End" << endl;
Class to change color to the output stream.
Class to implement Molten Salt Reactor multiphysics problem.
PtrList< volScalarField > NSFREC
PtrList< volScalarField > Tsnapshots
PtrList< volScalarField > DEC3REC
msrProblem * problem
Pointer to the FOM problem.
PtrList< volScalarField > vREC
void reconstruct_n(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct neutronics
newton_msr_n newton_object_n
PtrList< volScalarField > DEC2REC
PtrList< volScalarField > Dec1modes
PtrList< volScalarField > Pmodes
PtrList< volVectorField > Usnapshots
List of pointers to store the snapshots for each field.
PtrList< volScalarField > Psnapshots
PtrList< volScalarField > FLUXREC
PtrList< volVectorField > Umodes
List of pointers to store the modes for each field.
PtrList< volScalarField > Prec8snapshots
PtrList< volScalarField > Dec2modes
scalar nu
characteristic constants of the problem
PtrList< volScalarField > TXSREC
bool recall
boolean variable to check if the user wants to reconstruct all the three physics of the problem
PtrList< volScalarField > PREC4REC
PtrList< volScalarField > Prec2modes
void reconstruct_C(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct temperature dependent constants
PtrList< volScalarField > Prec6snapshots
PtrList< volScalarField > PREC3REC
PtrList< volScalarField > vsnapshots
PtrList< volVectorField > UREC
Recontructed fields.
List< Eigen::MatrixXd > online_solution_n
PtrList< volScalarField > Prec7modes
PtrList< volScalarField > Prec6modes
PtrList< volScalarField > SPmodes
void reconstruct_t(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct thermal fields
PtrList< volScalarField > Prec5snapshots
PtrList< volScalarField > TXSsnapshots
PtrList< volScalarField > Dec3snapshots
PtrList< volScalarField > Prec3modes
PtrList< volScalarField > Prec3snapshots
PtrList< volScalarField > Prec7snapshots
newton_msr_fd newton_object_fd
Newton object used to solve the non linear problem.
PtrList< volScalarField > Amodes
void loadConstants(msrProblem *problem)
Method to load all the constants needed in the ROM from ///the FOM.
PtrList< volScalarField > DEC1REC
PtrList< volScalarField > NSFmodes
PtrList< volScalarField > Fluxsnapshots
void reconstruct_fd(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct fluid-dynamics
int N_BC
Number of parametrized boundary conditions.
PtrList< volScalarField > Dec1snapshots
PtrList< volScalarField > SPsnapshots
PtrList< volScalarField > Asnapshots
PtrList< volScalarField > Prec4modes
PtrList< volScalarField > PREC1REC
PtrList< volScalarField > Prec1snapshots
PtrList< volScalarField > Dsnapshots
PtrList< volScalarField > Prec8modes
List< Eigen::MatrixXd > online_solution_fd
List of Eigen matrices to store the online solution.
PtrList< volScalarField > POWERDENSREC
PtrList< volScalarField > PREC6REC
PtrList< volScalarField > TREC
PtrList< volScalarField > NSFsnapshots
PtrList< volScalarField > TXSmodes
PtrList< volScalarField > Dec3modes
int count_online_solve
Counter to count the online solutions.
PtrList< volScalarField > AREC
PtrList< volScalarField > vmodes
PtrList< volScalarField > Prec4snapshots
PtrList< volScalarField > PREC2REC
PtrList< volScalarField > Tmodes
PtrList< volScalarField > Dmodes
int Nphi_u
Number of modes for each field.
List< Eigen::MatrixXd > online_solution_t
void reconstructAP(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Methods to reconstruct a solution from an online solve with a PPE stabilisation technique.
PtrList< volScalarField > Prec5modes
PtrList< volScalarField > Prec2snapshots
newton_msr_t newton_object_t
PtrList< volScalarField > PREC
PtrList< volScalarField > SPREC
void clearFields()
Method to clear all the fields of MSR (sets the size to zero)
PtrList< volScalarField > Dec2snapshots
PtrList< volScalarField > Prec1modes
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
PtrList< volScalarField > PREC5REC
PtrList< volScalarField > Fluxmodes
PtrList< volScalarField > PREC8REC
PtrList< volScalarField > DREC
List< Eigen::MatrixXd > online_solution_C
PtrList< volScalarField > PREC7REC
virtual void solveOnline()
Virtual Method to perform and online Solve.
dimensionedScalar & decLam1
dimensionedScalar & decLam3
dimensionedScalar & decLam2
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
Structure to implement a newton object for a stationary MSR problem.
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &n, Eigen::MatrixXd &fjacn) const
int operator()(const Eigen::VectorXd &n, Eigen::VectorXd &fvecn) const
int df(const Eigen::VectorXd &n, Eigen::MatrixXd &fjacn) const
int operator()(const Eigen::VectorXd &n, Eigen::VectorXd &fvecn) const