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 newtonSteadyNSTurbSUP(int Nx, int Ny,
56 SteadyNSTurb& 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 Eigen::VectorXd gNut;
74 scalar nu;
75 Eigen::MatrixXd tauU;
76 Eigen::VectorXd bc;
77 std::vector<SPLINTER::RBFSpline*> SPLINES;
78};
79
81{
82 public:
84 newtonSteadyNSTurbPPE(int Nx, int Ny,
85 SteadyNSTurb& problem): newton_argument<double>(Nx, Ny),
87 Nphi_u(problem.NUmodes + problem.liftfield.size()),
88 nphiNut(problem.nNutModes),
89 Nphi_p(problem.NPmodes),
90 N_BC(problem.inletIndex.rows()),
91 gNut(problem.nNutModes)
92 {}
93
94 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
95 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
96
98 int Nphi_u;
101 int N_BC;
102 Eigen::VectorXd gNut;
103 scalar nu;
104 Eigen::MatrixXd tauU;
105 Eigen::VectorXd bc;
106 std::vector<SPLINTER::RBFSpline*> SPLINES;
107};
108
109
110/*---------------------------------------------------------------------------*\
111 Class reducedProblem Declaration
112\*---------------------------------------------------------------------------*/
113
114
116
119{
120 private:
121
122 public:
123 // Constructors
131
133
134
136 PtrList<volScalarField> nutFields;
137
139 PtrList<volScalarField> nutModes;
140
142 Eigen::MatrixXd rbfCoeffMat;
143
145 Eigen::VectorXd rbfCoeff;
146
149
152
155
158
160 PtrList<volScalarField> nutRecFields;
161
170 void solveOnlineSUP(Eigen::MatrixXd velNow);
171
178 void solveOnlinePPE(Eigen::MatrixXd velNow);
179
186 void reconstruct(bool exportFields = false,
187 fileName folder = "./ITHACAoutput/online_rec",
188 int printevery = 1);
189
197 Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel);
198
199
200};
201
202// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203
204
205
206#endif
207
208
209
210
211
212
Header file of the ITHACAutilities namespace.
Header file of the reducedProblem class.
Header file of the reducedSteadyNS class.
Header file of the SteadyNSTurb class.
Class where it is implemented a reduced problem for the steady turbulent Navier-stokes problem.
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 ...
Template object created to solve non linear problems using the Eigen library.
Class where it is implemented a reduced problem for the steady Navier-stokes problem.
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)