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 {}
58 volScalarField& T;
60 volScalarField& nu;
62 volScalarField& S;
63
65 void offlineSolve(word folder = "./ITHACAoutput/Offline/")
66 {
67 if (offline)
68 {
71 ITHACAstream::readMatrix(folder + "/mu_samples_mat.txt");
72 }
73 else
74 {
75 List<scalar> mu_now(9);
76 scalar IF = 0;
77
78 for (label i = 0; i < mu.rows(); i++)
79 {
80 for (label j = 0; j < mu.cols() ; j++)
81 {
82 mu_now[j] = mu(i, j);
83 theta[j] = mu(i, j);
84 }
85
86 assignIF(T, IF);
87 Info << i << endl;
88 truthSolve(mu_now, folder);
89 }
90 }
91 }
92
94 void SetSource()
95 {
96 volScalarField yPos = T.mesh().C().component(vector::Y).ref();
97 volScalarField xPos = T.mesh().C().component(vector::X).ref();
99 {
100 S[counter] = Foam::sin(xPos[counter] / 0.9 * M_PI) + Foam::sin(
101 yPos[counter] / 0.9 * M_PI);
102 }
103 }
104
107 {
108 nu_list.resize(9);
109 volScalarField nu1(nu);
110 volScalarField nu2(nu);
111 volScalarField nu3(nu);
112 volScalarField nu4(nu);
113 volScalarField nu5(nu);
114 volScalarField nu6(nu);
115 volScalarField nu7(nu);
116 volScalarField nu8(nu);
117 volScalarField nu9(nu);
118 Eigen::MatrixXd Box1(2, 3);
119 Box1 << 0, 0, 0, 0.3, 0.3, 0.1;
120 Eigen::MatrixXd Box2(2, 3);
121 Box2 << 0.3, 0, 0, 0.6, 0.3, 0.1;
122 Eigen::MatrixXd Box3(2, 3);
123 Box3 << 0.6, 0, 0, 0.91, 0.3, 0.1;
124 Eigen::MatrixXd Box4(2, 3);
125 Box4 << 0, 0.3, 0, 0.3, 0.6, 0.1;
126 Eigen::MatrixXd Box5(2, 3);
127 Box5 << 0.3, 0.3, 0, 0.6, 0.6, 0.1;
128 Eigen::MatrixXd Box6(2, 3);
129 Box6 << 0.6, 0.3, 0, 0.91, 0.6, 0.1;
130 Eigen::MatrixXd Box7(2, 3);
131 Box7 << 0, 0.6, 0, 0.3, 0.91, 0.1;
132 Eigen::MatrixXd Box8(2, 3);
133 Box8 << 0.3, 0.61, 0, 0.6, 0.91, 0.1;
134 Eigen::MatrixXd Box9(2, 3);
135 Box9 << 0.6, 0.6, 0, 0.9, 0.91, 0.1;
136 ITHACAutilities::setBoxToValue(nu1, Box1, 1.0);
137 ITHACAutilities::setBoxToValue(nu2, Box2, 1.0);
138 ITHACAutilities::setBoxToValue(nu3, Box3, 1.0);
139 ITHACAutilities::setBoxToValue(nu4, Box4, 1.0);
140 ITHACAutilities::setBoxToValue(nu5, Box5, 1.0);
141 ITHACAutilities::setBoxToValue(nu6, Box6, 1.0);
142 ITHACAutilities::setBoxToValue(nu7, Box7, 1.0);
143 ITHACAutilities::setBoxToValue(nu8, Box8, 1.0);
144 ITHACAutilities::setBoxToValue(nu9, Box9, 1.0);
145 nu_list.set(0, (nu1).clone());
146 nu_list.set(1, (nu2).clone());
147 nu_list.set(2, (nu3).clone());
148 nu_list.set(3, (nu4).clone());
149 nu_list.set(4, (nu5).clone());
150 nu_list.set(5, (nu6).clone());
151 nu_list.set(6, (nu7).clone());
152 nu_list.set(7, (nu8).clone());
153 nu_list.set(8, (nu9).clone());
154 }
155
158 {
159 for (int i = 0; i < nu_list.size(); i++)
160 {
161 operator_list.append(fvm::laplacian(nu_list[i], T));
162 }
163 }
164
166 volScalarField solveFull(double _mu)
167 {
168 word folder = "./ITHACAoutput/test/";
169 scalar IF = 0;
170 List<scalar> mu_now(9);
171 volScalarField& T = _T();
172
173 for (label j = 0; j < mu.cols() ; j++)
174 {
175 mu_now[j] = _mu;
176 theta[j] = _mu;
177 }
178
179 assignIF(T, IF);
180 truthSolve(mu_now, folder);
181 return T;
182 }
183
184};
185
186
187int main(int argc, char* argv[])
188{
189 // Create the example object of the tutorialIPOD type
190 tutorialIPOD example(argc, argv);
191 // Read some parameters from file
193 example._runTime());
194 int NmodesTout = para->ITHACAdict->lookupOrDefault<int>("NmodesTout", 15);
195 int NmodesTproj = para->ITHACAdict->lookupOrDefault<int>("NmodesTproj", 10);
196 double tolleranceSVD =
197 para->ITHACAdict->lookupOrDefault<double>("tolleranceSVD", 1);
198 // Set the number of parameters
199 example.Pnumber = 9;
200 example.Tnumber = NmodesTout;
201 // Set the parameters
202 example.setParameters();
203 // Set the parameter ranges, in all the subdomains the diffusivity varies between
204 // 0.001 and 0.1
205 example.mu_range.col(0) = Eigen::MatrixXd::Ones(9, 1) * 0.001;
206 example.mu_range.col(1) = Eigen::MatrixXd::Ones(9, 1) * 0.1;
207 // Generate the Parameters
208 example.genRandPar(example.Tnumber);
209 // Set the size of the list of values that are multiplying the affine forms
210 example.theta.resize(9);
211 // Set the source term
212 example.SetSource();
213 // Compute the diffusivity field for each subdomain
214 example.compute_nu();
215 // Assemble all the operators of the affine decomposition
216 example.assemble_operator();
217 // Perform an Offline Solve
218 example.offlineSolve();
219 // Perform a POD decomposition and get the modes
220 ITHACAPOD::getModes(example.Tfield, example.Tmodes, example._T().name(),
221 example.podex, 0, 0,
222 NmodesTout);
223 // Set up the incremental POD space
224 scalarIncrementalPOD IPOD(example.Tfield[0], tolleranceSVD, "L2");
225
226 // Fill the incremental POD space
227 for (int fieldI = 1; fieldI < example.Tfield.size(); fieldI++)
228 {
229 IPOD.addSnapshot(example.Tfield[fieldI]);
230 }
231
232 IPOD.writeModes();
233 // Post processing
234 word folder = "./ITHACAoutput/testReconstruction";
235 // Compute new full order solution
236 volScalarField Tfull(example.solveFull(0.05));
237 PtrList<volScalarField> TfullList;
238 TfullList.append(Tfull.clone());
239 PtrList<volScalarField> Tproj;
240 // Project the full order solution onto the POD space
241 example.Tmodes.projectSnapshots(TfullList, Tproj, NmodesTproj);
242 ITHACAstream::exportSolution(Tfull, "1", folder, "Tfull");
243 ITHACAstream::exportSolution(Tproj[0], "1", folder, "Tpod");
244 // Compute the relative error between POD projected field and full order snapshot
245 double EPS = 1e-16;
246 volScalarField relativeErrorField(Tproj[0]);
247
248 for (label i = 0; i < relativeErrorField.internalField().size(); i++)
249 {
250 if (std::abs(Tfull.ref()[i]) < EPS)
251 {
252 relativeErrorField.ref()[i] = (std::abs(Tfull.ref()[i] - Tproj[0].ref()[i])) /
253 EPS;
254 }
255 else
256 {
257 relativeErrorField.ref()[i] = (std::abs(Tfull.ref()[i] - Tproj[0].ref()[i])) /
258 Tfull.ref()[i];
259 }
260 }
261
262 ITHACAstream::exportSolution(relativeErrorField,
263 "1", folder,
264 "relativeErrorField_POD");
265 Info << "Relative error L2 norm POD = " << ITHACAutilities::L2Norm(
266 relativeErrorField) << endl;
267 // Project the full order solution onto the incremental POD space
268 IPOD.projectSnapshots(TfullList, Tproj);
269 ITHACAstream::exportSolution(Tproj[0], "1", folder, "Tipod");
270 volScalarField Tipod = Tproj[0];
271
272 // Compute the relative error between incremental POD projected field and full order snapshot
273 for (label i = 0; i < relativeErrorField.internalField().size(); i++)
274 {
275 if (std::abs(Tfull.ref()[i]) < EPS)
276 {
277 relativeErrorField.ref()[i] = (std::abs(Tfull.ref()[i] - Tipod.ref()[i])) / EPS;
278 }
279 else
280 {
281 relativeErrorField.ref()[i] = (std::abs(Tfull.ref()[i] - Tipod.ref()[i])) /
282 Tfull.ref()[i];
283 }
284 }
285
286 ITHACAstream::exportSolution(relativeErrorField,
287 "1", folder,
288 "relativeErrorField_IPOD");
289 Info << "\n\nRelative error L2 norm incrementalPOD = " <<
290 ITHACAutilities::L2Norm(relativeErrorField) << endl;
291}
292//--------
296
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.
Header file of the reducedLaplacian class.
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.
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:408
Implementation of a incremental POD algorithm according to Oxberry et al.
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.
Class to implement a full order laplacian parametrized problem.
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.
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:93
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