Loading...
Searching...
No Matches
04unsteadyNS.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 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{
41 public:
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
49 void offlineSolve()
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);
66 change_viscosity(mu(0, i));
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 // Read parameters from ITHACAdict file
81 ITHACAparameters* para = ITHACAparameters::getInstance(example._mesh(),
82 example._runTime());
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 // Set the number of parameters
90 example.Pnumber = 1;
91 // Set samples
92 example.Tnumber = 1;
93 // Set the parameters infos
94 example.setParameters();
95 // Set the parameter ranges
96 example.mu_range(0, 0) = 0.005;
97 example.mu_range(0, 1) = 0.005;
98 // Generate equispaced samples inside the parameter range
99 example.genEquiPar();
100 // Set the inlet boundaries where we have non homogeneous boundary conditions
101 example.inletIndex.resize(1, 2);
102 example.inletIndex(0, 0) = 0;
103 example.inletIndex(0, 1) = 0;
104 // Time parameters
105 example.startTime = 60;
106 example.finalTime = 70;
107 example.timeStep = 0.01;
108 example.writeEvery = 0.1;
109 // Perform The Offline Solve
110 example.offlineSolve();
111
112 // Check if lift or penalty method should be used
113 if (example.bcMethod == "lift")
114 {
115 // Search the lift function
116 example.liftSolve();
117 // Normalize the lifting function
118 ITHACAutilities::normalizeFields(example.liftfield);
119 // Create homogeneous basis functions for velocity
120 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
121 // Perform a POD decomposition for velocity and pressure
122 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
123 example.podex, 0, 0, NmodesUout);
124 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
125 example.podex, 0, 0, NmodesPout);
126 }
127 else
128 {
129 // Perform a POD decomposition for velocity and pressure without lifting function
130 ITHACAPOD::getModes(example.Ufield, example.Umodes, example._U().name(),
131 example.podex, 0, 0, NmodesUout);
132 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
133 example.podex, 0, 0, NmodesPout);
134 }
135
136 // Check the method and perform the appropriate projection
137 if (example.method == "supremizer")
138 {
139 example.solvesupremizer();
140 ITHACAPOD::getModes(example.supfield, example.supmodes, example._U().name(),
141 example.podex, example.supex, 1, NmodesSUPout);
142 example.projectSUP("./Matrices", NmodesUproj, NmodesPproj, NmodesSUPproj);
143 }
144 else if (example.method == "PPE")
145 {
146 example.projectPPE("./Matrices", NmodesUproj, NmodesPproj);
147 }
148
149 reducedUnsteadyNS reduced(example);
150 // Set values of the reduced stuff
151 reduced.nu = 0.005;
152 reduced.tstart = 60;
153 reduced.finalTime = 70;
154 reduced.dt = 0.005;
155 reduced.storeEvery = 0.005;
156 reduced.exportEvery = 0.1;
157 // Set the online velocity
158 Eigen::MatrixXd vel_now(1, 1);
159 vel_now(0, 0) = 1;
160
161 // If using the penalty method, set tauU
162 if (example.bcMethod == "penalty")
163 {
164 reduced.tauU = Eigen::MatrixXd::Zero(1, 1);
165 reduced.tauU(0, 0) = 1e-1; // Example penalty coefficient
166 }
167
168 // Perform the online solve based on the method
169 if (example.method == "supremizer")
170 {
171 reduced.solveOnline_sup(vel_now, 1);
172 }
173 else if (example.method == "PPE")
174 {
175 reduced.solveOnline_PPE(vel_now, 1);
176 }
177
178 // Reconstruct the solution and export it
179 reduced.reconstruct(true,
180 "./ITHACAoutput/Reconstruction" + example.method + "/");
181 return 0;
182}
183
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...
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 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 truthSolve()
Perform a TruthSolve.
void change_viscosity(double mu)
Function to change the viscosity.
Definition steadyNS.C:2066
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< volVectorField > _U
Velocity field.
Definition steadyNS.H:272
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:269
unsteadyNS()
Construct Null.
Definition unsteadyNS.C:39
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.
Header file of the unsteadyNS class.