Loading...
Searching...
No Matches
inverseLaplacianProblem_paramBC Class Reference

Implementation of a parameterization method to solve inverse Laplacian problems. More...

#include <inverseLaplacianProblem_paramBC.H>

Inheritance diagram for inverseLaplacianProblem_paramBC:
inverseLaplacianProblem laplacianProblem reductionProblem IHTP01inverseLaplacian_paramBC inverseLaplacianProblemTotalHeatMeasure_paramBC

Public Member Functions

 inverseLaplacianProblem_paramBC ()
 Null constructor.
 
 inverseLaplacianProblem_paramBC (int argc, char *argv[])
 Construct with argc and argv. Reads the thermocouples dictionary.
 
virtual ~inverseLaplacianProblem_paramBC ()
 Destructor.
 
virtual void set_gBaseFunctions ()
 Define the base functions used for the parametrization of the heat flux g The center of each base function is the projection of each thermocouple on the boundary hotSide.
 
void set_gBaseFunctionsPOD (label Nmodes)
 Performs POD on the RBF basis functions.
 
void set_gParametrized (word baseFuncType, scalar _shapeParameter=1)
 Set initial heat flux for the parameterized BC method.
 
void update_gParametrized (List< scalar > weigths)
 Update the boundary heat flux.
 
void parameterizedBCoffline (bool force=0)
 Performs offline computation for the parameterized BC method, if the offline directory "./ITHACAoutputs/offlineParamBC" exists, it reads the solution from there.
 
Eigen::VectorXd parameterizedBC (word linSys_solver="fullPivLU", double regPar=0)
 Solve the online phase.
 
void parameterizedBCpostProcess (Eigen::VectorXd weigths)
 Reconstruct the temperature field and compute the cost function J.
 
virtual void solveAdditional ()
 Set BC the additional problem and solves it.
 
void reconstructT ()
 Reconstuct the field T using the offline computed fields.
 
- Public Member Functions inherited from inverseLaplacianProblem
 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

word folderOffline = "./ITHACAoutput/offlineParamBC/"
 Folder where the offline solutions and matrices are saved.
 
PtrList< volScalarField > Tbasis
 Solutions of the direct problem with the basis of the parameterization as boundary heat flux.
 
PtrList< volScalarField > Tad_base
 Solution of the additional problem.
 
Eigen::VectorXd residual
 Residual of the linear system.
 
Eigen::MatrixXd Theta
 Theta matrix of the lynear system.
 
Eigen::MatrixXd gPODmodes
 Modes of the POD.
 
Eigen::VectorXd addSol
 Solution of the additional problem at the thermocouples positions.
 
bool offlineReady = 0
 Flag to know if the offline phase was computed.
 
List< List< scalar > > gBaseFunctions
 Basis of the heat flux parameterization.
 
List< scalar > gWeights
 Weights of the heat flux parameterization.
 
double shapeParameter
 RBF shape parameter.
 
word baseFuncType
 Type of basis functions used for the parameterization of the heat flux.
 
- Public Attributes inherited from inverseLaplacianProblem
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 parameterization method to solve inverse Laplacian problems.

In this class, we implement the parameterization method for solving the inverse problem of estimating the boundary heat flux, given pointwise temperature measurements, in a Laplacian problem

Definition at line 72 of file inverseLaplacianProblem_paramBC.H.

Constructor & Destructor Documentation

◆ inverseLaplacianProblem_paramBC() [1/2]

inverseLaplacianProblem_paramBC::inverseLaplacianProblem_paramBC ( )

Null constructor.

Definition at line 39 of file inverseLaplacianProblem_paramBC.C.

◆ inverseLaplacianProblem_paramBC() [2/2]

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

Construct with argc and argv. Reads the thermocouples dictionary.

Definition at line 41 of file inverseLaplacianProblem_paramBC.C.

◆ ~inverseLaplacianProblem_paramBC()

virtual inverseLaplacianProblem_paramBC::~inverseLaplacianProblem_paramBC ( )
inlinevirtual

Destructor.

Definition at line 82 of file inverseLaplacianProblem_paramBC.H.

Member Function Documentation

◆ parameterizedBC()

Eigen::VectorXd inverseLaplacianProblem_paramBC::parameterizedBC ( word linSys_solver = "fullPivLU",
double regPar = 0 )

Solve the online phase.

Parameters
[in]linSys_solverType of linear solver are fullPivLU, jacobiSvd, householderQr, ldlt, inverse, TSVD, Tikhonov
[in]regParValue of the regularization parameter

Definition at line 309 of file inverseLaplacianProblem_paramBC.C.

◆ parameterizedBCoffline()

void inverseLaplacianProblem_paramBC::parameterizedBCoffline ( bool force = 0)

Performs offline computation for the parameterized BC method, if the offline directory "./ITHACAoutputs/offlineParamBC" exists, it reads the solution from there.

Parameters
[in]forceIf 1, force the offline phase to be computed

Definition at line 198 of file inverseLaplacianProblem_paramBC.C.

◆ parameterizedBCpostProcess()

void inverseLaplacianProblem_paramBC::parameterizedBCpostProcess ( Eigen::VectorXd weigths)

Reconstruct the temperature field and compute the cost function J.

Parameters
[in]weigthsWeights of the parameterization

Definition at line 367 of file inverseLaplacianProblem_paramBC.C.

◆ reconstructT()

void inverseLaplacianProblem_paramBC::reconstructT ( )

Reconstuct the field T using the offline computed fields.

Definition at line 448 of file inverseLaplacianProblem_paramBC.C.

◆ set_gBaseFunctions()

void inverseLaplacianProblem_paramBC::set_gBaseFunctions ( )
virtual

Define the base functions used for the parametrization of the heat flux g The center of each base function is the projection of each thermocouple on the boundary hotSide.

The type of basis function is defined by the value of baseFuncType (rbf = radial basis functions, pod = proper orthogonal decomposition).

Definition at line 50 of file inverseLaplacianProblem_paramBC.C.

◆ set_gBaseFunctionsPOD()

void inverseLaplacianProblem_paramBC::set_gBaseFunctionsPOD ( label Nmodes)

Performs POD on the RBF basis functions.

Parameters
[in]NmodesNumber of modes used for the POD

Definition at line 113 of file inverseLaplacianProblem_paramBC.C.

◆ set_gParametrized()

void inverseLaplacianProblem_paramBC::set_gParametrized ( word baseFuncType,
scalar _shapeParameter = 1 )

Set initial heat flux for the parameterized BC method.

Parameters
[in]baseFuncTypeType of basis function (rbf or pod)
[in]_shapeParameterValue of the RBF shape parameter

Definition at line 160 of file inverseLaplacianProblem_paramBC.C.

◆ solveAdditional()

void inverseLaplacianProblem_paramBC::solveAdditional ( )
virtual

Set BC the additional problem and solves it.

Reimplemented in IHTP01inverseLaplacian_paramBC.

Definition at line 399 of file inverseLaplacianProblem_paramBC.C.

◆ update_gParametrized()

void inverseLaplacianProblem_paramBC::update_gParametrized ( List< scalar > weigths)

Update the boundary heat flux.

Parameters
[in]weigthsNew values of the weights of the basis functions

Definition at line 183 of file inverseLaplacianProblem_paramBC.C.

Member Data Documentation

◆ addSol

Eigen::VectorXd inverseLaplacianProblem_paramBC::addSol

Solution of the additional problem at the thermocouples positions.

Definition at line 104 of file inverseLaplacianProblem_paramBC.H.

◆ baseFuncType

word inverseLaplacianProblem_paramBC::baseFuncType

Type of basis functions used for the parameterization of the heat flux.

Definition at line 119 of file inverseLaplacianProblem_paramBC.H.

◆ folderOffline

word inverseLaplacianProblem_paramBC::folderOffline = "./ITHACAoutput/offlineParamBC/"

Folder where the offline solutions and matrices are saved.

Definition at line 85 of file inverseLaplacianProblem_paramBC.H.

◆ gBaseFunctions

List<List<scalar> > inverseLaplacianProblem_paramBC::gBaseFunctions

Basis of the heat flux parameterization.

Definition at line 110 of file inverseLaplacianProblem_paramBC.H.

◆ gPODmodes

Eigen::MatrixXd inverseLaplacianProblem_paramBC::gPODmodes

Modes of the POD.

Definition at line 101 of file inverseLaplacianProblem_paramBC.H.

◆ gWeights

List<scalar> inverseLaplacianProblem_paramBC::gWeights

Weights of the heat flux parameterization.

Definition at line 113 of file inverseLaplacianProblem_paramBC.H.

◆ offlineReady

bool inverseLaplacianProblem_paramBC::offlineReady = 0

Flag to know if the offline phase was computed.

Definition at line 107 of file inverseLaplacianProblem_paramBC.H.

◆ residual

Eigen::VectorXd inverseLaplacianProblem_paramBC::residual

Residual of the linear system.

Definition at line 95 of file inverseLaplacianProblem_paramBC.H.

◆ shapeParameter

double inverseLaplacianProblem_paramBC::shapeParameter

RBF shape parameter.

Definition at line 116 of file inverseLaplacianProblem_paramBC.H.

◆ Tad_base

PtrList<volScalarField> inverseLaplacianProblem_paramBC::Tad_base

Solution of the additional problem.

Definition at line 92 of file inverseLaplacianProblem_paramBC.H.

◆ Tbasis

PtrList<volScalarField> inverseLaplacianProblem_paramBC::Tbasis

Solutions of the direct problem with the basis of the parameterization as boundary heat flux.

Definition at line 88 of file inverseLaplacianProblem_paramBC.H.

◆ Theta

Eigen::MatrixXd inverseLaplacianProblem_paramBC::Theta

Theta matrix of the lynear system.

Definition at line 98 of file inverseLaplacianProblem_paramBC.H.


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