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:
57 Nphi_u(problem.nModesOnline),
58 N_BC(problem.inletIndex.rows())
59 {}
60
61 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
62 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
63
65 int Nphi_u;
66 int N_BC;
67 scalar nu;
68 scalar dt;
69 Eigen::VectorXd y_old;
70 Eigen::VectorXd yOldOld;
71 Eigen::VectorXd bc;
72 Eigen::MatrixXd tauU;
73};
74
76{
77 public:
82 Nphi_u(problem.NUmodes),
83 Nphi_p(problem.NPmodes),
84 N_BC(problem.inletIndex.rows())
85 {}
86
87 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
88 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
89
91 int Nphi_u;
92 int Nphi_p;
93 int N_BC;
94 scalar nu;
95 scalar dt;
96 Eigen::VectorXd y_old;
97 Eigen::VectorXd yOldOld;
98 Eigen::VectorXd bc;
99 Eigen::MatrixXd tauU;
100};
101
102/*---------------------------------------------------------------------------*\
103 Class reducedProblem Declaration
104\*---------------------------------------------------------------------------*/
105
107
110{
111 private:
112
113 public:
114 // Constructors
122
124
127
129 PtrList<volScalarField> nutModes;
130
133
136
140 Eigen::MatrixXd initCond;
141
143 PtrList<volScalarField> nutRecFields;
144
145 // Functions
146
147 //--------------------------------------------------------------------------
152 void solveOnline(Eigen::MatrixXd vel);
153
154 //--------------------------------------------------------------------------
161 void solveOnlinePPE(Eigen::MatrixXd vel);
162
163 //--------------------------------------------------------------------------
169 void reconstruct(bool exportFields = false,
170 fileName folder = "./online_rec");
171
172 //--------------------------------------------------------------------------
178 void reconstructPPE(bool exportFields = false,
179 fileName folder = "./online_rec");
180
181 //--------------------------------------------------------------------------
188 Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel);
189
190 //--------------------------------------------------------------------------
197 fileName folder);
198
199};
200
201
202// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203
204
205
206#endif
207
208
209
210
211
212
Header file of the reducedSteadyNS class.
Header file of the reducedUnsteadyNS class.
Header file of the UnsteadyNSTurbIntrusive class.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
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 ...
Template object created to solve non linear problems using the Eigen library.
virtual void solveOnline()
Virtual Method to perform and online Solve.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
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