Loading...
Searching...
No Matches
ReducedUnsteadyBB.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 ReducedUnsteadyBB
26Description
27 A reduced problem for the unsteady Boussinesq equations
28SourceFiles
29 ReducedUnsteadyBB.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ReducedUnsteadyBB_H
38#define ReducedUnsteadyBB_H
39
40#include "fvCFD.H"
41#include "IOmanip.H"
42#include "ReducedUnsteadyNS.H"
43#include "UnsteadyBB.H"
44#include <Eigen/Dense>
45#include <unsupported/Eigen/NonLinearOptimization>
46#include <unsupported/Eigen/NumericalDiff>
47
50{
51 public:
53 newton_unsteadyBB_sup(int Nx, int Ny,
54 UnsteadyBB& problem): newton_argument<double>(Nx, Ny),
56 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
57 Nphi_t(problem.NTmodes + problem.liftfieldT.size()),
58 N_BC_t(problem.inletIndexT.rows()),
59 N_BC(problem.inletIndex.rows()),
60 Nphi_prgh(problem.NPrghmodes)
61 {}
62
63 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
64 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
66 int Nphi_u;
67 int Nphi_t;
68 int N_BC_t;
69 int N_BC;
71 scalar nu;
72 scalar Pr;
73 scalar dt;
74
75 Eigen::VectorXd y_old;
76 Eigen::VectorXd BC_t;
77 Eigen::VectorXd BC;
78
79};
80
82{
83 public:
85 newton_unsteadyBB_PPE(int Nx, int Ny,
86 UnsteadyBB& problem): newton_argument<double>(Nx, Ny),
88 Nphi_u(problem.NUmodes + problem.liftfield.size()),
89 Nphi_t(problem.NTmodes + problem.liftfieldT.size()),
90 N_BC_t(problem.inletIndexT.rows()),
91 N_BC(problem.inletIndex.rows()),
92 Nphi_p(problem.NPmodes),
93 Nphi_prgh(problem.NPrghmodes)
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
105 int N_BC;
106 scalar nu;
107 scalar Pr;
108 scalar dt;
111
112 Eigen::VectorXd y_old;
113 Eigen::VectorXd BC_t;
114 Eigen::VectorXd BC;
115};
116
117
118
119/*---------------------------------------------------------------------------*\
120 Class reducedProblem Declaration
121\*---------------------------------------------------------------------------*/
122
124
127{
128 private:
129
130 public:
131 // Constructors
140
142 Eigen::MatrixXd online_solutiont;
143
145 PtrList<volScalarField> LTmodes;
146
148 PtrList<volVectorField> LUmodes;
149
152
155
157 PtrList<volScalarField> TREC;
158
159 PtrList<volScalarField> PREC;
160
163
166
169
171 scalar Pr;
172
175
176 // Functions
177
190 Eigen::MatrixXd solveOnline_sup(Eigen::MatrixXd& temp_now_BC,
191 Eigen::MatrixXd& vel_now_BC, int NParaSet = 0, int startSnap = 0);
192
193
206 Eigen::MatrixXd solveOnline_PPE(Eigen::MatrixXd& temp_now_BC,
207 Eigen::MatrixXd& vel_now_BC, int NParaSet = 0,
208 int startSnap = 0);
209
216 void reconstruct_sup(fileName folder = "./ITHACAOutput/online_rec",
217 int printevery = 1);
218};
219
220
221// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222
223
224#endif
225
226
227
Header file of the reducedUnsteadyNS class.
Header file of the UnsteadyBB class.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
int Nphi_t
Number of temperature modes.
int N_BC_t
Number of parametrized boundary conditions related to temperature field.
PtrList< volScalarField > PREC
void reconstruct_sup(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique.
newton_unsteadyBB_PPE newton_object_PPE
Function object to call the non linear solver PPE approach.
scalar Pr
DimensionedScalar Pr;.
PtrList< volVectorField > LUmodes
List of pointers to store the modes for velocity including lift field modes.
Eigen::MatrixXd online_solutiont
List of Eigen matrices to store current online solution for temperature equation.
Eigen::MatrixXd solveOnline_sup(Eigen::MatrixXd &temp_now_BC, Eigen::MatrixXd &vel_now_BC, int NParaSet=0, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
PtrList< volScalarField > LTmodes
List of pointers to store the modes for temperature including lift field modes.
int Nphi_prgh
Number of shifted pressure modes.
PtrList< volScalarField > TREC
Reconstructed temperature field.
Eigen::MatrixXd solveOnline_PPE(Eigen::MatrixXd &temp_now_BC, Eigen::MatrixXd &vel_now_BC, int NParaSet=0, int startSnap=0)
Method to perform an online solve using a PPE stabilisation method.
newton_unsteadyBB_sup newton_object_sup
Function object to call the non linear solver sup approach.
ReducedUnsteadyBB()
Construct Null.
UnsteadyBB * problem
Pointer to the FOM problem.
Implementation of a parametrized full order unsteady Boussinesq problem and preparation of the the ...
Definition UnsteadyBB.H:60
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 df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
newton_unsteadyBB_PPE(int Nx, int Ny, UnsteadyBB &problem)
Newton object for the resolution of the reduced problem using a supremizer approach.
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
newton_unsteadyBB_sup(int Nx, int Ny, UnsteadyBB &problem)
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const