Loading...
Searching...
No Matches
laplacianProblem Class Reference

Class to implement a full order laplacian parametrized problem. More...

#include <laplacianProblem.H>

Inheritance diagram for laplacianProblem:
reductionProblem DEIMLaplacian EnKF_1DinverseHeatTransfer inverseLaplacianProblem tutorial02 tutorialIPOD inverseLaplacianProblem_CG inverseLaplacianProblem_paramBC IHTP01inverseLaplacian_CG inverseLaplacianProblemTotalHeatMeasure_CG IHTP01inverseLaplacian_paramBC inverseLaplacianProblemTotalHeatMeasure_paramBC

Public Member Functions

 laplacianProblem ()
 
 laplacianProblem (int argc, char *argv[])
 Construct with argc and argv.
 
 ~laplacianProblem ()
 
void truthSolve (List< scalar > mu_now, word folder="./ITHACAoutput/Offline/")
 Perform a truthsolve.
 
void project (label Nmodes)
 Perform a projection onto the POD modes.
 
- Public Member Functions inherited from reductionProblem
 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

PtrList< volScalarField > Tfield
 List of snapshots for the solution.
 
PtrList< volScalarField > Tonline
 List of snapshots for the solution.
 
volScalarModes Tmodes
 List of POD modes.
 
PtrList< fvScalarMatrix > operator_list
 List of operators.
 
List< scalar > theta
 Theta (coefficients of the affine expansion)
 
PtrList< volScalarField > nu_list
 Nu (diffusivity)
 
label NTmodes
 Number of modes reduced problem.
 
List< Eigen::MatrixXd > A_matrices
 A matrices.
 
Eigen::MatrixXd source
 Source vector.
 
autoPtr< volScalarField > _T
 Temperature field.
 
autoPtr< volScalarField > _S
 Source Term.
 
autoPtr< volScalarField > _nu
 Diffusivity.
 
autoPtr< fvMesh > _mesh
 Mesh.
 
autoPtr< Time > _runTime
 Time.
 
- Public Attributes inherited from reductionProblem
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

Class to implement a full order laplacian parametrized problem.

Examples
02thermalBlock.C, 09DEIM_ROM.C, and 20incrementalPOD.C.

Definition at line 50 of file laplacianProblem.H.

Constructor & Destructor Documentation

◆ laplacianProblem() [1/2]

laplacianProblem::laplacianProblem ( )
Examples
02thermalBlock.C, and 20incrementalPOD.C.

Definition at line 40 of file laplacianProblem.C.

◆ laplacianProblem() [2/2]

laplacianProblem::laplacianProblem ( int argc,
char * argv[] )

Construct with argc and argv.

Definition at line 41 of file laplacianProblem.C.

◆ ~laplacianProblem()

laplacianProblem::~laplacianProblem ( )
inline

Definition at line 58 of file laplacianProblem.H.

Member Function Documentation

◆ project()

void laplacianProblem::project ( label Nmodes)

Perform a projection onto the POD modes.

Parameters
[in]NmodesThe number of modes used for the projection

Export the A matrices

Export the source term

Examples
02thermalBlock.C.

Definition at line 101 of file laplacianProblem.C.

◆ truthSolve()

void laplacianProblem::truthSolve ( List< scalar > mu_now,
word folder = "./ITHACAoutput/Offline/" )

Perform a truthsolve.

Definition at line 64 of file laplacianProblem.C.

Member Data Documentation

◆ _mesh

autoPtr<fvMesh> laplacianProblem::_mesh
mutable

Mesh.

Examples
02enKF_1DinverseHeatTransfer.C, 02thermalBlock.C, and 09DEIM_ROM.C.

Definition at line 96 of file laplacianProblem.H.

◆ _nu

autoPtr<volScalarField> laplacianProblem::_nu

Diffusivity.

Examples
02thermalBlock.C, 09DEIM_ROM.C, and 20incrementalPOD.C.

Definition at line 94 of file laplacianProblem.H.

◆ _runTime

autoPtr<Time> laplacianProblem::_runTime

Time.

Examples
02enKF_1DinverseHeatTransfer.C, and 02thermalBlock.C.

Definition at line 98 of file laplacianProblem.H.

◆ _S

autoPtr<volScalarField> laplacianProblem::_S

Source Term.

Examples
02thermalBlock.C, 09DEIM_ROM.C, and 20incrementalPOD.C.

Definition at line 92 of file laplacianProblem.H.

◆ _T

autoPtr<volScalarField> laplacianProblem::_T

Temperature field.

Examples
02enKF_1DinverseHeatTransfer.C, 02thermalBlock.C, 09DEIM_ROM.C, and 20incrementalPOD.C.

Definition at line 90 of file laplacianProblem.H.

◆ A_matrices

List<Eigen::MatrixXd> laplacianProblem::A_matrices

A matrices.

Definition at line 84 of file laplacianProblem.H.

◆ NTmodes

label laplacianProblem::NTmodes

Number of modes reduced problem.

Definition at line 80 of file laplacianProblem.H.

◆ nu_list

PtrList<volScalarField> laplacianProblem::nu_list

Nu (diffusivity)

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

Definition at line 77 of file laplacianProblem.H.

◆ operator_list

PtrList<fvScalarMatrix> laplacianProblem::operator_list

List of operators.

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

Definition at line 71 of file laplacianProblem.H.

◆ source

Eigen::MatrixXd laplacianProblem::source

Source vector.

Definition at line 86 of file laplacianProblem.H.

◆ Tfield

PtrList<volScalarField> laplacianProblem::Tfield

List of snapshots for the solution.

Examples
02thermalBlock.C, 09DEIM_ROM.C, and 20incrementalPOD.C.

Definition at line 62 of file laplacianProblem.H.

◆ theta

List<scalar> laplacianProblem::theta

Theta (coefficients of the affine expansion)

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

Definition at line 74 of file laplacianProblem.H.

◆ Tmodes

volScalarModes laplacianProblem::Tmodes

List of POD modes.

Examples
02thermalBlock.C, and 09DEIM_ROM.C.

Definition at line 68 of file laplacianProblem.H.

◆ Tonline

PtrList<volScalarField> laplacianProblem::Tonline

List of snapshots for the solution.

Examples
09DEIM_ROM.C.

Definition at line 65 of file laplacianProblem.H.


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