Loading...
Searching...
No Matches
ReducedUnsteadyNST.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 reducedUnsteadyNST
26Description
27 A reduced problem for the unsteady NS plus the energy equations
28SourceFiles
29 reducedUnsteadyNST.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ReducedUnsteadyNST_H
38#define ReducedUnsteadyNST_H
39#include "fvCFD.H"
40#include "IOmanip.H"
41#include "ReducedUnsteadyNS.H"
42#include "unsteadyNST.H"
43#include <Eigen/Dense>
44#include <unsupported/Eigen/NonLinearOptimization>
45#include <unsupported/Eigen/NumericalDiff>
46
49{
50 public:
52
53 newton_unsteadyNST_sup(int Nx, int Ny,
54 unsteadyNST& problem): newton_argument<double>(Nx, Ny),
56 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
57 Nphi_p(problem.NPmodes),
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;
64 int Nphi_u;
65 int Nphi_p;
66 int N_BC;
67 scalar nu;
68 scalar dt;
69 Eigen::VectorXd y_old;
70 Eigen::VectorXd BC;
71
72};
73
75{
76 public:
78
80 unsteadyNST& problem): newton_argument<double>(Nx, Ny),
82 Nphi_t(problem.NTmodes + problem.liftfieldT.size()),
83 N_BC_t(problem.inletIndexT.rows())
84 {}
85
86 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
87 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
89 int Nphi_t;
90 int N_BC_t;
91 int Nphi_u;
92 scalar nu;
93 scalar DT;
94 scalar dt;
95 Eigen::VectorXd a_tmp;
96 Eigen::VectorXd z_old;
97 Eigen::VectorXd BC_t;
98
99};
100
101/*---------------------------------------------------------------------------*\
102 Class reducedProblem Declaration
103\*---------------------------------------------------------------------------*/
104
106
109{
110 private:
111
112 public:
113 // Constructors
122
124
126 Eigen::MatrixXd B_matrix;
127
129 Eigen::MatrixXd Y_matrix;
130
132 Eigen::MatrixXd MT_matrix;
133
135 Eigen::MatrixXd M_matrix;
136
138 Eigen::MatrixXd K_matrix;
139
141 List <Eigen::MatrixXd> S_matrix;
142
144 List <Eigen::MatrixXd> C_matrix;
145
147 List <Eigen::MatrixXd> Q_matrix;
148
150 Eigen::MatrixXd P_matrix;
151
153 Eigen::MatrixXd D_matrix;
154
156 List <Eigen::MatrixXd> G_matrix;
157
159 List < Eigen::MatrixXd> online_solutiont;
160
162 PtrList<volScalarField> Tmodes;
163
165 PtrList<volScalarField> Tsnapshots;
166
169
172
174 PtrList<volScalarField> T_rec;
175
178
180 PtrList<volScalarField> TREC;
181
184
186 scalar time;
187
189 scalar dt;
190
192 Eigen::VectorXd z_old;
193
195 Eigen::VectorXd z;
196
198 scalar DT;
199
201 scalar finalTime;
202
204 scalar tstart;
205
208
209 // Functions
210
222 void solveOnline_sup(Eigen::MatrixXd& vel_now, Eigen::MatrixXd& temp_now,
223 int startSnap = 0);
224
231 void reconstruct_sup(fileName folder = "./ITHACAOutput/online_rec",
232 int printevery = 1);
233
240 void reconstruct_supt(fileName folder = "./ITHACAOutput/online_rec",
241 int printevery = 1);
242};
243
244// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245
246#endif
247
248
249
Header file of the reducedUnsteadyNS class.
Eigen::MatrixXd vel_now
Online inlet velocity vector.
Eigen::MatrixXd B_matrix
Diffusion Term.
Eigen::MatrixXd MT_matrix
Mass Matrix Term for temperature.
PtrList< volScalarField > TREC
Reconstructed temperature field.
scalar finalTime
Scalar to store the final time if the online simulation.
int N_BC_t
Number of parametrized boundary conditions related to temperature field.
Eigen::VectorXd z_old
Vector to store the previous temperature solution during the Newton procedure.
Eigen::MatrixXd P_matrix
Divergence of velocity.
PtrList< volScalarField > Tmodes
List of pointers to store the modes for temperature.
Eigen::MatrixXd K_matrix
Gradient of pressure.
void reconstruct_supt(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique for t...
PtrList< volScalarField > T_rec
Reconstructed temperature field.
reducedUnsteadyNST()
Construct Null.
List< Eigen::MatrixXd > online_solutiont
List of Eigen matrices to store the online solution for temperature equation.
List< Eigen::MatrixXd > S_matrix
Turbulent thermal term.
newton_unsteadyNST_sup_t newton_object_sup_t
Functor object to call the non linear solver sup. approach for temperature case.
Eigen::MatrixXd D_matrix
Laplacian of pressure.
scalar tstart
Scalar to store the final time if the online simulation.
unsteadyNST * problem
pointer to the FOM problem
List< Eigen::MatrixXd > Q_matrix
Convective Term for temperature transport equation.
List< Eigen::MatrixXd > G_matrix
Divergence of momentum.
Eigen::MatrixXd M_matrix
Mass Matrix Term.
int Nphi_t
Number of temperature modes.
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.
scalar dt
Scalar to store the time increment.
Eigen::MatrixXd Y_matrix
Laplacian Term for temperature trasnport equation.
newton_unsteadyNST_sup newton_object_sup
Functor object to call the non linear solver sup. approach for velocity case.
Eigen::VectorXd z
Vector to store the temperature solution during the Newton procedure.
List< Eigen::MatrixXd > C_matrix
Convective Term.
PtrList< volScalarField > Tsnapshots
List of pointers to store the snapshots for temperature.
scalar time
Scalar to store the current time.
void solveOnline_sup(Eigen::MatrixXd &vel_now, Eigen::MatrixXd &temp_now, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
reducedUnsteadyNS()
Construct Null.
Implementation of a parametrized full order unsteady NS problem weakly coupled with the energy equat...
Definition unsteadyNST.H:61
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newton_unsteadyNST_sup_t(int Nx, int Ny, unsteadyNST &problem)
Newton object for the resolution of the reduced problem using a supremizer approach.
newton_unsteadyNST_sup(int Nx, int Ny, unsteadyNST &problem)
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
Header file of the unsteadyNST class.