Loading...
Searching...
No Matches
21unsteadyNSTurb_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 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 an unsteady NS Reduction Problem
26SourceFiles
27 04unsteadyNS.C
28\*---------------------------------------------------------------------------*/
29
30#include "unsteadyNS.H"
31#include "UnsteadyNSTurb.H"
32#include "ITHACAPOD.H"
33#include "ReducedUnsteadyNS.H"
35#include "ITHACAstream.H"
36#include <chrono>
37#include<math.h>
38#include<iomanip>
39
40#include "fvCFD.H"
41#include "fvOptions.H"
42#include "IOmanip.H"
43#include "ITHACAutilities.H"
44#include "Foam2Eigen.H"
45#include <Eigen/Dense>
46#include <Eigen/SVD>
47#include <Eigen/SparseLU>
48#include "EigenFunctions.H"
49#include "reductionProblem.H"
50#include "steadyNS.H"
51
53{
54 public:
55 explicit tutorial21(int argc, char* argv[])
56 :
57 UnsteadyNSTurb(argc, argv),
58 U(_U()),
59 p(_p()),
60 nut(_nut())
61 {}
62
63 // Fields To Perform
64 volVectorField& U;
65 volScalarField& p;
66 volScalarField& nut;
67
68 void offlineSolve(std::string offlinepath)
69 {
70 Vector<double> inl(1, 0, 0);
71 List<scalar> mu_now(1);
72 label BCind = 0;
73
74 if ((offline) && (ITHACAutilities::check_folder(offlinepath) == true))
75 {
76 ITHACAstream::read_fields(Ufield, U, offlinepath);
77 ITHACAstream::read_fields(Pfield, p, offlinepath);
79 }
80 else
81 {
82 for (label i = 0; i < mu.cols(); i++)
83 {
84 inl[0] = mu(0, i);
85 mu_now[0] = mu(0, i);
86 assignBC(U, BCind, inl);
87 truthSolve(mu_now, offlinepath);
88 }
89 }
90 }
91};
92
93/*---------------------------------------------------------------------------*\
94 Starting the MAIN
95\*---------------------------------------------------------------------------*/
96int main(int argc, char* argv[])
97{
98 // Construct the tutorial21 object
99 tutorial21 example(argc, argv);
100 // Read parameters from ITHACAdict file
102 example._runTime());
103 int NmodesU = para->ITHACAdict->lookupOrDefault<int>("NmodesU", 10);
104 int NmodesP = para->ITHACAdict->lookupOrDefault<int>("NmodesP", 10);
105 int NmodesSUP = para->ITHACAdict->lookupOrDefault<int>("NmodesSUP", 10);
106 int NmodesNUT = para->ITHACAdict->lookupOrDefault<int>("NmodesNUT", 10);
107 int NmodesProject = para->ITHACAdict->lookupOrDefault<int>("NmodesProject", 10);
109 example.Pnumber = 1;
111 example.Tnumber = 2;
113 example.setParameters();
114 // Set the parameter ranges
115 example.mu_range(0, 0) = 1.0;
116 example.mu_range(0, 1) = 1.1;
117 // Generate equispaced samples inside the parameter range
118 example.genEquiPar();
119 // Set the inlet boundaries where we have non homogeneous boundary conditions
120 example.inletIndex.resize(1, 2);
121 example.inletIndex(0, 0) = 0;
122 example.inletIndex(0, 1) = 0;
123 // Time parameters
124 example.startTime = 0;
125 example.finalTime = 5;
126 example.timeStep = 0.001;
127 example.writeEvery = 0.1;
128 // Perform The Offline Solve;
129 example.offlineSolve("./ITHACAoutput/Offline/");
130 // Define velRBF
131 auto mu_mat =
132 ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt");
133 example.velRBF = mu_mat.col(1);
134 // Solve the supremizer problem
135 example.solvesupremizer();
136 // Search the lift function
137 example.liftSolve();
138 // Normalize the lifting function
140 // Create homogeneous basis functions for velocity
141 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
142 // Export the homogeneous velocity snapshots
143 ITHACAstream::exportFields(example.Uomfield, "./ITHACAoutput/Offline",
144 "Uofield");
145 // Perform a POD decomposition for the velocity, the pressure and the eddy viscosity
146 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
147 example.podex,
148 example.supex, 0, NmodesProject);
149 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example.p().name(),
150 example.podex,
151 example.supex, 0, NmodesProject);
152 ITHACAPOD::getModes(example.nutFields, example.nutModes, example._nut().name(),
153 example.podex,
154 example.supex, 0, NmodesProject);
155 // Solve the supremizer problem based on the pressure modes
156 example.solvesupremizer("modes");
157 // Compute the reduced order matrices
158 example.projectSUP("./Matrices", NmodesU, NmodesP, NmodesSUP,
159 NmodesNUT, true);
160 // Create an object of the turbulent class
161 ReducedUnsteadyNSTurb pod_rbf(
162 example);
163 // Set value of the reduced viscosity and the penalty factor
164 pod_rbf.nu = 1e-05;
165 pod_rbf.tauU.resize(1, 1);
166 pod_rbf.tstart = 0;
167 pod_rbf.finalTime = 10.1;
168 pod_rbf.dt = 0.005;
169 pod_rbf.storeEvery = 0.005;
170 pod_rbf.exportEvery = 0.1;
171 // We create the matrix rbfCoeff which will store the values of the interpolation results for the eddy viscosity field
172 Eigen::MatrixXd rbfCoeff;
173 rbfCoeff.resize(NmodesNUT + 1, 1);
174
175 // Perform an online solve for the new values of inlet velocities
176 for (label k = 0; k < 1; k++)
177 {
178 Eigen::MatrixXd velNow(1, 1);
179 velNow(0, 0) = 1.05;
180 pod_rbf.tauU(0, 0) = 0;
181 pod_rbf.solveOnlineSUP(velNow);
182 rbfCoeff.col(k) = pod_rbf.rbfCoeffMat.col(k);
183 Eigen::MatrixXd tmp_sol(pod_rbf.y.rows() + 1, 1);
184 tmp_sol(0) = k + 1;
185 tmp_sol.col(0).tail(pod_rbf.y.rows()) = pod_rbf.y;
186 pod_rbf.online_solution.append(tmp_sol);
187 }
188
189 // Save the matrix of interpolated eddy viscosity coefficients
190 ITHACAstream::exportMatrix(rbfCoeff, "rbfCoeff", "python",
191 "./ITHACAoutput/Matrices/");
192 ITHACAstream::exportMatrix(rbfCoeff, "rbfCoeff", "matlab",
193 "./ITHACAoutput/Matrices/");
194 // Save the online solution
195 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "python",
196 "./ITHACAoutput/red_coeff");
197 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "matlab",
198 "./ITHACAoutput/red_coeff");
199 ITHACAstream::exportMatrix(pod_rbf.online_solution, "red_coeff", "eigen",
200 "./ITHACAoutput/red_coeff");
201 // Reconstruct and export the solution
202 pod_rbf.reconstruct(true, "./ITHACAoutput/Reconstruction/");
203 // Solve the full order problem for the online velocity values for the purpose of comparison
204 tutorial21 example2(argc, argv);
206 example2.Pnumber = 1;
208 example2.Tnumber = 2;
210 example2.setParameters();
211 // Set the parameter ranges
212 example2.mu_range(0, 0) = 1.05;
213 example2.mu_range(0, 1) = 1.05;
214 // Generate equispaced samples inside the parameter range
215 example2.genEquiPar();
216 // Set the inlet boundaries where we have non homogeneous boundary conditions
217 example2.inletIndex.resize(1, 2);
218 example2.inletIndex(0, 0) = 0;
219 example2.inletIndex(0, 1) = 0;
220 // Time parameters
221 example2.startTime = 0;
222 example2.finalTime = 5;
223 example2.timeStep = 0.001;
224 example2.writeEvery = 0.1;
225 example2.offlineSolve("./ITHACAoutput/Offline_check/");
226 Eigen::MatrixXd errFrobU = ITHACAutilities::errorFrobRel(example2.Ufield,
227 pod_rbf.uRecFields);
228 Eigen::MatrixXd errFrobP = ITHACAutilities::errorFrobRel(example2.Pfield,
229 pod_rbf.pRecFields);
230 Eigen::MatrixXd errFrobNut = ITHACAutilities::errorFrobRel(example2.nutFields,
231 pod_rbf.nutRecFields);
232 ITHACAstream::exportMatrix(errFrobU, "errFrobU", "matlab",
233 "./ITHACAoutput/ErrorsFrob/");
234 ITHACAstream::exportMatrix(errFrobP, "errFrobP", "matlab",
235 "./ITHACAoutput/ErrorsFrob/");
236 ITHACAstream::exportMatrix(errFrobNut, "errFrobNut", "matlab",
237 "./ITHACAoutput/ErrorsFrob/");
238 Eigen::MatrixXd errL2U = ITHACAutilities::errorL2Rel(example2.Ufield,
239 pod_rbf.uRecFields);
240 Eigen::MatrixXd errL2P = ITHACAutilities::errorL2Rel(example2.Pfield,
241 pod_rbf.pRecFields);
242 Eigen::MatrixXd errL2Nut = ITHACAutilities::errorL2Rel(example2.nutFields,
243 pod_rbf.nutRecFields);
244 ITHACAstream::exportMatrix(errL2U, "errL2U", "matlab",
245 "./ITHACAoutput/ErrorsL2/");
246 ITHACAstream::exportMatrix(errL2P, "errL2P", "matlab",
247 "./ITHACAoutput/ErrorsL2/");
248 ITHACAstream::exportMatrix(errL2Nut, "errL2Nut", "matlab",
249 "./ITHACAoutput/ErrorsL2/");
250 exit(0);
251}
int main(int argc, char *argv[])
Header file of the EigenFunctions class.
Header file of the Foam2Eigen class.
Header file of the ITHACAPOD class.
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 ReducedUnsteadyNSTurb class.
Header file of the reducedUnsteadyNS class.
Header file of the UnsteadyNSTurb 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.
PtrList< volScalarField > nutRecFields
Reconstructed eddy viscosity fields list.
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="./online_rec")
Method to reconstruct the solutions from an online solve with any of the two techniques SUP or the PP...
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes, bool rbfInterp=true)
Project using a supremizer approach.
autoPtr< volScalarField > _nut
Eddy viscosity field.
volScalarModes nutModes
List of POD modes for eddy viscosity.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
Eigen::MatrixXd velRBF
Velocity coefficients for RBF interpolation.
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)
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
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.
double exportEvery
A variable for exporting the fields.
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.
double storeEvery
A variable for storing the reduced coefficients.
label Pnumber
Number of parameters.
void assignBC(volVectorField &s, label BC_ind, Vector< double > &value)
Assign Boundary Condition to a volVectorField.
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.
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.
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:256
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
void solvesupremizer(word type="snapshots")
solve the supremizer either with the use of the pressure snaphots or the pressure modes
Definition steadyNS.C:136
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
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:281
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
void liftSolve()
Perform a lift solve.
Definition steadyNS.C:252
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 & p
volVectorField & U
void offlineSolve(std::string offlinepath)
volScalarField & nut
tutorial21(int argc, char *argv[])
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 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.
void normalizeFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &fields)
Normalize list of Geometric fields.
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
bool check_folder(word folder)
Checks if a folder exists.
Header file of the reductionProblem class.
label i
Definition pEqn.H:46
Header file of the steadyNS class.
Header file of the unsteadyNS class.