Loading...
Searching...
No Matches
UnsteadyNSTurb.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 UnsteadyNSTurb
26Description
27 Reduction class for a NON-Stationary turbulent NS problem
28SourceFiles
29 UnsteadyNSTurb.C
30
31\*---------------------------------------------------------------------------*/
32
37
38#ifndef UnsteadyNSTurb_H
39#define UnsteadyNSTurb_H
40#include "fvCFD.H"
41#include "singlePhaseTransportModel.H"
42#include "turbulentTransportModel.H"
43#include "pimpleControl.H"
44#include "fvOptions.H"
45#include "IOporosityModelList.H"
46#include "IOMRFZoneList.H"
47#include "fixedFluxPressureFvPatchScalarField.H"
48#include "steadyNS.H"
49#include "unsteadyNS.H"
50#include <iostream>
51#include <datatable.h>
52#include <bspline.h>
53#include <bsplinebuilder.h>
54#include <rbfspline.h>
55#include <spline.h>
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59/*---------------------------------------------------------------------------*\
60 Class SteadyNS Declaration
61\*---------------------------------------------------------------------------*/
62
64
68{
69 public:
70 // Constructors
73
75 UnsteadyNSTurb(int argc, char* argv[]);
76
77 // Specific variable for the unstationary case
78 // Member Functions
80 PtrList<volScalarField> nutFields;
81
84
86 PtrList<volScalarField> nutAve;
87
89 std::vector<SPLINTER::DataTable*> samples;
90
92 std::vector<SPLINTER::RBFSpline*> rbfSplines;
93
95 Eigen::MatrixXd btMatrix;
96
98 Eigen::MatrixXd z;
99
101 Eigen::Tensor<double, 3 > ct1Tensor;
102
104 Eigen::Tensor<double, 3 > ct1AveTensor;
105
107 Eigen::Tensor<double, 3 > ct1PPETensor;
108
110 Eigen::Tensor<double, 3 > ct1PPEAveTensor;
111
113 Eigen::Tensor<double, 3 > ct2Tensor;
114
116 Eigen::Tensor<double, 3 > ct2AveTensor;
117
119 Eigen::Tensor<double, 3 > ct2PPETensor;
120
122 Eigen::Tensor<double, 3 > ct2PPEAveTensor;
123
125 Eigen::Tensor<double, 3 > cTotalTensor;
126
128 Eigen::Tensor<double, 3 > cTotalAveTensor;
129
131 Eigen::Tensor<double, 3 > cTotalPPETensor;
132
134 Eigen::Tensor<double, 3 > cTotalPPEAveTensor;
135
137 Eigen::MatrixXd bTotalMatrix;
138
140 Eigen::MatrixXd coeffL2;
141
143 Eigen::MatrixXd velRBF;
144
146 Eigen::VectorXd radii;
147
149 double e = 1;
150
153
155 label interChoice = 1;
156
158 autoPtr<volScalarField> _nut;
159
161 autoPtr<surfaceVectorField> _Uf;
162
165
168
169 //--------------------------------------------------------------------------
177 Eigen::MatrixXd btTurbulence(label NUmodes, label NSUPmodes);
178
179 //--------------------------------------------------------------------------
188 Eigen::Tensor<double, 3 > turbulenceTensor1(label NUmodes,
189 label NSUPmodes, label nNutModes);
190
191 //--------------------------------------------------------------------------
200 Eigen::Tensor<double, 3 > turbulenceTensor2(label NUmodes,
201 label NSUPmodes, label nNutModes);
202
203 //--------------------------------------------------------------------------
213 Eigen::Tensor<double, 3 > turbulencePPETensor1(label NUmodes,
214 label NSUPmodes, label NPmodes, label nNutModes);
215
216 //--------------------------------------------------------------------------
226 Eigen::Tensor<double, 3 > turbulencePPETensor2(label NUmodes,
227 label NSUPmodes, label NPmodes, label nNutModes);
228
229 //--------------------------------------------------------------------------
238 Eigen::Tensor<double, 3> turbulenceAveTensor1(label NUmodes,
239 label NSUPmodes);
240
241 //--------------------------------------------------------------------------
250 Eigen::Tensor<double, 3> turbulenceAveTensor2(label NUmodes,
251 label NSUPmodes);
252 //--------------------------------------------------------------------------
263 Eigen::Tensor<double, 3> turbulencePPEAveTensor1(label NUmodes,
264 label NSUPmodes, label NPmodes);
265 //--------------------------------------------------------------------------
276 Eigen::Tensor<double, 3> turbulencePPEAveTensor2(label NUmodes,
277 label NSUPmodes, label NPmodes);
278 //--------------------------------------------------------------------------
285 void truthSolve(List<scalar> mu_now, std::string& offlinepath);
286
287 //----------------------------------------------------------------------
297 void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes,
298 label nNutModes, bool rbfInterp = true);
299
300 //----------------------------------------------------------------------
311 void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes,
312 label nNutModes, bool rbfInterp = true);
313
314 //--------------------------------------------------------------------------
327 List < Eigen::MatrixXd > velDerivativeCoeff(Eigen::MatrixXd A,
328 Eigen::MatrixXd G,
329 Eigen::VectorXd initSnapInd, Eigen::VectorXd timeSnap);
330
331 //--------------------------------------------------------------------------
341 List < Eigen::MatrixXd > velParCoeff(Eigen::MatrixXd A, Eigen::MatrixXd G);
342
343
344 //--------------------------------------------------------------------------
358 List < Eigen::MatrixXd > velParDerivativeCoeff(Eigen::MatrixXd A,
359 Eigen::MatrixXd G,
360 Eigen::VectorXd initSnapInd, Eigen::VectorXd timeSnap);
361
362 //--------------------------------------------------------------------------
375 Eigen::MatrixXd velParDerivativeCoeff(Eigen::MatrixXd A,
376 Eigen::VectorXd par, double timeSnap);
377
378};
379
380// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
381
382#endif
383
384
385
386
387
388
389
390
391
392
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
Eigen::Tensor< double, 3 > ct1AveTensor
Turbulent average viscosity tensor for the splitting approach.
autoPtr< surfaceVectorField > _Uf
Face velocity field.
Eigen::Tensor< double, 3 > turbulenceAveTensor2(label NUmodes, label NSUPmodes)
ct2Ave added tensor for approximation of the averaged part of the eddy viscosity
Eigen::Tensor< double, 3 > cTotalTensor
Turbulent total viscosity tensor.
PtrList< volScalarField > nutAve
List of for eddy viscosity time-averaged fields.
scalar _pRefValue
Pressure reference value.
Eigen::MatrixXd btTurbulence(label NUmodes, label NSUPmodes)
bt added matrix for the turbulence treatement
List< Eigen::MatrixXd > velParDerivativeCoeff(Eigen::MatrixXd A, Eigen::MatrixXd G, Eigen::VectorXd initSnapInd, Eigen::VectorXd timeSnap)
A method to compute the two matrices needed for the RBF interpolation by combining the parameter valu...
Eigen::Tensor< double, 3 > turbulencePPETensor1(label NUmodes, label NSUPmodes, label NPmodes, label nNutModes)
ct1PPE added tensor for the turbulence treatement in the PPE method
std::vector< SPLINTER::RBFSpline * > rbfSplines
Create a samples for interpolation.
Eigen::Tensor< double, 3 > turbulenceAveTensor1(label NUmodes, label NSUPmodes)
ct1Ave added tensor for approximation of the averaged part of the eddy viscosity
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes, bool rbfInterp=true)
Project using a supremizer approach.
Eigen::Tensor< double, 3 > turbulencePPEAveTensor1(label NUmodes, label NSUPmodes, label NPmodes)
ct1PPEAve added tensor for approximation of the averaged part of the eddy viscosity with the usage of...
UnsteadyNSTurb()
Construct Null.
Eigen::Tensor< double, 3 > turbulenceTensor2(label NUmodes, label NSUPmodes, label nNutModes)
ct2 added tensor for the turbulence treatement
autoPtr< volScalarField > _nut
Eddy viscosity field.
Eigen::Tensor< double, 3 > ct2PPEAveTensor
Turbulent average viscosity tensor for the splitting approach in the PPE equation.
label interChoice
Interpolation independent variable choice.
Eigen::MatrixXd coeffL2
The matrix of L2 projection coefficients for the eddy viscosity.
Eigen::MatrixXd btMatrix
Turbulent viscosity term.
volScalarModes nutModes
List of POD modes for eddy viscosity.
Eigen::Tensor< double, 3 > turbulencePPEAveTensor2(label NUmodes, label NSUPmodes, label NPmodes)
ct2PPEAve added tensor for approximation of the averaged part of the eddy viscosity with the usage of...
Eigen::Tensor< double, 3 > cTotalPPEAveTensor
Turbulent total average viscosity tensor for the splitting approach in the PPE equation.
double e
RBF functions radius.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
Eigen::Tensor< double, 3 > ct1PPEAveTensor
Turbulent average viscosity tensor for the splitting approach in the PPE equation.
Eigen::Tensor< double, 3 > ct2PPETensor
Turbulent viscosity tensor in the PPE equation.
Eigen::Tensor< double, 3 > ct1PPETensor
Turbulent viscosity tensor in the PPE equation.
Eigen::Tensor< double, 3 > ct2Tensor
Turbulent viscosity tensor.
Eigen::Tensor< double, 3 > ct2AveTensor
Turbulent average viscosity tensor for the splitting approach.
label _pRefCell
Pressure reference cell.
Eigen::Tensor< double, 3 > ct1Tensor
Turbulent viscosity tensor.
Eigen::Tensor< double, 3 > cTotalPPETensor
Turbulent total viscosity tensor in the PPE equation.
List< Eigen::MatrixXd > velParCoeff(Eigen::MatrixXd A, Eigen::MatrixXd G)
A method to compute the two matrices needed for the RBF interpolation by combining the parameter samp...
Eigen::MatrixXd bTotalMatrix
Total B Matrix.
List< Eigen::MatrixXd > velDerivativeCoeff(Eigen::MatrixXd A, Eigen::MatrixXd G, Eigen::VectorXd initSnapInd, Eigen::VectorXd timeSnap)
A method to compute the two matrices needed for the RBF interpolation by combining the velocity L2 pr...
label nNutModes
Number of viscoisty modes used for the projection.
Eigen::Tensor< double, 3 > cTotalAveTensor
Turbulent total average viscosity tensor for the splitting approach.
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes, bool rbfInterp=true)
Project using the Poisson Equation for pressure.
Eigen::VectorXd radii
RBF shape parameters vector.
Eigen::Tensor< double, 3 > turbulencePPETensor2(label NUmodes, label NSUPmodes, label NPmodes, label nNutModes)
ct2PPE added tensor for the turbulence treatement in the PPE method
Eigen::MatrixXd z
Time-parameter combined matrix.
std::vector< SPLINTER::DataTable * > samples
Create a Rbf splines for interpolation.
Eigen::MatrixXd velRBF
Velocity coefficients for RBF interpolation.
Eigen::Tensor< double, 3 > turbulenceTensor1(label NUmodes, label NSUPmodes, label nNutModes)
ct1 added tensor for the turbulence treatement
void truthSolve()
Perform a TruthSolve.
label NPmodes
Number of pressure modes used for the projection.
Definition steadyNS.H:143
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 and preparation of the the reduced ...
Definition unsteadyNS.H:62
volScalarField & A
Header file of the steadyNS class.
Header file of the unsteadyNS class.