Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
reductionProblem Class Reference

A general class for the implementation of a full order parametrized problem. More...

#include <reductionProblem.H>

Inheritance diagram for reductionProblem:
Burgers laplacianProblem msrProblem steadyNS tutorial23 DEIMLaplacian EnKF_1DinverseHeatTransfer inverseLaplacianProblem tutorial02 tutorialIPOD usmsrProblem SteadyNSSimple SteadyNSTurb SteadyNSTurbIntrusive tutorial03 tutorial05 unsteadyNS

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
 
ITHACAparallelparal
 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
 

Detailed Description

A general class for the implementation of a full order parametrized problem.

Definition at line 72 of file reductionProblem.H.

Constructor & Destructor Documentation

◆ reductionProblem()

reductionProblem::reductionProblem ( )

Construct Null.

Definition at line 41 of file reductionProblem.C.

◆ ~reductionProblem()

reductionProblem::~reductionProblem ( )
inline

Definition at line 81 of file reductionProblem.H.

Member Function Documentation

◆ assignBC() [1/2]

void reductionProblem::assignBC ( volScalarField & s,
label BC_ind,
double & value )

Assign Boundary Condition to a volScalarField.

Parameters
[in]sfield where you want to assign the BC.
[in]BC_indThe BC index.
[in]valueThe value you want to assign (it must be a double).

Definition at line 105 of file reductionProblem.C.

◆ assignBC() [2/2]

void reductionProblem::assignBC ( volVectorField & s,
label BC_ind,
Vector< double > & value )

Assign Boundary Condition to a volVectorField.

Parameters
[in]sfield where you want to assign the BC.
[in]BC_indThe BC index.
[in]valueThe value you want to assign (it must be an OpenFOAM vector).
Examples
06POD_RBF.C, 10UnsteadyBBEnclosed.C, 11UnsteadyBBOpen.C, and 17YJunction.C.

Definition at line 160 of file reductionProblem.C.

◆ assignIF()

template<typename T , typename G >
void reductionProblem::assignIF ( T & s,
G & value )

Assign internal field condition.

Parameters
[in,out]sfield where you want to assign the internal field condition
valueThe value you want to assign
Template Parameters
Ttype of field (volVectorField or volScalarField)
Gtype of value you want to assign (OpenFOAM vector or scalar)
Examples
02thermalBlock.C, 04unsteadyNS.C, 10UnsteadyBBEnclosed.C, 11UnsteadyBBOpen.C, 17YJunction.C, and 20incrementalPOD.C.

Definition at line 34 of file reductionProblemTemplates.C.

◆ computeLift()

template<typename T >
void reductionProblem::computeLift ( T & Lfield,
T & liftfield,
T & omfield )

Homogenize the snapshot matrix, it works with PtrList of volVectorField and volScalarField.

Parameters
[in]LfieldThe list of snapshots to be homogenized.
[in]liftfieldThe list containing.
[out]omfieldThe homogenized snapshot matrix.
Template Parameters
T{ description }
Examples
03steadyNS.C.

Definition at line 43 of file reductionProblemTemplates.C.

◆ computeLiftT()

template<typename T >
void reductionProblem::computeLiftT ( T & Lfield,
T & liftfield,
T & omfield )

Virtual function to compute the lifting function.

Definition at line 80 of file reductionProblemTemplates.C.

◆ genEquiPar()

void reductionProblem::genEquiPar ( )

Generate Equidistributed Numbers.

Definition at line 86 of file reductionProblem.C.

◆ genRandPar() [1/2]

void reductionProblem::genRandPar ( )

Generate Random Numbers.

Examples
02thermalBlock.C.

Definition at line 56 of file reductionProblem.C.

◆ genRandPar() [2/2]

void reductionProblem::genRandPar ( label tsize)

Generate Random Numbers given the dimension of the training set.

Parameters
[in]tsizeDimension of the training set.

Definition at line 68 of file reductionProblem.C.

◆ getCoeffManifoldRBF() [1/2]

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.

Parameters
[in]snapshotsSnapshots scalar fields, used to compute the coefficient matrix
[in]modesPOD modes scalar fields, used to compute the coefficient matrix
[in]rbfBasisThe RBF basis type. Implemented bases are "GAUSSIAN", "THIN_PLATE", "MULTI_QUADRIC", "INVERSE_QUADRIC", and "INVERSE_MULTI_QUADRIC". Default basis is "Gaussian"
Returns
Vector of objects to the RBF splines corresponding to each mode

Definition at line 354 of file reductionProblem.C.

◆ getCoeffManifoldRBF() [2/2]

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.

Parameters
[in]snapshotsSnapshots vector fields, used to compute the coefficient matrix
[in]modesPOD modes vector fields, used to compute the coefficient matrix
[in]rbfBasisThe RBF basis type. Implemented bases are "GAUSSIAN", "THIN_PLATE", "MULTI_QUADRIC", "INVERSE_QUADRIC", and "INVERSE_MULTI_QUADRIC". Default basis is "Gaussian"
Returns
Vector of objects to the RBF splines corresponding to each mode

Definition at line 270 of file reductionProblem.C.

◆ getCoeffManifoldSPL() [1/2]

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.

Parameters
[in]snapshotsSnapshots scalar fields, used to compute the coefficient matrix
[in]modesPOD modes scalar fields, used to compute the coefficient matrix
[in]splDegThe B-spline degree. Default value is 3
Returns
Vector of objects to the B-splines corresponding to each mode

Definition at line 492 of file reductionProblem.C.

◆ getCoeffManifoldSPL() [2/2]

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.

Parameters
[in]snapshotsSnapshots vector fields, used to compute the coefficient matrix
[in]modesPOD modes vector fields, used to compute the coefficient matrix
[in]splDegThe B-spline degree. Default value is 3
Returns
Vector of objects to the B-splines corresponding to each mode

Definition at line 438 of file reductionProblem.C.

◆ liftSolve()

void reductionProblem::liftSolve ( )

Virtual function to compute the lifting function for scalar field.

Definition at line 256 of file reductionProblem.C.

◆ liftSolveT()

void reductionProblem::liftSolveT ( )

Definition at line 263 of file reductionProblem.C.

◆ project()

void reductionProblem::project ( )

General projection operation.

Definition at line 235 of file reductionProblem.C.

◆ reconstructFromMatrix() [1/2]

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)

Parameters
[out]rec_field2The reconstructed field as PtrList of volScalarField.
[in]modesThe modes used for reconstruction as PtrList of volScalarField.
[in]NmodesThe number of modes you want to use.
[in]coeff_matrixThe matrix of coefficients.

Exact reconstruction using a certain number of modes for scalar list of fields

Definition at line 215 of file reductionProblem.C.

◆ reconstructFromMatrix() [2/2]

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)

Parameters
[out]rec_field2The reconstructed field as PtrList of volVectorField.
[in]modesThe modes used for reconstruction as PtrList of volVectorField.
[in]NmodesThe number of modes you want to use.
[in]coeff_matrixThe matrix of coefficients.

Definition at line 193 of file reductionProblem.C.

◆ setParameters()

void reductionProblem::setParameters ( )

Set Parameters Problems.

Examples
02thermalBlock.C.

Definition at line 49 of file reductionProblem.C.

◆ truthSolve()

void reductionProblem::truthSolve ( )

◆ writeMu()

void reductionProblem::writeMu ( List< scalar > mu_now)

Write out a list of scalar corresponding to the parameters used in the truthSolve.

Parameters
[in]mu_nowThe list of scalars.

Definition at line 242 of file reductionProblem.C.

Member Data Documentation

◆ _args

autoPtr<argList> reductionProblem::_args

argList

Examples
03steadyNS.C.

Definition at line 103 of file reductionProblem.H.

◆ counter

label reductionProblem::counter = 1

Counter used for the output of the full order solutions.

Examples
02thermalBlock.C, 06POD_RBF.C, and 20incrementalPOD.C.

Definition at line 110 of file reductionProblem.H.

◆ folderN

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.

◆ inletIndex

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} \]

Examples
03steadyNS.C.

Definition at line 137 of file reductionProblem.H.

◆ inletIndexT

Eigen::MatrixXi reductionProblem::inletIndexT
Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 187 of file reductionProblem.H.

◆ inletPatch

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} \]

Examples
17YJunction.C.

Definition at line 158 of file reductionProblem.H.

◆ ITHACAdict

IOdictionary* reductionProblem::ITHACAdict

dictionary to store input output infos

Examples
09DEIM_ROM.C.

Definition at line 101 of file reductionProblem.H.

◆ mu

Eigen::MatrixXd reductionProblem::mu

◆ mu_cur

double reductionProblem::mu_cur

Current value of the parameter.

Definition at line 95 of file reductionProblem.H.

◆ mu_range

Eigen::MatrixXd reductionProblem::mu_range

Range of the parameter spaces.

Examples
02thermalBlock.C.

Definition at line 91 of file reductionProblem.H.

◆ mu_samples

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.

Examples
02thermalBlock.C, 03steadyNS.C, and 20incrementalPOD.C.

Definition at line 93 of file reductionProblem.H.

◆ offline

bool reductionProblem::offline

Boolean variable, it is 1 if the Offline phase has already been computed, else 0.

Examples
02thermalBlock.C, 03steadyNS.C, 04unsteadyNS.C, 06POD_RBF.C, 09DEIM_ROM.C, 10UnsteadyBBEnclosed.C, 11UnsteadyBBOpen.C, 17YJunction.C, 19UnsteadyNSExplicit.C, and 20incrementalPOD.C.

Definition at line 99 of file reductionProblem.H.

◆ paral

ITHACAparallel* reductionProblem::paral

parallel handling

Definition at line 105 of file reductionProblem.H.

◆ Pnumber

label reductionProblem::Pnumber

Number of parameters.

Examples
02thermalBlock.C.

Definition at line 85 of file reductionProblem.H.

◆ podex

bool reductionProblem::podex

Boolean variable, it is 1 if the POD has already been computed, else 0.

Examples
02thermalBlock.C, and 03steadyNS.C.

Definition at line 97 of file reductionProblem.H.

◆ Tnumber

label reductionProblem::Tnumber

Dimension of the training set (used only when gerating parameters without input)

Examples
02thermalBlock.C.

Definition at line 87 of file reductionProblem.H.


The documentation for this class was generated from the following files: