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 {}
66 // Relevant Fields
67 volVectorField& U;
68 volScalarField& p;
69 volScalarField& nut;
72 {
73 Vector<double> inl(0, 0, 0);
74 List<scalar> mu_now(2);
75
76 // if the offline solution is already performed read the fields
77 if (offline)
78 {
79 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
80 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
81 ITHACAstream::read_fields(nutFields, nut, "./ITHACAoutput/Offline/");
82
83 // if the offline stage is not completed then resume it
84 if (Ufield.size() < mu.rows())
85 {
86 Vector<double> Uinl(0, 0, 0);
87 label BCind = 0;
88
89 for (label i = Ufield.size(); i < mu.rows(); i++)
90 {
91 Uinl[0] = mu(i, 0);
92 Uinl[1] = mu(i, 1);
93 assignBC(_U(), BCind, Uinl);
94 counter = Ufield.size() + 1;
95 truthSolve("./ITHACAoutput/Offline/");
96 }
97 }
98 }
99 else
100 {
101 Vector<double> Uinl(0, 0, 0);
102 label BCind = 0;
103
104 for (label i = 0; i < mu.rows(); i++)
105 {
106 Uinl[0] = mu(i, 0);
107 Uinl[1] = mu(i, 1);
108 assignBC(_U(), BCind, Uinl);
109 truthSolve("./ITHACAoutput/Offline/");
110 }
111 }
112 }
114 void offlineSolve(Eigen::MatrixXd par, fileName folder)
115 {
116 Vector<double> inl(0, 0, 0);
117 List<scalar> mu_now(2);
118 Vector<double> Uinl(0, 0, 0);
119 label BCind = 0;
120
121 for (label i = 0; i < par.rows(); i++)
122 {
123 Uinl[0] = par(i, 0);
124 Uinl[1] = par(i, 1);
125 assignBC(_U(), BCind, Uinl);
126 truthSolve(folder);
127 }
128 }
129 void truthSolve(fileName folder)
130 {
131 Time& runTime = _runTime();
132 fvMesh& mesh = _mesh();
133 volScalarField& p = _p();
134 volVectorField& U = _U();
135 surfaceScalarField& phi = _phi();
136 fv::options& fvOptions = _fvOptions();
137 simpleControl& simple = _simple();
138 IOMRFZoneList& MRF = _MRF();
139 singlePhaseTransportModel& laminarTransport = _laminarTransport();
140#include "NLsolveSteadyNSTurb.H"
141 ITHACAstream::exportSolution(U, name(counter), folder);
142 ITHACAstream::exportSolution(p, name(counter), folder);
143 volScalarField _nut(turbulence->nut());
145 Ufield.append((U).clone());
146 Pfield.append((p).clone());
147 nutFields.append((_nut).clone());
148 counter++;
149 }
150};
151
152int main(int argc, char* argv[])
153{
154 // Construct the tutorial06 object
155 tutorial06 example(argc, argv);
156 // Read parameters samples for the offline stage and the ones for the online stages
157 word par_offline("./par_offline");
158 word par_new("./par_online");
159 Eigen::MatrixXd par_online = ITHACAstream::readMatrix(par_new);
160 example.mu = ITHACAstream::readMatrix(par_offline);
161 // Set the inlet boundaries where we have parameterized boundary conditions
162 example.inletIndex.resize(2, 2);
163 example.inletIndex(0, 0) = 0;
164 example.inletIndex(0, 1) = 0;
165 example.inletIndex(1, 0) = 0;
166 example.inletIndex(1, 1) = 1;
167 ITHACAparameters* para = ITHACAparameters::getInstance(example._mesh(),
168 example._runTime());
169 // Read parameters from ITHACAdict file
170 int NmodesU = para->ITHACAdict->lookupOrDefault<int>("NmodesU", 5);
171 int NmodesP = para->ITHACAdict->lookupOrDefault<int>("NmodesP", 5);
172 int NmodesSUP = para->ITHACAdict->lookupOrDefault<int>("NmodesSUP", 5);
173 int NmodesNUT = para->ITHACAdict->lookupOrDefault<int>("NmodesNUT", 5);
174 int NmodesProject = para->ITHACAdict->lookupOrDefault<int>("NmodesProject", 5);
175 word stabilization = para->ITHACAdict->lookupOrDefault<word>("Stabilization",
176 "supremizer");
177 // Perform The Offline Solve;
178 example.offlineSolve();
179 // Read the lift functions
180 ITHACAstream::read_fields(example.liftfield, example.U, "./lift/");
181 // Create the homogeneous set of snapshots for the velocity field
182 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
183 // Export the homogeneous velocity snapshots
184 ITHACAstream::exportFields(example.Uomfield, "./ITHACAoutput/Offline",
185 "Uofield");
186 // Perform a POD decomposition for the velocity, the pressure and the eddy viscosity
187 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
188 example.podex,
189 example.supex, 0, NmodesProject);
190 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example.p().name(),
191 example.podex,
192 example.supex, 0, NmodesProject);
193 ITHACAPOD::getModes(example.nutFields, example.nutModes, example._nut().name(),
194 example.podex,
195 example.supex, 0, NmodesProject);
196
197 // Solve the supremizer problem based on the pressure modes
198 if (stabilization == "supremizer")
199 {
200 example.solvesupremizer("modes");
201 }
202
203 // Compute the reduced order matrices
204 // Get reduced matrices
205 if (stabilization == "supremizer")
206 {
207 example.projectSUP("./Matrices", NmodesU, NmodesP, NmodesSUP,
208 NmodesNUT);
209 }
210 else if (stabilization == "PPE")
211 {
212 example.projectPPE("./Matrices", NmodesU, NmodesP, NmodesSUP,
213 NmodesNUT);
214 }
215
216 // Create an object of the turbulent class
217 ReducedSteadyNSTurb pod_rbf(
218 example);
219 // Set value of the reduced viscosity and the penalty factor
220 pod_rbf.nu = 1e-3;
221 pod_rbf.tauU.resize(2, 1);
222 // We create the matrix rbfCoeff which will store the values of the interpolation results for the eddy viscosity field
223 Eigen::MatrixXd rbfCoeff;
224 rbfCoeff.resize(NmodesNUT, par_online.rows());
225
226 // Perform an online solve for the new values of inlet velocities
227 for (label k = 0; k < par_online.rows(); k++)
228 {
229 Eigen::MatrixXd velNow(2, 1);
230 velNow(0, 0) = par_online(k, 0);
231 velNow(1, 0) = par_online(k, 1);
232 pod_rbf.tauU(0, 0) = 0;
233 pod_rbf.tauU(1, 0) = 0;
234
235 if (stabilization == "supremizer")
236 {
237 pod_rbf.solveOnlineSUP(velNow);
238 }
239 else if (stabilization == "PPE")
240 {
241 pod_rbf.solveOnlinePPE(velNow);
242 }
243
244 rbfCoeff.col(k) = pod_rbf.rbfCoeff;
245 Eigen::MatrixXd tmp_sol(pod_rbf.y.rows() + 1, 1);
246 tmp_sol(0) = k + 1;
247 tmp_sol.col(0).tail(pod_rbf.y.rows()) = pod_rbf.y;
248 pod_rbf.online_solution.append(tmp_sol);
249 }
250
251 // Save the matrix of interpolated eddy viscosity coefficients
252 ITHACAstream::exportMatrix(rbfCoeff, "rbfCoeff", "python",
253 "./ITHACAoutput/Matrices/");
254 ITHACAstream::exportMatrix(rbfCoeff, "rbfCoeff", "matlab",
255 "./ITHACAoutput/Matrices/");
256 // Save the online solution
257 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "python",
258 "./ITHACAoutput/red_coeff");
259 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "matlab",
260 "./ITHACAoutput/red_coeff");
261 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "eigen",
262 "./ITHACAoutput/red_coeff");
263 pod_rbf.rbfCoeffMat = rbfCoeff;
264 // Reconstruct and export the solution
265 pod_rbf.reconstruct(true, "./ITHACAoutput/Reconstruction/");
266 // Create an object of the laminar class
267 reducedSteadyNS pod_normal(
268 example);
269 // Set value of the reduced viscosity and the penalty factor
270 pod_normal.nu = 1e-3;
271 pod_normal.tauU.resize(2, 1);
272
273 // Perform an online solve for the new values of inlet velocities
274 for (label k = 0; k < par_online.rows(); k++)
275 {
276 Eigen::MatrixXd vel_now(2, 1);
277 vel_now(0, 0) = par_online(k, 0);
278 vel_now(1, 0) = par_online(k, 1);
279 pod_normal.tauU(0, 0) = 0;
280 pod_normal.tauU(1, 0) = 0;
281 pod_normal.solveOnline_sup(vel_now);
282 Eigen::MatrixXd tmp_sol(pod_normal.y.rows() + 1, 1);
283 tmp_sol(0) = k + 1;
284 tmp_sol.col(0).tail(pod_normal.y.rows()) = pod_normal.y;
285 pod_normal.online_solution.append(tmp_sol);
286 }
287
288 // Save the online solution
289 ITHACAstream::exportMatrix(pod_normal.online_solution, "red_coeffnew", "python",
290 "./ITHACAoutput/red_coeffnew");
291 ITHACAstream::exportMatrix(pod_normal.online_solution, "red_coeffnew", "matlab",
292 "./ITHACAoutput/red_coeffnew");
293 ITHACAstream::exportMatrix(pod_normal.online_solution, "red_coeffnew", "eigen",
294 "./ITHACAoutput/red_coeffnew");
295 // Reconstruct and export the solution
296 pod_normal.reconstruct(true, "./ITHACAoutput/Lam_Rec/");
297 // Solve the full order problem for the online velocity values for the purpose of comparison
298 // if (ITHACAutilities::check_folder("./ITHACAoutput/Offline_check") == false)
299 // {
300 // example.offlineSolve(par_online, "./ITHACAoutput/Offline_check/");
301 // ITHACAutilities::createSymLink("./ITHACAoutput/Offline_check");
302 // }
303 Eigen::MatrixXd errFrobU = ITHACAutilities::errorFrobRel(example.Ufield,
304 pod_rbf.uRecFields);
305 Eigen::MatrixXd errFrobP = ITHACAutilities::errorFrobRel(example.Pfield,
306 pod_rbf.pRecFields);
307 Eigen::MatrixXd errFrobNut = ITHACAutilities::errorFrobRel(example.nutFields,
308 pod_rbf.nutRecFields);
309 ITHACAstream::exportMatrix(errFrobU, "errFrobU", "matlab",
310 "./ITHACAoutput/ErrorsFrob/");
311 ITHACAstream::exportMatrix(errFrobP, "errFrobP", "matlab",
312 "./ITHACAoutput/ErrorsFrob/");
313 ITHACAstream::exportMatrix(errFrobNut, "errFrobNut", "matlab",
314 "./ITHACAoutput/ErrorsFrob/");
315 Eigen::MatrixXd errL2U = ITHACAutilities::errorL2Rel(example.Ufield,
316 pod_rbf.uRecFields);
317 Eigen::MatrixXd errL2P = ITHACAutilities::errorL2Rel(example.Pfield,
318 pod_rbf.pRecFields);
319 Eigen::MatrixXd errL2Nut = ITHACAutilities::errorL2Rel(example.nutFields,
320 pod_rbf.nutRecFields);
321 ITHACAstream::exportMatrix(errL2U, "errL2U", "matlab",
322 "./ITHACAoutput/ErrorsL2/");
323 ITHACAstream::exportMatrix(errL2P, "errL2P", "matlab",
324 "./ITHACAoutput/ErrorsL2/");
325 ITHACAstream::exportMatrix(errL2Nut, "errL2Nut", "matlab",
326 "./ITHACAoutput/ErrorsL2/");
327 exit(0);
328}
329
330//--------
334
495
496
497
498
int main(int argc, char *argv[])
Definition 06POD_RBF.C:152
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.
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...
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 steady turbulent Navier-stokes problem.
Implementation of a parametrized full order steady turbulent Navier Stokes problem and preparation ...
autoPtr< volScalarField > _nut
Eddy viscosity field.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
Class where it is implemented a reduced problem for the steady Navier-stokes problem.
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.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
Eigen::MatrixXd mu
Row matrix of parameters.
void truthSolve()
Perform a TruthSolve.
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:284
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:287
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:290
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition steadyNS.H:89
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:281
autoPtr< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model)
Definition steadyNS.H:302
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:305
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition steadyNS.H:299
autoPtr< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
Definition steadyNS.H:272
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
Class where the tutorial number 6 is implemented.
Definition 06POD_RBF.C:56
volVectorField & U
[tutorial06]
Definition 06POD_RBF.C:67
volScalarField & p
Definition 06POD_RBF.C:68
tutorial06(int argc, char *argv[])
Definition 06POD_RBF.C:58
volScalarField & nut
Definition 06POD_RBF.C:69
void truthSolve(fileName folder)
Definition 06POD_RBF.C:129
void offlineSolve()
Perform an Offline solve.
Definition 06POD_RBF.C:71
void offlineSolve(Eigen::MatrixXd par, fileName folder)
Perform an Offline solve for a special set of parameter samples called par.
Definition 06POD_RBF.C:114
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 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.