1#ifndef IHTP01inverseLaplacian_H
2#define IHTP01inverseLaplacian_H
81 volScalarField Tad(
_T());
84 List<scalar> RobinBC = -
Tf;
87 if (patchI ==
mesh.boundaryMesh().findPatchID(
"coldSide"))
91 else if (patchI ==
mesh.boundaryMesh().findPatchID(
"hotSide"))
95 else if (patchI ==
mesh.boundaryMesh().findPatchID(
"gammaEx1"))
99 else if (patchI ==
mesh.boundaryMesh().findPatchID(
"gammaEx2"))
103 else if (patchI ==
mesh.boundaryMesh().findPatchID(
"gammaEx3"))
107 else if (patchI ==
mesh.boundaryMesh().findPatchID(
"gammaEx4"))
112#if defined(OFVER) && (OFVER == 6)
119 while (
simple.correctNonOrthogonal())
123 fvm::laplacian(
DT, Tad)
131 instantList Times =
runTime.times();
148 void postProcess(word folder, word heatFluxFieldName, scalar innerField = 0.0)
167 word outputFolder =
"./ITHACAoutput/true/";
180 Eigen::MatrixXd mass(Nbasis, Nbasis);
181 Eigen::VectorXd
source(Nbasis);
192 return mass.fullPivLu().solve(
source);
204 Eigen::VectorXd coeff(Nbasis);
205 Eigen::MatrixXd Phi = Eigen::MatrixXd::Zero(Nbasis, Nbasis);
214 double radius = Foam::sqrt((tcX - centerX) * (tcX - centerX) +
215 (tcZ - centerZ) * (tcZ - centerZ));
222 return Phi.fullPivLu().solve(coeff);
316 void postProcess(word folder, word heatFluxFieldName, scalar innerField = 0.0)
335 word outputFolder =
"./ITHACAoutput/true/";
forAll(example_CG.gList, solutionI)
Class where the first inverse heat transfer problem tutorial is solved using Alifanov's regularizatio...
label gammaEx3_ind
gammaEx3 boundary patches label
label gammaEx4_ind
gammaEx4 boundary patches label
IHTP01inverseLaplacian_CG(int argc, char *argv[])
double a
Constant to define boundary conditions.
volScalarField & deltaT
Sensitivity field.
PtrList< volScalarField > gField
List of estimated heat fluxes used for post processing.
volScalarField & T
Temperature field.
double b
Constant to define boundary conditions.
List< scalar > heatFlux_gammaEx1
Boundary heat flux at the gammaEx1 patch.
void assignDirectBC()
Assign the correct boundary condition to direct problem and assign -g/k to the hotSide patch.
void postProcess(word folder, word heatFluxFieldName, scalar innerField=0.0)
Post process all the estimated heat fluxes heatFluxFieldName inside folder computing the relative err...
label gammaEx2_ind
gammaEx2 boundary patches label
double d
Constant to define boundary conditions.
double c
Constant to define boundary conditions.
List< scalar > heatFlux_gammaEx4
Boundary heat flux at the gammaEx4 patch.
label gammaEx1_ind
gammaEx1 boundary patches label
volScalarField & lambda
Adjoint field.
List< scalar > heatFlux_gammaEx2
Boundary heat flux at the gammaEx2 patch.
void solveTrue()
Perform a solution of the direct problem with the correct boundary conditions.
List< scalar > heatFlux_gammaEx3
Boundary heat flux at the gammaEx3 patch.
Class where the first inverse heat transfer problem tutorial is solved using the Parameterization of ...
double a
Constant to define boundary conditions.
double c
Constant to define boundary conditions.
volScalarField & T
Temperature field.
IHTP01inverseLaplacian_paramBC(int argc, char *argv[])
label gammaEx3_ind
gammaEx3 boundary patches label
label gammaEx1_ind
gammaEx1 boundary patches label
void solveTrue()
Perform a solution of the direct problem with the correct boundary conditions.
label gammaEx2_ind
gammaEx2 boundary patches label
void assignDirectBC()
Assign the correct boundary condition to direct problem and assign -g/k to the hotSide patch.
Eigen::VectorXd bestInterpolation()
Compute the best interpolation of the true heat flux gTrue in the parameterized space.
double d
Constant to define boundary conditions.
void postProcess(word folder, word heatFluxFieldName, scalar innerField=0.0)
Post process all the estimated heat fluxes heatFluxFieldName inside folder computing the relative err...
List< scalar > heatFlux_gammaEx2
Boundary heat flux at the gammaEx2 patch.
List< scalar > heatFlux_gammaEx4
Boundary heat flux at the gammaEx4 patch.
Eigen::VectorXd bestApproximator()
Compute the best approximation of the true heat flux gTrue in the parameterized space.
label gammaEx4_ind
gammaEx4 boundary patches label
List< scalar > heatFlux_gammaEx3
Boundary heat flux at the gammaEx3 patch.
PtrList< volScalarField > gField
List of estimated heat fluxes used for post processing.
List< scalar > heatFlux_gammaEx1
Boundary heat flux at the gammaEx1 patch.
double b
Constant to define boundary conditions.
void solveAdditional()
Solve the additional problem.
Implementation of a Alifanov's regularization to solve inverse Laplacian problems.
bool interpolationPlaneDefined
Flag for definition of interpolation plane.
int cgIter
Conjugate Gradient (CG) interations counter.
autoPtr< volScalarField > _deltaT
Sensibility temperature field.
autoPtr< volScalarField > _lambda
Adjoint field.
Implementation of a parameterization method to solve inverse Laplacian problems.
Eigen::VectorXd addSol
Solution of the additional problem at the thermocouples positions.
List< List< scalar > > gBaseFunctions
Basis of the heat flux parameterization.
double shapeParameter
RBF shape parameter.
word folderOffline
Folder where the offline solutions and matrices are saved.
PtrList< volScalarField > Tad_base
Solution of the additional problem.
label hotSide_ind
Index of the hotSide patch.
dimensionedScalar DT
Dummy thermal conductivity with unitary value.
List< vector > thermocouplesPos
List containing the positions of the thermocouples.
List< scalar > Tf
Temperature at coldSide [K].
autoPtr< simpleControl > _simple
simpleControl
bool thermocouplesRead
Flag to know if thermocouples file was read.
void restart()
Restart fields.
List< scalar > valueFraction
Value fraction for the Robin BC.
autoPtr< Time > _runTime
Time.
scalar homogeneousBC
Homogenerous BC.
void set_valueFraction()
Set valueFraction list values for Robin condition.
label coldSide_ind
Index of the coldSide patch.
autoPtr< fvMesh > _mesh
Mesh.
Eigen::VectorXd fieldValueAtThermocouples(volScalarField &field)
Interpolates the field value at the thermocouples points.
List< scalar > gTrue
True heat flux at hotSide [w/m2].
List< scalar > refGrad
Reference gradient for the Robin BC.
double k
Thermal diffusivity [W/(m K)].
autoPtr< volScalarField > _T
Temperature field.
Eigen::MatrixXd source
Source vector.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void assignMixedBC(GeometricField< Type, fvPatchField, volMesh > &field, label BC_ind, List< Type > &value, List< Type > &grad, List< scalar > &valueFrac)
Assign value of a boundary condition of type "mixed".
void assignBC(GeometricField< scalar, fvPatchField, volMesh > &s, label BC_ind, double value)
Assign uniform Boundary Condition to a volScalarField.
double L2productOnPatch(fvMesh &mesh, List< scalar > &field1, List< scalar > &field2, word patch)
Evaluate the L2 inner product between two scalarLists.
simpleControl simple(mesh)