Loading...
Searching...
No Matches
11UnsteadyBBOpen.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2Copyright (C) 2017 by the ITHACA-FV authors
3
4License
5 This file is part of ITHACA-FV
6
7 ITHACA-FV is free software: you can redistribute it and/or modify
8 it under the terms of the GNU Lesser General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 ITHACA-FV is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU Lesser General Public License for more details.
16
17 You should have received a copy of the GNU Lesser General Public License
18 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
19
20Description
21 Example of Boussinesq approximation for two way coupling NS-momentum equation
22 and heat transport equation for open flows where pressure connot be neglected.
23SourceFiles
24 11UnsteadyBBOpen.C
25\*---------------------------------------------------------------------------*/
26
27#include "UnsteadyBB.H"
28#include "ITHACAPOD.H"
29#include "ReducedUnsteadyBB.H"
30#include "ITHACAstream.H"
31#include <chrono>
32#include <math.h>
33#include <iomanip>
34
36{
37 public:
38 explicit tutorial11(int argc, char* argv[])
39 :
40 UnsteadyBB(argc, argv),
41 U(_U()),
42 p(_p()),
43 p_rgh(_p_rgh()),
44 T(_T())
45 {}
46
47 // Fields To Perform
48 volVectorField& U;
49 volScalarField& p;
50 volScalarField& p_rgh;
51 volScalarField& T;
52
53 void offlineSolve(Eigen::MatrixXd par_BC)
54 {
55 List<scalar> mu_now(1);
56
57 if (offline)
58 {
59 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
60 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
61 ITHACAstream::read_fields(Prghfield, p_rgh, "./ITHACAoutput/Offline/");
62 ITHACAstream::read_fields(Tfield, T, "./ITHACAoutput/Offline/");
63 }
64 else
65 {
66 for (label k = 0; k < par_BC.rows(); k++)
67 {
68 for (label j = 0; j < par_BC.cols(); j++)
69 {
70 for (label i = 0; i < mu.cols(); i++)
71 {
72 mu_now[0] = mu(0, i);
73 }
74
75 assignBC(T, inletIndexT(j, 0), par_BC(k, j));
76 }
77
78 truthSolve(mu_now);
79 }
80 }
81 }
82
83
84 void onlineSolveFull(Eigen::MatrixXd par_BC, label para_set_BC, fileName folder)
85 {
87 {
88 }
89 else
90 {
92 label i = para_set_BC;
93
94 for (label j = 0; j < par_BC.cols(); j++)
95 {
96 assignBC(T, inletIndexT(j, 0), par_BC(i, j));
97 }
98
99 truthSolve(folder);
100 }
101 }
102
103 void onlineSolveRead(fileName folder)
104 {
106 {
109 }
110 else
111 {
112 }
113 }
114
115
116 // Method to compute the lifting function for temperature
118 {
119 for (label k = 0; k < inletIndexT.rows(); k++)
120 {
121 Time& runTime = _runTime();
122 fvMesh& mesh = _mesh();
123 volScalarField& T = _T();
124 volVectorField& U = _U();
125 surfaceScalarField& phi = _phi();
126 phi = linearInterpolate(U) & mesh.Sf();
127 simpleControl simple(mesh);
128 // IOMRFZoneList& MRF = _MRF();
129 // singlePhaseTransportModel& laminarTransport = _laminarTransport();
130 // volScalarField& nut = _nut();
131 volScalarField& alphat = _alphat();
132 // dimensionedScalar& nu = _nu();
133 dimensionedScalar& Pr = _Pr();
134 dimensionedScalar& Prt = _Prt();
135 label BCind = inletIndexT(k, 0);
136 volScalarField Tlift("Tlift" + name(k), T);
137 instantList Times = runTime.times();
138 runTime.setTime(Times[1], 1);
139 Info << "Solving a lifting Problem" << endl;
140 scalar t1 = 1;
141 scalar t0 = 0;
142 alphat = turbulence->nut() / Prt;
143 alphat.correctBoundaryConditions();
144 volScalarField alphaEff("alphaEff", turbulence->nu() / Pr + alphat);
145
146 for (label j = 0; j < T.boundaryField().size(); j++)
147 {
148 if (j == BCind)
149 {
150 assignBC(Tlift, j, t1);
151 assignIF(Tlift, t0);
152 }
153 else if (T.boundaryField()[BCind].type() == "fixedValue")
154 {
155 assignBC(Tlift, j, t0);
156 assignIF(Tlift, t0);
157 }
158 else
159 {
160 }
161 }
162
163 while (simple.correctNonOrthogonal())
164 {
165 fvScalarMatrix TEqn
166 (
167 fvm::div(phi, Tlift)
168 - fvm::laplacian(alphaEff, Tlift)
169 );
170 TEqn.solve();
171 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
172 << " ClockTime = " << runTime.elapsedClockTime() << " s"
173 << nl << endl;
174 }
175
176 Tlift.write();
177 liftfieldT.append((Tlift).clone());
178 }
179 }
180
181};
182
183
184/*---------------------------------------------------------------------------*\
185 Starting the MAIN
186\*---------------------------------------------------------------------------*/
187int main(int argc, char* argv[])
188{
189 // Construct the tutorial object
190 tutorial11 example(argc, argv);
191 // the offline samples for the boundary conditions
192 word par_offline_BC("./par_offline_BC");
193 // the samples which will be used for setting the boundary condition in the online stage
194 word par_online_BC("./par_online_BC");
195 Eigen::MatrixXd par_off_BC = ITHACAstream::readMatrix(par_offline_BC);
196 Eigen::MatrixXd par_on_BC = ITHACAstream::readMatrix(par_online_BC);
197 // Read some parameters from file
199 example._runTime());
200 word stabilization = para->ITHACAdict->lookupOrDefault<word>("Stabilization",
201 "supremizer");
202 int NmodesUproj = para->ITHACAdict->lookupOrDefault<int>("NmodesUproj", 5);
203 int NmodesPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesPproj", 5);
204 int NmodesPrghproj = para->ITHACAdict->lookupOrDefault<int>("NmodesPrghproj",
205 5);
206 int NmodesTproj = para->ITHACAdict->lookupOrDefault<int>("NmodesTproj", 5);
207 int NmodesSUPproj = 0;
208
209 if (stabilization == "supremizer")
210 {
211 NmodesSUPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesSUPproj", 5);
212 }
213
214 int NmodesOut = para->ITHACAdict->lookupOrDefault<int>("NmodesOut", 15);
215 // Set the number of parameters
216 example.Pnumber = 1;
217 // Set samples
218 example.Tnumber = 1;
219 // Set the parameters infos
220 example.setParameters();
221 // Set the parameter ranges
222 example.mu_range(0, 0) = 0.00001;
223 example.mu_range(0, 1) = 0.00001;
224 // Generate equispaced samples inside the parameter range
225 example.genEquiPar();
226 // Set the inlet Temperature boundaries where we have non homogeneous boundary conditions
227 example.inletIndexT.resize(3, 1);
228 example.inletIndexT << 1, 2, 3;
229 // Set the inlet Velocity boundaries where we have non homogeneous boundary conditions
230 example.inletIndex.resize(1, 2);
231 example.inletIndex << 3, 0;
232 // Time parameters
233 example.startTime = 0.0;
234 example.finalTime = 5.0;
235 example.timeStep = 0.002;
236 example.writeEvery = 0.01;
237 // Perform the Offline Solve;
238 example.offlineSolve(par_off_BC);
239 // Search the lift function for the velocity
240 example.liftSolve();
241 // Search the lift function for the temperature
242 example.liftSolveT();
243 // Create homogeneous basis functions for velocity
244 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
245 // Create homogeneous basis functions for temperature
246 example.computeLiftT(example.Tfield, example.liftfieldT, example.Tomfield);
247 // Perform a POD decomposition for velocity ,temperature and pressure fields
249 example.Uomfield, example.Umodes, example._U().name(),
250 example.podex, 0, 0, NmodesOut, false);
252 example.Pfield, example.Pmodes, example._p().name(),
253 example.podex, 0, 0, NmodesOut, false);
255 example.Prghfield, example.Prghmodes, example._p_rgh().name(),
256 example.podex, 0, 0, NmodesOut, false);
258 example.Tomfield, example.Tmodes, example._T().name(),
259 example.podex, 0, 0, NmodesOut, false);
260
261 // Solve the supremizer problem
262 if (stabilization == "supremizer")
263 {
264 example.solvesupremizer("modes");
265 }
266
267 // Create a list with number of modes for which the projection needs to be performed
268 Eigen::MatrixXd List_of_modes(NmodesOut - 5, 1);
269
270 for (int i = 0; i < List_of_modes.rows(); i++)
271 {
272 List_of_modes(i, 0) = i + 1;
273 }
274
275 // Export with number of modes for which the projection needs to be performed
276 ITHACAstream::exportMatrix(List_of_modes, "List_of_modes", "eigen",
277 "./ITHACAoutput/l2error");
278 // Create locally the temperature modes
279 PtrList<volScalarField> TLmodes;
280
281 for (label k = 0; k < example.liftfieldT.size(); k++)
282 {
283 TLmodes.append((example.liftfieldT[k]).clone());
284 }
285
286 for (label k = 0; k < List_of_modes.size(); k++)
287 {
288 TLmodes.append((example.Tmodes[k]).clone());
289 }
290
291 // Create locally the velocity modes
292 PtrList<volVectorField> ULmodes;
293
294 for (label k = 0; k < example.liftfield.size(); k++)
295 {
296 ULmodes.append((example.liftfield[k]).clone());
297 }
298
299 for (label k = 0; k < List_of_modes.size(); k++)
300 {
301 ULmodes.append((example.Umodes[k]).clone());
302 }
303
304 if (stabilization == "supremizer")
305 {
306 for (label k = 0; k < NmodesSUPproj; k++)
307 {
308 ULmodes.append((example.supmodes[k]).clone());
309 }
310 }
311
312 // Perform the projection for all number of modes in List_of_modes for temperature and velocity
313 Eigen::MatrixXd L2errorProjMatrixU(example.Ufield.size(), List_of_modes.rows());
314 Eigen::MatrixXd L2errorProjMatrixT(example.Tfield.size(), List_of_modes.rows());
315
316 // Calculate the coefficients and L2 error and store the error in a matrix for each number of modes
317 for (int i = 0; i < List_of_modes.rows(); i++)
318 {
319 Eigen::MatrixXd coeffU = ITHACAutilities::getCoeffs(example.Ufield, ULmodes,
320 List_of_modes(i, 0) + example.liftfield.size() + NmodesSUPproj);
321 Eigen::MatrixXd coeffT = ITHACAutilities::getCoeffs(example.Tfield, TLmodes,
322 List_of_modes(i, 0) + example.liftfieldT.size());
323 PtrList<volVectorField> rec_fieldU = ITHACAutilities::reconstructFromCoeff(
324 ULmodes, coeffU, List_of_modes(i,
325 0) + example.liftfield.size() + NmodesSUPproj);
326 PtrList<volScalarField> rec_fieldT = ITHACAutilities::reconstructFromCoeff(
327 TLmodes, coeffT, List_of_modes(i, 0) + example.liftfieldT.size());
328 Eigen::MatrixXd L2errorProjU = ITHACAutilities::errorL2Rel(example.Ufield,
329 rec_fieldU);
330 Eigen::MatrixXd L2errorProjT = ITHACAutilities::errorL2Rel(example.Tfield,
331 rec_fieldT);
332 L2errorProjMatrixU.col(i) = L2errorProjU;
333 L2errorProjMatrixT.col(i) = L2errorProjT;
334 }
335
336 // Export the matrix containing the error
337 ITHACAstream::exportMatrix(L2errorProjMatrixU, "L2errorProjMatrixU", "eigen",
338 "./ITHACAoutput/l2error");
339 ITHACAstream::exportMatrix(L2errorProjMatrixT, "L2errorProjMatrixT", "eigen",
340 "./ITHACAoutput/l2error");
341
342 // Get reduced matrices
343 if (stabilization == "supremizer")
344 {
345 example.projectSUP("./Matrices", NmodesUproj, NmodesPrghproj, NmodesTproj,
346 NmodesSUPproj);
347 }
348 else if (stabilization == "PPE")
349 {
350 example.projectPPE("./Matrices", NmodesUproj, NmodesPrghproj, NmodesTproj,
351 NmodesSUPproj);
352 }
353 else
354 {
355 // TODO: warning message
356 }
357
358 // Resize the modes for projection
359 example.Tmodes.resize(NmodesTproj);
360 example.Umodes.resize(NmodesUproj);
361 example.Pmodes.resize(NmodesPproj);
362 example.Prghmodes.resize(NmodesPrghproj);
363 // Online part
364 ReducedUnsteadyBB reduced(example);
365 // Set values of the online solve
366 reduced.nu = 0.00001;
367 reduced.Pr = 0.71;
368 reduced.tstart = 0.0;
369 reduced.finalTime = 5.0;
370 reduced.dt = 0.002;
371 // Set the online velocity
372 Eigen::MatrixXd vel_now_BC(1, 1);
373 vel_now_BC(0, 0) = 0.0157;
374
375 // Set the online temperature BC and solve reduced model
376 for (label k = 0; k < par_on_BC.rows(); k++)
377 {
378 Eigen::MatrixXd temp_now_BC(3, 1);
379 temp_now_BC(0, 0) = par_on_BC(k, 0);
380 temp_now_BC(1, 0) = par_on_BC(k, 1);
381 temp_now_BC(2, 0) = par_on_BC(k, 2);
382
383
384 if (stabilization == "supremizer")
385 {
386 reduced.solveOnline_sup(temp_now_BC, vel_now_BC, k, par_on_BC.rows());
387 reduced.reconstruct_sup("./ITHACAoutput/ReconstructionSUP/", 5);
388 }
389 else if (stabilization == "PPE")
390 {
391 reduced.solveOnline_PPE(temp_now_BC, vel_now_BC, k, par_on_BC.rows());
392 reduced.reconstruct_PPE("./ITHACAoutput/ReconstructionPPE/", 5);
393 }
394 }
395
396 // Performing full order simulation for second parameter set of temp_now_BC
397 tutorial11 HFonline2(argc, argv);
398 HFonline2.Pnumber = 1;
399 HFonline2.Tnumber = 1;
400 HFonline2.setParameters();
401 HFonline2.mu_range(0, 0) = 0.00001;
402 HFonline2.mu_range(0, 1) = 0.00001;
403 HFonline2.genEquiPar();
404 HFonline2.inletIndexT.resize(3, 1);
405 HFonline2.inletIndexT << 1, 2, 3;
406 HFonline2.inletIndex.resize(1, 2);
407 HFonline2.inletIndex << 3, 0;
408 HFonline2.startTime = 0.0;
409 HFonline2.finalTime = 5.0;
410 HFonline2.timeStep = 0.002;
411 HFonline2.writeEvery = 0.01;
412 // Reconstruct the online solution
413 HFonline2.onlineSolveFull(par_on_BC, 1,
414 "./ITHACAoutput/high_fidelity_online2");
415 // Reading in the high-fidelity solutions for the parameter set
416 // for which the offline solve has been performed
417 example.onlineSolveRead("./ITHACAoutput/Offline/");
418 // Reading in the high-fidelity solutions for the second parameter set
419 example.onlineSolveRead("./ITHACAoutput/high_fidelity_online2/");
420 // Calculate error between online- and corresponding full order solution
421 Eigen::MatrixXd L2errorMatrixU = ITHACAutilities::errorL2Rel(
422 example.Ufield_on, reduced.UREC);
423 Eigen::MatrixXd L2errorMatrixT = ITHACAutilities::errorL2Rel(
424 example.Tfield_on, reduced.TREC);
425 //Export the matrix containing the error
426 ITHACAstream::exportMatrix(L2errorMatrixU, "L2errorMatrixU", "eigen",
427 "./ITHACAoutput/l2error");
428 ITHACAstream::exportMatrix(L2errorMatrixT, "L2errorMatrixT", "eigen",
429 "./ITHACAoutput/l2error");
430 exit(0);
431}
432
436
int main(int argc, char *argv[])
Header file of the ITHACAPOD class.
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ReducedUnsteadyBB class.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Header file of the UnsteadyBB 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 where it is implemented a reduced problem for the unsteady Navier-stokes problem.
void reconstruct_sup(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique.
scalar Pr
DimensionedScalar Pr;.
Eigen::MatrixXd solveOnline_sup(Eigen::MatrixXd &temp_now_BC, Eigen::MatrixXd &vel_now_BC, int NParaSet=0, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
PtrList< volScalarField > TREC
Reconstructed temperature field.
Eigen::MatrixXd solveOnline_PPE(Eigen::MatrixXd &temp_now_BC, Eigen::MatrixXd &vel_now_BC, int NParaSet=0, int startSnap=0)
Method to perform an online solve using a PPE stabilisation method.
Implementation of a parametrized full order unsteady Boussinesq problem and preparation of the the ...
Definition UnsteadyBB.H:60
autoPtr< volScalarField > _p_rgh
Shifted Pressure field.
Definition UnsteadyBB.H:131
void solvesupremizer(word type="snapshots")
solve the supremizer either with the use of the pressure snaphots or the pressure modes
Definition UnsteadyBB.C:281
PtrList< volScalarField > Prghmodes
List of pointers used to form the shifted pressure modes.
Definition UnsteadyBB.H:79
autoPtr< volScalarField > _T
Temperature field.
Definition UnsteadyBB.H:134
autoPtr< dimensionedScalar > _Prt
dimensionedScalar Prt;
Definition UnsteadyBB.H:146
autoPtr< dimensionedScalar > _Pr
dimensionedScalar Pr;
Definition UnsteadyBB.H:143
PtrList< volScalarField > Tomfield
List of pointers used to form the homogeneous velocity snapshots.
Definition UnsteadyBB.H:97
PtrList< volScalarField > Prghfield
List of pointers used to form the shifted pressure snapshots matrix.
Definition UnsteadyBB.H:73
PtrList< volVectorField > Ufield_on
List of pointers used to form the temperature snapshots matrix.
Definition UnsteadyBB.H:91
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NTmodes, label NSUPmodes)
Project using a supremizer approach.
Definition UnsteadyBB.C:488
void liftSolve()
Perform a lift solve for velocity field.
autoPtr< fvMesh > _mesh
Mesh.
Definition UnsteadyBB.H:128
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
Definition UnsteadyBB.H:76
autoPtr< volScalarField > _alphat
dimensionedScalar alphat;
Definition UnsteadyBB.H:149
PtrList< volScalarField > Tfield_on
List of pointers used to form the temperature snapshots matrix.
Definition UnsteadyBB.H:85
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of lifting functions.
Definition UnsteadyBB.H:94
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
Definition UnsteadyBB.H:82
void projectPPE(fileName folder, label NUmodes, label NPrghmodes, label NTmodes, label NSUPmodes)
Project using a PPE approach.
Definition UnsteadyBB.C:652
scalar startTime
Start Time (initial time to start storing the snapshots)
scalar writeEvery
Time step of the writing procedure.
scalar timeStep
Time step of the simulation.
scalar finalTime
Final time (final time of the simulation and consequently of the acquisition of the snapshots)
scalar nu
Reduced viscosity in case of parametrized viscosity.
PtrList< volVectorField > UREC
Recontructed velocity field.
void reconstruct_PPE(fileName folder="./ITHACAoutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a PPE stabilisation technique.
scalar finalTime
Scalar to store the final time if the online simulation.
scalar tstart
Scalar to store the initial time if the online simulation.
double dt
Scalar to store the time increment.
Eigen::MatrixXi inletIndexT
label Pnumber
Number of parameters.
void assignBC(volVectorField &s, label BC_ind, Vector< double > &value)
Assign Boundary Condition to a volVectorField.
void assignIF(T &s, G &value)
Assign internal field condition.
label Tnumber
Dimension of the training set (used only when gerating parameters without input)
void computeLift(T &Lfield, T &liftfield, T &omfield)
Homogenize the snapshot matrix, it works with PtrList of volVectorField and volScalarField.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
void computeLiftT(T &Lfield, T &liftfield, T &omfield)
Virtual function to compute the lifting function.
Eigen::MatrixXd mu
Row matrix of parameters.
Eigen::MatrixXd mu_range
Range of the parameter spaces.
void setParameters()
Set Parameters Problems.
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
void genEquiPar()
Generate Equidistributed Numbers.
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
volVectorModes supmodes
List of pointers used to form the supremizer modes.
Definition steadyNS.H:104
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:290
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition steadyNS.H:101
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition steadyNS.H:89
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
PtrList< volVectorField > Uomfield
List of pointers used to form the homogeneous velocity snapshots.
Definition steadyNS.H:113
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
volScalarField & T
void onlineSolveFull(Eigen::MatrixXd par_BC, label para_set_BC, fileName folder)
volScalarField & p_rgh
void onlineSolveRead(fileName folder)
volVectorField & U
void offlineSolve(Eigen::MatrixXd par_BC)
tutorial11(int argc, char *argv[])
volScalarField & p
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition unsteadyNS.H:76
fvScalarMatrix & TEqn
Definition TEqn.H:15
volScalarField & alphat
dimensionedScalar & Prt
dimensionedScalar & Pr
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 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.
Eigen::VectorXd getCoeffs(GeometricField< Type, PatchField, GeoMesh > &snapshot, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Projects a snapshot on a basis function and gets the coefficients of the projection.
PtrList< GeometricField< Type, PatchField, GeoMesh > > reconstructFromCoeff(PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, Eigen::MatrixXd &coeff_matrix, label Nmodes)
Exact reconstruction using a certain number of modes for a list of fields and the projection coeffici...
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.
bool check_folder(word folder)
Checks if a folder exists.
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
simpleControl simple(mesh)
surfaceScalarField & phi
label i
Definition pEqn.H:46