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:
55 UnsteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
57 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
58 nphiNut(problem.nNutModes),
59 Nphi_p(problem.NPmodes),
60 N_BC(problem.inletIndex.rows()),
61 gNut(problem.nNutModes)
62 {}
63
64 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
65 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
66
68 int Nphi_u;
70 int Nphi_p;
71 int N_BC;
72 scalar nu;
73 scalar dt;
74 Eigen::VectorXd y_old;
75 Eigen::VectorXd yOldOld;
76 Eigen::VectorXd bc;
77 Eigen::MatrixXd tauU;
78 Eigen::VectorXd gNut;
79 std::vector<SPLINTER::RBFSpline*> SPLINES;
80};
81
82
83
85{
86 public:
89 UnsteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
91 Nphi_u(problem.NUmodes + problem.liftfield.size()),
92 nphiNut(problem.nNutModes),
93 Nphi_p(problem.NPmodes),
94 N_BC(problem.inletIndex.rows()),
95 gNut(problem.nNutModes)
96 {}
97
98 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
99 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
100
105 int N_BC;
106 scalar nu;
107 scalar dt;
108 Eigen::VectorXd y_old;
109 Eigen::VectorXd yOldOld;
110 Eigen::VectorXd bc;
111 Eigen::MatrixXd tauU;
112 Eigen::VectorXd gNut;
113 std::vector<SPLINTER::RBFSpline*> SPLINES;
114};
115
117{
118 public:
121 UnsteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
123 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
124 nphiNut(problem.nNutModes),
125 Nphi_p(problem.NPmodes),
126 N_BC(problem.inletIndex.rows()),
127 gNut(problem.nNutModes),
128 gNutAve(problem.nutAve.size())
129 {}
130
131 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
132 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
133
138 int N_BC;
139 scalar nu;
140 scalar dt;
141 Eigen::VectorXd y_old;
142 Eigen::VectorXd yOldOld;
143 Eigen::VectorXd bc;
144 Eigen::MatrixXd tauU;
145 Eigen::VectorXd gNut;
146 Eigen::VectorXd gNutAve;
147 std::vector<SPLINTER::RBFSpline*> SPLINES;
148};
149
151{
152 public:
155 UnsteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
157 Nphi_u(problem.NUmodes + problem.liftfield.size()),
158 nphiNut(problem.nNutModes),
159 Nphi_p(problem.NPmodes),
160 N_BC(problem.inletIndex.rows()),
161 gNut(problem.nNutModes),
162 gNutAve(problem.nutAve.size())
163
164 {}
165
166 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
167 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
168
173 int N_BC;
174 scalar nu;
175 scalar dt;
176 Eigen::VectorXd y_old;
177 Eigen::VectorXd yOldOld;
178 Eigen::VectorXd bc;
179 Eigen::MatrixXd tauU;
180 Eigen::VectorXd gNut;
181 Eigen::VectorXd gNutAve;
182 std::vector<SPLINTER::RBFSpline*> SPLINES;
183};
184
185
186
187/*---------------------------------------------------------------------------*\
188 Class reducedProblem Declaration
189\*---------------------------------------------------------------------------*/
190
192
195{
196 private:
197
198 public:
199 // Constructors
207
209
212
215
217 int dimA;
218
220 PtrList<volScalarField> nutModes;
221
223 Eigen::MatrixXd rbfCoeffMat;
224
226 Eigen::MatrixXd initCond;
227
229 Eigen::VectorXd nut0;
230
233
236
239
242
244 Eigen::VectorXd muStar;
245
247 Eigen::VectorXd gNutAve;
248
250 int interChoice = 1;
251
253 bool skipLift = false;
254
256 PtrList<volScalarField> nutRecFields;
257
258 // Functions
259
267 void solveOnlinePPE(Eigen::MatrixXd velNow);
268
276 void solveOnlineSUP(Eigen::MatrixXd velNow);
277
284 void solveOnlineSUPAve(Eigen::MatrixXd velNow);
285
292 void solveOnlinePPEAve(Eigen::MatrixXd velNow);
293
300 void reconstruct(bool exportFields = false,
301 fileName folder = "./online_rec");
302
310 Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel);
311
312};
313
314
315// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316
317
318
319#endif
320
321
322
323
324
325
Header file of the reducedSteadyNS class.
Header file of the reducedUnsteadyNS class.
Header file of the UnsteadyNSTurb class.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
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 ...
Template object created to solve non linear problems using the Eigen library.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
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