Loading...
Searching...
No Matches
03steadyNS.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 steady NS Reduction Problem
26SourceFiles
27 03steadyNS.C
28\*---------------------------------------------------------------------------*/
29
30#include "IOmanip.H"
31#include "ITHACAPOD.H"
32#include "ITHACAstream.H"
33#include "ReducedSteadyNS.H"
34#include "forces.H"
35#include "steadyNS.H"
36
37class tutorial03 : public steadyNS
38{
39 public:
41 explicit tutorial03(int argc, char* argv[])
42 : steadyNS(argc, argv), U(_U()), p(_p()), args(_args()) {}
43
45 volVectorField& U;
47 volScalarField& p;
49 argList& args;
50
53 {
54 Vector<double> inl(0, 0, 0);
55 List<scalar> mu_now(1);
56
57 // if the offline solution is already performed read the fields
58 if (offline)
59 {
60 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
61 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
63 ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt");
64 }
65 else
66 {
67 Vector<double> Uinl(0, 0, 0);
68
69 for (label i = 0; i < mu.cols(); i++)
70 {
71 mu_now[0] = mu(0, i);
73 truthSolve(mu_now);
74 restart();
75 }
76 }
77 }
78};
79
80void offline_stage(tutorial03& example);
81void online_stage(tutorial03& example);
82
83int main(int argc, char* argv[])
84{
85#if OPENFOAM > 1812
86 // load stage to perform
87 argList::addOption("stage", "offline", "Perform offline stage");
88 argList::addOption("stage", "online", "Perform online stage");
89 // add options for parallel run
90 HashTable<string> validParOptions;
91 validParOptions.set
92 (
93 "stage",
94 "offline"
95 );
96 validParOptions.set
97 (
98 "stage",
99 "online"
100 );
101 Pstream::addValidParOptions(validParOptions);
102 // Construct the tutorial object
103 tutorial03 example(argc, argv);
104
105 if (example._args().get("stage").match("offline"))
106 {
107 // perform the offline stage, extracting the modes from the snapshots'
108 // dataset corresponding to parOffline
109 offline_stage(example);
110 }
111 else if (example._args().get("stage").match("online"))
112 {
113 // load precomputed modes and reduced matrices
114 offline_stage(example);
115 // perform online solve with respect to the parameters in parOnline
116 online_stage(example);
117 }
118 else
119 {
120 Info << "Pass '-stage offline', '-stage online'" << endl;
121 }
122
123#else
124
125 if (argc == 1)
126 {
127 Info << "Pass 'offline' or 'online' as first arguments."
128 << endl;
129 exit(0);
130 }
131
132 // process arguments removing "offline" or "online" keywords
133 int argc_proc = argc - 1;
134 char* argv_proc[argc_proc];
135 argv_proc[0] = argv[0];
136
137 if (argc > 2)
138 {
139 std::copy(argv + 2, argv + argc, argv_proc + 1);
140 }
141
142 argc--;
143 // Construct the tutorial object
144 tutorial03 example(argc, argv);
145
146 if (std::strcmp(argv[1], "offline") == 0)
147 {
148 // perform the offline stage, extracting the modes from the snapshots' dataset corresponding to parOffline
149 offline_stage(example);
150 }
151 else if (std::strcmp(argv[1], "online") == 0)
152 {
153 // load precomputed modes and reduced matrices
154 offline_stage(example);
155 // perform online solve with respect to the parameters in parOnline
156 online_stage(example);
157 }
158 else
159 {
160 Info << "Pass offline, online" << endl;
161 }
162
163#endif
164 exit(0);
165}
166
168{
169 // Read some parameters from file
170 ITHACAparameters* para =
171 ITHACAparameters::getInstance(example._mesh(), example._runTime());
172 int NmodesUout = para->ITHACAdict->lookupOrDefault<int>("NmodesUout", 15);
173 int NmodesPout = para->ITHACAdict->lookupOrDefault<int>("NmodesPout", 15);
174 int NmodesSUPout = para->ITHACAdict->lookupOrDefault<int>("NmodesSUPout", 15);
175 int NmodesUproj = para->ITHACAdict->lookupOrDefault<int>("NmodesUproj", 10);
176 int NmodesPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesPproj", 10);
177 int NmodesSUPproj =
178 para->ITHACAdict->lookupOrDefault<int>("NmodesSUPproj", 10);
179 // Read the par file where the training parameters are stored
180 word filename("./parOffline");
181 example.mu = ITHACAstream::readMatrix(filename);
182 // Set the inlet boundaries patch 0 directions x and y
183 example.inletIndex.resize(1, 2);
184 example.inletIndex(0, 0) = 0;
185 example.inletIndex(0, 1) = 0;
186 // Perform the offline solve
187 example.offlineSolve();
188 // Solve the supremizer problem
189 example.solvesupremizer();
190
191 if (example.bcMethod == "lift")
192 {
193 ITHACAstream::read_fields(example.liftfield, example._U(), "./lift/");
195 // Homogenize the snapshots
196 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
197 // Perform POD on the velocity snapshots
198 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
199 example.podex, 0, 0, NmodesUout);
200 }
201 else
202 {
203 ITHACAPOD::getModes(example.Ufield, example.Umodes, example._U().name(),
204 example.podex, 0, 0, NmodesUout);
205 }
206
207 // Perform POD on velocity pressure and supremizers and store the first 10
208 // modes
209 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
210 example.podex, 0, 0, NmodesPout);
211 ITHACAPOD::getModes(example.supfield, example.supmodes, example._U().name(),
212 example.podex, example.supex, 1, NmodesSUPout);
213 // Perform the Galerkin Projection
214 example.projectSUP("./Matrices", NmodesUproj, NmodesPproj, NmodesSUPproj);
215}
216
218{
219 // Create the reduced object
220 reducedSteadyNS reduced(example);
221 // Read the par file where the test parameters are stored
222 word filename("./parOnline");
223 example.mu = ITHACAstream::readMatrix(filename);
224 // Set the inlet velocity
225 Eigen::MatrixXd vel_now(1, 1);
226 vel_now(0, 0) = 1;
227 // used only for penalty approach
228 reduced.tauU = Eigen::MatrixXd::Zero(1, 1);
229 reduced.tauU(0, 0) = 1e-1;
230
231 // Perform an online solve for the new values of inlet velocities
232 for (label k = 0; k < example.mu.size(); k++)
233 {
234 Info << "Evaluation of the reduced order model on the test set" << endl;
235 Info << "Inlet Ux = " << vel_now(0, 0) << " nu = " << example.mu(0, k)
236 << endl;
237 // Set the reduced viscosity
238 reduced.nu = example.mu(0, k);
239 reduced.solveOnline_sup(vel_now);
240 Eigen::MatrixXd tmp_sol(reduced.y.rows() + 1, 1);
241 tmp_sol(0) = k + 1;
242 tmp_sol.col(0).tail(reduced.y.rows()) = reduced.y;
243 reduced.online_solution.append(tmp_sol);
244 }
245
246 // Save the online solution
247 ITHACAstream::exportMatrix(reduced.online_solution, "red_coeff", "python",
248 "./ITHACAoutput/red_coeff");
249 ITHACAstream::exportMatrix(reduced.online_solution, "red_coeff", "matlab",
250 "./ITHACAoutput/red_coeff");
251 ITHACAstream::exportMatrix(reduced.online_solution, "red_coeff", "eigen",
252 "./ITHACAoutput/red_coeff");
253 // Reconstruct and export the solution
254 reduced.reconstruct(true, "./ITHACAoutput/Reconstruction/");
255
256 if (Pstream::parRun())
257 {
258 bool endedOnline = true;
259 reduce(endedOnline, sumOp<label>());
260 }
261}
262
263//--------
267
int main(int argc, char *argv[])
Definition 03steadyNS.C:83
void offline_stage(tutorial03 &example)
Definition 03steadyNS.C:167
void online_stage(tutorial03 &example)
Definition 03steadyNS.C:217
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.
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 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.
void reconstruct(bool exportFields=false, fileName folder="./ITHACAoutput/online_rec", int printevery=1)
Method to reconstruct the solutions from an online solve.
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
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.
autoPtr< argList > _args
argList
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.
Implementation of a parametrized full order steady NS problem and preparation of the the reduced ma...
Definition steadyNS.H:70
void change_viscosity(double mu)
Function to change the viscosity.
Definition steadyNS.C:2014
void restart()
set U and P back to the values into the 0 folder
Definition steadyNS.C:2253
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
volVectorModes supmodes
List of pointers used to form the supremizer modes.
Definition steadyNS.H:104
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
PtrList< volVectorField > supfield
List of pointers used to form the supremizer snapshots matrix.
Definition steadyNS.H:92
autoPtr< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
Definition steadyNS.H:272
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes)
Project using a supremizer approach.
Definition steadyNS.C:543
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
word bcMethod
Boundary Method.
Definition steadyNS.H:317
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
void offlineSolve()
Perform an Offline solve.
Definition 03steadyNS.C:52
volScalarField & p
Pressure field.
Definition 03steadyNS.C:47
argList & args
Arg List.
Definition 03steadyNS.C:49
volVectorField & U
Velocity field.
Definition 03steadyNS.C:45
tutorial03(int argc, char *argv[])
Constructor.
Definition 03steadyNS.C:41
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 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.
label i
Definition pEqn.H:46
Header file of the steadyNS class.