Loading...
Searching...
No Matches
18simpleTurbNS.C
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 turbulent steady NS Reduction Problem solved by the use of the SIMPLE algorithm
26SourceFiles
27 18simpleTurbNS.C
28\*---------------------------------------------------------------------------*/
29
30#include "SteadyNSSimple.H"
31#include "ITHACAstream.H"
32#include "ITHACAPOD.H"
34#include "forces.H"
35#include "IOmanip.H"
36
37
39{
40 public:
42 explicit tutorial18(int argc, char* argv[])
43 :
44 SteadyNSSimple(argc, argv),
45 U(_U()),
46 p(_p()),
47 phi(_phi())
48 {}
49
51 volVectorField& U;
53 volScalarField& p;
55 surfaceScalarField& phi;
56
59 {
60 Vector<double> inl(0, 0, 0);
61 List<scalar> mu_now(1);
62
63 // if the offline solution is already performed read the fields
65 {
66 ITHACAstream::readMiddleFields(Ufield, U, "./ITHACAoutput/Offline/");
67 ITHACAstream::readMiddleFields(Pfield, p, "./ITHACAoutput/Offline/");
68 auto nut = _mesh().lookupObject<volScalarField>("nut");
69 ITHACAstream::readLastFields(nutFields, nut, "./ITHACAoutput/Offline/");
71 ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt");
72 }
73 else if (offline)
74 {
75 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
76 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
78 ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt");
79 }
80 else
81 {
82 Vector<double> Uinl(1, 0, 0);
83 label BCind = 0;
84
85 for (label i = 0; i < mu.rows(); i++)
86 {
87 mu_now[0] = mu(i, 0);
88 change_viscosity(mu_now[0]);
89 assignIF(U, Uinl);
90 truthSolve2(mu_now);
91 nutFields.clear();
92 auto nut = _mesh().lookupObject<volScalarField>("nut");
93 ITHACAstream::readConvergedFields(nutFields, nut, "./ITHACAoutput/Offline/");
94 }
95 }
96 }
97
98};
99
100int main(int argc, char* argv[])
101{
102 // Construct the tutorial object
103 tutorial18 example(argc, argv);
104 // Read some parameters from file
105 ITHACAparameters* para = ITHACAparameters::getInstance(example._mesh(),
106 example._runTime());
107 // Read the par file where the parameters are stored
108 std::ifstream exFileOff("./parsOff_mat.txt");
109
110 if (exFileOff)
111 {
112 example.mu = ITHACAstream::readMatrix("./parsOff_mat.txt");
113 }
114 else
115 {
116 example.mu = Eigen::VectorXd::LinSpaced(50, 1.00e-04, 1.00e-05);
117 ITHACAstream::exportMatrix(example.mu, "parsOff", "eigen", "./");
118 }
119
120 Eigen::MatrixXd parOn;
121 std::ifstream exFileOn("./parsOn_mat.txt");
122
123 if (exFileOn)
124 {
125 parOn = ITHACAstream::readMatrix("./parsOn_mat.txt");
126 }
127 else
128 {
129 parOn = ITHACAutilities::rand(20, 1, 1.00e-04, 1.00e-05);
130 ITHACAstream::exportMatrix(parOn, "parsOn", "eigen", "./");
131 }
132
133 // Set the inlet boundaries patch 0 directions x and y
134 example.inletIndex.resize(1, 2);
135 example.inletIndex(0, 0) = 0;
136 example.inletIndex(0, 1) = 0;
137 // Set the maximum iterations number for the offline phase
138 example.maxIter = para->ITHACAdict->lookupOrDefault<int>("maxIter", 2000);
139 // Perform the offline solve
140 example.offlineSolve();
141 ITHACAstream::read_fields(example.liftfield, example.U, "./lift/");
142 // Homogenize the snapshots
143 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
144 // Perform POD on velocity and pressure and store the first 10 modes
145 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
146 example.podex, 0, 0,
147 example.NUmodesOut);
148 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
149 example.podex, 0, 0,
150 example.NPmodesOut);
151
153 {
154 ITHACAPOD::getModes(example.nutFields, example.nutModes, "nut",
155 example.podex, 0, 0, example.NNutModesOut);
156 // Create the RBF for turbulence
157 example.getTurbRBF(example.NNutModes);
158 }
159
160 // Create the reduced object
161 reducedSimpleSteadyNS reduced(example);
162 PtrList<volVectorField> U_rec_list;
163 PtrList<volScalarField> P_rec_list;
164 // Reads inlet velocities boundary conditions.
165 word vel_file(para->ITHACAdict->lookup("online_velocities"));
166 Eigen::MatrixXd vel = ITHACAstream::readMatrix(vel_file);
167 // Set the maximum iterations number for the online phase
168 reduced.maxIterOn = para->ITHACAdict->lookupOrDefault<int>("maxIterOn", 2000);
169
170 //Perform the online solutions
171 for (label k = 0; k < parOn.size(); k++)
172 {
173 scalar mu_now = parOn(k, 0);
174 example.restart();
175 example.change_viscosity(mu_now);
176 reduced.setOnlineVelocity(vel);
177
179 {
180 reduced.solveOnline_Simple(mu_now, example.NUmodes, example.NPmodes,
181 example.NNutModes);
182 }
183 else
184 {
185 reduced.solveOnline_Simple(mu_now, example.NUmodes, example.NPmodes,
186 example.NNutModes);
187 }
188 }
189
190 return 0;
191}
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 reducedSteadyNS class.
Header file of the steadyNS 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.
SteadyNSSimple()
Null constructor.
void truthSolve2(List< scalar > mu_now, word Folder="./ITHACAoutput/Offline/")
Offline solver for the whole Navier-Stokes problem.
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.
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
void assignIF(T &s, G &value)
Assign internal field condition.
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 change_viscosity(double mu)
Function to change the viscosity.
Definition steadyNS.C:2066
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:302
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
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< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
Definition steadyNS.H:281
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:272
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:269
void offlineSolve()
Perform an Offline solve.
volVectorField & U
Velocity field.
volScalarField & p
Pressure field.
tutorial18(int argc, char *argv[])
Constructor.
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 readConvergedFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, GeometricField< Type, PatchField, GeoMesh > &field, fileName casename)
Function to read a list of volVectorField from name of the field including only converged snapshots.
void readMiddleFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, GeometricField< Type, PatchField, GeoMesh > &field, fileName casename)
Funtion to read a list of volVectorField from name of the field including all the intermediate snapsh...
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 readLastFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, const GeometricField< Type, PatchField, GeoMesh > &field, const fileName casename)
Funtion to read a list of volVectorField from name of the field including only the last snapshots.
Eigen::MatrixXd rand(label rows, label cols, double min, double max)
Generates random matrix with random values in an interval.
bool isTurbulent()
This function checks if the case is turbulent.