39#define _USE_MATH_DEFINES
74 List<scalar> mu_now(9);
77 for (label
i = 0;
i <
mu.rows();
i++)
79 for (label j = 0; j <
mu.cols() ; j++)
95 volScalarField
yPos =
T.mesh().C().component(vector::Y).ref();
96 volScalarField
xPos =
T.mesh().C().component(vector::X).ref();
108 volScalarField nu1(
nu);
109 volScalarField nu2(
nu);
110 volScalarField nu3(
nu);
111 volScalarField nu4(
nu);
112 volScalarField nu5(
nu);
113 volScalarField nu6(
nu);
114 volScalarField nu7(
nu);
115 volScalarField nu8(
nu);
116 volScalarField nu9(
nu);
117 Eigen::MatrixXd Box1(2, 3);
118 Box1 << 0, 0, 0, 0.3, 0.3, 0.1;
119 Eigen::MatrixXd Box2(2, 3);
120 Box2 << 0.3, 0, 0, 0.6, 0.3, 0.1;
121 Eigen::MatrixXd Box3(2, 3);
122 Box3 << 0.6, 0, 0, 0.91, 0.3, 0.1;
123 Eigen::MatrixXd Box4(2, 3);
124 Box4 << 0, 0.3, 0, 0.3, 0.6, 0.1;
125 Eigen::MatrixXd Box5(2, 3);
126 Box5 << 0.3, 0.3, 0, 0.6, 0.6, 0.1;
127 Eigen::MatrixXd Box6(2, 3);
128 Box6 << 0.6, 0.3, 0, 0.91, 0.6, 0.1;
129 Eigen::MatrixXd Box7(2, 3);
130 Box7 << 0, 0.6, 0, 0.3, 0.91, 0.1;
131 Eigen::MatrixXd Box8(2, 3);
132 Box8 << 0.3, 0.61, 0, 0.6, 0.91, 0.1;
133 Eigen::MatrixXd Box9(2, 3);
134 Box9 << 0.6, 0.6, 0, 0.9, 0.91, 0.1;
169int main(
int argc,
char* argv[])
173 argList::addOption(
"stage",
"offline",
"Perform offline stage");
174 argList::addOption(
"stage",
"online",
"Perform online stage");
176 HashTable<string> validParOptions;
187 Pstream::addValidParOptions(validParOptions);
193 if (example.
_args().get(
"stage").match(
"offline"))
199 else if (example.
_args().get(
"stage").match(
"online"))
208 Info <<
"Pass '-stage offline', '-stage online'" << endl;
215 Info <<
"Pass 'offline' or 'online' as first arguments."
221 int argc_proc = argc - 1;
222 char* argv_proc[argc_proc];
223 argv_proc[0] = argv[0];
227 std::copy(argv + 2, argv + argc, argv_proc + 1);
236 if (std::strcmp(argv[1],
"offline") == 0)
242 else if (std::strcmp(argv[1],
"online") == 0)
251 Info <<
"Pass offline, online" << endl;
263 int NmodesTout = para->
ITHACAdict->lookupOrDefault<
int>(
"NmodesTout", 15);
264 int NmodesTproj = para->
ITHACAdict->lookupOrDefault<
int>(
"NmodesTproj", 10);
272 example.
mu_range.col(0) = Eigen::MatrixXd::Ones(9, 1) * 0.001;
273 example.
mu_range.col(1) = Eigen::MatrixXd::Ones(9, 1) * 0.1;
277 example.
theta.resize(9);
299 FOM_test.
mu_range.col(0) = Eigen::MatrixXd::Ones(9, 1) * 0.001;
300 FOM_test.
mu_range.col(1) = Eigen::MatrixXd::Ones(9, 1) * 0.1;
304 FOM_test.
theta.resize(9);
321 for (
int i = 0;
i < FOM_test.
mu.rows();
i++)
327 reduced.
reconstruct(
"./ITHACAoutput/Reconstruction");
int main(int argc, char *argv[])
void online_stage(tutorial02 &example, tutorial02 &FOM_test)
void offline_stage(tutorial02 &example, tutorial02 &FOM_test)
forAll(example_CG.gList, solutionI)
Header file of the ITHACAPOD class.
Header file of the ITHACAutilities namespace.
Header file of the reducedLaplacian class.
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.
Class to implement a full order laplacian parametrized problem.
PtrList< volScalarField > nu_list
Nu (diffusivity)
List< scalar > theta
Theta (coefficients of the affine expansion)
PtrList< fvScalarMatrix > operator_list
List of operators.
autoPtr< volScalarField > _T
Temperature field.
void project(label Nmodes)
Perform a projection onto the POD modes.
autoPtr< volScalarField > _nu
Diffusivity.
autoPtr< fvMesh > _mesh
Mesh.
volScalarModes Tmodes
List of POD modes.
autoPtr< volScalarField > _S
Source Term.
PtrList< volScalarField > Tfield
List of snapshots for the solution.
autoPtr< Time > _runTime
Time.
Class to solve the online reduced problem associated with a Laplace problem.
void reconstruct(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Function to recover the solution given the online solution.
void solveOnline(Eigen::MatrixXd mu)
Function to perform an online solve given a certain mu.
PtrList< volScalarField > Trec
Reduced reconstructed Solution.
label Pnumber
Number of parameters.
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 assignIF(T &s, G &value)
Assign internal field condition.
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.
autoPtr< argList > _args
argList
void setParameters()
Set Parameters Problems.
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
Class where the tutorial number 2 is implemented.
volScalarField & S
Source term field.
void offlineSolve(word folder="./ITHACAoutput/Offline/")
It perform an offline Solve.
volScalarField & nu
Diffusivity field.
void compute_nu()
Compute the diffusivity in each subdomain.
void assemble_operator()
Construct the operator_list where each term of the affine decomposition is stored.
tutorial02(int argc, char *argv[])
volScalarField & T
[tutorial02] Temperature field
void SetSource()
Define the source term function.
Header file of the laplacianProblem class.
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.
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.
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.
void setBoxToValue(GeometricField< Type, fvPatchField, volMesh > &field, Eigen::MatrixXd Box, Type value)
Set value of a volScalarField to a constant inside a given box.