Loading...
Searching...
No Matches
20incrementalPOD.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 of a heat transfer Reduction Problem
26SourceFiles
27 20incrementalPOD.C
28\*---------------------------------------------------------------------------*/
29
30#include <iostream>
31#include "fvCFD.H"
32#include "IOmanip.H"
33#include "Time.H"
34#include "laplacianProblem.H"
35#include "ReducedLaplacian.H"
36#include "ITHACAPOD.H"
37#include "ITHACAutilities.H"
38#include "incrementalPOD.H"
39#include <Eigen/Dense>
40#define _USE_MATH_DEFINES
41#include <cmath>
42
47{
48 public:
49 explicit tutorialIPOD(int argc, char* argv[])
50 :
51 laplacianProblem(argc, argv),
52 T(_T()),
53 nu(_nu()),
54 S(_S())
55 {}
56
59 volScalarField& T;
61 volScalarField& nu;
63 volScalarField& S;
64
66 void offlineSolve(word folder = "./ITHACAoutput/Offline/")
67 {
68 if (offline)
69 {
72 ITHACAstream::readMatrix(folder + "/mu_samples_mat.txt");
73 }
74 else
75 {
76 List<scalar> mu_now(9);
77 scalar IF = 0;
78
79 for (label i = 0; i < mu.rows(); i++)
80 {
81 for (label j = 0; j < mu.cols() ; j++)
82 {
83 mu_now[j] = mu(i, j);
84 theta[j] = mu(i, j);
85 }
86
87 assignIF(T, IF);
88 Info << i << endl;
89 truthSolve(mu_now, folder);
90 }
91 }
92 }
93
95 void SetSource()
96 {
97 volScalarField yPos = T.mesh().C().component(vector::Y).ref();
98 volScalarField xPos = T.mesh().C().component(vector::X).ref();
100 {
101 S[counter] = Foam::sin(xPos[counter] / 0.9 * M_PI) + Foam::sin(
102 yPos[counter] / 0.9 * M_PI);
103 }
104 }
105
108 {
109 nu_list.resize(9);
110 volScalarField nu1(nu);
111 volScalarField nu2(nu);
112 volScalarField nu3(nu);
113 volScalarField nu4(nu);
114 volScalarField nu5(nu);
115 volScalarField nu6(nu);
116 volScalarField nu7(nu);
117 volScalarField nu8(nu);
118 volScalarField nu9(nu);
119 Eigen::MatrixXd Box1(2, 3);
120 Box1 << 0, 0, 0, 0.3, 0.3, 0.1;
121 Eigen::MatrixXd Box2(2, 3);
122 Box2 << 0.3, 0, 0, 0.6, 0.3, 0.1;
123 Eigen::MatrixXd Box3(2, 3);
124 Box3 << 0.6, 0, 0, 0.91, 0.3, 0.1;
125 Eigen::MatrixXd Box4(2, 3);
126 Box4 << 0, 0.3, 0, 0.3, 0.6, 0.1;
127 Eigen::MatrixXd Box5(2, 3);
128 Box5 << 0.3, 0.3, 0, 0.6, 0.6, 0.1;
129 Eigen::MatrixXd Box6(2, 3);
130 Box6 << 0.6, 0.3, 0, 0.91, 0.6, 0.1;
131 Eigen::MatrixXd Box7(2, 3);
132 Box7 << 0, 0.6, 0, 0.3, 0.91, 0.1;
133 Eigen::MatrixXd Box8(2, 3);
134 Box8 << 0.3, 0.61, 0, 0.6, 0.91, 0.1;
135 Eigen::MatrixXd Box9(2, 3);
136 Box9 << 0.6, 0.6, 0, 0.9, 0.91, 0.1;
137 ITHACAutilities::setBoxToValue(nu1, Box1, 1.0);
138 ITHACAutilities::setBoxToValue(nu2, Box2, 1.0);
139 ITHACAutilities::setBoxToValue(nu3, Box3, 1.0);
140 ITHACAutilities::setBoxToValue(nu4, Box4, 1.0);
141 ITHACAutilities::setBoxToValue(nu5, Box5, 1.0);
142 ITHACAutilities::setBoxToValue(nu6, Box6, 1.0);
143 ITHACAutilities::setBoxToValue(nu7, Box7, 1.0);
144 ITHACAutilities::setBoxToValue(nu8, Box8, 1.0);
145 ITHACAutilities::setBoxToValue(nu9, Box9, 1.0);
146 nu_list.set(0, (nu1).clone());
147 nu_list.set(1, (nu2).clone());
148 nu_list.set(2, (nu3).clone());
149 nu_list.set(3, (nu4).clone());
150 nu_list.set(4, (nu5).clone());
151 nu_list.set(5, (nu6).clone());
152 nu_list.set(6, (nu7).clone());
153 nu_list.set(7, (nu8).clone());
154 nu_list.set(8, (nu9).clone());
155 }
156
159 {
160 for (int i = 0; i < nu_list.size(); i++)
161 {
162 operator_list.append(fvm::laplacian(nu_list[i], T));
163 }
164 }
165
167 volScalarField solveFull(double _mu)
168 {
169 word folder = "./ITHACAoutput/test/";
170 scalar IF = 0;
171 List<scalar> mu_now(9);
172 volScalarField& T = _T();
173
174 for (label j = 0; j < mu.cols() ; j++)
175 {
176 mu_now[j] = _mu;
177 theta[j] = _mu;
178 }
179
180 assignIF(T, IF);
181 truthSolve(mu_now, folder);
182 return T;
183 }
184
185};
186
187
188int main(int argc, char* argv[])
189{
190 // Create the example object of the tutorialIPOD type
191 tutorialIPOD example(argc, argv);
192 // Read some parameters from file
194 example._runTime());
195 int NmodesTout = para->ITHACAdict->lookupOrDefault<int>("NmodesTout", 15);
196 int NmodesTproj = para->ITHACAdict->lookupOrDefault<int>("NmodesTproj", 10);
197 double tolleranceSVD =
198 para->ITHACAdict->lookupOrDefault<double>("tolleranceSVD", 1);
199 // Set the number of parameters
200 example.Pnumber = 9;
201 example.Tnumber = NmodesTout;
202 // Set the parameters
203 example.setParameters();
204 // Set the parameter ranges, in all the subdomains the diffusivity varies between
205 // 0.001 and 0.1
206 example.mu_range.col(0) = Eigen::MatrixXd::Ones(9, 1) * 0.001;
207 example.mu_range.col(1) = Eigen::MatrixXd::Ones(9, 1) * 0.1;
208 // Generate the Parameters
209 example.genRandPar(example.Tnumber);
210 // Set the size of the list of values that are multiplying the affine forms
211 example.theta.resize(9);
212 // Set the source term
213 example.SetSource();
214 // Compute the diffusivity field for each subdomain
215 example.compute_nu();
216 // Assemble all the operators of the affine decomposition
217 example.assemble_operator();
218 // Perform an Offline Solve
219 example.offlineSolve();
220 // Perform a POD decomposition and get the modes
221 ITHACAPOD::getModes(example.Tfield, example.Tmodes, example._T().name(),
222 example.podex, 0, 0,
223 NmodesTout);
224 // Set up the incremental POD space
225 scalarIncrementalPOD IPOD(example.Tfield[0], tolleranceSVD, "L2");
226
227 // Fill the incremental POD space
228 for (int fieldI = 1; fieldI < example.Tfield.size(); fieldI++)
229 {
230 IPOD.addSnapshot(example.Tfield[fieldI]);
231 }
232
233 IPOD.writeModes();
234 // Post processing
235 word folder = "./ITHACAoutput/testReconstruction";
236 // Compute new full order solution
237 volScalarField Tfull(example.solveFull(0.05));
238 PtrList<volScalarField> TfullList;
239 TfullList.append(Tfull.clone());
240 PtrList<volScalarField> Tproj;
241 // Project the full order solution onto the POD space
242 example.Tmodes.projectSnapshots(TfullList, Tproj, NmodesTproj);
243 ITHACAstream::exportSolution(Tfull, "1", folder, "Tfull");
244 ITHACAstream::exportSolution(Tproj[0], "1", folder, "Tpod");
245 // Compute the relative error between POD projected field and full order snapshot
246 double EPS = 1e-16;
247 volScalarField relativeErrorField(Tproj[0]);
248
249 for (label i = 0; i < relativeErrorField.internalField().size(); i++)
250 {
251 if (std::abs(Tfull.ref()[i]) < EPS)
252 {
253 relativeErrorField.ref()[i] = (std::abs(Tfull.ref()[i] - Tproj[0].ref()[i])) /
254 EPS;
255 }
256 else
257 {
258 relativeErrorField.ref()[i] = (std::abs(Tfull.ref()[i] - Tproj[0].ref()[i])) /
259 Tfull.ref()[i];
260 }
261 }
262
263 ITHACAstream::exportSolution(relativeErrorField,
264 "1", folder,
265 "relativeErrorField_POD");
266 Info << "Relative error L2 norm POD = " << ITHACAutilities::L2Norm(
267 relativeErrorField) << endl;
268 // Project the full order solution onto the incremental POD space
269 IPOD.projectSnapshots(TfullList, Tproj);
270 ITHACAstream::exportSolution(Tproj[0], "1", folder, "Tipod");
271 volScalarField Tipod = Tproj[0];
272
273 // Compute the relative error between incremental POD projected field and full order snapshot
274 for (label i = 0; i < relativeErrorField.internalField().size(); i++)
275 {
276 if (std::abs(Tfull.ref()[i]) < EPS)
277 {
278 relativeErrorField.ref()[i] = (std::abs(Tfull.ref()[i] - Tipod.ref()[i])) / EPS;
279 }
280 else
281 {
282 relativeErrorField.ref()[i] = (std::abs(Tfull.ref()[i] - Tipod.ref()[i])) /
283 Tfull.ref()[i];
284 }
285 }
286
287 ITHACAstream::exportSolution(relativeErrorField,
288 "1", folder,
289 "relativeErrorField_IPOD");
290 Info << "\n\nRelative error L2 norm incrementalPOD = " <<
291 ITHACAutilities::L2Norm(relativeErrorField) << endl;
292}
293
294//--------
298
int main(int argc, char *argv[])
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Header file of the ITHACAPOD class.
Header file of the ITHACAutilities namespace.
ITHACAparameters * para
Definition NLsolve.H:40
Header file of the reducedLaplacian class.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
void projectSnapshots(PtrList< GeometricField< Type, PatchField, GeoMesh > > snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &projSnapshots, label numberOfModes=0, word innerProduct="L2")
Function to project a list of fields into the modes manifold.
Definition Modes.C:410
void addSnapshot(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Add a snapshot to the POD space.
void writeModes()
Write to modes to file.
void projectSnapshots(PtrList< GeometricField< Type, PatchField, GeoMesh > > snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &projSnapshots, label numberOfModes=0)
Project a list of snapshots onto the reduced space.
PtrList< volScalarField > nu_list
Nu (diffusivity)
List< scalar > theta
Theta (coefficients of the affine expansion)
PtrList< fvScalarMatrix > operator_list
List of operators.
autoPtr< volScalarField > _T
Temperature field.
autoPtr< volScalarField > _nu
Diffusivity.
autoPtr< fvMesh > _mesh
Mesh.
volScalarModes Tmodes
List of POD modes.
autoPtr< volScalarField > _S
Source Term.
PtrList< volScalarField > Tfield
List of snapshots for the solution.
autoPtr< Time > _runTime
Time.
label Pnumber
Number of parameters.
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
label counter
Counter used for the output of the full order solutions.
void assignIF(T &s, G &value)
Assign internal field condition.
label Tnumber
Dimension of the training set (used only when gerating parameters without input)
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
void genRandPar()
Generate Random Numbers.
Eigen::MatrixXd mu
Row matrix of parameters.
Eigen::MatrixXd mu_range
Range of the parameter spaces.
void setParameters()
Set Parameters Problems.
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
Class where the tutorial number 20 is implemented.
tutorialIPOD(int argc, char *argv[])
void compute_nu()
Compute the diffusivity in each subdomain.
volScalarField & S
Source term field.
volScalarField & nu
Diffusivity field.
void SetSource()
Define the source term function.
volScalarField & T
[tutorialIPOD] Temperature field
volScalarField solveFull(double _mu)
Performs a full order solution for a uniform value of mu.
void assemble_operator()
Construct the operator_list where each term of the affine decomposition is stored.
void offlineSolve(word folder="./ITHACAoutput/Offline/")
It perform an offline Solve.
Header file of the incrementalPOD class.
incrementalPOD< scalar, fvPatchField, volMesh > scalarIncrementalPOD
Header file of the laplacianProblem class.
void getModes(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Computes the bases or reads them for a field.
Definition ITHACAPOD.C:90
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
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 L2Norm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:41
void setBoxToValue(GeometricField< Type, fvPatchField, volMesh > &field, Eigen::MatrixXd Box, Type value)
Set value of a volScalarField to a constant inside a given box.
label i
Definition pEqn.H:46
volScalarField xPos
volScalarField yPos