Loading...
Searching...
No Matches
04unsteadyNS.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 "ITHACAPOD.H"
32#include "ReducedUnsteadyNS.H"
33#include "ITHACAstream.H"
34#include <chrono>
35#include <math.h>
36#include <iomanip>
37#include <iostream>
38
39class tutorial04 : public unsteadyNS
40{
41public:
42 explicit tutorial04(int argc, char* argv[])
43 : unsteadyNS(argc, argv), U(_U()), p(_p()) {}
44
45 // Fields To Perform
46 volVectorField& U;
47 volScalarField& p;
48
50 {
51 Vector<double> inl(1, 0, 0);
52 List<scalar> mu_now(1);
53 label BCind = 0;
54
55 if (offline)
56 {
57 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
58 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
59 }
60 else
61 {
62 for (label i = 0; i < mu.cols(); i++)
63 {
64 mu_now[0] = mu(0, i);
65 assignIF(U, inl);
67 truthSolve(mu_now);
68 }
69 }
70 }
71};
72
73/*---------------------------------------------------------------------------*\
74 Starting the MAIN
75\*---------------------------------------------------------------------------*/
76int main(int argc, char* argv[])
77{
78 // Construct the tutorial04 object
79 tutorial04 example(argc, argv);
80
81 // Read parameters from ITHACAdict file
83 int NmodesUout = para->ITHACAdict->lookupOrDefault<int>("NmodesUout", 15);
84 int NmodesPout = para->ITHACAdict->lookupOrDefault<int>("NmodesPout", 15);
85 int NmodesSUPout = para->ITHACAdict->lookupOrDefault<int>("NmodesSUPout", 15);
86 int NmodesUproj = para->ITHACAdict->lookupOrDefault<int>("NmodesUproj", 10);
87 int NmodesPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesPproj", 10);
88 int NmodesSUPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesSUPproj", 10);
89
90 // Set the number of parameters
91 example.Pnumber = 1;
92 // Set samples
93 example.Tnumber = 1;
94 // Set the parameters infos
95 example.setParameters();
96 // Set the parameter ranges
97 example.mu_range(0, 0) = 0.005;
98 example.mu_range(0, 1) = 0.005;
99 // Generate equispaced samples inside the parameter range
100 example.genEquiPar();
101 // Set the inlet boundaries where we have non homogeneous boundary conditions
102 example.inletIndex.resize(1, 2);
103 example.inletIndex(0, 0) = 0;
104 example.inletIndex(0, 1) = 0;
105 // Time parameters
106 example.startTime = 60;
107 example.finalTime = 70;
108 example.timeStep = 0.01;
109 example.writeEvery = 0.1;
110 // Perform The Offline Solve
111 example.offlineSolve();
112
113 // Check if lift or penalty method should be used
114 if (example.bcMethod == "lift")
115 {
116 // Search the lift function
117 example.liftSolve();
118 // Normalize the lifting function
120 // Create homogeneous basis functions for velocity
121 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
122 // Perform a POD decomposition for velocity and pressure
123 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(), example.podex, 0, 0, NmodesUout);
124 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(), example.podex, 0, 0, NmodesPout);
125 }
126 else
127 {
128 // Perform a POD decomposition for velocity and pressure without lifting function
129 ITHACAPOD::getModes(example.Ufield, example.Umodes, example._U().name(), example.podex, 0, 0, NmodesUout);
130 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(), example.podex, 0, 0, NmodesPout);
131 }
132 // Check the method and perform the appropriate projection
133 if (example.method == "supremizer")
134 {
135 example.solvesupremizer();
136 ITHACAPOD::getModes(example.supfield, example.supmodes, example._U().name(), example.podex, example.supex, 1, NmodesSUPout);
137 example.projectSUP("./Matrices", NmodesUproj, NmodesPproj, NmodesSUPproj);
138 }
139 else if (example.method == "PPE")
140 {
141 example.projectPPE("./Matrices", NmodesUproj, NmodesPproj);
142 }
143
144 reducedUnsteadyNS reduced(example);
145
146 // Set values of the reduced stuff
147 reduced.nu = 0.005;
148 reduced.tstart = 60;
149 reduced.finalTime = 70;
150 reduced.dt = 0.005;
151 reduced.storeEvery = 0.005;
152 reduced.exportEvery = 0.1;
153
154 // Set the online velocity
155 Eigen::MatrixXd vel_now(1, 1);
156 vel_now(0, 0) = 1;
157
158 // If using the penalty method, set tauU
159 if (example.bcMethod == "penalty")
160 {
161 reduced.tauU = Eigen::MatrixXd::Zero(1, 1);
162 reduced.tauU(0, 0) = 1e-1; // Example penalty coefficient
163 }
164
165 // Perform the online solve based on the method
166 if (example.method == "supremizer")
167 {
168 reduced.solveOnline_sup(vel_now, 1);
169 }
170 else if (example.method == "PPE")
171 {
172 reduced.solveOnline_PPE(vel_now, 1);
173 }
174
175 // Reconstruct the solution and export it
176 reduced.reconstruct(true, "./ITHACAoutput/Reconstruction" + example.method + "/");
177
178 exit(0);
179}
180
181
182
183
184
188
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 reducedUnsteadyNS 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.
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.
Eigen::MatrixXd tauU
Penalty Factor.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
double exportEvery
A variable for exporting the fields.
scalar finalTime
Scalar to store the final time if the online simulation.
void solveOnline_PPE(Eigen::MatrixXd vel_now, int startSnap=0)
Method to perform an online solve using a PPE stabilisation method.
void solveOnline_sup(Eigen::MatrixXd vel_now, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
scalar tstart
Scalar to store the initial time if the online simulation.
void reconstruct(bool exportFields=false, fileName folder="./online_rec")
Method to reconstruct the solutions from an online solve with a supremizer stabilisation technique.
double dt
Scalar to store the time increment.
double storeEvery
A variable for storing the reduced coefficients.
label Pnumber
Number of parameters.
void assignIF(T &s, G &value)
Assign internal field condition.
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.
void change_viscosity(double mu)
Function to change the viscosity.
Definition steadyNS.C:2014
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:256
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Poisson Equation for pressure.
Definition steadyNS.C:353
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
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
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
word bcMethod
Boundary Method.
Definition steadyNS.H:317
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
tutorial04(int argc, char *argv[])
void offlineSolve()
volVectorField & U
volScalarField & p
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
Definition unsteadyNS.H:62
word method
Definition unsteadyNS.H:90
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.
void normalizeFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &fields)
Normalize list of Geometric fields.
label i
Definition pEqn.H:46
Header file of the unsteadyNS class.