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:
62 UnsteadyNSTTurb& problem): newton_argument<double>(Nx, Ny),
64 Nphi_u(problem.NUmodes + problem.liftfield.size() + problem.NSUPmodes),
65 Nphi_nut(problem.Nnutmodes),
66 Nphi_p(problem.NPmodes),
67 N_BC(problem.inletIndex.rows()),
68 nu_c(problem.Nnutmodes)
69 {}
70
71 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
72 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
74 int Nphi_u;
76 int Nphi_p;
77 int N_BC;
78 Eigen::VectorXd nu_c;
80 Eigen::VectorXd mu_now;
81 scalar nu;
82 scalar dt;
83 Eigen::VectorXd y_old;
84 Eigen::VectorXd BC;
85 std::vector<SPLINTER::RBFSpline*> SPLINES;
86};
87
89{
90 public:
93 UnsteadyNSTTurb& problem): newton_argument<double>(Nx, Ny),
95 Nphi_t(problem.NTmodes + problem.liftfieldT.size()),
96 N_BC_t(problem.inletIndexT.rows()),
97 nu_c(problem.Nnutmodes)
98 {}
99
100 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
101 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
107 Eigen::VectorXd nu_c;
108 scalar nu;
109 scalar DT;
110 scalar Prt;
111 scalar Pr;
112 scalar dt;
113 Eigen::VectorXd a_tmp;
114 Eigen::MatrixXd Y_matrix;
115 Eigen::MatrixXd MT_matrix;
116 Eigen::VectorXd mu_now;
117 Eigen::VectorXd z;
118 Eigen::VectorXd z_old;
119 Eigen::VectorXd BC_t;
120 std::vector<SPLINTER::RBFSpline*> SPLINES;
121};
122
123/*---------------------------------------------------------------------------*\
124 Class reducedProblem Declaration
125\*---------------------------------------------------------------------------*/
126
128
131{
132 private:
133
134 public:
135 // Constructors
137 //
139
146
148
150 PtrList<volScalarField> nutFields;
151
153 PtrList<volScalarField> nuTmodes;
154
156 PtrList<volScalarField> nutREC;
157
160
163
165 scalar Prt;
166
168 scalar Pr;
169
171 double time;
172
174 double dt;
175
177 double finalTime;
178
180 double tstart;
181
184
187
188 // Functions
200 //void solveOnline_sup(Eigen::MatrixXd, int startSnap=0);
201 void solveOnlineSup(Eigen::MatrixXd& vel_now, Eigen::MatrixXd& temp_now,
202 int startSnap = 0);
203
210 void reconstructSup(fileName folder = "./ITHACAOutput/online_rec",
211 int printevery = 1);
212
219 void reconstructSupt(fileName folder = "./ITHACAOutput/online_rec",
220 int printevery = 1);
221};
222
223// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224
225#endif
Header file of the reducedSteadyNS class.
Header file of the reducedUnsteadyNST class.
Header file of the unsteadyNSTTurb class.
Class where it is implemented a reduced problem for the unsteady Navier-stokes weakly coupled with t...
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...
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...
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