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
54 newton_unsteadyBB_sup(int Nx, int Ny,
55 UnsteadyBB& problem): newton_argument<double>(Nx, Ny),
57 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
58 Nphi_t(problem.NTmodes + problem.liftfieldT.size()),
59 N_BC_t(problem.inletIndexT.rows()),
60 N_BC(problem.inletIndex.rows()),
61 Nphi_prgh(problem.NPrghmodes)
62 {}
63
64 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
65 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
67 int Nphi_u;
68 int Nphi_t;
69 int N_BC_t;
70 int N_BC;
72 scalar nu;
73 scalar Pr;
74 scalar dt;
75
76 Eigen::VectorXd y_old;
77 Eigen::VectorXd BC_t;
78 Eigen::VectorXd BC;
79
80};
81
83{
84 public:
86
87 newton_unsteadyBB_PPE(int Nx, int Ny,
88 UnsteadyBB& problem): newton_argument<double>(Nx, Ny),
90 Nphi_u(problem.NUmodes + problem.liftfield.size()),
91 Nphi_t(problem.NTmodes + problem.liftfieldT.size()),
92 N_BC_t(problem.inletIndexT.rows()),
93 N_BC(problem.inletIndex.rows()),
94 Nphi_p(problem.NPmodes),
95 Nphi_prgh(problem.NPrghmodes)
96 {}
97
98 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
99 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
100
107 int N_BC;
108 scalar nu;
109 scalar Pr;
110 scalar dt;
113
114 Eigen::VectorXd y_old;
115 Eigen::VectorXd BC_t;
116 Eigen::VectorXd BC;
117};
118
119
120
121/*---------------------------------------------------------------------------*\
122 Class reducedProblem Declaration
123\*---------------------------------------------------------------------------*/
124
126
129{
130 private:
131
132 public:
133 // Constructors
142
144 Eigen::MatrixXd online_solutiont;
145
147 PtrList<volScalarField> LTmodes;
148
150 PtrList<volVectorField> LUmodes;
151
154
157
159 PtrList<volScalarField> TREC;
160
161 PtrList<volScalarField> PREC;
162
165
168
171
173 scalar Pr;
174
177
178 // Functions
179
192 Eigen::MatrixXd solveOnline_sup(Eigen::MatrixXd& temp_now_BC,
193 Eigen::MatrixXd& vel_now_BC, int NParaSet = 0, int startSnap = 0);
194
195
208 Eigen::MatrixXd solveOnline_PPE(Eigen::MatrixXd& temp_now_BC,
209 Eigen::MatrixXd& vel_now_BC, int NParaSet = 0,
210 int startSnap = 0);
211
218 void reconstruct_sup(fileName folder = "./ITHACAOutput/online_rec",
219 int printevery = 1);
220};
221
222
223// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224
225
226#endif
227
228
229
Header file of the reducedUnsteadyNS class.
Header file of the UnsteadyBB class.
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
reducedUnsteadyNS()
Construct Null.
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