110 Eigen::VectorXd& fvec)
const
112 Eigen::VectorXd a_dot(
Nphi_u);
113 Eigen::VectorXd a_tmp(
Nphi_u);
114 Eigen::VectorXd b_tmp(
Nphi_p);
119 Eigen::MatrixXd cc(1, 1);
132 fvec(
i) = - M5(
i) + M1(
i) - cc(0, 0) - M2(
i);
135 for (
int j = 0; j <
Nphi_p; j++)
141 for (
int j = 0; j <
N_BC; j++)
143 fvec(j) = x(j) -
BC(j);
151 Eigen::MatrixXd& fjac)
const
153 Eigen::NumericalDiff<newton_unsteadyNST_sup> numDiff(*
this);
159 Eigen::VectorXd& fvect)
const
161 Eigen::VectorXd c_dot(
Nphi_t);
162 Eigen::VectorXd c_tmp(
Nphi_t);
166 Eigen::MatrixXd qq(1, 1);
175 fvect(
i) = -M8(
i) + M6(
i) - qq(0, 0);
178 for (
int j = 0; j <
N_BC_t; j++)
180 fvect(j) = t(j) -
BC_t(j);
186 Eigen::MatrixXd& fjact)
const
188 Eigen::NumericalDiff<newton_unsteadyNST_sup_t> numDiff(*
this);
189 numDiff.df(t, fjact);
195 Eigen::MatrixXd& temp_now,
int startSnap)
208 for (
int j = 0; j <
N_BC; j++)
213 for (
int j = 0; j <
N_BC_t; j++)
215 z(j) = temp_now(j, 0);
228 for (
int j = 0; j <
N_BC_t; j++)
233 for (
int j = 0; j <
N_BC; j++)
249 tmp_sol.col(0).tail(
y.rows()) =
y;
251 Eigen::MatrixXd tmp_solt(
Nphi_t + 1, 1);
253 tmp_solt.col(0).tail(
z.rows()) =
z;
258 Eigen::HybridNonLinearSolver<newton_unsteadyNST_sup_t> hnlst(
269 Eigen::VectorXd res(
y);
270 Eigen::VectorXd rest(
z);
275 for (
int j = 0; j <
N_BC; j++)
285 for (
int j = 0; j <
N_BC_t; j++)
287 z(j) = temp_now(j, 0);
295 " ##################" << std::endl;
296 Info <<
"Time = " <<
time << endl;
297 std::cout <<
"Solving for the parameter: " <<
vel_now << std::endl;
298 std::cout <<
"Solving for the parameter: " << temp_now << std::endl;
300 if (res.norm() < 1e-5)
302 std::cout << green <<
"|F(x)| = " << res.norm() <<
" - Minimun reached in " <<
303 hnls.iter <<
" iterations " << def << std::endl << std::endl;
307 std::cout << red <<
"|F(x)| = " << res.norm() <<
" - Minimun reached in " <<
308 hnls.iter <<
" iterations " << def << std::endl << std::endl;
311 if (rest.norm() < 1e-5)
313 std::cout << green <<
"|F(x)| = " << rest.norm() <<
" - Minimun reached in " <<
314 hnlst.iter <<
" iterations " << def << std::endl << std::endl;
318 std::cout << red <<
"|F(x)| = " << rest.norm() <<
" - Minimun reached in " <<
319 hnlst.iter <<
" iterations " << def << std::endl << std::endl;
324 tmp_sol.col(0).tail(
y.rows()) =
y;
336 tmp_solt.col(0).tail(
z.rows()) =
z;
352 "./ITHACAoutput/red_coeff");
354 "./ITHACAoutput/red_coeff");
356 "./ITHACAoutput/red_coeff_t");
358 "./ITHACAoutput/red_coeff_t");
373 if (counter == nextwrite)
375 volVectorField U_rec(
"U_rec",
Umodes[0] * 0);
377 for (
int j = 0; j <
Nphi_u; j++)
384 volScalarField P_rec(
"P_rec",
Pmodes[0] * 0);
386 for (
int j = 0; j <
Nphi_p; j++)
393 nextwrite += printevery;
395 UREC.append((U_rec).clone());
396 PREC.append((P_rec).clone());
415 if (counter == nextwrite)
419 for (
int j = 0; j <
Nphi_t; j++)
425 nextwrite += printevery;
Header file of the reducedUnsteadyNST class.
Class to change color to the output stream.
PtrList< volScalarField > PREC
Reconstructed pressure field.
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
PtrList< volScalarField > Pmodes
List of pointers to store the modes for pressure.
List< Eigen::MatrixXd > online_solution
List of Eigen matrices to store the online solution.
int Nphi_p
Number of pressure modes.
scalar nu
Reduced viscosity in case of parametrized viscosity.
int count_online_solve
Counter to count the online solutions.
Eigen::MatrixXd vel_now
Online inlet velocity vector.
PtrList< volScalarField > Psnapshots
List of pointers to store the snapshots for pressure.
PtrList< volVectorField > Usnapshots
List of pointers to store the snapshots for velocity.
PtrList< volVectorField > Umodes
List of pointers to store the modes for velocity.
int N_BC
Number of parametrized boundary conditions.
PtrList< volVectorField > UREC
Recontructed velocity field.
int Nphi_u
Number of velocity modes.
PtrList< volScalarField > TREC
Reconstructed temperature field.
scalar finalTime
Scalar to store the final time if the online simulation.
int N_BC_t
Number of parametrized boundary conditions related to temperature field.
PtrList< volScalarField > Tmodes
List of pointers to store the modes for temperature.
void reconstruct_supt(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique for t...
PtrList< volScalarField > T_rec
Reconstructed temperature field.
reducedUnsteadyNST()
Construct Null.
List< Eigen::MatrixXd > online_solutiont
List of Eigen matrices to store the online solution for temperature equation.
newton_unsteadyNST_sup_t newton_object_sup_t
Functor object to call the non linear solver sup. approach for temperature case.
scalar tstart
Scalar to store the final time if the online simulation.
unsteadyNST * problem
pointer to the FOM problem
int Nphi_t
Number of temperature modes.
void reconstruct_sup(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique.
scalar dt
Scalar to store the time increment.
newton_unsteadyNST_sup newton_object_sup
Functor object to call the non linear solver sup. approach for velocity case.
Eigen::VectorXd z
Vector to store the temperature solution during the Newton procedure.
PtrList< volScalarField > Tsnapshots
List of pointers to store the snapshots for temperature.
scalar time
Scalar to store the current time.
void solveOnline_sup(Eigen::MatrixXd &vel_now, Eigen::MatrixXd &temp_now, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
Eigen::MatrixXi inletIndexT
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
label NPmodes
Number of pressure modes used for the projection.
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
volVectorModes supmodes
List of pointers used to form the supremizer modes.
volVectorModes Umodes
List of pointers used to form the velocity modes.
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
label NUmodes
Number of velocity modes used for the projection.
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Eigen::MatrixXd B_matrix
Diffusion term.
label NSUPmodes
Number of supremizer modes used for the projection.
Eigen::MatrixXd K_matrix
Gradient of pressure matrix.
List< Eigen::MatrixXd > C_matrix
Non linear term.
Eigen::MatrixXd P_matrix
Div of velocity.
Eigen::MatrixXd M_matrix
Mass Matrix.
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Implementation of a parametrized full order unsteady NS problem weakly coupled with the energy equat...
Eigen::MatrixXd Y_matrix
Gradient of pressure matrix.
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of the temperature lifting functions.
label NTmodes
Number of temperature modes used for the projection.
Eigen::MatrixXd MT_matrix
Mass Matrix T.
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
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 ...
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.
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
Newton object for the resolution of the reduced problem using a supremizer approach.
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const