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 newton_unsteadyNST_sup(int Nx, int Ny,
53 unsteadyNST& problem): newton_argument<double>(Nx, Ny),
55 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
56 Nphi_p(problem.NPmodes),
57 N_BC(problem.inletIndex.rows())
58 {}
59
60 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
61 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
63 int Nphi_u;
64 int Nphi_p;
65 int N_BC;
66 scalar nu;
67 scalar dt;
68 Eigen::VectorXd y_old;
69 Eigen::VectorXd BC;
70
71};
72
74{
75 public:
78 unsteadyNST& problem): newton_argument<double>(Nx, Ny),
80 Nphi_t(problem.NTmodes + problem.liftfieldT.size()),
81 N_BC_t(problem.inletIndexT.rows())
82 {}
83
84 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
85 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
87 int Nphi_t;
88 int N_BC_t;
89 int Nphi_u;
90 scalar nu;
91 scalar DT;
92 scalar dt;
93 Eigen::VectorXd a_tmp;
94 Eigen::VectorXd z_old;
95 Eigen::VectorXd BC_t;
96
97};
98
99/*---------------------------------------------------------------------------*\
100 Class reducedProblem Declaration
101\*---------------------------------------------------------------------------*/
102
104
107{
108 private:
109
110 public:
111 // Constructors
120
122
124 Eigen::MatrixXd B_matrix;
125
127 Eigen::MatrixXd Y_matrix;
128
130 Eigen::MatrixXd MT_matrix;
131
133 Eigen::MatrixXd M_matrix;
134
136 Eigen::MatrixXd K_matrix;
137
139 List <Eigen::MatrixXd> S_matrix;
140
142 List <Eigen::MatrixXd> C_matrix;
143
145 List <Eigen::MatrixXd> Q_matrix;
146
148 Eigen::MatrixXd P_matrix;
149
151 Eigen::MatrixXd D_matrix;
152
154 List <Eigen::MatrixXd> G_matrix;
155
157 List < Eigen::MatrixXd> online_solutiont;
158
160 PtrList<volScalarField> Tmodes;
161
163 PtrList<volScalarField> Tsnapshots;
164
167
170
172 PtrList<volScalarField> T_rec;
173
176
178 PtrList<volScalarField> TREC;
179
182
184 scalar time;
185
187 scalar dt;
188
190 Eigen::VectorXd z_old;
191
193 Eigen::VectorXd z;
194
196 scalar DT;
197
199 scalar finalTime;
200
202 scalar tstart;
203
206
207 // Functions
208
220 void solveOnline_sup(Eigen::MatrixXd& vel_now, Eigen::MatrixXd& temp_now,
221 int startSnap = 0);
222
229 void reconstruct_sup(fileName folder = "./ITHACAOutput/online_rec",
230 int printevery = 1);
231
238 void reconstruct_supt(fileName folder = "./ITHACAOutput/online_rec",
239 int printevery = 1);
240};
241// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242
243#endif
244
245
246
Header file of the reducedUnsteadyNS class.
Template object created to solve non linear problems using the Eigen library.
Eigen::MatrixXd vel_now
Online inlet velocity vector.
Class where it is implemented a reduced problem for the unsteady Navier-stokes weakly coupled with t...
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.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
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.