Loading...
Searching...
No Matches
ReducedUnsteadyNSTTurb.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-------------------------------------------------------------------------------
12
13License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29Class
30 reducedProblem
31
32Description
33 A reduced problem for the unsteady NS plus the energy equations for turbulence modelling
34
35SourceFiles
36 reducedUnsteadyNSTTurb.C
37
38\*---------------------------------------------------------------------------*/
39
44
45#ifndef ReducedUnsteadyNSTTurb_H
46#define ReducedUnsteadyNSTTurb_H
47
48#include "fvCFD.H"
49#include "IOmanip.H"
50#include "ReducedSteadyNS.H"
51#include "ReducedUnsteadyNST.H"
52#include "UnsteadyNSTTurb.H"
53#include <Eigen/Dense>
54#include <unsupported/Eigen/NonLinearOptimization>
55#include <unsupported/Eigen/NumericalDiff>
56
58{
59 public:
61
63 UnsteadyNSTTurb& problem): newton_argument<double>(Nx, Ny),
65 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
66 Nphi_nut(problem.Nnutmodes),
67 Nphi_p(problem.NPmodes),
68 N_BC(problem.inletIndex.rows()),
69 nu_c(problem.Nnutmodes)
70 {}
71
72 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
73 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
75 int Nphi_u;
77 int Nphi_p;
78 int N_BC;
79 Eigen::VectorXd nu_c;
81 Eigen::VectorXd mu_now;
82 scalar nu;
83 scalar dt;
84 Eigen::VectorXd y_old;
85 Eigen::VectorXd BC;
86 std::vector<SPLINTER::RBFSpline*> SPLINES;
87};
88
90{
91 public:
93
95 UnsteadyNSTTurb& problem): newton_argument<double>(Nx, Ny),
97 Nphi_t(problem.NTmodes + problem.liftfieldT.size()),
98 N_BC_t(problem.inletIndexT.rows()),
99 nu_c(problem.Nnutmodes)
100 {}
101
102 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
103 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
109 Eigen::VectorXd nu_c;
110 scalar nu;
111 scalar DT;
112 scalar Prt;
113 scalar Pr;
114 scalar dt;
115 Eigen::VectorXd a_tmp;
116 Eigen::MatrixXd Y_matrix;
117 Eigen::MatrixXd MT_matrix;
118 Eigen::VectorXd mu_now;
119 Eigen::VectorXd z;
120 Eigen::VectorXd z_old;
121 Eigen::VectorXd BC_t;
122 std::vector<SPLINTER::RBFSpline*> SPLINES;
123};
124
125/*---------------------------------------------------------------------------*\
126 Class reducedProblem Declaration
127\*---------------------------------------------------------------------------*/
128
130
133{
134 private:
135
136 public:
137 // Constructors
139 //
141
148
150
152 PtrList<volScalarField> nutFields;
153
155 PtrList<volScalarField> nuTmodes;
156
158 PtrList<volScalarField> nutREC;
159
162
165
167 scalar Prt;
168
170 scalar Pr;
171
173 double time;
174
176 double dt;
177
179 double finalTime;
180
182 double tstart;
183
186
189
190 // Functions
202 //void solveOnline_sup(Eigen::MatrixXd, int startSnap=0);
203 void solveOnlineSup(Eigen::MatrixXd& vel_now, Eigen::MatrixXd& temp_now,
204 int startSnap = 0);
205
212 void reconstructSup(fileName folder = "./ITHACAOutput/online_rec",
213 int printevery = 1);
214
221 void reconstructSupt(fileName folder = "./ITHACAOutput/online_rec",
222 int printevery = 1);
223};
224
225// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226
227#endif
Header file of the reducedSteadyNS class.
Header file of the reducedUnsteadyNST class.
Header file of the unsteadyNSTTurb class.
ReducedUnsteadyNSTTurb()
Construct Null.
void solveOnlineSup(Eigen::MatrixXd &vel_now, Eigen::MatrixXd &temp_now, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
double time
Scalar to store the current time.
newton_unsteadyNSTTurb_sup_t newton_object_sup_t
Functor object to call the non linear solver sup. approach for temperature field.
double dt
Scalar to store the time increment.
newton_unsteadyNSTTurb_sup newton_object_sup
Functor object to call the non linear solver sup. approach.
scalar Prt
Scalar to store the turbulent Prandtl number.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
PtrList< volScalarField > nuTmodes
List of POD modes for eddy viscosity.
int Nphi_nut
Number of nut field modes.
UnsteadyNSTTurb * problem
Pointer to the FOM problem.
void reconstructSup(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution for velocity and pressure from an online solve with a supremizer sta...
PtrList< volScalarField > nutREC
Reconstructed eddy viscosity field.
double finalTime
Scalar to store the final time if the online simulation.
void reconstructSupt(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution for temperature from an online solve with a supremizer stabilisation...
scalar Pr
Scalar to store the Prandtl number.
double tstart
Scalar to store the final time if the online simulation.
Implementation of a parametrized full order unsteady NST problem weakly coupled with the energy equa...
Eigen::MatrixXd vel_now
Online inlet velocity vector.
reducedUnsteadyNST()
Construct Null.
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
std::vector< SPLINTER::RBFSpline * > SPLINES
newton_unsteadyNSTTurb_sup_t(int Nx, int Ny, UnsteadyNSTTurb &problem)
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newton_unsteadyNSTTurb_sup(int Nx, int Ny, UnsteadyNSTTurb &problem)
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
std::vector< SPLINTER::RBFSpline * > SPLINES