Loading...
Searching...
No Matches
02thermalBlock.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 02thermalBlock.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 <Eigen/Dense>
39#define _USE_MATH_DEFINES
40#include <cmath>
41
46{
47 public:
48 explicit tutorial02(int argc, char* argv[])
49 :
50 laplacianProblem(argc, argv),
51 T(_T()),
52 nu(_nu()),
53 S(_S())
54 {}
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[i] =
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
165};
166
167void offline_stage(tutorial02& example, tutorial02& FOM_test);
168void online_stage(tutorial02& example, tutorial02& FOM_test);
169
170int main(int argc, char* argv[])
171{
172#if OPENFOAM > 1812
173 // load stage to perform
174 argList::addOption("stage", "offline", "Perform offline stage");
175 argList::addOption("stage", "online", "Perform online stage");
176 // add options for parallel run
177 HashTable<string> validParOptions;
178 validParOptions.set
179 (
180 "stage",
181 "offline"
182 );
183 validParOptions.set
184 (
185 "stage",
186 "online"
187 );
188 Pstream::addValidParOptions(validParOptions);
189 // Construct the tutorial object
190 tutorial02 example(argc, argv);
192 tutorial02 FOM_test(argc, argv);
193
194 if (example._args().get("stage").match("offline"))
195 {
196 // perform the offline stage, extracting the modes from the snapshots'
197 // dataset corresponding to parOffline
198 offline_stage(example, FOM_test);
199 }
200 else if (example._args().get("stage").match("online"))
201 {
202 // load precomputed modes and reduced matrices
203 offline_stage(example, FOM_test);
204 // perform online solve with respect to the parameters in parOnline
205 online_stage(example, FOM_test);
206 }
207 else
208 {
209 Info << "Pass '-stage offline', '-stage online'" << endl;
210 }
211
212#else
213
214 if (argc == 1)
215 {
216 Info << "Pass 'offline' or 'online' as first arguments."
217 << endl;
218 exit(0);
219 }
220
221 // process arguments removing "offline" or "online" keywords
222 int argc_proc = argc - 1;
223 char* argv_proc[argc_proc];
224 argv_proc[0] = argv[0];
225
226 if (argc > 2)
227 {
228 std::copy(argv + 2, argv + argc, argv_proc + 1);
229 }
230
231 argc--;
232 // Construct the tutorial object
233 tutorial02 example(argc, argv);
235 tutorial02 FOM_test(argc, argv);
236
237 if (std::strcmp(argv[1], "offline") == 0)
238 {
239 // perform the offline stage, extracting the modes from the snapshots'
240 // dataset corresponding to parOffline
241 offline_stage(example, FOM_test);
242 }
243 else if (std::strcmp(argv[1], "online") == 0)
244 {
245 // load precomputed modes and reduced matrices
246 offline_stage(example, FOM_test);
247 // perform online solve with respect to the parameters in parOnline
248 online_stage(example, FOM_test);
249 }
250 else
251 {
252 Info << "Pass offline, online" << endl;
253 }
254
255#endif
256 exit(0);
257}
258
259void offline_stage(tutorial02& example, tutorial02& FOM_test)
260{
261 // Read some parameters from file
263 example._runTime());
264 int NmodesTout = para->ITHACAdict->lookupOrDefault<int>("NmodesTout", 15);
265 int NmodesTproj = para->ITHACAdict->lookupOrDefault<int>("NmodesTproj", 10);
266 // Set the number of parameters
267 example.Pnumber = 9;
268 example.Tnumber = 500;
269 // Set the parameters
270 example.setParameters();
271 // Set the parameter ranges, in all the subdomains the diffusivity varies between
272 // 0.001 and 0.1
273 example.mu_range.col(0) = Eigen::MatrixXd::Ones(9, 1) * 0.001;
274 example.mu_range.col(1) = Eigen::MatrixXd::Ones(9, 1) * 0.1;
275 // Generate the Parameters
276 example.genRandPar(500);
277 // Set the size of the list of values that are multiplying the affine forms
278 example.theta.resize(9);
279 // Set the source term
280 example.SetSource();
281 // Compute the diffusivity field for each subdomain
282 example.compute_nu();
283 // Assemble all the operators of the affine decomposition
284 example.assemble_operator();
285 // Perform an Offline Solve
286 example.offlineSolve();
287 // Perform a POD decomposition and get the modes
288 ITHACAPOD::getModes(example.Tfield, example.Tmodes, example._T().name(),
289 example.podex, 0, 0,
290 NmodesTout);
291 // Perform the Galerkin projection onto the space spanned by the POD modes
292 example.project(NmodesTproj);
293 FOM_test.offline = false;
294 FOM_test.Pnumber = 9;
295 FOM_test.Tnumber = 50;
296 // Set the parameters
297 FOM_test.setParameters();
298 // Set the parameter ranges, in all the subdomains the diffusivity varies between
299 // 0.001 and 0.1
300 FOM_test.mu_range.col(0) = Eigen::MatrixXd::Ones(9, 1) * 0.001;
301 FOM_test.mu_range.col(1) = Eigen::MatrixXd::Ones(9, 1) * 0.1;
302 // Generate the Parameters
303 FOM_test.genRandPar(50);
304 // Set the size of the list of values that are multiplying the affine forms
305 FOM_test.theta.resize(9);
306 // Set the source term
307 FOM_test.SetSource();
308 // Compute the diffusivity field for each subdomain
309 FOM_test.compute_nu();
310 // Assemble all the operators of the affine decomposition
311 FOM_test.assemble_operator();
312 // Perform an Offline Solve
313 FOM_test.offlineSolve("./ITHACAoutput/FOMtest");
314}
315
316void online_stage(tutorial02& example, tutorial02& FOM_test)
317{
318 // Create a reduced object
319 reducedLaplacian reduced(example);
320
321 // Solve the online reduced problem some new values of the parameters
322 for (int i = 0; i < FOM_test.mu.rows(); i++)
323 {
324 reduced.solveOnline(FOM_test.mu.row(i));
325 }
326
327 // Reconstruct the solution and store it into Reconstruction folder
328 reduced.reconstruct("./ITHACAoutput/Reconstruction");
329 // Compute the error on the testing set
330 Eigen::MatrixXd error = ITHACAutilities::errorL2Rel(FOM_test.Tfield,
331 reduced.Trec);
332}
333
334//--------
338
int main(int argc, char *argv[])
void online_stage(tutorial02 &example, tutorial02 &FOM_test)
void offline_stage(tutorial02 &example, tutorial02 &FOM_test)
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.
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.
void project(label Nmodes)
Perform a projection onto the POD modes.
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.
Class to solve the online reduced problem associated with a Laplace problem.
void reconstruct(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Function to recover the solution given the online solution.
void solveOnline(Eigen::MatrixXd mu)
Function to perform an online solve given a certain mu.
PtrList< volScalarField > Trec
Reduced reconstructed Solution.
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.
autoPtr< argList > _args
argList
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 2 is implemented.
volScalarField & S
Source term field.
void offlineSolve(word folder="./ITHACAoutput/Offline/")
It perform an offline Solve.
volScalarField & nu
Diffusivity field.
void compute_nu()
Compute the diffusivity in each subdomain.
void assemble_operator()
Construct the operator_list where each term of the affine decomposition is stored.
tutorial02(int argc, char *argv[])
volScalarField & T
[tutorial02] Temperature field
void SetSource()
Define the source term function.
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
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 errorL2Rel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in L2 norm.
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