Loading...
Searching...
No Matches
offline.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8 * In real Time Highly Advanced Computational Applications for Finite Volumes
9 * Copyright (C) 2017 by the ITHACA-FV authors
10-------------------------------------------------------------------------------
11License
12 This file is part of ITHACA-FV
13
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
19 ITHACA-FV is distributed in the hope that it will be useful,
20 but WITHOUT ANY WARRANTY; without even the implied warranty of
21 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22 GNU Lesser General Public License for more details.
23
24 You should have received a copy of the GNU Lesser General Public License
25 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
26
27Description
28 Example of NS-Stokes Reduction Problem
29
30\*---------------------------------------------------------------------------*/
31#include "steadyNS.H"
32#include "UnsteadyNSTurb.H"
33#include <math.h>
34
35#include "ITHACAutilities.H"
36#include "reductionProblem.H"
37
38#include "ReducedSteadyNS.H"
40#include "ReducedUnsteadyNS.H"
41
42#include <Eigen/Dense>
43#include <Eigen/Core>
44#include "forces.H"
45#include "IOmanip.H"
46#include "forces.H"
47#include "forceCoeffs.H"
48#include "fvCFD.H"
49#include "fvOptions.H"
50#include "simpleControl.H"
51#include "Foam2Eigen.H"
52#include <chrono>
53#include "ITHACAstream.H"
54#include "ITHACAPOD.H"
55#include "EigenFunctions.H"
56#include <chrono>
57#include <Eigen/SVD>
58#include <Eigen/SparseLU>
59
60#include <stdio.h>
61#include <sys/types.h>
62#include <dirent.h>
63#include <algorithm>
64#include <fstream>
65#include <string>
66#include <stdexcept>
67#include <sstream>
68#include <vector>
69#include <cstdio>
70#include <typeinfo>
71#include <iostream>
72#include <cassert>
73#include <zlib.h>
74
75
77{
78 public:
79 explicit tutorial22(int argc, char* argv[])
80 :
81 UnsteadyNSTurb(argc, argv),
82 U(_U()),
83 p(_p()),
84 nut(_nut())
85 {}
86
87 // Relevant Fields
88 volVectorField& U;
89 volScalarField& p;
90 volScalarField& nut;
91
92 void offlineSolve(std::string offlinepath)
93 {
94 Vector<double> inl(0, 0, 0);
95 List<scalar> mu_now(1);
96 mu_now[0] = 1e-05;
97
98 if ((offline) && (ITHACAutilities::check_folder(offlinepath) == true))
99 {
100 ITHACAstream::read_fields(Ufield, U, offlinepath);
101 ITHACAstream::read_fields(Pfield, p, offlinepath);
103 }
104 else
105 {
106 truthSolve(mu_now, offlinepath);
107 }
108 }
109 Eigen::MatrixXd vectorTensorMult(Eigen::VectorXd g, Eigen::Tensor<double, 3> c,
110 Eigen::VectorXd a)
111 {
112 int prodDim = c.dimension(0);
113 Eigen::MatrixXd prod;
114 prod.resize(prodDim, 1);
115
116 for (int i = 0; i < prodDim; i++)
117 {
118 prod(i, 0) = g.transpose() *
119 Eigen::SliceFromTensor(c, 0, i) * a;
120 }
121
122 return prod;
123 }
124};
125
126
127void supremizer_approach(tutorial22& example);
128void poisson_approach(tutorial22& example);
129
130/*---------------------------------------------------------------------------*\
131 Starting the MAIN
132\*---------------------------------------------------------------------------*/
133
134int main(int argc, char* argv[])
135{
136 if (argc == 1)
137 {
138 std::cout << "Pass 'supremizer' or 'poisson' as first arguments."
139 << std::endl;
140 exit(0);
141 }
142
143 // process arguments removing "offline" or "online" keywords
144 int argc_proc = argc - 1;
145 char* argv_proc[argc_proc];
146 argv_proc[0] = argv[0];
147
148 if (argc > 2)
149 {
150 std::copy(argv + 2, argv + argc, argv_proc + 1);
151 }
152
153 argc--;
154 tutorial22 example(argc, argv);
155
156 if (std::strcmp(argv[1], "supremizer") == 0)
157 {
158 // perform the offline stage, extracting the modes from the snapshots' dataset corresponding to parOffline
159 supremizer_approach(example);
160 }
161 else if (std::strcmp(argv[1], "poisson") == 0)
162 {
163 poisson_approach(example);
164 }
165 else
166 {
167 std::cout << "Pass supremizer, poisson" << std::endl;
168 }
169
170 exit(0);
171}
172
174
175{
176 word filename("./par");
177 Eigen::VectorXd par;
178 example.inletIndex.resize(1, 2);
179 example.inletIndex << 0, 0;
180 example.inletIndexT.resize(1, 1);
181 example.inletIndexT << 1;
183 example._runTime());
184 int NmodesU = para->ITHACAdict->lookupOrDefault<int>("NmodesU", 5);
185 int NmodesP = para->ITHACAdict->lookupOrDefault<int>("NmodesP", 5);
186 int NmodesSUP = para->ITHACAdict->lookupOrDefault<int>("NmodesSUP", 5);
187 int NmodesNUT = para->ITHACAdict->lookupOrDefault<int>("NmodesNUT", 5);
188 int NmodesProject = para->ITHACAdict->lookupOrDefault<int>("NmodesProject", 5);
189 int NmodesMatrixRec = para->ITHACAdict->lookupOrDefault<int>("NmodesMatrixRec",
190 5);
191 double penaltyFactor =
192 para->ITHACAdict->lookupOrDefault<double>("penaltyFactor", 5);
193 double U_BC = para->ITHACAdict->lookupOrDefault<double>("U_BC", 0.001);
194 double romStartTime = para->ITHACAdict->lookupOrDefault<double>("romStartTime",
195 0);
196 double romEndTime = para->ITHACAdict->lookupOrDefault<double>("romEndTime", 3);
197 double romTimeStep = para->ITHACAdict->lookupOrDefault<double>("romTimeStep",
198 0.001);
199 double e = para->ITHACAdict->lookupOrDefault<double>("RBFradius", 1);
200 example.startTime = 79.992;
201 example.finalTime = 99.996;
202 example.timeStep = 0.0002;
203 example.writeEvery = 0.004;
204 example.offlineSolve("./ITHACAoutput/Offline/");
205 example.solvesupremizer();
206 ITHACAPOD::getModes(example.nutFields, example.nutModes, example._nut().name(),
207 example.podex, 0, 0, NmodesProject);
208 ITHACAPOD::getModes(example.Ufield, example.Umodes, example._U().name(),
209 example.podex, 0, 0, NmodesProject);
210 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
211 example.podex, 0, 0, NmodesProject);
212 ITHACAPOD::getModes(example.supfield, example.supmodes, example._U().name(),
213 example.podex,
214 example.supex, 1, NmodesProject);
215 example.projectSUP("./Matrices", NmodesU, NmodesP, NmodesSUP, NmodesNUT);
216 Eigen::MatrixXd coeefs = ITHACAutilities::getCoeffs(example.Ufield,
217 example.L_U_SUPmodes);
218 Eigen::MatrixXd coeefsNut = ITHACAutilities::getCoeffs(example.nutFields,
219 example.nutModes);
220 Eigen::MatrixXd coeefsP = ITHACAutilities::getCoeffs(example.Pfield,
221 example.Pmodes);
222 cnpy::save(coeefs, "./ITHACAoutput/Matrices/coeefs.npy");
223 cnpy::save(coeefsNut, "./ITHACAoutput/Matrices/coeefsNut.npy");
224 cnpy::save(coeefsP, "./ITHACAoutput/Matrices/coeefsP.npy");
225}
226
227
229{
230 word filename("./par");
231 Eigen::VectorXd par;
232 example.inletIndex.resize(1, 2);
233 example.inletIndex << 0, 0;
234 example.inletIndexT.resize(1, 1);
235 example.inletIndexT << 1;
237 example._runTime());
238 int NmodesU = para->ITHACAdict->lookupOrDefault<int>("NmodesU", 5);
239 int NmodesP = para->ITHACAdict->lookupOrDefault<int>("NmodesP", 5);
240 int NmodesSUP = para->ITHACAdict->lookupOrDefault<int>("NmodesSUP", 5);
241 int NmodesNUT = para->ITHACAdict->lookupOrDefault<int>("NmodesNUT", 5);
242 int NmodesProject = para->ITHACAdict->lookupOrDefault<int>("NmodesProject", 5);
243 int NmodesMatrixRec = para->ITHACAdict->lookupOrDefault<int>("NmodesMatrixRec",
244 5);
245 double penaltyFactor =
246 para->ITHACAdict->lookupOrDefault<double>("penaltyFactor", 5);
247 double U_BC = para->ITHACAdict->lookupOrDefault<double>("U_BC", 0.001);
248 double romStartTime = para->ITHACAdict->lookupOrDefault<double>("romStartTime",
249 0);
250 double romEndTime = para->ITHACAdict->lookupOrDefault<double>("romEndTime", 3);
251 double romTimeStep = para->ITHACAdict->lookupOrDefault<double>("romTimeStep",
252 0.001);
253 double e = para->ITHACAdict->lookupOrDefault<double>("RBFradius", 1);
254 example.startTime = 79.992;
255 example.finalTime = 99.996;
256 example.timeStep = 0.0002;
257 example.writeEvery = 0.004;
258 example.offlineSolve("./ITHACAoutput/Offline");
259 example.solvesupremizer();
260 //modNut = ITHACAPOD::getModes(example.nutFields, example.nutModes, example._nut().name(),
261 // example.podex, 0, 0, NmodesProject);
262 ITHACAPOD::getModes(example.Ufield, example.Umodes, example._U().name(),
263 example.podex, 0, 0, NmodesProject);
264 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
265 example.podex, 0, 0, NmodesProject);
266 ITHACAPOD::getModes(example.supfield, example.supmodes, example._U().name(),
267 example.podex,
268 example.supex, 1, NmodesProject);
269 example.projectPPE("./Matrices", NmodesU, NmodesP, NmodesSUP, NmodesNUT);
270 Eigen::MatrixXd coeefs = ITHACAutilities::getCoeffs(example.Ufield,
271 example.L_U_SUPmodes);
272 Eigen::MatrixXd coeefsNut = ITHACAutilities::getCoeffs(example.nutFields,
273 example.nutModes);
274 Eigen::MatrixXd coeefsP = ITHACAutilities::getCoeffs(example.Pfield,
275 example.Pmodes);
276 cnpy::save(coeefs, "./ITHACAoutput/Matrices/coeefs.npy");
277 cnpy::save(coeefsNut, "./ITHACAoutput/Matrices/coeefsNut.npy");
278 cnpy::save(coeefsP, "./ITHACAoutput/Matrices/coeefsP.npy");
279}
280
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 reducedSteadyNS class.
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.
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.
double e
RBF functions radius.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes, bool rbfInterp=true)
Project using the Poisson Equation for pressure.
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::MatrixXi inletIndexT
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
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: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
volVectorModes L_U_SUPmodes
List of pointers containing the total number of lift, supremizer and velocity modes.
Definition steadyNS.H:116
PtrList< volVectorField > supfield
List of pointers used to form the supremizer snapshots matrix.
Definition steadyNS.H:92
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
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 & nut
Definition offline.C:90
void offlineSolve(std::string offlinepath)
Definition offline.C:92
volScalarField & p
Definition offline.C:89
tutorial22(int argc, char *argv[])
Definition offline.C:79
Eigen::MatrixXd vectorTensorMult(Eigen::VectorXd g, Eigen::Tensor< double, 3 > c, Eigen::VectorXd a)
Definition offline.C:109
volVectorField & U
Definition offline.C:88
Matrix< VectorType, Dynamic, Dynamic > SliceFromTensor(Eigen::Tensor< VectorType, 3 > &tensor, label dim, label index1)
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 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.
Eigen::VectorXd getCoeffs(GeometricField< Type, PatchField, GeoMesh > &snapshot, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Projects a snapshot on a basis function and gets the coefficients of the projection.
bool check_folder(word folder)
Checks if a folder exists.
void save(const Eigen::Matrix< T, -1, dim > &mat, const std::string fname)
void supremizer_approach(tutorial22 &example)
Definition offline.C:173
int main(int argc, char *argv[])
Definition offline.C:134
void poisson_approach(tutorial22 &example)
Definition offline.C:228
Header file of the reductionProblem class.
dimensionedVector & g
label i
Definition pEqn.H:46
Header file of the steadyNS class.