Loading...
Searching...
No Matches
UnsteadyNSTTurb.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 unsteadyNS
31
32Description
33 Reduction class for a NON-Stationary NS plus Energy Equation turbulent problem
34
35SourceFiles
36 unsteadyNSTTurb.C
37
38\*---------------------------------------------------------------------------*/
39
44
45#ifndef unsteadyNSTTurb_H
46#define unsteadyNSTTurb_H
47#include "fvCFD.H"
48#include "singlePhaseTransportModel.H"
49#include "turbulentTransportModel.H"
50#include "pimpleControl.H"
51#include "fvOptions.H"
52#include "IOporosityModelList.H"
53#include "IOMRFZoneList.H"
54#include "fixedFluxPressureFvPatchScalarField.H"
55#include "steadyNS.H"
56#include "unsteadyNST.H"
57#include <iostream>
58#include <datatable.h>
59#include <bspline.h>
60#include <bsplinebuilder.h>
61#include <rbfspline.h>
62#include <spline.h>
63
64// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65
66/*---------------------------------------------------------------------------*\
67 Class Declaration
68\*---------------------------------------------------------------------------*/
69
71
75{
76 public:
77 // Constructors
80
82 UnsteadyNSTTurb(int argc, char* argv[]);
84
85 // Specific variable for the unstationary case
86
87 // Member Functions
88
90 PtrList<volScalarField> nutFields;
91
93 PtrList<volScalarField> liftfieldnut;
94
96 PtrList<volScalarField> nutomFields;
97
99 PtrList<volScalarField> nuTmodes;
100
102 PtrList<volScalarField> alphatfield;
103
105 scalar Prt;
106
108 scalar Pr;
109
111 label Nphi_nut;
112
114 std::vector<SPLINTER::DataTable*> SAMPLES;
115
117 std::vector<SPLINTER::RBFSpline*> rbfsplines;
118
121
126 //
128 Eigen::MatrixXd BT_matrix;
129
131 List <Eigen::MatrixXd> CT2_matrix;
132
134 List <Eigen::MatrixXd> CT1_matrix;
135
137 List <Eigen::MatrixXd> S_matrix;
138
140 List <Eigen::MatrixXd> C_total_matrix;
141
143 Eigen::MatrixXd B_total_matrix;
145
147 autoPtr<volScalarField> _nut;
148
150 autoPtr<volScalarField> _nuTilda;
151
153 autoPtr<dimensionedScalar> _Prt;
154
156 autoPtr<dimensionedScalar> _Pr;
157
159 autoPtr<volScalarField> _alphat;
160
168 Eigen::MatrixXd BTturbulence(label NUmodes, label NSUPmodes);
169
179 List < Eigen::MatrixXd > turbulenceTerm1(label NUmodes, label NSUPmodes,
180 label Nnutmodes);
181
190 List < Eigen::MatrixXd > temperatureTurbulenceTerm(label NTmodes,
191 label Nnutmodes);
192
202 List < Eigen::MatrixXd > turbulenceTerm2(label NUmodes, label NSUPmodes,
203 label Nnutmodes);
204
205 //--------------------------------------------------------------------------
211 void truthSolve(List<scalar> mu_now);
212
219 bool checkWrite(Time& timeObject);
220
222 void liftSolveT();
223
233 void projectSUP(fileName folder, label NU, label NP, label NSUP, label Nnut,
234 label NT);
235};
236
237// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238
239#endif
Implementation of a parametrized full order unsteady NST problem weakly coupled with the energy equa...
label Nnutmodes
Number of viscoisty modes used for the projection.
scalar Prt
Scalar to store the turbulent Prandtl number.
PtrList< volScalarField > nutomFields
List of snapshots to form the homogeneous eddy viscosity fields.
Eigen::MatrixXd BT_matrix
Turbulent viscosity term.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
label Nphi_nut
Number of nut field modes.
List< Eigen::MatrixXd > CT1_matrix
Turbulent viscosity term.
bool checkWrite(Time &timeObject)
Function to check if the solution must be exported.
List< Eigen::MatrixXd > S_matrix
Turbulent diffusivity term.
scalar Pr
Scalar to store the Prandtl number.
List< Eigen::MatrixXd > C_total_matrix
Total C Matrix.
void projectSUP(fileName folder, label NU, label NP, label NSUP, label Nnut, label NT)
Project using a supremizer approach.
List< Eigen::MatrixXd > CT2_matrix
Turbulent viscosity term.
autoPtr< volScalarField > _nuTilda
Eddy viscosity for the Spalart-Allmaras turbulence model.
Eigen::MatrixXd B_total_matrix
Total B Matrix.
autoPtr< volScalarField > _alphat
Turbulent thermal diffusivity.
List< Eigen::MatrixXd > turbulenceTerm1(label NUmodes, label NSUPmodes, label Nnutmodes)
CT1 added matrix for the turbulence treatement.
PtrList< volScalarField > alphatfield
Eddy diffusivity field.
List< Eigen::MatrixXd > turbulenceTerm2(label NUmodes, label NSUPmodes, label Nnutmodes)
CT2 added matrix for the turbulence treatement.
PtrList< volScalarField > liftfieldnut
List of snapshots for the solution for eddy viscosity.
std::vector< SPLINTER::RBFSpline * > rbfsplines
Create a SAMPLES for interpolation.
Eigen::MatrixXd BTturbulence(label NUmodes, label NSUPmodes)
BT added matrix for the turbulence treatement.
std::vector< SPLINTER::DataTable * > SAMPLES
Create a Rbf splines for interpolation.
PtrList< volScalarField > nuTmodes
List of POD modes for eddy viscosity.
autoPtr< dimensionedScalar > _Prt
Turbulent Prandtl number.
autoPtr< volScalarField > _nut
Eddy viscosity field.
List< Eigen::MatrixXd > temperatureTurbulenceTerm(label NTmodes, label Nnutmodes)
S added matrix for the thermal turbulence treatement.
autoPtr< dimensionedScalar > _Pr
Prandtl number.
void liftSolveT()
Perform a lift solve for temperature field.
UnsteadyNSTTurb()
Construct Null.
void truthSolve()
Perform a TruthSolve.
label NUmodes
Number of velocity modes used for the projection.
Definition steadyNS.H:140
label NSUPmodes
Number of supremizer modes used for the projection.
Definition steadyNS.H:146
Implementation of a parametrized full order unsteady NS problem weakly coupled with the energy equat...
Definition unsteadyNST.H:61
label NTmodes
Number of temperature modes used for the projection.
Header file of the steadyNS class.
Header file of the unsteadyNST class.