Loading...
Searching...
No Matches
inverseLaplacianProblemTotalHeatMeasure_paramBC.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-------------------------------------------------------------------------------
12
13 License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
33
34
36
37// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
38
39// Constructors
41
48
49// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
50
51
52
54 bool force)
55{
56 fvMesh& mesh = _mesh();
57 Tbasis.resize(0);
58 Tad_base.resize(0);
59 char recomputeOffline;
60
61 if (ITHACAutilities::check_file(folderOffline + "Theta_mat.txt") && force == 0)
62 {
63 do
64 {
65 metaData_offline metaData;
66 std::ifstream fin(folderOffline + "metaData.txt");
67 fin >> metaData.numberTC >> metaData.numberBasis >>
68 metaData.basisType >> metaData.shapeParameter;
69 fin.close();
70 std::cout << "\nOffline FOUND with parameter:\n" <<
71 "Number of thermocouples = " << metaData.numberTC <<
72 "\nNumber of basis functions = " << metaData.numberBasis <<
73 "\nType of basis functions = " << metaData.basisType <<
74 "\nRBF shape parameters = " << metaData.shapeParameter <<
75 "\n\nShould I recompute it? [y/n]" << std::endl;
76 std::cin >> recomputeOffline;
77 }
78 while ( !cin.fail() && recomputeOffline != 'y' && recomputeOffline != 'n' );
79 }
80
81 if (recomputeOffline == 'y')
82 {
83 force = 1;
84 }
85
86 if (ITHACAutilities::check_file(folderOffline + "Theta_mat.txt") && force == 0)
87 {
88 Info << "\nOffline already computed." << endl;
89 Info << "Check that the basis used for the parameterized BC are correct (RBF, POD, etc.)\n";
90 Theta = ITHACAstream::readMatrix(folderOffline + "Theta_mat.txt");
93 addSol = ITHACAstream::readMatrix(folderOffline + "addSol_mat.txt");
94 volScalarField& T(_T());
98 }
99 else
100 {
101 Info << "\nComputing offline" << endl;
104 M_Assert(Tmeas.size() > 0, "Initialize Tmeas");
105 M_Assert(gWeights.size() > 0, "Initialize gWeights");
106 Theta.resize(Tmeas.size(), gWeights.size());
107 metaData_offline metaData(Tmeas.size(), gWeights.size(), baseFuncType,
109 std::ofstream fout(folderOffline + "metaData.txt");
110 fout << metaData.numberTC << ' ' <<
111 metaData.numberBasis << ' ' <<
112 metaData.basisType << ' ' <<
113 metaData.shapeParameter << ' ';
114 fout.close();
115 int Nbasis = Theta.cols();
116 Phi.resize(gWeights.size(), gWeights.size());
117 phi.resize(Nbasis);
118
119 for (label j = 0; j < Nbasis; j++)
120 {
121 gWeights = Foam::zero();
122 gWeights[j] = 1;
124 Info << "Solving for j = " << j << endl;
125 solveDirect();
126 phi(j) = ITHACAutilities::integralOnPatch(mesh, g, "hotSide");
127 volScalarField& T = _T();
128 Tbasis.append(T.clone());
130
131 for (label i = 0; i < Theta.rows(); i++)
132 {
133 Theta(i, j) = Tdirect(i) + addSol(i);
134 }
135
136 volScalarField gParametrizedField = list2Field(g);
137 ITHACAstream::exportSolution(gParametrizedField, std::to_string(j + 1),
139 "gParametrized");
140 }
141
142 for (label baseI = 0; baseI < Nbasis; baseI++)
143 {
144 for (label baseJ = 0; baseJ < Nbasis; baseJ++)
145 {
146 Phi(baseI, baseJ) = phi(baseI) * phi(baseJ);
147 }
148 }
149
154 Info << "\nOffline part ENDED\n" << endl;
155 }
156
157 Eigen::MatrixXd A = Theta.transpose() * Theta + gIntegralWeight * Phi;
158 Eigen::JacobiSVD<Eigen::MatrixXd> svd(A,
159 Eigen::ComputeThinU | Eigen::ComputeThinV);
160 Eigen::MatrixXd singularValues = svd.singularValues();
161
162 if (singularValues.minCoeff() > 0)
163 {
164 double conditionNumber = singularValues.maxCoeff() / singularValues.minCoeff();
165 Info << "Condition number = " << conditionNumber << endl;
166 }
167
168 ITHACAstream::exportMatrix(singularValues, "singularValues", "eigen",
170}
171
172Eigen::VectorXd
174 word linSys_solver,
175 double regPar)
176{
177 Info << endl << "Using quasilinearity of direct problem" << endl;
178 List<Eigen::MatrixXd> linSys;
179 linSys.resize(2);
180 Info << "debug: Theta size = " << Theta.rows() << ", " << Theta.cols() << endl;
181 linSys[0] = Theta.transpose() * Theta + gIntegralWeight * Phi;
182 linSys[1] = gIntegralWeight * gIntegral_meas * phi + Theta.transpose() *
183 (Tmeas + addSol);
184 Eigen::VectorXd weigths;
185 M_Assert(std::abs(gIntegral_meas) > 1e-16, "First set up gIntegral_meas");
186 M_Assert(std::abs(gIntegralWeight) > 1e-16, "First set up gIntegralWeight");
187
188 if (linSys_solver == "fullPivLU")
189 {
190 weigths = linSys[0].fullPivLu().solve(linSys[1]);
191 }
192 else if (linSys_solver == "jacobiSvd")
193 {
194 Eigen::JacobiSVD<Eigen::MatrixXd> svd(linSys[0],
195 Eigen::ComputeThinU | Eigen::ComputeThinV);
196 weigths = svd.solve(linSys[1]);
197 }
198 else if (linSys_solver == "householderQr")
199 {
200 weigths = linSys[0].householderQr().solve(linSys[1]);
201 }
202 else if (linSys_solver == "ldlt")
203 {
204 weigths = linSys[0].ldlt().solve(linSys[1]);
205 }
206 else if (linSys_solver == "inverse")
207 {
208 weigths = linSys[0].inverse() * linSys[1];
209 }
210 else if (linSys_solver == "TSVD")
211 {
212 weigths = ITHACAregularization::TSVD(linSys[0], linSys[1], int(regPar));
213 }
214 else if (linSys_solver == "Tikhonov")
215 {
216 weigths = ITHACAregularization::Tikhonov(linSys[0], linSys[1], regPar);
217 }
218 else
219 {
220 Info << "Select a linear system solver in this list:" << endl
221 << "fullPivLU, jacobiSvd, householderQr, ldlt, inverse, TSVD" << endl;
222 exit(1);
223 }
224
226 return weigths;
227}
Foam::fvMesh & mesh
Definition createMesh.H:47
#define M_Assert(Expr, Msg)
void parameterizedBCoffline(bool force=0)
Performs offline computation for the parameterized BC method, if the offline directory "".
Eigen::VectorXd parameterizedBC(word linSys_solver="fullPivLU", double regPar=0)
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::MatrixXd Theta
Theta matrix of the lynear system.
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.
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.
word folderOffline
Folder where the offline solutions and matrices are saved.
PtrList< volScalarField > Tad_base
Solution of the additional problem.
List< scalar > g
Heat flux at hotSide [W/m2].
Eigen::VectorXd Tdirect
Vector of computed temperatures at the thermocouples locations [K].
volScalarField list2Field(List< scalar > list, scalar innerField=0.0)
Create a field with the hotSide boundary heat flux at the hotSide bounday cells for visualization.
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.
void solveDirect()
Solve direct problem.
autoPtr< volScalarField > _T
Temperature field.
volScalarField & T
Definition createT.H:46
Header file of the inverseLaplacianProblemTotalHeatMeasure_paramBC class.
volScalarField & A
Eigen::VectorXd Tikhonov(Eigen::MatrixXd A, Eigen::MatrixXd b, double regularizationParameter)
Truncated Singular Value regularization.
Eigen::VectorXd TSVD(Eigen::MatrixXd A, Eigen::MatrixXd b, int filter)
Truncated Singular Value regularization.
void exportFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &field, word folder, word fieldname)
Function to export a scalar of vector field.
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 exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
List< Eigen::MatrixXd > readMatrix(word folder, word mat_name)
Read a three dimensional matrix from a txt file in Eigen format.
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
double integralOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the integral on a patch.
bool check_file(std::string fileName)
Function that returns true if a file exists.
Info<< endl;Info<< "*********************************************************"<< endl;Info<< "Performing test for the parameterized BC inverse solver"<< endl;Info<< endl;word outputFolder="./ITHACAoutput/parameterizedBCtest_RBFparameter/";volScalarField gTrueField=example_paramBC.list2Field(example_paramBC.gTrue);ITHACAstream::exportSolution(gTrueField, "1", outputFolder, "gTrue");Eigen::VectorXd rbfWidth=EigenFunctions::ExpSpaced(0.001, 100, rbfWidthTest_size);ITHACAstream::exportMatrix(rbfWidth, "rbfShapeParameters", "eigen", outputFolder);Eigen::VectorXd residualNorms;residualNorms.resize(rbfWidthTest_size);Eigen::VectorXd heatFluxL2norm(rbfWidthTest_size);Eigen::VectorXd heatFluxLinfNorm=heatFluxL2norm;Eigen::VectorXd condNumber=heatFluxL2norm;Eigen::MatrixXd singVal;for(int i=0;i< rbfWidthTest_size;i++){ Info<< "*********************************************************"<< endl;Info<< "RBF parameter "<< rbfWidth(i)<< endl;Info<< "Test "<< i+1<< endl;Info<< endl;example_paramBC.set_gParametrized("rbf", rbfWidth(i));example_paramBC.parameterizedBCoffline(1);example_paramBC.parameterizedBC("fullPivLU");volScalarField gParametrizedField=example_paramBC.list2Field(example_paramBC.g);ITHACAstream::exportSolution(gParametrizedField, std::to_string(i+1), outputFolder, "gParametrized");volScalarField &T(example_paramBC._T());ITHACAstream::exportSolution(T, std::to_string(i+1), outputFolder, "T");residualNorms(i)=Foam::sqrt(example_paramBC.residual.squaredNorm());Eigen::MatrixXd A=example_paramBC.Theta.transpose() *example_paramBC.Theta;Eigen::JacobiSVD< Eigen::MatrixXd > svd(A, Eigen::ComputeThinU|Eigen::ComputeThinV)
Eigen::MatrixXd singularValues
double conditionNumber
label i
Definition pEqn.H:46