40#ifndef inverseLaplacianProblem_CG_H
41#define inverseLaplacianProblem_CG_H
44#define _USE_MATH_DEFINES
188 void solve(
const char* problemID);
245 int isInPlane(
double cx,
double cy,
double cz,
246 Foam::vector thermocoupleCellDim);
254 void writeFields(label folderNumber,
const char* folder);
286 void restart(word fieldName =
"all");
Implementation of a Alifanov's regularization to solve inverse Laplacian problems.
int cgIterMax
Maximum CG iterations.
List< Eigen::MatrixXd > ArbT
Temperature reduced matrix.
Eigen::VectorXd Tsens
Vector of solutions of the sensitivity problem at the thermocouples points.
Eigen::MatrixXd Jlist
Vector to store the const function J.
void set_valueFraction()
Set valueFraction list values for Robin condition.
PtrList< volScalarField > lambdaField
List of adjoint solutions.
double gamma
Conjugate coefficient.
bool interpolationPlaneDefined
Flag for definition of interpolation plane.
virtual ~inverseLaplacianProblem_CG()
Destructor.
bool saveSolInLists
Flag to save solutions in lists.
List< Eigen::MatrixXd > ArbLambda
Adjoint reduced matrix.
double Jtol
Absolute stopping criterion for the CG.
void defineThermocouplesPlane()
Identifies the plane defined by the thermocouples.
void assignSensitivityBC()
Set BC and IF of the sensitivity problem.
thermocouplesPlane interpolationPlane
Interpolation plane.
void thermocouplesInterpolation()
Interpolates the thermocouples measurements in the plane defined in readThermocouples() using radial ...
int isInPlane(double cx, double cy, double cz, Foam::vector thermocoupleCellDim)
Checks if a cell crosses the interpolation plane.
int cgIter
Conjugate Gradient (CG) interations counter.
double gamma_den
Denoinator of the conjugate coefficient.
List< Eigen::MatrixXd > ArbDeltaT
Sensitivity reduced matrix.
void sensibilitySolAtThermocouplesLocations()
Fill the Foam::vector containing the values of the sensitivity solution at the thermocouples location...
void assignAdjointBC()
Set BC of the adjoint problem.
double beta
CG search step size.
void updateHeatFlux()
Updates the heat flux in the conjugate gradient iterations.
int conjugateGradient()
Conjugate gradient method.
inverseLaplacianProblem_CG()
Null constructor.
void writeFields(label folderNumber, const char *folder)
Writes fields to file.
List< scalar > gradJ
Gradient of the cost function.
bool interpolation
Flag for interpolation of the temperature measurements.
List< scalar > P
Search direction.
void solveSensitivity()
Solve sensibility problem.
double gradJ_L2norm
L2 norm of the gradient of J.
PtrList< volScalarField > Tfield
List of temperature solutions.
double JtolRel
Relative stopping criterion for the CG.
autoPtr< volScalarField > _deltaT
Sensibility temperature field.
void computeGradJ()
Computes the gradient of cost function J and its L2 norm.
void solveAdjoint()
Solve adjoint problem.
void computeSearchStep()
Compute the search step beta.
Eigen::VectorXd cellsInPlane
IDs of the cells in the interpolation plane.
PtrList< volScalarField > deltaTfield
List of sensitivity solutions.
volScalarField assignAdjointBCandSource()
Assign the BC for the adjoint problem and returs the source field.
List< scalar > valueFractionAdj
Value fraction for the adjoint Robin BC.
autoPtr< volScalarField > _lambda
Adjoint field.
void searchDirection()
Computes the search direction P.
void differenceBetweenDirectAndMeasure()
Computes the difference between direct problem solution and measure.
int conjugateGradientConvergenceCheck()
Convergence cher for the conjugate gradient method.
Implementation of a inverse Laplacian problem .
void restart()
Restart fields.
Header file of the inverseLaplacianProblem class.
Class defining a plane for thermocouples measurements interpolation.