Loading...
Searching...
No Matches
reductionProblem.H
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/>.
24Class
25 reductionProblem
26Description
27 A general reduction problem class is implemented
28SourceFiles
29 reductionProblem.C
30 reductionProblemTemplates.C
31\*---------------------------------------------------------------------------*/
32
37
38
39#ifndef reductionProblem_H
40#define reductionProblem_H
41
42#include <random>
43#include "fvCFD.H"
44#include "IOmanip.H"
45#include "fixedFluxPressureFvPatchScalarField.H"
46#include "freestreamFvPatchField.H"
47#include <sys/stat.h>
48#include "ITHACAutilities.H"
49#include "ITHACAparallel.H"
50#include "ITHACAstream.H"
51#pragma GCC diagnostic push
52#pragma GCC diagnostic ignored "-Wold-style-cast"
53#pragma GCC diagnostic ignored "-Wnon-virtual-dtor"
54#pragma GCC diagnostic ignored "-Wreturn-type"
55#include <datatable.h>
56#include <bspline.h>
57#include <bsplinebuilder.h>
58#include <rbfspline.h>
59#pragma GCC diagnostic pop
60
61// #include <spline.h>
62
63
64// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65
66
67/*---------------------------------------------------------------------------*\
68 Class reductionProblem Declaration
69\*---------------------------------------------------------------------------*/
70
73{
74 private:
75
76 public:
77 // Constructors
78
82
83 // Members
85 label Pnumber;
87 label Tnumber;
89 Eigen::MatrixXd mu;
91 Eigen::MatrixXd mu_range;
93 Eigen::MatrixXd mu_samples;
95 double mu_cur;
97 bool podex;
99 bool offline;
101 IOdictionary* ITHACAdict;
103 autoPtr<argList> _args;
107 label folderN = 1;
108
110 label counter = 1;
111
112#pragma GCC diagnostic push
113#pragma GCC diagnostic ignored "-Wcomment"
137 Eigen::MatrixXi inletIndex;
138
158 Eigen::MatrixXi inletPatch;
159
160#pragma GCC diagnostic pop
161
162 // Functions
163
164 //--------------------------------------------------------------------------
166 void setParameters();
167
168 //--------------------------------------------------------------------------
170 void genRandPar();
171
172 //--------------------------------------------------------------------------
177 void genRandPar(label tsize);
178
179 //--------------------------------------------------------------------------
181 void genEquiPar();
182
183 //--------------------------------------------------------------------------
185 void truthSolve();
186
187 Eigen::MatrixXi inletIndexT;
188
189 //--------------------------------------------------------------------------
196 void assignBC(volVectorField& s, label BC_ind, Vector<double>& value);
197
198 //--------------------------------------------------------------------------
205 void assignBC(volScalarField& s, label BC_ind, double& value);
206
207 //--------------------------------------------------------------------------
216 void reconstructFromMatrix(PtrList<volVectorField>& rec_field2,
217 PtrList<volVectorField>& modes, label Nmodes, Eigen::MatrixXd coeff_matrix);
218
219 //--------------------------------------------------------------------------
229 void reconstructFromMatrix(PtrList<volScalarField>& rec_field2,
230 PtrList<volScalarField>& modes, label Nmodes, Eigen::MatrixXd coeff_matrix);
231
232 //--------------------------------------------------------------------------
241 template<typename T, typename G>
242 void assignIF(T& s, G& value);
243
244 //--------------------------------------------------------------------------
253 template<typename T>
254 void computeLift(T& Lfield, T& liftfield, T& omfield);
255
256 //--------------------------------------------------------------------------
259
260 template<typename T>
261 void computeLiftT(T& Lfield, T& liftfield, T& omfield);
262 //--------------------------------------------------------------------------
265 void liftSolve();
266 void liftSolveT();
267 //--------------------------------------------------------------------------
270 void project();
271
272 //--------------------------------------------------------------------------
277 void writeMu(List<scalar> mu_now);
278
279 //--------------------------------------------------------------------------
288 std::vector<SPLINTER::RBFSpline> getCoeffManifoldRBF(PtrList<volVectorField>
289 snapshots, PtrList<volVectorField>& modes, word rbfBasis = "GAUSSIAN");
290
291 //--------------------------------------------------------------------------
300 std::vector<SPLINTER::RBFSpline> getCoeffManifoldRBF(PtrList<volScalarField>
301 snapshots, PtrList<volScalarField>& modes, word rbfBasis = "GAUSSIAN");
302
303 //--------------------------------------------------------------------------
311 // std::vector<SPLINTER::BSpline*> getCoeffManifoldSPL(PtrList<volVectorField> snapshots, PtrList<volVectorField>& modes, label splDeg=3);
312 std::vector<SPLINTER::BSpline> getCoeffManifoldSPL(PtrList<volVectorField>
313 snapshots, PtrList<volVectorField>& modes, label splDeg = 3);
314
315 //--------------------------------------------------------------------------
323 std::vector<SPLINTER::BSpline> getCoeffManifoldSPL(PtrList<volScalarField>
324 snapshots, PtrList<volScalarField>& modes, label splDeg = 3);
325
326};
327
328#ifdef NoRepository
330#endif
331
332
333// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334
335#endif
336
337// Additional Doxygen Documentation
341
342// Doxigen Documentation (Main Page)
343
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
Class for parallel handling, it has several functions to deal with parallel problems,...
A general class for the implementation of a full order parametrized problem.
Eigen::MatrixXi inletIndexT
Eigen::MatrixXi inletPatch
Matrix that contains informations about the inlet boundaries without specifing the direction Rows = N...
void project()
General projection operation.
void writeMu(List< scalar > mu_now)
Write out a list of scalar corresponding to the parameters used in the truthSolve.
label Pnumber
Number of parameters.
void reconstructFromMatrix(PtrList< volVectorField > &rec_field2, PtrList< volVectorField > &modes, label Nmodes, Eigen::MatrixXd coeff_matrix)
Exact reconstruction using a certain number of modes for vector list of fields and the projection coe...
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
reductionProblem()
Construct Null.
label counter
Counter used for the output of the full order solutions.
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.
ITHACAparallel * paral
parallel handling
void computeLiftT(T &Lfield, T &liftfield, T &omfield)
Virtual function to compute the lifting function.
label folderN
Counter to save intermediate steps in the correct folder, for unsteady and some stationary cases.
void genRandPar()
Generate Random Numbers.
IOdictionary * ITHACAdict
dictionary to store input output infos
double mu_cur
Current value of the parameter.
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.
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
void liftSolve()
Virtual function to compute the lifting function for scalar field.
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
std::vector< SPLINTER::BSpline > getCoeffManifoldSPL(PtrList< volVectorField > snapshots, PtrList< volVectorField > &modes, label splDeg=3)
Constructs the parameters-coefficients manifold for vector fields, based on the B-spline model.
std::vector< SPLINTER::RBFSpline > getCoeffManifoldRBF(PtrList< volVectorField > snapshots, PtrList< volVectorField > &modes, word rbfBasis="GAUSSIAN")
Constructs the parameters-coefficients manifold for vector fields, based on RBF-spline model.
void genEquiPar()
Generate Equidistributed Numbers.
volScalarField & T
Definition createT.H:46
Template file of the reductionProblem class.