Loading...
Searching...
No Matches
19UnsteadyNSExplicit.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) 2021 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 with the Discretize-Then-
26 Project approach for the lid driven cavity problem
27SourceFiles
28 19UnsteadyNSExplicit.C
29\*---------------------------------------------------------------------------*/
30
31#include "UnsteadyNSExplicit.H"
32#include "ITHACAPOD.H"
34#include "ITHACAstream.H"
35#include <chrono>
36#include<math.h>
37#include<iomanip>
38
40{
41 public:
42 explicit tutorial19(int argc, char* argv[])
43 :
44 UnsteadyNSExplicit(argc, argv),
45 U(_U()),
46 p(_p()),
47 phi(_phi())
48 {}
49
50 // Fields To Perform
51 volVectorField& U;
52 volScalarField& p;
53 surfaceScalarField& phi;
54
56 {
57 Vector<double> inl(1, 0, 0);
58 List<scalar> mu_now(1);
59
60 if (offline)
61 {
62 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
63 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
64
65 if (fluxMethod == "consistent")
66 {
67 ITHACAstream::read_fields(Phifield, phi, "./ITHACAoutput/Offline/");
68 }
69 }
70 else
71 {
72 for (label i = 0; i < mu.cols(); i++)
73 {
74 mu_now[0] = mu(0, i);
75 truthSolve(mu_now);
76 }
77 }
78 }
79};
80
81/*---------------------------------------------------------------------------*\
82 Starting the MAIN
83\*---------------------------------------------------------------------------*/
84int main(int argc, char* argv[])
85{
86 // Construct the tutorial19 object
87 tutorial19 example(argc, argv);
88 // Read parameters from ITHACAdict file
90 example._runTime());
91 int NmodesUout = para->ITHACAdict->lookupOrDefault<int>("NmodesUout", 15);
92 int NmodesPout = para->ITHACAdict->lookupOrDefault<int>("NmodesPout", 15);
93 int NmodesSUPout = 0;
94 int NmodesUproj = para->ITHACAdict->lookupOrDefault<int>("NmodesUproj", 10);
95 int NmodesPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesPproj", 10);
96 int NmodesSUPproj = 0;
98 example.Pnumber = 1;
100 example.Tnumber = 1;
102 example.setParameters();
103 // Set the parameter ranges
104 example.mu_range(0, 0) = 0.01;
105 example.mu_range(0, 1) = 0.01;
106 // Generate equispaced samples inside the parameter range
107 example.genEquiPar();
108 // Set the inlet boundaries where we have non homogeneous boundary conditions
109 example.inletIndex.resize(1, 2);
110 example.inletIndex(0, 0) = 0;
111 example.inletIndex(0, 1) = 0;
112 // Time parameters
113 example.startTime = 0.0;
114 example.finalTime = 1.0;
115 example.timeStep = 0.005;
116 example.writeEvery = 0.005;
117 // Perform The Offline Solve;
118 example.offlineSolve();
119 // Perform a POD decomposition for velocity and pressure
120 ITHACAPOD::getModes(example.Ufield, example.Umodes, example._U().name(),
121 example.podex, 0, 0,
122 NmodesUout);
123 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
124 example.podex, 0, 0,
125 NmodesPout);
126
127 if (example.fluxMethod == "consistent")
128 {
129 ITHACAPOD::getModes(example.Phifield, example.Phimodes, example._phi().name(),
130 example.podex, 0, 0,
131 NmodesUout);
132 }
133
134 // Galerkin Projection
135 example.discretizeThenProject("./Matrices", NmodesUproj, NmodesPproj,
136 NmodesSUPproj);
137 ReducedUnsteadyNSExplicit reduced(example);
138 // Set values of the reduced order model
139 reduced.nu = 0.01;
140 reduced.tstart = 0.0;
141 reduced.finalTime = 1.0;
142 reduced.dt = 0.005;
143 reduced.storeEvery = 0.005;
144 reduced.exportEvery = 0.005;
145 // Set the online velocity
146 Eigen::MatrixXd vel_now(1, 1);
147 vel_now(0, 0) = 1;
148 reduced.solveOnline(vel_now, 1);
149 // Reconstruct the solution and export it
150 reduced.reconstruct(false, "./ITHACAoutput/Reconstruction/");
151 exit(0);
152}
153
int main(int argc, char *argv[])
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 ReducedUnsteadyNSExplicit class.
Header file of the UnsteadyNSExplicit 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.
void reconstruct(bool exportFields=false, fileName folder="./online_rec")
Method to reconstruct the solutions from an online solve with a supremizer stabilisation technique.
void solveOnline(Eigen::MatrixXd vel_now, label startSnap=0)
Method to perform an online solve without a pressure stabilisation method.
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
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)
scalar nu
Reduced viscosity in case of parametrized viscosity.
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.
label Tnumber
Dimension of the training set (used only when gerating parameters without input)
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.
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
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< surfaceScalarField > Phifield
List of pointers used to form the flux snapshots matrix.
Definition steadyNS.H:95
word fluxMethod
Flux Method.
Definition steadyNS.H:320
void discretizeThenProject(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Discretize-then-project approach.
Definition steadyNS.C:687
surfaceScalarModes Phimodes
List of pointers used to form the flux modes.
Definition steadyNS.H:107
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 & p
volVectorField & U
tutorial19(int argc, char *argv[])
surfaceScalarField & phi
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.
label i
Definition pEqn.H:46