Loading...
Searching...
No Matches
ReducedSteadyNSTurb.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 ReducedSteadyNSTurb
26Description
27 A reduced problem for the stationary turbulent NS equations
28SourceFiles
29 ReducedSteadyNSTurb.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ReducedSteadyNSTurb_H
38#define ReducedSteadyNSTurb_H
39
40#include "fvCFD.H"
41#include "IOmanip.H"
42#include "ReducedProblem.H"
43#include "ReducedSteadyNS.H"
44#include "SteadyNSTurb.H"
45#include "ITHACAutilities.H"
46#include <Eigen/Dense>
47#include <unsupported/Eigen/NonLinearOptimization>
48#include <unsupported/Eigen/NumericalDiff>
49
52{
53 public:
55
56 newtonSteadyNSTurbSUP(int Nx, int Ny,
57 SteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
59 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
60 nphiNut(problem.nNutModes),
61 Nphi_p(problem.NPmodes),
62 N_BC(problem.inletIndex.rows()),
63 gNut(problem.nNutModes)
64 {}
65
66 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
67 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
68
70 int Nphi_u;
72 int Nphi_p;
73 int N_BC;
74 Eigen::VectorXd gNut;
75 scalar nu;
76 Eigen::MatrixXd tauU;
77 Eigen::VectorXd bc;
78 std::vector<SPLINTER::RBFSpline*> SPLINES;
79};
80
82{
83 public:
85
86 newtonSteadyNSTurbPPE(int Nx, int Ny,
87 SteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
89 Nphi_u(problem.NUmodes + problem.liftfield.size()),
90 nphiNut(problem.nNutModes),
91 Nphi_p(problem.NPmodes),
92 N_BC(problem.inletIndex.rows()),
93 gNut(problem.nNutModes)
94 {}
95
96 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
97 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
98
103 int N_BC;
104 Eigen::VectorXd gNut;
105 scalar nu;
106 Eigen::MatrixXd tauU;
107 Eigen::VectorXd bc;
108 std::vector<SPLINTER::RBFSpline*> SPLINES;
109};
110
111
112/*---------------------------------------------------------------------------*\
113 Class reducedProblem Declaration
114\*---------------------------------------------------------------------------*/
115
116
118
121{
122 private:
123
124 public:
125 // Constructors
133
135
136
138 PtrList<volScalarField> nutFields;
139
141 PtrList<volScalarField> nutModes;
142
144 Eigen::MatrixXd rbfCoeffMat;
145
147 Eigen::VectorXd rbfCoeff;
148
151
154
157
160
162 PtrList<volScalarField> nutRecFields;
163
172 void solveOnlineSUP(Eigen::MatrixXd velNow);
173
180 void solveOnlinePPE(Eigen::MatrixXd velNow);
181
188 void reconstruct(bool exportFields = false,
189 fileName folder = "./ITHACAoutput/online_rec",
190 int printevery = 1);
191
199 Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel);
200
201
202};
203
204// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205
206
207
208#endif
209
210
211
212
213
214
Header file of the ITHACAutilities namespace.
Header file of the reducedProblem class.
Header file of the reducedSteadyNS class.
Header file of the SteadyNSTurb class.
PtrList< volScalarField > nutRecFields
Reconstructed eddy viscosity fields list.
newtonSteadyNSTurbSUP newtonObjectSUP
Newton Object to solve the nonlinear problem sup approach.
newtonSteadyNSTurbPPE newtonObjectPPE
Newton Object to solve the nonlinear problem PPE approach.
void solveOnlinePPE(Eigen::MatrixXd velNow)
Method to perform an online solve using a PPE stabilisation method.
ReducedSteadyNSTurb()
Construct Null.
int nphiNut
Number of viscosity modes.
Eigen::MatrixXd rbfCoeffMat
The matrix of the eddy viscosity RBF interoplated coefficients.
void solveOnlineSUP(Eigen::MatrixXd velNow)
Method to perform an online solve using a supremizer stabilisation method.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
SteadyNSTurb * problem
Pointer to the FOM problem.
void reconstruct(bool exportFields=false, fileName folder="./ITHACAoutput/online_rec", int printevery=1)
Method to reconstruct the solutions from an online solve.
Eigen::VectorXd rbfCoeff
Vector of eddy viscosity RBF interoplated coefficients.
PtrList< volScalarField > nutModes
List of POD modes for eddy viscosity.
Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel)
Sets the online velocity.
Implementation of a parametrized full order steady turbulent Navier Stokes problem and preparation ...
reducedSteadyNS()
Construct Null.
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
newtonSteadyNSTurbPPE(int Nx, int Ny, SteadyNSTurb &problem)
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
std::vector< SPLINTER::RBFSpline * > SPLINES
Structure to implement a newton object for a stationary NS problem.
std::vector< SPLINTER::RBFSpline * > SPLINES
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newtonSteadyNSTurbSUP(int Nx, int Ny, SteadyNSTurb &problem)