Loading...
Searching...
No Matches
inverseLaplacianProblem Class Reference

Implementation of a inverse Laplacian problem . More...

#include <inverseLaplacianProblem.H>

Inheritance diagram for inverseLaplacianProblem:
laplacianProblem reductionProblem inverseLaplacianProblem_CG inverseLaplacianProblem_paramBC IHTP01inverseLaplacian_CG inverseLaplacianProblemTotalHeatMeasure_CG IHTP01inverseLaplacian_paramBC inverseLaplacianProblemTotalHeatMeasure_paramBC

Public Member Functions

 inverseLaplacianProblem ()
 Null constructor.
 
 inverseLaplacianProblem (int argc, char *argv[])
 Construct with argc and argv. Reads the thermocouples dictionary.
 
virtual ~inverseLaplacianProblem ()
 Destructor.
 
void set_g ()
 Set the right g size and fills it with zeros.
 
volScalarField list2Field (List< scalar > list, scalar innerField=0.0)
 Create a field with the hotSide boundary heat flux at the hotSide bounday cells for visualization.
 
void set_valueFraction ()
 Set valueFraction list values for Robin condition.
 
virtual void solveTrue ()
 Solve the direct problem with the true heat flux as boundary condition and fills the measured temeprature vector.
 
virtual void assignDirectBC ()
 Set boundary condition of the direct problem.
 
void solveDirect ()
 Solve direct problem.
 
void solve (const char *problem)
 Solve Laplacian problem without source term.
 
virtual void readThermocouples ()
 Identifies in the mesh the cells corresponding to the termocouples locations.
 
Eigen::VectorXd fieldValueAtThermocouples (volScalarField &field)
 Interpolates the field value at the thermocouples points.
 
void differenceBetweenDirectAndMeasure ()
 Computes the difference between direct problem solution and measures Saves the difference vector in Tdiff.
 
Foam::vector cellDim (const faceList &ff, const pointField &pp, const cell &cc, labelList pLabels, pointField pLocal)
 Compute maximum cell dimension in x, y and z.
 
void restart ()
 Restart fields.
 
- Public Member Functions inherited from laplacianProblem
 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

ITHACAparameterspara
 
autoPtr< volScalarField > _T
 Temperature field.
 
autoPtr< fvMesh > _mesh
 Mesh.
 
autoPtr< simpleControl > _simple
 simpleControl
 
autoPtr< fv::options > _fvOptions
 fvOptions
 
autoPtr< Time > _runTime
 Time.
 
dimensionedScalar DT
 Dummy thermal conductivity with unitary value.
 
double k
 Thermal diffusivity [W/(m K)].
 
double H
 Heat transfer coefficient [W/(m2 K)].
 
bool thermocouplesRead = 0
 Flag to know if thermocouples file was read.
 
int thermocouplesNum
 Number of thermocouples.
 
double J
 Cost funtion [K^2].
 
double L2norm
 
double LinfNorm
 
scalar homogeneousBC = 0.0
 Homogenerous BC.
 
List< scalar > homogeneousBCcoldSide
 List of zeros of the size of coldSide patch.
 
List< scalar > Tf
 Temperature at coldSide [K].
 
List< scalar > refGrad
 Reference gradient for the Robin BC.
 
List< scalar > valueFraction
 Value fraction for the Robin BC.
 
label hotSide_ind
 Index of the hotSide patch.
 
label coldSide_ind
 Index of the coldSide patch.
 
List< scalar > g
 Heat flux at hotSide [W/m2].
 
List< List< scalar > > gList
 List of boundary heat fluxes.
 
List< scalar > gTrue
 True heat flux at hotSide [w/m2].
 
List< scalar > faceCellArea
 List of patch faces areas [m2].
 
List< vector > thermocouplesPos
 List containing the positions of the thermocouples.
 
List< int > thermocouplesCellID
 List of cells indices containing a thermocouple.
 
List< int > thermocouplesCellProc
 List of incedes of the processors contining each thermocouple.
 
Eigen::VectorXd Tmeas
 Vector of measured temperatures at the thermocouples locations [K].
 
Eigen::VectorXd Tdirect
 Vector of computed temperatures at the thermocouples locations [K].
 
Eigen::VectorXd Tdiff
 Difference between computed and measured temperatures at the thermocouples.
 
label nProcs
 Number of processors.
 
- Public Attributes inherited from laplacianProblem
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

Implementation of a inverse Laplacian problem .

In this class, a general inverse Laplacian problem is set up. In this problem, the objective is to estimate the boundary heat flux at the hotSide patch given pointwise temeprature measurements in the interior of the domain. The different solution methodologies are implemented as child of this class

Definition at line 74 of file inverseLaplacianProblem.H.

Constructor & Destructor Documentation

◆ inverseLaplacianProblem() [1/2]

inverseLaplacianProblem::inverseLaplacianProblem ( )

Null constructor.

Definition at line 40 of file inverseLaplacianProblem.C.

◆ inverseLaplacianProblem() [2/2]

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

Construct with argc and argv. Reads the thermocouples dictionary.

Definition at line 42 of file inverseLaplacianProblem.C.

◆ ~inverseLaplacianProblem()

virtual inverseLaplacianProblem::~inverseLaplacianProblem ( )
inlinevirtual

Destructor.

Definition at line 84 of file inverseLaplacianProblem.H.

Member Function Documentation

◆ assignDirectBC()

void inverseLaplacianProblem::assignDirectBC ( )
virtual

Set boundary condition of the direct problem.

Reimplemented in IHTP01inverseLaplacian_CG, and IHTP01inverseLaplacian_paramBC.

Definition at line 139 of file inverseLaplacianProblem.C.

◆ cellDim()

Foam::vector inverseLaplacianProblem::cellDim ( const faceList & ff,
const pointField & pp,
const cell & cc,
labelList pLabels,
pointField pLocal )

Compute maximum cell dimension in x, y and z.

Parameters
[in]ff
[in]pp
[in]cc
[in]pLabels
[in]pLocal
[out]dimMaximum cell dimensions

Definition at line 283 of file inverseLaplacianProblem.C.

◆ differenceBetweenDirectAndMeasure()

void inverseLaplacianProblem::differenceBetweenDirectAndMeasure ( )

Computes the difference between direct problem solution and measures Saves the difference vector in Tdiff.

Definition at line 275 of file inverseLaplacianProblem.C.

◆ fieldValueAtThermocouples()

Eigen::VectorXd inverseLaplacianProblem::fieldValueAtThermocouples ( volScalarField & field)

Interpolates the field value at the thermocouples points.

Parameters
[in]fieldField to read the values
Returns
Vector of field values at thermocouples points
Examples
IHTP01inverseLaplacian.C.

Definition at line 252 of file inverseLaplacianProblem.C.

◆ list2Field()

volScalarField inverseLaplacianProblem::list2Field ( List< scalar > list,
scalar innerField = 0.0 )

Create a field with the hotSide boundary heat flux at the hotSide bounday cells for visualization.

Parameters
[in]listList of boundary heat flux
[in]innerFieldValue of the inner field
Returns
Field with the heat flux value at the boundary cells

Definition at line 100 of file inverseLaplacianProblem.C.

◆ readThermocouples()

void inverseLaplacianProblem::readThermocouples ( )
virtual

Identifies in the mesh the cells corresponding to the termocouples locations.

Examples
IHTP01inverseLaplacian.C.

Definition at line 205 of file inverseLaplacianProblem.C.

◆ restart()

void inverseLaplacianProblem::restart ( )

Restart fields.

Definition at line 300 of file inverseLaplacianProblem.C.

◆ set_g()

void inverseLaplacianProblem::set_g ( )

Set the right g size and fills it with zeros.

Definition at line 90 of file inverseLaplacianProblem.C.

◆ set_valueFraction()

void inverseLaplacianProblem::set_valueFraction ( )

Set valueFraction list values for Robin condition.

Definition at line 121 of file inverseLaplacianProblem.C.

◆ solve()

void inverseLaplacianProblem::solve ( const char * problem)

Solve Laplacian problem without source term.

Definition at line 168 of file inverseLaplacianProblem.C.

◆ solveDirect()

void inverseLaplacianProblem::solveDirect ( )

Solve direct problem.

Definition at line 161 of file inverseLaplacianProblem.C.

◆ solveTrue()

virtual void inverseLaplacianProblem::solveTrue ( )
inlinevirtual

Solve the direct problem with the true heat flux as boundary condition and fills the measured temeprature vector.

Reimplemented in IHTP01inverseLaplacian_CG, and IHTP01inverseLaplacian_paramBC.

Definition at line 205 of file inverseLaplacianProblem.H.

Member Data Documentation

◆ _fvOptions

autoPtr<fv::options> inverseLaplacianProblem::_fvOptions

fvOptions

Definition at line 98 of file inverseLaplacianProblem.H.

◆ _mesh

autoPtr<fvMesh> inverseLaplacianProblem::_mesh
mutable

Mesh.

Examples
IHTP01inverseLaplacian.C.

Definition at line 92 of file inverseLaplacianProblem.H.

◆ _runTime

autoPtr<Time> inverseLaplacianProblem::_runTime

Time.

Examples
IHTP01inverseLaplacian.C.

Definition at line 101 of file inverseLaplacianProblem.H.

◆ _simple

autoPtr<simpleControl> inverseLaplacianProblem::_simple

simpleControl

Definition at line 95 of file inverseLaplacianProblem.H.

◆ _T

autoPtr<volScalarField> inverseLaplacianProblem::_T

Temperature field.

Examples
IHTP01inverseLaplacian.C.

Definition at line 89 of file inverseLaplacianProblem.H.

◆ coldSide_ind

label inverseLaplacianProblem::coldSide_ind

Index of the coldSide patch.

Definition at line 142 of file inverseLaplacianProblem.H.

◆ DT

dimensionedScalar inverseLaplacianProblem::DT

Dummy thermal conductivity with unitary value.

Definition at line 104 of file inverseLaplacianProblem.H.

◆ faceCellArea

List<scalar> inverseLaplacianProblem::faceCellArea

List of patch faces areas [m2].

Definition at line 154 of file inverseLaplacianProblem.H.

◆ g

List<scalar> inverseLaplacianProblem::g

Heat flux at hotSide [W/m2].

Examples
IHTP01inverseLaplacian.C.

Definition at line 145 of file inverseLaplacianProblem.H.

◆ gList

List<List<scalar> > inverseLaplacianProblem::gList

List of boundary heat fluxes.

Definition at line 148 of file inverseLaplacianProblem.H.

◆ gTrue

List<scalar> inverseLaplacianProblem::gTrue

True heat flux at hotSide [w/m2].

Examples
IHTP01inverseLaplacian.C.

Definition at line 151 of file inverseLaplacianProblem.H.

◆ H

double inverseLaplacianProblem::H

Heat transfer coefficient [W/(m2 K)].

Examples
IHTP01inverseLaplacian.C.

Definition at line 110 of file inverseLaplacianProblem.H.

◆ homogeneousBC

scalar inverseLaplacianProblem::homogeneousBC = 0.0

Homogenerous BC.

Definition at line 124 of file inverseLaplacianProblem.H.

◆ homogeneousBCcoldSide

List<scalar> inverseLaplacianProblem::homogeneousBCcoldSide

List of zeros of the size of coldSide patch.

Definition at line 127 of file inverseLaplacianProblem.H.

◆ hotSide_ind

label inverseLaplacianProblem::hotSide_ind

Index of the hotSide patch.

Examples
IHTP01inverseLaplacian.C.

Definition at line 139 of file inverseLaplacianProblem.H.

◆ J

double inverseLaplacianProblem::J

Cost funtion [K^2].

Definition at line 119 of file inverseLaplacianProblem.H.

◆ k

double inverseLaplacianProblem::k

Thermal diffusivity [W/(m K)].

Examples
IHTP01inverseLaplacian.C.

Definition at line 107 of file inverseLaplacianProblem.H.

◆ L2norm

double inverseLaplacianProblem::L2norm

Definition at line 120 of file inverseLaplacianProblem.H.

◆ LinfNorm

double inverseLaplacianProblem::LinfNorm

Definition at line 121 of file inverseLaplacianProblem.H.

◆ nProcs

label inverseLaplacianProblem::nProcs

Number of processors.

Definition at line 175 of file inverseLaplacianProblem.H.

◆ para

ITHACAparameters* inverseLaplacianProblem::para

Definition at line 86 of file inverseLaplacianProblem.H.

◆ refGrad

List<scalar> inverseLaplacianProblem::refGrad

Reference gradient for the Robin BC.

Definition at line 133 of file inverseLaplacianProblem.H.

◆ Tdiff

Eigen::VectorXd inverseLaplacianProblem::Tdiff

Difference between computed and measured temperatures at the thermocouples.

Definition at line 172 of file inverseLaplacianProblem.H.

◆ Tdirect

Eigen::VectorXd inverseLaplacianProblem::Tdirect

Vector of computed temperatures at the thermocouples locations [K].

Definition at line 169 of file inverseLaplacianProblem.H.

◆ Tf

List<scalar> inverseLaplacianProblem::Tf

Temperature at coldSide [K].

Definition at line 130 of file inverseLaplacianProblem.H.

◆ thermocouplesCellID

List<int> inverseLaplacianProblem::thermocouplesCellID

List of cells indices containing a thermocouple.

Definition at line 160 of file inverseLaplacianProblem.H.

◆ thermocouplesCellProc

List<int> inverseLaplacianProblem::thermocouplesCellProc

List of incedes of the processors contining each thermocouple.

Definition at line 163 of file inverseLaplacianProblem.H.

◆ thermocouplesNum

int inverseLaplacianProblem::thermocouplesNum

Number of thermocouples.

Definition at line 116 of file inverseLaplacianProblem.H.

◆ thermocouplesPos

List<vector> inverseLaplacianProblem::thermocouplesPos

List containing the positions of the thermocouples.

Definition at line 157 of file inverseLaplacianProblem.H.

◆ thermocouplesRead

bool inverseLaplacianProblem::thermocouplesRead = 0

Flag to know if thermocouples file was read.

Definition at line 113 of file inverseLaplacianProblem.H.

◆ Tmeas

Eigen::VectorXd inverseLaplacianProblem::Tmeas

Vector of measured temperatures at the thermocouples locations [K].

Examples
IHTP01inverseLaplacian.C.

Definition at line 166 of file inverseLaplacianProblem.H.

◆ valueFraction

List<scalar> inverseLaplacianProblem::valueFraction

Value fraction for the Robin BC.

Definition at line 136 of file inverseLaplacianProblem.H.


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