A general class for the implementation of a full order parametrized problem. More...
#include <reductionProblem.H>
Public Member Functions | |
reductionProblem () | |
Construct Null. | |
~reductionProblem () | |
void | setParameters () |
Set Parameters Problems. | |
void | genRandPar () |
Generate Random Numbers. | |
void | genRandPar (label tsize) |
Generate Random Numbers given the dimension of the training set. | |
void | genEquiPar () |
Generate Equidistributed Numbers. | |
void | truthSolve () |
Perform a TruthSolve. | |
void | assignBC (volVectorField &s, label BC_ind, Vector< double > &value) |
Assign Boundary Condition to a volVectorField. | |
void | assignBC (volScalarField &s, label BC_ind, double &value) |
Assign Boundary Condition to a volScalarField. | |
void | reconstructFromMatrix (PtrList< volVectorField > &rec_field2, PtrList< volVectorField > &modes, label Nmodes, Eigen::MatrixXd coeff_matrix) |
Exact reconstruction using a certain number of modes for vector list of fields and the projection coefficients (volVectorField) | |
void | reconstructFromMatrix (PtrList< volScalarField > &rec_field2, PtrList< volScalarField > &modes, label Nmodes, Eigen::MatrixXd coeff_matrix) |
Exact reconstruction using a certain number of modes for vector list of fields and the projection coefficients (volScalarField) | |
template<typename T , typename G > | |
void | assignIF (T &s, G &value) |
Assign internal field condition. | |
template<typename T > | |
void | computeLift (T &Lfield, T &liftfield, T &omfield) |
Homogenize the snapshot matrix, it works with PtrList of volVectorField and volScalarField. | |
template<typename T > | |
void | computeLiftT (T &Lfield, T &liftfield, T &omfield) |
Virtual function to compute the lifting function. | |
void | liftSolve () |
Virtual function to compute the lifting function for scalar field. | |
void | liftSolveT () |
void | project () |
General projection operation. | |
void | writeMu (List< scalar > mu_now) |
Write out a list of scalar corresponding to the parameters used in the truthSolve. | |
std::vector< SPLINTER::RBFSpline > | getCoeffManifoldRBF (PtrList< volVectorField > snapshots, PtrList< volVectorField > &modes, word rbfBasis="GAUSSIAN") |
Constructs the parameters-coefficients manifold for vector fields, based on RBF-spline model. | |
std::vector< SPLINTER::RBFSpline > | getCoeffManifoldRBF (PtrList< volScalarField > snapshots, PtrList< volScalarField > &modes, word rbfBasis="GAUSSIAN") |
Constructs the parameters-coefficients manifold for scalar fields, based on RBF-spline model. | |
std::vector< SPLINTER::BSpline > | getCoeffManifoldSPL (PtrList< volVectorField > snapshots, PtrList< volVectorField > &modes, label splDeg=3) |
Constructs the parameters-coefficients manifold for vector fields, based on the B-spline model. | |
std::vector< SPLINTER::BSpline > | getCoeffManifoldSPL (PtrList< volScalarField > snapshots, PtrList< volScalarField > &modes, label splDeg=3) |
Constructs the parameters-coefficients manifold for scalar fields, based on the B-spline model. | |
Public Attributes | |
label | Pnumber |
Number of parameters. | |
label | Tnumber |
Dimension of the training set (used only when gerating parameters without input) | |
Eigen::MatrixXd | mu |
Row matrix of parameters. | |
Eigen::MatrixXd | mu_range |
Range of the parameter spaces. | |
Eigen::MatrixXd | mu_samples |
Matrix of parameters to be used for PODI, where each row corresponds to a sample point. In this matrix the time dimension is regarded as a parameter for unsteady problems. | |
double | mu_cur |
Current value of the parameter. | |
bool | podex |
Boolean variable, it is 1 if the POD has already been computed, else 0. | |
bool | offline |
Boolean variable, it is 1 if the Offline phase has already been computed, else 0. | |
IOdictionary * | ITHACAdict |
dictionary to store input output infos | |
autoPtr< argList > | _args |
argList | |
ITHACAparallel * | paral |
parallel handling | |
label | folderN = 1 |
Counter to save intermediate steps in the correct folder, for unsteady and some stationary cases. | |
label | counter = 1 |
Counter used for the output of the full order solutions. | |
Eigen::MatrixXi | inletIndex |
Matrix that contains informations about the inlet boundaries. | |
Eigen::MatrixXi | inletPatch |
Matrix that contains informations about the inlet boundaries without specifing the direction Rows = Number of parametrized boundary conditions Cols = 1 Example: example.inletIndex.resize(2, 1); example.inletIndex(0, 0) = 0; example.inletIndex(1, 0) = 1; Means that there are two parametrized boundary conditions of which the first row is of patch 0 and the second row of patch 1. | |
Eigen::MatrixXi | inletIndexT |
A general class for the implementation of a full order parametrized problem.
Definition at line 72 of file reductionProblem.H.
reductionProblem::reductionProblem | ( | ) |
Construct Null.
Definition at line 41 of file reductionProblem.C.
|
inline |
Definition at line 81 of file reductionProblem.H.
void reductionProblem::assignBC | ( | volScalarField & | s, |
label | BC_ind, | ||
double & | value ) |
Assign Boundary Condition to a volScalarField.
[in] | s | field where you want to assign the BC. |
[in] | BC_ind | The BC index. |
[in] | value | The value you want to assign (it must be a double). |
Definition at line 105 of file reductionProblem.C.
void reductionProblem::assignBC | ( | volVectorField & | s, |
label | BC_ind, | ||
Vector< double > & | value ) |
Assign Boundary Condition to a volVectorField.
[in] | s | field where you want to assign the BC. |
[in] | BC_ind | The BC index. |
[in] | value | The value you want to assign (it must be an OpenFOAM vector). |
Definition at line 160 of file reductionProblem.C.
Assign internal field condition.
[in,out] | s | field where you want to assign the internal field condition |
value | The value you want to assign |
T | type of field (volVectorField or volScalarField) |
G | type of value you want to assign (OpenFOAM vector or scalar) |
Definition at line 34 of file reductionProblemTemplates.C.
Homogenize the snapshot matrix, it works with PtrList of volVectorField and volScalarField.
[in] | Lfield | The list of snapshots to be homogenized. |
[in] | liftfield | The list containing. |
[out] | omfield | The homogenized snapshot matrix. |
T | { description } |
Definition at line 43 of file reductionProblemTemplates.C.
void reductionProblem::computeLiftT | ( | T & | Lfield, |
T & | liftfield, | ||
T & | omfield ) |
Virtual function to compute the lifting function.
Definition at line 80 of file reductionProblemTemplates.C.
void reductionProblem::genEquiPar | ( | ) |
Generate Equidistributed Numbers.
Definition at line 86 of file reductionProblem.C.
void reductionProblem::genRandPar | ( | ) |
Generate Random Numbers.
Definition at line 56 of file reductionProblem.C.
void reductionProblem::genRandPar | ( | label | tsize | ) |
Generate Random Numbers given the dimension of the training set.
[in] | tsize | Dimension of the training set. |
Definition at line 68 of file reductionProblem.C.
std::vector< SPLINTER::RBFSpline > reductionProblem::getCoeffManifoldRBF | ( | PtrList< volScalarField > | snapshots, |
PtrList< volScalarField > & | modes, | ||
word | rbfBasis = "GAUSSIAN" ) |
Constructs the parameters-coefficients manifold for scalar fields, based on RBF-spline model.
[in] | snapshots | Snapshots scalar fields, used to compute the coefficient matrix |
[in] | modes | POD modes scalar fields, used to compute the coefficient matrix |
[in] | rbfBasis | The RBF basis type. Implemented bases are "GAUSSIAN", "THIN_PLATE", "MULTI_QUADRIC", "INVERSE_QUADRIC", and "INVERSE_MULTI_QUADRIC". Default basis is "Gaussian" |
Definition at line 354 of file reductionProblem.C.
std::vector< SPLINTER::RBFSpline > reductionProblem::getCoeffManifoldRBF | ( | PtrList< volVectorField > | snapshots, |
PtrList< volVectorField > & | modes, | ||
word | rbfBasis = "GAUSSIAN" ) |
Constructs the parameters-coefficients manifold for vector fields, based on RBF-spline model.
[in] | snapshots | Snapshots vector fields, used to compute the coefficient matrix |
[in] | modes | POD modes vector fields, used to compute the coefficient matrix |
[in] | rbfBasis | The RBF basis type. Implemented bases are "GAUSSIAN", "THIN_PLATE", "MULTI_QUADRIC", "INVERSE_QUADRIC", and "INVERSE_MULTI_QUADRIC". Default basis is "Gaussian" |
Definition at line 270 of file reductionProblem.C.
std::vector< SPLINTER::BSpline > reductionProblem::getCoeffManifoldSPL | ( | PtrList< volScalarField > | snapshots, |
PtrList< volScalarField > & | modes, | ||
label | splDeg = 3 ) |
Constructs the parameters-coefficients manifold for scalar fields, based on the B-spline model.
[in] | snapshots | Snapshots scalar fields, used to compute the coefficient matrix |
[in] | modes | POD modes scalar fields, used to compute the coefficient matrix |
[in] | splDeg | The B-spline degree. Default value is 3 |
Definition at line 492 of file reductionProblem.C.
std::vector< SPLINTER::BSpline > reductionProblem::getCoeffManifoldSPL | ( | PtrList< volVectorField > | snapshots, |
PtrList< volVectorField > & | modes, | ||
label | splDeg = 3 ) |
Constructs the parameters-coefficients manifold for vector fields, based on the B-spline model.
[in] | snapshots | Snapshots vector fields, used to compute the coefficient matrix |
[in] | modes | POD modes vector fields, used to compute the coefficient matrix |
[in] | splDeg | The B-spline degree. Default value is 3 |
Definition at line 438 of file reductionProblem.C.
void reductionProblem::liftSolve | ( | ) |
Virtual function to compute the lifting function for scalar field.
Definition at line 256 of file reductionProblem.C.
void reductionProblem::liftSolveT | ( | ) |
Definition at line 263 of file reductionProblem.C.
void reductionProblem::project | ( | ) |
General projection operation.
Definition at line 235 of file reductionProblem.C.
void reductionProblem::reconstructFromMatrix | ( | PtrList< volScalarField > & | rec_field2, |
PtrList< volScalarField > & | modes, | ||
label | Nmodes, | ||
Eigen::MatrixXd | coeff_matrix ) |
Exact reconstruction using a certain number of modes for vector list of fields and the projection coefficients (volScalarField)
[out] | rec_field2 | The reconstructed field as PtrList of volScalarField. |
[in] | modes | The modes used for reconstruction as PtrList of volScalarField. |
[in] | Nmodes | The number of modes you want to use. |
[in] | coeff_matrix | The matrix of coefficients. |
Exact reconstruction using a certain number of modes for scalar list of fields
Definition at line 215 of file reductionProblem.C.
void reductionProblem::reconstructFromMatrix | ( | PtrList< volVectorField > & | rec_field2, |
PtrList< volVectorField > & | modes, | ||
label | Nmodes, | ||
Eigen::MatrixXd | coeff_matrix ) |
Exact reconstruction using a certain number of modes for vector list of fields and the projection coefficients (volVectorField)
[out] | rec_field2 | The reconstructed field as PtrList of volVectorField. |
[in] | modes | The modes used for reconstruction as PtrList of volVectorField. |
[in] | Nmodes | The number of modes you want to use. |
[in] | coeff_matrix | The matrix of coefficients. |
Definition at line 193 of file reductionProblem.C.
void reductionProblem::setParameters | ( | ) |
Set Parameters Problems.
Definition at line 49 of file reductionProblem.C.
void reductionProblem::truthSolve | ( | ) |
Perform a TruthSolve.
Definition at line 97 of file reductionProblem.C.
void reductionProblem::writeMu | ( | List< scalar > | mu_now | ) |
Write out a list of scalar corresponding to the parameters used in the truthSolve.
[in] | mu_now | The list of scalars. |
Definition at line 242 of file reductionProblem.C.
autoPtr<argList> reductionProblem::_args |
label reductionProblem::counter = 1 |
Counter used for the output of the full order solutions.
Definition at line 110 of file reductionProblem.H.
label reductionProblem::folderN = 1 |
Counter to save intermediate steps in the correct folder, for unsteady and some stationary cases.
Definition at line 107 of file reductionProblem.H.
Eigen::MatrixXi reductionProblem::inletIndex |
Matrix that contains informations about the inlet boundaries.
The dimension is:
Rows = Number of parametrized boundary conditions
Cols = 2
Example:
example.inletIndex.resize(1, 2);
example.inletIndex(0, 0) = 0;
example.inletIndex(0, 1) = 0;
Means that there is one parametrized boundary conditions, the first col contains the index of the parametrized boundary, the second col contains the direction along which the BC is parametrized:
0 for x
1 for y
2 for z
The Matrix should be implemented in the following way: \[ \mbox{inletIndex}=\begin{bmatrix} \mbox{BC1_index} & \mbox{Direction of BC1} \\ \ \mbox{BC2_index} & \mbox{Direction of BC3} \\ \ \vdots & \vdots \\ \ \mbox{BCN_index} & \mbox{Direction of BCN} \ \ \end{bmatrix} \]
Definition at line 137 of file reductionProblem.H.
Eigen::MatrixXi reductionProblem::inletIndexT |
Definition at line 187 of file reductionProblem.H.
Eigen::MatrixXi reductionProblem::inletPatch |
Matrix that contains informations about the inlet boundaries without specifing the direction Rows = Number of parametrized boundary conditions
Cols = 1
Example:
example.inletIndex.resize(2, 1);
example.inletIndex(0, 0) = 0;
example.inletIndex(1, 0) = 1;
Means that there are two parametrized boundary conditions of which the first row is of patch 0 and the second row of patch 1.
The Matrix should be implemented in the following way: \[ \mbox{inletIndex}=\begin{bmatrix} \mbox{BC1_index} \\ \ \mbox{BC2_index} \\ \ \vdots & \vdots \\ \ \mbox{BCN_index} \ \ \end{bmatrix} \]
Definition at line 158 of file reductionProblem.H.
IOdictionary* reductionProblem::ITHACAdict |
dictionary to store input output infos
Definition at line 101 of file reductionProblem.H.
Eigen::MatrixXd reductionProblem::mu |
Row matrix of parameters.
Definition at line 89 of file reductionProblem.H.
double reductionProblem::mu_cur |
Current value of the parameter.
Definition at line 95 of file reductionProblem.H.
Eigen::MatrixXd reductionProblem::mu_range |
Range of the parameter spaces.
Definition at line 91 of file reductionProblem.H.
Eigen::MatrixXd reductionProblem::mu_samples |
Matrix of parameters to be used for PODI, where each row corresponds to a sample point. In this matrix the time dimension is regarded as a parameter for unsteady problems.
Definition at line 93 of file reductionProblem.H.
bool reductionProblem::offline |
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
Definition at line 99 of file reductionProblem.H.
ITHACAparallel* reductionProblem::paral |
parallel handling
Definition at line 105 of file reductionProblem.H.
label reductionProblem::Pnumber |
bool reductionProblem::podex |
Boolean variable, it is 1 if the POD has already been computed, else 0.
Definition at line 97 of file reductionProblem.H.
label reductionProblem::Tnumber |
Dimension of the training set (used only when gerating parameters without input)
Definition at line 87 of file reductionProblem.H.