Loading...
Searching...
No Matches
IHTP01inverseLaplacian.H
Go to the documentation of this file.
1#ifndef IHTP01inverseLaplacian_H
2#define IHTP01inverseLaplacian_H
3
10{
11 public:
12 explicit IHTP01inverseLaplacian_paramBC(int argc, char* argv[])
13 :
15 T(_T()),
16 mesh(_mesh()),
18 {
19 hotSide_ind = mesh.boundaryMesh().findPatchID("hotSide");
20 coldSide_ind = mesh.boundaryMesh().findPatchID("coldSide");
22 }
23
25 volScalarField& T;
26
28 fvMesh& mesh;
29
31 Time& runTime;
32
34 double a;
35
37 double b;
38
40 double c;
41
43 double d;
44
46 List<scalar> heatFlux_gammaEx1;
47
49 List<scalar> heatFlux_gammaEx2;
50
52 List<scalar> heatFlux_gammaEx3;
53
55 List<scalar> heatFlux_gammaEx4;
56
59
62
65
68
70 PtrList<volScalarField> gField;
71
72
73 //--------------------------------------------------------------------------
77 {
78 restart();
79 fvMesh& mesh = _mesh();
80 simpleControl& simple = _simple();
81 volScalarField Tad(_T());
82 Foam::Time& runTime = _runTime();
84 List<scalar> RobinBC = - Tf;
85 forAll(mesh.boundaryMesh(), patchI)
86 {
87 if (patchI == mesh.boundaryMesh().findPatchID("coldSide"))
88 {
90 }
91 else if (patchI == mesh.boundaryMesh().findPatchID("hotSide"))
92 {
94 }
95 else if (patchI == mesh.boundaryMesh().findPatchID("gammaEx1"))
96 {
98 }
99 else if (patchI == mesh.boundaryMesh().findPatchID("gammaEx2"))
100 {
102 }
103 else if (patchI == mesh.boundaryMesh().findPatchID("gammaEx3"))
104 {
106 }
107 else if (patchI == mesh.boundaryMesh().findPatchID("gammaEx4"))
108 {
110 }
111 }
112#if defined(OFVER) && (OFVER == 6)
113
114 while (simple.loop(runTime))
115#else
116 while (simple.loop())
117#endif
118 {
119 while (simple.correctNonOrthogonal())
120 {
121 fvScalarMatrix TEqn
122 (
123 fvm::laplacian(DT, Tad)
124 );
125 TEqn.solve();
126 }
127 }
128
129 Tad_base.append(Tad.clone());
130 //Reinitializing runTime
131 instantList Times = runTime.times();
132 runTime.setTime(Times[1], 1);
136 "Tad");
137 }
138
139 //--------------------------------------------------------------------------
148 void postProcess(word folder, word heatFluxFieldName, scalar innerField = 0.0)
149 {
150#include"postProcess.H"
151 }
152
153
154 //--------------------------------------------------------------------------
158 {
159#include"directBC.H"
160 }
161
162 //--------------------------------------------------------------------------
166 {
167 word outputFolder = "./ITHACAoutput/true/";
168#include"solveTrue.H"
169 }
170
171 //--------------------------------------------------------------------------
176 Eigen::VectorXd bestApproximator()
177 {
178 fvMesh& mesh = _mesh();
179 int Nbasis = gBaseFunctions.size();
180 Eigen::MatrixXd mass(Nbasis, Nbasis);
181 Eigen::VectorXd source(Nbasis);
182 forAll(gBaseFunctions, baseI)
183 {
184 forAll(gBaseFunctions, baseJ)
185 {
186 mass(baseI, baseJ) = ITHACAutilities::L2productOnPatch(mesh,
187 gBaseFunctions[baseI], gBaseFunctions[baseJ], "hotSide");
188 }
190 gTrue, "hotSide");
191 }
192 return mass.fullPivLu().solve(source); //Vector of coefficients
193 }
194
195
196 //--------------------------------------------------------------------------
201 Eigen::VectorXd bestInterpolation()
202 {
203 int Nbasis = gBaseFunctions.size();
204 Eigen::VectorXd coeff(Nbasis);
205 Eigen::MatrixXd Phi = Eigen::MatrixXd::Zero(Nbasis, Nbasis);
206 forAll(gBaseFunctions, baseI)
207 {
208 scalar centerX = thermocouplesPos[baseI].x();
209 scalar centerZ = thermocouplesPos[baseI].z();
211 {
212 scalar tcX = thermocouplesPos[TCi].x();
213 scalar tcZ = thermocouplesPos[TCi].z();
214 double radius = Foam::sqrt((tcX - centerX) * (tcX - centerX) +
215 (tcZ - centerZ) * (tcZ - centerZ));
216 Phi(TCi, baseI) = Foam::exp(- (shapeParameter *
218 * radius * radius));
219 }
220 coeff(baseI) = k * (b * thermocouplesPos[baseI].x() + c);
221 }
222 return Phi.fullPivLu().solve(coeff);
223 }
224
225};
226
227
228
235{
236 public:
237 explicit IHTP01inverseLaplacian_CG(int argc, char* argv[])
238 :
239 inverseLaplacianProblem_CG(argc, argv),
240 T(_T()),
241 lambda(_lambda()),
242 deltaT(_deltaT()),
243 mesh(_mesh()),
245 {
246 hotSide_ind = mesh.boundaryMesh().findPatchID("hotSide");
247 coldSide_ind = mesh.boundaryMesh().findPatchID("coldSide");
249 cgIter = 0;
251 }
252
254 volScalarField& T;
255
257 volScalarField& lambda;
258
260 volScalarField& deltaT;
261
263 fvMesh& mesh;
264
266 Time& runTime;
267
269 double a;
270
272 double b;
273
275 double c;
276
278 double d;
279
281 List<scalar> heatFlux_gammaEx1;
282
284 List<scalar> heatFlux_gammaEx2;
285
287 List<scalar> heatFlux_gammaEx3;
288
290 List<scalar> heatFlux_gammaEx4;
291
294
297
300
303
305 PtrList<volScalarField> gField;
306
307 //--------------------------------------------------------------------------
316 void postProcess(word folder, word heatFluxFieldName, scalar innerField = 0.0)
317 {
318#include"postProcess.H"
319 }
320
321
322 //--------------------------------------------------------------------------
326 {
327#include"directBC.H"
328 }
329
330 //--------------------------------------------------------------------------
334 {
335 word outputFolder = "./ITHACAoutput/true/";
336#include"solveTrue.H"
337 }
338};
339#endif
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
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.
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.
List< scalar > valueFraction
Value fraction for the Robin BC.
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.
fvScalarMatrix & TEqn
Definition TEqn.H:15
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)