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 {}
57 volScalarField& T;
59 volScalarField& nu;
61 volScalarField& S;
62
64 void offlineSolve(word folder = "./ITHACAoutput/Offline/")
65 {
66 if (offline)
67 {
70 ITHACAstream::readMatrix(folder + "/mu_samples_mat.txt");
71 }
72 else
73 {
74 List<scalar> mu_now(9);
75 scalar IF = 0;
76
77 for (label i = 0; i < mu.rows(); i++)
78 {
79 for (label j = 0; j < mu.cols() ; j++)
80 {
81 // mu_now[i] =
82 theta[j] = mu(i, j);
83 }
84
85 assignIF(T, IF);
86 Info << i << endl;
87 truthSolve(mu_now, folder);
88 }
89 }
90 }
91
93 void SetSource()
94 {
95 volScalarField yPos = T.mesh().C().component(vector::Y).ref();
96 volScalarField xPos = T.mesh().C().component(vector::X).ref();
98 {
99 S[counter] = Foam::sin(xPos[counter] / 0.9 * M_PI) + Foam::sin(
100 yPos[counter] / 0.9 * M_PI);
101 }
102 }
103
106 {
107 nu_list.resize(9);
108 volScalarField nu1(nu);
109 volScalarField nu2(nu);
110 volScalarField nu3(nu);
111 volScalarField nu4(nu);
112 volScalarField nu5(nu);
113 volScalarField nu6(nu);
114 volScalarField nu7(nu);
115 volScalarField nu8(nu);
116 volScalarField nu9(nu);
117 Eigen::MatrixXd Box1(2, 3);
118 Box1 << 0, 0, 0, 0.3, 0.3, 0.1;
119 Eigen::MatrixXd Box2(2, 3);
120 Box2 << 0.3, 0, 0, 0.6, 0.3, 0.1;
121 Eigen::MatrixXd Box3(2, 3);
122 Box3 << 0.6, 0, 0, 0.91, 0.3, 0.1;
123 Eigen::MatrixXd Box4(2, 3);
124 Box4 << 0, 0.3, 0, 0.3, 0.6, 0.1;
125 Eigen::MatrixXd Box5(2, 3);
126 Box5 << 0.3, 0.3, 0, 0.6, 0.6, 0.1;
127 Eigen::MatrixXd Box6(2, 3);
128 Box6 << 0.6, 0.3, 0, 0.91, 0.6, 0.1;
129 Eigen::MatrixXd Box7(2, 3);
130 Box7 << 0, 0.6, 0, 0.3, 0.91, 0.1;
131 Eigen::MatrixXd Box8(2, 3);
132 Box8 << 0.3, 0.61, 0, 0.6, 0.91, 0.1;
133 Eigen::MatrixXd Box9(2, 3);
134 Box9 << 0.6, 0.6, 0, 0.9, 0.91, 0.1;
135 ITHACAutilities::setBoxToValue(nu1, Box1, 1.0);
136 ITHACAutilities::setBoxToValue(nu2, Box2, 1.0);
137 ITHACAutilities::setBoxToValue(nu3, Box3, 1.0);
138 ITHACAutilities::setBoxToValue(nu4, Box4, 1.0);
139 ITHACAutilities::setBoxToValue(nu5, Box5, 1.0);
140 ITHACAutilities::setBoxToValue(nu6, Box6, 1.0);
141 ITHACAutilities::setBoxToValue(nu7, Box7, 1.0);
142 ITHACAutilities::setBoxToValue(nu8, Box8, 1.0);
143 ITHACAutilities::setBoxToValue(nu9, Box9, 1.0);
144 nu_list.set(0, (nu1).clone());
145 nu_list.set(1, (nu2).clone());
146 nu_list.set(2, (nu3).clone());
147 nu_list.set(3, (nu4).clone());
148 nu_list.set(4, (nu5).clone());
149 nu_list.set(5, (nu6).clone());
150 nu_list.set(6, (nu7).clone());
151 nu_list.set(7, (nu8).clone());
152 nu_list.set(8, (nu9).clone());
153 }
154
157 {
158 for (int i = 0; i < nu_list.size(); i++)
159 {
160 operator_list.append(fvm::laplacian(nu_list[i], T));
161 }
162 }
163
164};
165
166void offline_stage(tutorial02& example, tutorial02& FOM_test);
167void online_stage(tutorial02& example, tutorial02& FOM_test);
168
169int main(int argc, char* argv[])
170{
171#if OPENFOAM > 1812
172 // load stage to perform
173 argList::addOption("stage", "offline", "Perform offline stage");
174 argList::addOption("stage", "online", "Perform online stage");
175 // add options for parallel run
176 HashTable<string> validParOptions;
177 validParOptions.set
178 (
179 "stage",
180 "offline"
181 );
182 validParOptions.set
183 (
184 "stage",
185 "online"
186 );
187 Pstream::addValidParOptions(validParOptions);
188 // Construct the tutorial object
189 tutorial02 example(argc, argv);
191 tutorial02 FOM_test(argc, argv);
192
193 if (example._args().get("stage").match("offline"))
194 {
195 // perform the offline stage, extracting the modes from the snapshots'
196 // dataset corresponding to parOffline
197 offline_stage(example, FOM_test);
198 }
199 else if (example._args().get("stage").match("online"))
200 {
201 // load precomputed modes and reduced matrices
202 offline_stage(example, FOM_test);
203 // perform online solve with respect to the parameters in parOnline
204 online_stage(example, FOM_test);
205 }
206 else
207 {
208 Info << "Pass '-stage offline', '-stage online'" << endl;
209 }
210
211#else
212
213 if (argc == 1)
214 {
215 Info << "Pass 'offline' or 'online' as first arguments."
216 << endl;
217 exit(0);
218 }
219
220 // process arguments removing "offline" or "online" keywords
221 int argc_proc = argc - 1;
222 char* argv_proc[argc_proc];
223 argv_proc[0] = argv[0];
224
225 if (argc > 2)
226 {
227 std::copy(argv + 2, argv + argc, argv_proc + 1);
228 }
229
230 argc--;
231 // Construct the tutorial object
232 tutorial02 example(argc, argv);
234 tutorial02 FOM_test(argc, argv);
235
236 if (std::strcmp(argv[1], "offline") == 0)
237 {
238 // perform the offline stage, extracting the modes from the snapshots'
239 // dataset corresponding to parOffline
240 offline_stage(example, FOM_test);
241 }
242 else if (std::strcmp(argv[1], "online") == 0)
243 {
244 // load precomputed modes and reduced matrices
245 offline_stage(example, FOM_test);
246 // perform online solve with respect to the parameters in parOnline
247 online_stage(example, FOM_test);
248 }
249 else
250 {
251 Info << "Pass offline, online" << endl;
252 }
253
254#endif
255 exit(0);
256}
257
258void offline_stage(tutorial02& example, tutorial02& FOM_test)
259{
260 // Read some parameters from file
262 example._runTime());
263 int NmodesTout = para->ITHACAdict->lookupOrDefault<int>("NmodesTout", 15);
264 int NmodesTproj = para->ITHACAdict->lookupOrDefault<int>("NmodesTproj", 10);
265 // Set the number of parameters
266 example.Pnumber = 9;
267 example.Tnumber = 500;
268 // Set the parameters
269 example.setParameters();
270 // Set the parameter ranges, in all the subdomains the diffusivity varies between
271 // 0.001 and 0.1
272 example.mu_range.col(0) = Eigen::MatrixXd::Ones(9, 1) * 0.001;
273 example.mu_range.col(1) = Eigen::MatrixXd::Ones(9, 1) * 0.1;
274 // Generate the Parameters
275 example.genRandPar(500);
276 // Set the size of the list of values that are multiplying the affine forms
277 example.theta.resize(9);
278 // Set the source term
279 example.SetSource();
280 // Compute the diffusivity field for each subdomain
281 example.compute_nu();
282 // Assemble all the operators of the affine decomposition
283 example.assemble_operator();
284 // Perform an Offline Solve
285 example.offlineSolve();
286 // Perform a POD decomposition and get the modes
287 ITHACAPOD::getModes(example.Tfield, example.Tmodes, example._T().name(),
288 example.podex, 0, 0,
289 NmodesTout);
290 // Perform the Galerkin projection onto the space spanned by the POD modes
291 example.project(NmodesTproj);
292 FOM_test.offline = false;
293 FOM_test.Pnumber = 9;
294 FOM_test.Tnumber = 50;
295 // Set the parameters
296 FOM_test.setParameters();
297 // Set the parameter ranges, in all the subdomains the diffusivity varies between
298 // 0.001 and 0.1
299 FOM_test.mu_range.col(0) = Eigen::MatrixXd::Ones(9, 1) * 0.001;
300 FOM_test.mu_range.col(1) = Eigen::MatrixXd::Ones(9, 1) * 0.1;
301 // Generate the Parameters
302 FOM_test.genRandPar(50);
303 // Set the size of the list of values that are multiplying the affine forms
304 FOM_test.theta.resize(9);
305 // Set the source term
306 FOM_test.SetSource();
307 // Compute the diffusivity field for each subdomain
308 FOM_test.compute_nu();
309 // Assemble all the operators of the affine decomposition
310 FOM_test.assemble_operator();
311 // Perform an Offline Solve
312 FOM_test.offlineSolve("./ITHACAoutput/FOMtest");
313}
314
315void online_stage(tutorial02& example, tutorial02& FOM_test)
316{
317 // Create a reduced object
318 reducedLaplacian reduced(example);
319
320 // Solve the online reduced problem some new values of the parameters
321 for (int i = 0; i < FOM_test.mu.rows(); i++)
322 {
323 reduced.solveOnline(FOM_test.mu.row(i));
324 }
325
326 // Reconstruct the solution and store it into Reconstruction folder
327 reduced.reconstruct("./ITHACAoutput/Reconstruction");
328 // Compute the error on the testing set
329 Eigen::MatrixXd error = ITHACAutilities::errorL2Rel(FOM_test.Tfield,
330 reduced.Trec);
331}
332//--------
336
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.
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.
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.
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:93
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