Loading...
Searching...
No Matches
inverseLaplacianProblem_CG.H
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/>.
24Class
25 inverseLaplacianProblem_CG
26Description
27 Implementation of the Alifanov's regularization for the inverse problem of
28 estimating the boundary heat flux, given pointwise temperature measurements,
29 in a Laplacian problem
30SourceFiles
31 inverseLaplacianProblem_CG.C
32\*---------------------------------------------------------------------------*/
33
38
39
40#ifndef inverseLaplacianProblem_CG_H
41#define inverseLaplacianProblem_CG_H
42#include "thermocouplesPlane.H"
44#define _USE_MATH_DEFINES
45
46using namespace SPLINTER;
47
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51/*---------------------------------------------------------------------------*\
52 Class inverseLaplacianProblem_CG Declaration
53\*---------------------------------------------------------------------------*/
54
56
60{
61
62 public:
63 // Constructors
67 inverseLaplacianProblem_CG(int argc, char* argv[]);
70
71
73 autoPtr<volScalarField> _lambda;
74
76 autoPtr<volScalarField> _deltaT;
77
79 List<scalar> valueFractionAdj;
80
82 PtrList<volScalarField> Tfield;
83
85 PtrList<volScalarField> lambdaField;
86
88 PtrList<volScalarField> deltaTfield;
89
92
95
97 int cgIter;
98
101
104
106 Eigen::MatrixXd Jlist;
107
110
112 double Jtol;
113
115 double JtolRel;
116
118 double gamma;
119
121 double gamma_den;
122
124 double beta;
125
126
128 List<scalar> P;
129
131 List<scalar> gradJ;
132
134 Eigen::VectorXd Tsens;
135
137 List<Eigen::MatrixXd> ArbT;
138
140 List<Eigen::MatrixXd> ArbLambda;
141
143 List<Eigen::MatrixXd> ArbDeltaT;
144
147
149 Eigen::VectorXd cellsInPlane;
150
151 // Functions
152
153 //--------------------------------------------------------------------------
156 void set_valueFraction();
157
158 //--------------------------------------------------------------------------
161 void assignAdjointBC();
162
163 //--------------------------------------------------------------------------
168 volScalarField assignAdjointBCandSource();
169
170 //--------------------------------------------------------------------------
173 void assignSensitivityBC();
174
175 //--------------------------------------------------------------------------
178 void solveAdjoint();
179
180 //--------------------------------------------------------------------------
183 void solveSensitivity();
184
185 //--------------------------------------------------------------------------
188 void solve(const char* problemID);
189
190 //--------------------------------------------------------------------------
194
195 //--------------------------------------------------------------------------
200
201 //--------------------------------------------------------------------------
206 int conjugateGradient();
207
208 //--------------------------------------------------------------------------
211 void computeGradJ();
212
213 //--------------------------------------------------------------------------
216 void searchDirection();
217
218 //--------------------------------------------------------------------------
221 void computeSearchStep();
222
223 //--------------------------------------------------------------------------
226 void updateHeatFlux();
227
228 //--------------------------------------------------------------------------
234
235 //--------------------------------------------------------------------------
245 int isInPlane(double cx, double cy, double cz,
246 Foam::vector thermocoupleCellDim);
247
248 //--------------------------------------------------------------------------
254 void writeFields(label folderNumber, const char* folder);
255
256 //--------------------------------------------------------------------------
263
264 //--------------------------------------------------------------------------
273 void thermocouplesInterpolation(DenseMatrix& RBFweights, DenseMatrix& RBFbasis);
274
275 //--------------------------------------------------------------------------
279
280 //--------------------------------------------------------------------------
286 void restart(word fieldName = "all");
287
288
289 //--------------------------------------------------------------------------
292 //void offlineSolve();
293
294};
295
296#endif
TEqn solve()
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.
void updateHeatFlux()
Updates the heat flux in the conjugate gradient iterations.
int conjugateGradient()
Conjugate gradient method.
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 .
Header file of the inverseLaplacianProblem class.
Class defining a plane for thermocouples measurements interpolation.