Loading...
Searching...
No Matches
ReducedUnsteadyNSTurb.H
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/>.
24Class
25 ReducedUnsteadyNSTurb
26Description
27 A reduced problem for the unsteady turbulent NS equations
28SourceFiles
29 ReducedUnsteadyNSTurb.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ReducedUnsteadyNSTurb_H
38#define ReducedUnsteadyNSTurb_H
39
40#include "fvCFD.H"
41#include "IOmanip.H"
42#include "ReducedSteadyNS.H"
43#include "ReducedUnsteadyNS.H"
44#include "UnsteadyNSTurb.H"
45#include <Eigen/Dense>
46#include <unsupported/Eigen/NonLinearOptimization>
47#include <unsupported/Eigen/NumericalDiff>
48
49
51{
52 public:
54
56 UnsteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
58 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
59 nphiNut(problem.nNutModes),
60 Nphi_p(problem.NPmodes),
61 N_BC(problem.inletIndex.rows()),
62 gNut(problem.nNutModes)
63 {}
64
65 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
66 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
67
69 int Nphi_u;
71 int Nphi_p;
72 int N_BC;
73 scalar nu;
74 scalar dt;
75 Eigen::VectorXd y_old;
76 Eigen::VectorXd yOldOld;
77 Eigen::VectorXd bc;
78 Eigen::MatrixXd tauU;
79 Eigen::VectorXd gNut;
80 std::vector<SPLINTER::RBFSpline*> SPLINES;
81};
82
83
84
86{
87 public:
89
91 UnsteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
93 Nphi_u(problem.NUmodes + problem.liftfield.size()),
94 nphiNut(problem.nNutModes),
95 Nphi_p(problem.NPmodes),
96 N_BC(problem.inletIndex.rows()),
97 gNut(problem.nNutModes)
98 {}
99
100 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
101 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
102
107 int N_BC;
108 scalar nu;
109 scalar dt;
110 Eigen::VectorXd y_old;
111 Eigen::VectorXd yOldOld;
112 Eigen::VectorXd bc;
113 Eigen::MatrixXd tauU;
114 Eigen::VectorXd gNut;
115 std::vector<SPLINTER::RBFSpline*> SPLINES;
116};
117
119{
120 public:
122
124 UnsteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
125 problem(& problem),
126 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
127 nphiNut(problem.nNutModes),
128 Nphi_p(problem.NPmodes),
129 N_BC(problem.inletIndex.rows()),
130 gNut(problem.nNutModes),
131 gNutAve(problem.nutAve.size())
132 {}
133
134 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
135 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
136
141 int N_BC;
142 scalar nu;
143 scalar dt;
144 Eigen::VectorXd y_old;
145 Eigen::VectorXd yOldOld;
146 Eigen::VectorXd bc;
147 Eigen::MatrixXd tauU;
148 Eigen::VectorXd gNut;
149 Eigen::VectorXd gNutAve;
150 std::vector<SPLINTER::RBFSpline*> SPLINES;
151};
152
154{
155 public:
157
159 UnsteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
160 problem(& problem),
161 Nphi_u(problem.NUmodes + problem.liftfield.size()),
162 nphiNut(problem.nNutModes),
163 Nphi_p(problem.NPmodes),
164 N_BC(problem.inletIndex.rows()),
165 gNut(problem.nNutModes),
166 gNutAve(problem.nutAve.size())
167
168 {}
169
170 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
171 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
172
177 int N_BC;
178 scalar nu;
179 scalar dt;
180 Eigen::VectorXd y_old;
181 Eigen::VectorXd yOldOld;
182 Eigen::VectorXd bc;
183 Eigen::MatrixXd tauU;
184 Eigen::VectorXd gNut;
185 Eigen::VectorXd gNutAve;
186 std::vector<SPLINTER::RBFSpline*> SPLINES;
187};
188
189
190
191/*---------------------------------------------------------------------------*\
192 Class reducedProblem Declaration
193\*---------------------------------------------------------------------------*/
194
196
199{
200 private:
201
202 public:
203 // Constructors
211
213
216
219
221 int dimA;
222
224 PtrList<volScalarField> nutModes;
225
227 Eigen::MatrixXd rbfCoeffMat;
228
230 Eigen::MatrixXd initCond;
231
233 Eigen::VectorXd nut0;
234
237
240
243
246
248 Eigen::VectorXd muStar;
249
251 Eigen::VectorXd gNutAve;
252
254 int interChoice = 1;
255
257 bool skipLift = false;
258
260 PtrList<volScalarField> nutRecFields;
261
262 // Functions
263
271 void solveOnlinePPE(Eigen::MatrixXd velNow);
272
280 void solveOnlineSUP(Eigen::MatrixXd velNow);
281
288 void solveOnlineSUPAve(Eigen::MatrixXd velNow);
289
296 void solveOnlinePPEAve(Eigen::MatrixXd velNow);
297
304 void reconstruct(bool exportFields = false,
305 fileName folder = "./online_rec");
306
314 Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel);
315
316};
317
318
319// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320
321
322
323#endif
324
325
326
327
328
329
Header file of the reducedSteadyNS class.
Header file of the reducedUnsteadyNS class.
Header file of the UnsteadyNSTurb class.
void solveOnlinePPEAve(Eigen::MatrixXd velNow)
Method to perform an online solve using a PPE stabilisation method with the use of the average splitt...
ReducedUnsteadyNSTurb()
Construct Null.
int nphiNut
Number of viscosity modes.
newtonUnsteadyNSTurbSUP newtonObjectSUP
Function object to call the non linear solver sup approach.
newtonUnsteadyNSTurbSUPAve newtonObjectSUPAve
Function object to call the non linear solver sup approach with the splitting of the eddy viscosity.
PtrList< volScalarField > nutRecFields
Reconstructed eddy viscosity fields list.
void solveOnlineSUPAve(Eigen::MatrixXd velNow)
Method to perform an online solve using a supremizer stabilisation method with the use of the average...
Eigen::MatrixXd rbfCoeffMat
The matrix of the eddy viscosity RBF interoplated coefficients.
bool skipLift
Interpolation boolean variable to skip lifting functions.
int dimA
Dimension of the interpolation independent variable.
Eigen::VectorXd muStar
Online parameter value.
Eigen::MatrixXd initCond
The matrix of the initial velocity and pressure reduced coefficients.
Eigen::VectorXd gNutAve
The reduced average vector for viscosity coefficients.
void solveOnlineSUP(Eigen::MatrixXd velNow)
Method to perform an online solve using a supremizer stabilisation method.
int interChoice
Interpolation independent variable choice.
newtonUnsteadyNSTurbPPE newtonObjectPPE
Function object to call the non linear solver PPE approach.
void reconstruct(bool exportFields=false, fileName folder="./online_rec")
Method to reconstruct the solutions from an online solve with any of the two techniques SUP or the PP...
PtrList< volScalarField > nutModes
List of pointers to store the modes for the eddy viscosity.
void solveOnlinePPE(Eigen::MatrixXd velNow)
Method to perform an online solve using a PPE stabilisation method.
UnsteadyNSTurb * problem
Pointer to the FOM problem.
newtonUnsteadyNSTurbPPEAve newtonObjectPPEAve
Function object to call the non linear solver PPE approach with the splitting of the eddy viscosity.
Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel)
Sets the online velocity.
Eigen::VectorXd nut0
The initial eddy viscosity reduced coefficients.
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
reducedUnsteadyNS()
Construct Null.
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newtonUnsteadyNSTurbPPEAve(int Nx, int Ny, UnsteadyNSTurb &problem)
std::vector< SPLINTER::RBFSpline * > SPLINES
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
newtonUnsteadyNSTurbPPE(int Nx, int Ny, UnsteadyNSTurb &problem)
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
std::vector< SPLINTER::RBFSpline * > SPLINES
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
std::vector< SPLINTER::RBFSpline * > SPLINES
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
newtonUnsteadyNSTurbSUPAve(int Nx, int Ny, UnsteadyNSTurb &problem)
std::vector< SPLINTER::RBFSpline * > SPLINES
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newtonUnsteadyNSTurbSUP(int Nx, int Ny, UnsteadyNSTurb &problem)
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const