Loading...
Searching...
No Matches
inverseLaplacianProblem_paramBC.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_paramBC
26Description
27 Implementation of the parameterization method for the solution of the
28 inverse problem of estimating the boundary heat flux, given
29 pointwise temperature measurements, in a Laplacian problem
30SourceFiles
31 inverseLaplacianProblem_paramBC.C
32\*---------------------------------------------------------------------------*/
33
38
39#ifndef inverseLaplacianProblem_paramBC_H
40#define inverseLaplacianProblem_paramBC_H
41
43
44#define _USE_MATH_DEFINES
45
47{
52
53
55
56 metaData_offline(int _numberTC, int _numberBasis, word _basisType,
57 double _shapeParameter)
58 : numberTC(_numberTC), numberBasis(_numberBasis),
59 basisType(_basisType), shapeParameter(_shapeParameter) {}
60};
61
62// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63
64/*---------------------------------------------------------------------------*\
65 Class inverseLaplacianProblem_paramBC Declaration
66\*---------------------------------------------------------------------------*/
67
69
73{
74
75 public:
76 // Constructors
80 inverseLaplacianProblem_paramBC(int argc, char* argv[]);
83
85 word folderOffline = "./ITHACAoutput/offlineParamBC/";
86
88 PtrList<volScalarField> Tbasis;
89
90
92 PtrList<volScalarField> Tad_base;
93
95 Eigen::VectorXd residual;
96
98 Eigen::MatrixXd Theta;
99
101 Eigen::MatrixXd gPODmodes;
102
104 Eigen::VectorXd addSol;
105
107 bool offlineReady = 0;
108
110 List<List<scalar>> gBaseFunctions;
111
113 List<scalar> gWeights;
114
117
120
121
122 // Functions
123
124 //--------------------------------------------------------------------------
131 virtual void set_gBaseFunctions();
132
133 //--------------------------------------------------------------------------
138 void set_gBaseFunctionsPOD(label Nmodes);
139
140 //--------------------------------------------------------------------------
146 void set_gParametrized(word baseFuncType, scalar _shapeParameter = 1);
147
148 //--------------------------------------------------------------------------
153 void update_gParametrized(List<scalar> weigths);
154
155 //--------------------------------------------------------------------------
162 void parameterizedBCoffline(bool force = 0);
163
164 //--------------------------------------------------------------------------
170 Eigen::VectorXd parameterizedBC(word linSys_solver = "fullPivLU",
171 double regPar = 0);
172
173 //--------------------------------------------------------------------------
178 void parameterizedBCpostProcess(Eigen::VectorXd weigths);
179
180 //--------------------------------------------------------------------------
183 virtual void solveAdditional();
184
185 //--------------------------------------------------------------------------
188 void reconstructT();
189};
190
191#endif
Implementation of a parameterization method to solve inverse Laplacian problems.
void update_gParametrized(List< scalar > weigths)
Update the boundary heat flux.
Eigen::VectorXd addSol
Solution of the additional problem at the thermocouples positions.
Eigen::VectorXd parameterizedBC(word linSys_solver="fullPivLU", double regPar=0)
Solve the online phase.
virtual void set_gBaseFunctions()
Define the base functions used for the parametrization of the heat flux g The center of each base fun...
Eigen::MatrixXd Theta
Theta matrix of the lynear system.
List< List< scalar > > gBaseFunctions
Basis of the heat flux parameterization.
virtual void solveAdditional()
Set BC the additional problem and solves it.
void parameterizedBCpostProcess(Eigen::VectorXd weigths)
Reconstruct the temperature field and compute the cost function J.
word baseFuncType
Type of basis functions used for the parameterization of the heat flux.
void set_gBaseFunctionsPOD(label Nmodes)
Performs POD on the RBF basis functions.
PtrList< volScalarField > Tbasis
Solutions of the direct problem with the basis of the parameterization as boundary heat flux.
List< scalar > gWeights
Weights of the heat flux parameterization.
void set_gParametrized(word baseFuncType, scalar _shapeParameter=1)
Set initial heat flux for the parameterized BC method.
word folderOffline
Folder where the offline solutions and matrices are saved.
PtrList< volScalarField > Tad_base
Solution of the additional problem.
Eigen::VectorXd residual
Residual of the linear system.
void reconstructT()
Reconstuct the field T using the offline computed fields.
Eigen::MatrixXd gPODmodes
Modes of the POD.
bool offlineReady
Flag to know if the offline phase was computed.
Implementation of a inverse Laplacian problem .
Header file of the inverseLaplacianProblem class.
example_paramBC parameterizedBCoffline()
metaData_offline(int _numberTC, int _numberBasis, word _basisType, double _shapeParameter)