Loading...
Searching...
No Matches
ReducedUnsteadyNSTurbIntrusive.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 ReducedUnsteadyNSTurbIntrusive
26Description
27 A reduced problem for the unsteady turbulent NS equations using fully intrusive approach
28SourceFiles
29 ReducedUnsteadyNSTurbIntrusive.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ReducedUnsteadyNSTurbIntrusive_H
38#define ReducedUnsteadyNSTurbIntrusive_H
39
40#include "fvCFD.H"
41#include "IOmanip.H"
42#include "ReducedSteadyNS.H"
43#include "ReducedUnsteadyNS.H"
45#include <Eigen/Dense>
46#include <unsupported/Eigen/NonLinearOptimization>
47#include <unsupported/Eigen/NumericalDiff>
48
49
51{
52 public:
54
58 Nphi_u(problem.nModesOnline),
59 N_BC(problem.inletIndex.rows())
60 {}
61
62 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
63 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
64
66 int Nphi_u;
67 int N_BC;
68 scalar nu;
69 scalar dt;
70 Eigen::VectorXd y_old;
71 Eigen::VectorXd yOldOld;
72 Eigen::VectorXd bc;
73 Eigen::MatrixXd tauU;
74};
75
77{
78 public:
80
84 Nphi_u(problem.NUmodes),
85 Nphi_p(problem.NPmodes),
86 N_BC(problem.inletIndex.rows())
87 {}
88
89 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
90 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
91
93 int Nphi_u;
94 int Nphi_p;
95 int N_BC;
96 scalar nu;
97 scalar dt;
98 Eigen::VectorXd y_old;
99 Eigen::VectorXd yOldOld;
100 Eigen::VectorXd bc;
101 Eigen::MatrixXd tauU;
102};
103
104/*---------------------------------------------------------------------------*\
105 Class reducedProblem Declaration
106\*---------------------------------------------------------------------------*/
107
109
112{
113 private:
114
115 public:
116 // Constructors
124
126
129
131 PtrList<volScalarField> nutModes;
132
135
138
142 Eigen::MatrixXd initCond;
143
145 PtrList<volScalarField> nutRecFields;
146
147 // Functions
148
149 //--------------------------------------------------------------------------
154 void solveOnline(Eigen::MatrixXd vel);
155
156 //--------------------------------------------------------------------------
163 void solveOnlinePPE(Eigen::MatrixXd vel);
164
165 //--------------------------------------------------------------------------
171 void reconstruct(bool exportFields = false,
172 fileName folder = "./online_rec");
173
174 //--------------------------------------------------------------------------
180 void reconstructPPE(bool exportFields = false,
181 fileName folder = "./online_rec");
182
183 //--------------------------------------------------------------------------
190 Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel);
191
192 //--------------------------------------------------------------------------
199 fileName folder);
200
201};
202
203
204// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205
206
207
208#endif
209
210
211
212
213
214
Header file of the reducedSteadyNS class.
Header file of the reducedUnsteadyNS class.
Header file of the UnsteadyNSTurbIntrusive class.
void reconstructLiftAndDrag(UnsteadyNSTurbIntrusive &problem, fileName folder)
Method to compute the reduced order forces for same number of modes of velocity and pressure.
Eigen::MatrixXd initCond
Tha matrix containing the initial conditions for the reduced variables, in case of PPE approach it mu...
PtrList< volScalarField > nutRecFields
Reconstructed eddy viscosity fields list.
newtonUnsteadyNSTurbIntrusive newtonObject
Function object to call the non linear solver sup approach.
PtrList< volScalarField > nutModes
List of pointers to store the modes for the eddy viscosity.
UnsteadyNSTurbIntrusive * problem
Pointer to the FOM problem.
Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel)
Sets the online velocity.
void solveOnlinePPE(Eigen::MatrixXd vel)
Method to perform an online solve using an intrusive approach with the usage of PPE.
void reconstructPPE(bool exportFields=false, fileName folder="./online_rec")
Method to reconstruct the solutions from an online solve with a PPE semi-uniform approach.
void reconstruct(bool exportFields=false, fileName folder="./online_rec")
Method to reconstruct the solutions from an online solve with a fully uniform approach.
newtonUnsteadyNSTurbIntrusivePPE newtonObjectPPE
Function object to call the non linear solver sup approach.
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
virtual void solveOnline()
Virtual Method to perform and online Solve.
reducedUnsteadyNS()
Construct Null.
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newtonUnsteadyNSTurbIntrusivePPE(int Nx, int Ny, UnsteadyNSTurbIntrusive &problem)
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newtonUnsteadyNSTurbIntrusive(int Nx, int Ny, UnsteadyNSTurbIntrusive &problem)
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const