Loading...
Searching...
No Matches
06POD_RBF.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
15 ITHACA-FV is free software: you can redistribute it and/or modify
16 it under the terms of the GNU Lesser General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 ITHACA-FV is distributed in the hope that it will be useful,
21 but WITHOUT ANY WARRANTY; without even the implied warranty of
22 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23 GNU Lesser General Public License for more details.
24
25 You should have received a copy of the GNU Lesser General Public License
26 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
27
28Description
29 Example of NS-Stokes Reduction Problem for Turbulent Flow Case
30SourceFiles
31 06POD_RBF.C
32\*---------------------------------------------------------------------------*/
33#include "fvCFD.H"
34#include "fvOptions.H"
35#include "IOmanip.H"
36#include "simpleControl.H"
37#include "ITHACAutilities.H"
38#include "Foam2Eigen.H"
39#include "ITHACAstream.H"
40#include "ITHACAPOD.H"
41#include <Eigen/Dense>
42#include <Eigen/SVD>
43#include <Eigen/SparseLU>
44#include "EigenFunctions.H"
45#include <chrono>
46#include "reductionProblem.H"
47#include "steadyNS.H"
48#include "SteadyNSTurb.H"
49#include "ReducedSteadyNS.H"
50#include "ReducedSteadyNSTurb.H"
51
56{
57 public:
58 explicit tutorial06(int argc, char* argv[])
59 :
60 SteadyNSTurb(argc, argv),
61 U(_U()),
62 p(_p()),
63 nut(_nut())
64 {}
65
67 // Relevant Fields
68 volVectorField& U;
69 volScalarField& p;
70 volScalarField& nut;
73 {
74 Vector<double> inl(0, 0, 0);
75 List<scalar> mu_now(2);
76
77 // if the offline solution is already performed read the fields
78 if (offline)
79 {
80 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
81 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
82 ITHACAstream::read_fields(nutFields, nut, "./ITHACAoutput/Offline/");
83
84 // if the offline stage is not completed then resume it
85 if (Ufield.size() < mu.rows())
86 {
87 Vector<double> Uinl(0, 0, 0);
88 label BCind = 0;
89
90 for (label i = Ufield.size(); i < mu.rows(); i++)
91 {
92 Uinl[0] = mu(i, 0);
93 Uinl[1] = mu(i, 1);
94 assignBC(_U(), BCind, Uinl);
95 counter = Ufield.size() + 1;
96 truthSolve("./ITHACAoutput/Offline/");
97 }
98 }
99 }
100 else
101 {
102 Vector<double> Uinl(0, 0, 0);
103 label BCind = 0;
104
105 for (label i = 0; i < mu.rows(); i++)
106 {
107 Uinl[0] = mu(i, 0);
108 Uinl[1] = mu(i, 1);
109 assignBC(_U(), BCind, Uinl);
110 truthSolve("./ITHACAoutput/Offline/");
111 }
112 }
113 }
114
116 void offlineSolve(Eigen::MatrixXd par, fileName folder)
117 {
118 Vector<double> inl(0, 0, 0);
119 List<scalar> mu_now(2);
120 Vector<double> Uinl(0, 0, 0);
121 label BCind = 0;
122
123 for (label i = 0; i < par.rows(); i++)
124 {
125 Uinl[0] = par(i, 0);
126 Uinl[1] = par(i, 1);
127 assignBC(_U(), BCind, Uinl);
128 truthSolve(folder);
129 }
130 }
131
132 void truthSolve(fileName folder)
133 {
134 Time& runTime = _runTime();
135 fvMesh& mesh = _mesh();
136 volScalarField& p = _p();
137 volVectorField& U = _U();
138 surfaceScalarField& phi = _phi();
139 fv::options& fvOptions = _fvOptions();
140 simpleControl& simple = _simple();
141 IOMRFZoneList& MRF = _MRF();
142 singlePhaseTransportModel& laminarTransport = _laminarTransport();
143#include "NLsolveSteadyNSTurb.H"
144 ITHACAstream::exportSolution(U, name(counter), folder);
145 ITHACAstream::exportSolution(p, name(counter), folder);
146 volScalarField _nut(turbulence->nut());
148 Ufield.append((U).clone());
149 Pfield.append((p).clone());
150 nutFields.append((_nut).clone());
151 counter++;
152 }
153};
154
155int main(int argc, char* argv[])
156{
157 // Construct the tutorial06 object
158 tutorial06 example(argc, argv);
159 // Read parameters samples for the offline stage and the ones for the online stages
160 word par_offline("./par_offline");
161 word par_new("./par_online");
162 Eigen::MatrixXd par_online = ITHACAstream::readMatrix(par_new);
163 example.mu = ITHACAstream::readMatrix(par_offline);
164 // Set the inlet boundaries where we have parameterized boundary conditions
165 example.inletIndex.resize(2, 2);
166 example.inletIndex(0, 0) = 0;
167 example.inletIndex(0, 1) = 0;
168 example.inletIndex(1, 0) = 0;
169 example.inletIndex(1, 1) = 1;
171 example._runTime());
172 // Read parameters from ITHACAdict file
173 int NmodesU = para->ITHACAdict->lookupOrDefault<int>("NmodesU", 5);
174 int NmodesP = para->ITHACAdict->lookupOrDefault<int>("NmodesP", 5);
175 int NmodesSUP = para->ITHACAdict->lookupOrDefault<int>("NmodesSUP", 5);
176 int NmodesNUT = para->ITHACAdict->lookupOrDefault<int>("NmodesNUT", 5);
177 int NmodesProject = para->ITHACAdict->lookupOrDefault<int>("NmodesProject", 5);
178 word stabilization = para->ITHACAdict->lookupOrDefault<word>("Stabilization",
179 "supremizer");
180 // Perform The Offline Solve;
181 example.offlineSolve();
182 // Read the lift functions
183 ITHACAstream::read_fields(example.liftfield, example.U, "./lift/");
184 // Create the homogeneous set of snapshots for the velocity field
185 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
186 // Export the homogeneous velocity snapshots
187 ITHACAstream::exportFields(example.Uomfield, "./ITHACAoutput/Offline",
188 "Uofield");
189 // Perform a POD decomposition for the velocity, the pressure and the eddy viscosity
190 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
191 example.podex,
192 example.supex, 0, NmodesProject);
193 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example.p().name(),
194 example.podex,
195 example.supex, 0, NmodesProject);
196 ITHACAPOD::getModes(example.nutFields, example.nutModes, example._nut().name(),
197 example.podex,
198 example.supex, 0, NmodesProject);
199
200 // Solve the supremizer problem based on the pressure modes
201 if (stabilization == "supremizer")
202 {
203 example.solvesupremizer("modes");
204 }
205
206 // Compute the reduced order matrices
207 // Get reduced matrices
208 if (stabilization == "supremizer")
209 {
210 example.projectSUP("./Matrices", NmodesU, NmodesP, NmodesSUP,
211 NmodesNUT);
212 }
213 else if (stabilization == "PPE")
214 {
215 example.projectPPE("./Matrices", NmodesU, NmodesP, NmodesSUP,
216 NmodesNUT);
217 }
218
219 // Create an object of the turbulent class
220 ReducedSteadyNSTurb pod_rbf(
221 example);
222 // Set value of the reduced viscosity and the penalty factor
223 pod_rbf.nu = 1e-3;
224 pod_rbf.tauU.resize(2, 1);
225 // We create the matrix rbfCoeff which will store the values of the interpolation results for the eddy viscosity field
226 Eigen::MatrixXd rbfCoeff;
227 rbfCoeff.resize(NmodesNUT, par_online.rows());
228
229 // Perform an online solve for the new values of inlet velocities
230 for (label k = 0; k < par_online.rows(); k++)
231 {
232 Eigen::MatrixXd velNow(2, 1);
233 velNow(0, 0) = par_online(k, 0);
234 velNow(1, 0) = par_online(k, 1);
235 pod_rbf.tauU(0, 0) = 0;
236 pod_rbf.tauU(1, 0) = 0;
237
238 if (stabilization == "supremizer")
239 {
240 pod_rbf.solveOnlineSUP(velNow);
241 }
242 else if (stabilization == "PPE")
243 {
244 pod_rbf.solveOnlinePPE(velNow);
245 }
246
247 rbfCoeff.col(k) = pod_rbf.rbfCoeff;
248 Eigen::MatrixXd tmp_sol(pod_rbf.y.rows() + 1, 1);
249 tmp_sol(0) = k + 1;
250 tmp_sol.col(0).tail(pod_rbf.y.rows()) = pod_rbf.y;
251 pod_rbf.online_solution.append(tmp_sol);
252 }
253
254 // Save the matrix of interpolated eddy viscosity coefficients
255 ITHACAstream::exportMatrix(rbfCoeff, "rbfCoeff", "python",
256 "./ITHACAoutput/Matrices/");
257 ITHACAstream::exportMatrix(rbfCoeff, "rbfCoeff", "matlab",
258 "./ITHACAoutput/Matrices/");
259 // Save the online solution
260 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "python",
261 "./ITHACAoutput/red_coeff");
262 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "matlab",
263 "./ITHACAoutput/red_coeff");
264 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "eigen",
265 "./ITHACAoutput/red_coeff");
266 pod_rbf.rbfCoeffMat = rbfCoeff;
267 // Reconstruct and export the solution
268 pod_rbf.reconstruct(true, "./ITHACAoutput/Reconstruction/");
269 // Create an object of the laminar class
270 reducedSteadyNS pod_normal(
271 example);
272 // Set value of the reduced viscosity and the penalty factor
273 pod_normal.nu = 1e-3;
274 pod_normal.tauU.resize(2, 1);
275
276 // Perform an online solve for the new values of inlet velocities
277 for (label k = 0; k < par_online.rows(); k++)
278 {
279 Eigen::MatrixXd vel_now(2, 1);
280 vel_now(0, 0) = par_online(k, 0);
281 vel_now(1, 0) = par_online(k, 1);
282 pod_normal.tauU(0, 0) = 0;
283 pod_normal.tauU(1, 0) = 0;
284 pod_normal.solveOnline_sup(vel_now);
285 Eigen::MatrixXd tmp_sol(pod_normal.y.rows() + 1, 1);
286 tmp_sol(0) = k + 1;
287 tmp_sol.col(0).tail(pod_normal.y.rows()) = pod_normal.y;
288 pod_normal.online_solution.append(tmp_sol);
289 }
290
291 // Save the online solution
292 ITHACAstream::exportMatrix(pod_normal.online_solution, "red_coeffnew", "python",
293 "./ITHACAoutput/red_coeffnew");
294 ITHACAstream::exportMatrix(pod_normal.online_solution, "red_coeffnew", "matlab",
295 "./ITHACAoutput/red_coeffnew");
296 ITHACAstream::exportMatrix(pod_normal.online_solution, "red_coeffnew", "eigen",
297 "./ITHACAoutput/red_coeffnew");
298 // Reconstruct and export the solution
299 pod_normal.reconstruct(true, "./ITHACAoutput/Lam_Rec/");
300 // Solve the full order problem for the online velocity values for the purpose of comparison
301 // if (ITHACAutilities::check_folder("./ITHACAoutput/Offline_check") == false)
302 // {
303 // example.offlineSolve(par_online, "./ITHACAoutput/Offline_check/");
304 // ITHACAutilities::createSymLink("./ITHACAoutput/Offline_check");
305 // }
306 Eigen::MatrixXd errFrobU = ITHACAutilities::errorFrobRel(example.Ufield,
307 pod_rbf.uRecFields);
308 Eigen::MatrixXd errFrobP = ITHACAutilities::errorFrobRel(example.Pfield,
309 pod_rbf.pRecFields);
310 Eigen::MatrixXd errFrobNut = ITHACAutilities::errorFrobRel(example.nutFields,
311 pod_rbf.nutRecFields);
312 ITHACAstream::exportMatrix(errFrobU, "errFrobU", "matlab",
313 "./ITHACAoutput/ErrorsFrob/");
314 ITHACAstream::exportMatrix(errFrobP, "errFrobP", "matlab",
315 "./ITHACAoutput/ErrorsFrob/");
316 ITHACAstream::exportMatrix(errFrobNut, "errFrobNut", "matlab",
317 "./ITHACAoutput/ErrorsFrob/");
318 Eigen::MatrixXd errL2U = ITHACAutilities::errorL2Rel(example.Ufield,
319 pod_rbf.uRecFields);
320 Eigen::MatrixXd errL2P = ITHACAutilities::errorL2Rel(example.Pfield,
321 pod_rbf.pRecFields);
322 Eigen::MatrixXd errL2Nut = ITHACAutilities::errorL2Rel(example.nutFields,
323 pod_rbf.nutRecFields);
324 ITHACAstream::exportMatrix(errL2U, "errL2U", "matlab",
325 "./ITHACAoutput/ErrorsL2/");
326 ITHACAstream::exportMatrix(errL2P, "errL2P", "matlab",
327 "./ITHACAoutput/ErrorsL2/");
328 ITHACAstream::exportMatrix(errL2Nut, "errL2Nut", "matlab",
329 "./ITHACAoutput/ErrorsL2/");
330 exit(0);
331}
332
333//--------
337
498
499
500
501
int main(int argc, char *argv[])
Definition 06POD_RBF.C:155
Header file of the EigenFunctions class.
Header file of the Foam2Eigen class.
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 ITHACAutilities namespace.
fv::options & fvOptions
Definition NLsolve.H:25
ITHACAparameters * para
Definition NLsolve.H:40
Header file of the ReducedSteadyNSTurb class.
Header file of the reducedSteadyNS class.
Header file of the SteadyNSTurb 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.
Class where it is implemented a reduced problem for the steady turbulent Navier-stokes problem.
PtrList< volScalarField > nutRecFields
Reconstructed eddy viscosity fields list.
void solveOnlinePPE(Eigen::MatrixXd velNow)
Method to perform an online solve using a PPE stabilisation method.
Eigen::MatrixXd rbfCoeffMat
The matrix of the eddy viscosity RBF interoplated coefficients.
void solveOnlineSUP(Eigen::MatrixXd velNow)
Method to perform an online solve using a supremizer stabilisation method.
void reconstruct(bool exportFields=false, fileName folder="./ITHACAoutput/online_rec", int printevery=1)
Method to reconstruct the solutions from an online solve.
Eigen::VectorXd rbfCoeff
Vector of eddy viscosity RBF interoplated coefficients.
autoPtr< volScalarField > _nut
Eddy viscosity field.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes)
Project using a supremizer approach.
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes)
Project using a PPE approach.
volScalarModes nutModes
List of POD modes for eddy viscosity.
Class where it is implemented a reduced problem for the steady Navier-stokes problem.
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
void solveOnline_sup(Eigen::MatrixXd vel_now)
Method to perform an online solve using a supremizer stabilisation method.
List< Eigen::MatrixXd > online_solution
List of Eigen matrices to store the online solution.
scalar nu
Reduced viscosity in case of parametrized viscosity.
Eigen::MatrixXd tauU
Penalty Factor.
PtrList< volScalarField > pRecFields
Reconstructed pressure fields list.
PtrList< volVectorField > uRecFields
Recontructed velocity fields list.
void reconstruct(bool exportFields=false, fileName folder="./ITHACAoutput/online_rec", int printevery=1)
Method to reconstruct the solutions from an online solve.
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 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.
Eigen::MatrixXd mu
Row matrix of parameters.
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.
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:265
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:302
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:293
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
autoPtr< fv::options > _fvOptions
fvOptions
Definition steadyNS.H:296
void solvesupremizer(word type="snapshots")
solve the supremizer either with the use of the pressure snaphots or the pressure modes
Definition steadyNS.C:137
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:299
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
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:290
autoPtr< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model)
Definition steadyNS.H:311
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:116
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:314
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition steadyNS.H:308
autoPtr< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
Definition steadyNS.H:281
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:272
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:269
Class where the tutorial number 6 is implemented.
Definition 06POD_RBF.C:56
volVectorField & U
[tutorial06]
Definition 06POD_RBF.C:68
volScalarField & p
Definition 06POD_RBF.C:69
tutorial06(int argc, char *argv[])
Definition 06POD_RBF.C:58
volScalarField & nut
Definition 06POD_RBF.C:70
void truthSolve(fileName folder)
Definition 06POD_RBF.C:132
void offlineSolve()
Perform an Offline solve.
Definition 06POD_RBF.C:72
void offlineSolve(Eigen::MatrixXd par, fileName folder)
Perform an Offline solve for a special set of parameter samples called par.
Definition 06POD_RBF.C:116
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
void exportFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &field, word folder, word fieldname)
Function to export a scalar of vector field.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
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.
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.
double errorFrobRel(GeometricField< Type, PatchField, GeoMesh > &field1, GeometricField< Type, PatchField, GeoMesh > &field2, List< label > *labels)
Computes the relative error between two Fields in the Frobenius norm.
Definition ITHACAerror.C:78
Header file of the reductionProblem class.
simpleControl simple(mesh)
surfaceScalarField & phi
singlePhaseTransportModel & laminarTransport
label i
Definition pEqn.H:46
Header file of the steadyNS class.