Loading...
Searching...
No Matches
IHTP01inverseLaplacian.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Description
25 example_paramBC.of a heat transfer Reduction Problem
26SourceFiles
27 IHTP01inverseLaplacian.C
28\*---------------------------------------------------------------------------*/
29
30#include <iostream>
31#include <float.h>
32#include "interpolation.H"
33#include "fvCFD.H"
34#include "fvOptions.H"
35#include "simpleControl.H"
36#include "IOmanip.H"
37#include "Time.H"
40#include "ITHACAPOD.H"
41#include "ITHACAutilities.H"
42#include <Eigen/Dense>
43#define _USE_MATH_DEFINES
44#include <cmath>
45#include "Foam2Eigen.H"
46#include "mixedFvPatchFields.H"
47#include "cellDistFuncs.H"
48
50
51
52int main(int argc, char* argv[])
53{
54 solverPerformance::debug = 1; //No verbose output
55 IHTP01inverseLaplacian_paramBC example_paramBC(argc, argv);
56 IHTP01inverseLaplacian_CG example_CG(argc, argv);
57 //Setting parameters for the analytical benchmark
58 double a = 5;
59 double b = 10;
60 double c = 15;
61 double d = 20;
62 example_paramBC.a = a;
63 example_paramBC.b = b;
64 example_paramBC.c = c;
65 example_paramBC.d = d;
66 example_CG.a = a;
67 example_CG.b = b;
68 example_CG.c = c;
69 example_CG.d = d;
70 // Reading tests to perform
72 example_paramBC._runTime());
73 label CGtest = para->ITHACAdict->lookupOrDefault<int>("CGtest", 0);
74 label CGnoiseTest = para->ITHACAdict->lookupOrDefault<int>("CGnoiseTest", 0);
75 label CGnoiseLevelTest =
76 para->ITHACAdict->lookupOrDefault<int>("CGnoiseLevelTest", 0);
77 label parameterizedBCtest =
78 para->ITHACAdict->lookupOrDefault<int>("parameterizedBCtest", 0);
79 label parameterizedBCtest_RBFwidth =
80 para->ITHACAdict->lookupOrDefault<int>("parameterizedBCtest_RBFwidth", 0);
81 label thermocouplesLocationTest_CG =
82 para->ITHACAdict->lookupOrDefault<int>("thermocouplesLocationTest_CG", 0);
83 label thermocouplesLocationTest_paramBC =
84 para->ITHACAdict->lookupOrDefault<int>("thermocouplesLocationTest_paramBC", 0);
85 label thermocouplesNumberTest_CG =
86 para->ITHACAdict->lookupOrDefault<int>("thermocouplesNumberTest_CG", 0);
87 label thermocouplesNumberTest_paramBC =
88 para->ITHACAdict->lookupOrDefault<int>("thermocouplesNumberTest_paramBC", 0);
89 // Reading parameters from ITHACAdict
90 example_CG.cgIterMax = para->ITHACAdict->lookupOrDefault<int>("cgIterMax", 100);
91 example_CG.Jtol = para->ITHACAdict->lookupOrDefault<double>("Jtolerance",
92 0.000001);
93 example_CG.JtolRel =
94 para->ITHACAdict->lookupOrDefault<double>("JrelativeTolerance",
95 0.001);
96 double rbfShapePar = para->ITHACAdict->lookupOrDefault<double>("rbfShapePar",
97 0);
98 M_Assert(rbfShapePar > 0, "rbfShapePar not specified");
99 example_paramBC.k =
100 para->ITHACAdict->lookupOrDefault<double>("thermalConductivity", 0);
101 M_Assert(example_paramBC.k > 0, "thermalConductivity, k, not specified");
102 example_paramBC.H =
103 para->ITHACAdict->lookupOrDefault<double>("heatTranferCoeff", 0);
104 M_Assert(example_paramBC.H > 0, "Heat transfer coeff, H, not specified");
105 example_CG.k = example_paramBC.k;
106 example_CG.H = example_paramBC.H;
107 int rbfWidthTest_size =
108 para->ITHACAdict->lookupOrDefault<int>("rbfWidthTest_size", 0);
109 // Setting analytical solution
110 volScalarField T_true(example_paramBC._T());
111
112 for (label i = 0; i < T_true.internalField().size(); i++)
113 {
114 auto cx = T_true.mesh().C()[i].component(vector::X);
115 auto cy = T_true.mesh().C()[i].component(vector::Y);
116 auto cz = T_true.mesh().C()[i].component(vector::Z);
117 T_true.ref()[i] = a * cx * cx + b * cx * cy + c * cy - a * cz * cz + c;
118 }
119
120 // Perform true solution
121 fvMesh& mesh = example_paramBC._mesh();
122 example_paramBC.hotSide_ind = mesh.boundaryMesh().findPatchID("hotSide");
123 label hotSideSize = T_true.boundaryField()[example_paramBC.hotSide_ind].size();
124 example_paramBC.g.resize(hotSideSize);
125 example_paramBC.gTrue.resize(hotSideSize);
126 forAll(example_paramBC.g, faceI)
127 {
128 scalar faceX =
129 mesh.boundaryMesh()[example_paramBC.hotSide_ind].faceCentres()[faceI].x();
130 example_paramBC.g[faceI] = example_paramBC.k * (b * faceX + c) ;
131 }
132 example_paramBC.gTrue = example_paramBC.g;
133 example_paramBC.solveTrue();
134 example_CG.g = example_paramBC.g;
135 example_CG.gTrue = example_paramBC.g;
136 example_CG.solveTrue();
137 volScalarField& T(example_paramBC._T());
138 Info << "Exporting analytical solution" << endl;
139 ITHACAstream::exportSolution(T_true, "1", "./ITHACAoutput/true/",
140 "analyticalSol");
141 volScalarField error = (T_true - T).ref();
142 ITHACAstream::exportSolution(error, "1", "./ITHACAoutput/true/", "error");
143 Info << "L2 norm of the relative error = " << ITHACAutilities::errorL2Rel(
144 T_true, T) << endl;
145 Info << "Linf norm of the relative error = " <<
146 ITHACAutilities::errorLinfRel(T_true, T) << endl;
147 Info << "L2 norm of the error = " << ITHACAutilities::errorL2Abs(T_true,
148 T) << endl;
149 Info << "Linf norm of the error = " << ITHACAutilities::LinfNorm(error) << endl;
150 Info << "L2 norm of T = " << ITHACAutilities::L2Norm(T) << endl;
151 // Setting up The thermocouples
152 example_paramBC.readThermocouples();
153 example_paramBC.Tmeas = example_paramBC.fieldValueAtThermocouples(T_true);
154 example_CG.readThermocouples();
155 example_CG.Tmeas = example_CG.fieldValueAtThermocouples(T_true);
156
157 // Inverse problem tests
158
159 // Alifanov's regularization test
160 if (CGtest)
161 {
162#include "CGtest.H"
163 }
164
165 // Parameterized heat flux test
166 if (parameterizedBCtest)
167 {
168#include"parameterizedBCtest.H"
169 }
170
171 // Test that solves the inverse problem using the parameterization of the heat flux with different RBF shape parameters
172 if (parameterizedBCtest_RBFwidth)
173 {
175 }
176
177 // Test the Alifanov's regularization changing the distance of the thermocouples from the hotSide
178 if (thermocouplesLocationTest_CG)
179 {
181 }
182
183 // Test the parameterization method changing the distance of the thermocouples from the hotSide
184 if (thermocouplesLocationTest_paramBC)
185 {
187 }
188
189 // Test the Alifanov's regularization changing the number of thermocouples in the thermocouples plane
190 if (thermocouplesNumberTest_CG)
191 {
193 }
194
195 // Test the parameterization method changing the number of thermocouples in the thermocouples plane
196 if (thermocouplesNumberTest_paramBC)
197 {
199 }
200
201 return 0;
202}
203
204//--------
208
265
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Header file of the Foam2Eigen class.
int main(int argc, char *argv[])
Header file of the ITHACAPOD class.
Foam::fvMesh & mesh
Definition createMesh.H:47
#define M_Assert(Expr, Msg)
Header file of the ITHACAutilities namespace.
Class where the first inverse heat transfer problem tutorial is solved using Alifanov's regularizatio...
double a
Constant to define boundary conditions.
double b
Constant to define boundary conditions.
double d
Constant to define boundary conditions.
double c
Constant to define boundary conditions.
void solveTrue()
Perform a solution of the direct problem with the correct boundary conditions.
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.
void solveTrue()
Perform a solution of the direct problem with the correct boundary conditions.
double d
Constant to define boundary conditions.
double b
Constant to define boundary conditions.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
IOdictionary * ITHACAdict
Dictionary for input objects from file.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
int cgIterMax
Maximum CG iterations.
double Jtol
Absolute stopping criterion for the CG.
double JtolRel
Relative stopping criterion for the CG.
label hotSide_ind
Index of the hotSide patch.
double H
Heat transfer coefficient [W/(m2 K)].
List< scalar > g
Heat flux at hotSide [W/m2].
Eigen::VectorXd Tmeas
Vector of measured temperatures at the thermocouples locations [K].
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].
virtual void readThermocouples()
Identifies in the mesh the cells corresponding to the termocouples locations.
double k
Thermal diffusivity [W/(m K)].
autoPtr< volScalarField > _T
Temperature field.
volScalarField & T
Definition createT.H:46
Header file of the inverseLaplacianProblem_CG class.
Header file of the inverseLaplacianProblem_paramBC class.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
double LinfNorm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:57
double errorLinfRel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in Linf norm.
double L2Norm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:41
double errorL2Rel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in L2 norm.
double errorL2Abs(GeometricField< vector, fvPatchField, volMesh > &field1, GeometricField< vector, fvPatchField, volMesh > &field2, volScalarField &Volumes)
label i
Definition pEqn.H:46