Loading...
Searching...
No Matches
SteadyNSTurb.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 SteadyNSTurb
26Description
27 Reduction class for a Stationary turbulent Navier-Stokes problem.
28SourceFiles
29 SteadyNSTurb.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef SteadyNSTurb_H
38#define SteadyNSTurb_H
39#include "fvCFD.H"
40#include "singlePhaseTransportModel.H"
41#include "turbulentTransportModel.H"
42#include "simpleControl.H"
43#include "pisoControl.H"
44#include "fvOptions.H"
45#include "reductionProblem.H"
46#include "ITHACAstream.H"
47#include <iostream>
48#include "ithacaInterpolator.H"
49#include <datatable.h>
50#include <bspline.h>
51#include <bsplinebuilder.h>
52#include <rbfspline.h>
53#include <spline.h>
54#include <Eigen/Dense>
55#include <unsupported/Eigen/NonLinearOptimization>
56#include <unsupported/Eigen/NumericalDiff>
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59/*---------------------------------------------------------------------------*\
60 Class SteadyNSturb Declaration
61\*---------------------------------------------------------------------------*/
62
64
66class SteadyNSTurb: public steadyNS
67{
68
69
70 public:
71 // Constructors
72 SteadyNSTurb();
74 SteadyNSTurb(int argc, char* argv[]);
75
76 // Member Functions
78 PtrList<volScalarField> nutFields;
79
81 volScalarModes nutModes;
82
84 std::vector<std::shared_ptr<ithacaInterpolator>> rbfSplines;
85
87 Eigen::MatrixXd btMatrix;
88
90 List <Eigen::MatrixXd> ct2Matrix;
91 Eigen::Tensor<double, 3 > ct2Tensor;
92
94 List <Eigen::MatrixXd> ct1Matrix;
95 Eigen::Tensor<double, 3 > ct1Tensor;
96
98 Eigen::Tensor<double, 3 > ct1PPETensor;
99
101 Eigen::Tensor<double, 3 > ct2PPETensor;
102
104 List <Eigen::MatrixXd> cTotalMatrix;
105 Eigen::Tensor<double, 3 > cTotalTensor;
106
108 Eigen::Tensor<double, 3 > cTotalPPETensor;
109
111 Eigen::MatrixXd bTotalMatrix;
112
114 Eigen::MatrixXd coeffL2;
115
117 Eigen::VectorXd nutCoeff;
118
120 dictionary viscDict;
121
124
126 autoPtr<volScalarField> _nut;
127 autoPtr<volScalarField> _nuTilda;
128
129 //--------------------------------------------------------------------------
135 void truthSolve(List<scalar> mu_now);
136
137 //--------------------------------------------------------------------------
146 void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes,
147 label nNutModes);
148
149 //--------------------------------------------------------------------------
158 void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes,
159 label nNutModes);
160
161 //--------------------------------------------------------------------------
170 Eigen::MatrixXd btTurbulence(label NUmodes, label NSUPmodes);
171
172
173 //--------------------------------------------------------------------------
182 List < Eigen::MatrixXd > turbulenceTerm1(label NUmodes, label NSUPmodes,
183 label nNutModes);
184
185 //--------------------------------------------------------------------------
194 Eigen::Tensor<double, 3 > turbulenceTensor1(label NUmodes,
195 label NSUPmodes, label nNutModes);
196
197 //--------------------------------------------------------------------------
206 Eigen::Tensor<double, 3 > turbulenceTensor1_cache(label NUmodes,
207 label NSUPmodes, label nNutModes);
208
209 //--------------------------------------------------------------------------
218 List < Eigen::MatrixXd > turbulenceTerm2(label NUmodes, label NSUPmodes,
219 label nNutModes);
220
221 //--------------------------------------------------------------------------
230 Eigen::Tensor<double, 3 > turbulenceTensor2(label NUmodes,
231 label NSUPmodes, label nNutModes);
232
233 //--------------------------------------------------------------------------
242 Eigen::Tensor<double, 3 > turbulenceTensor2_cache(label NUmodes,
243 label NSUPmodes, label nNutModes);
244
245 //--------------------------------------------------------------------------
255 Eigen::Tensor<double, 3 > turbulencePPETensor1(label NUmodes,
256 label NSUPmodes, label NPmodes, label nNutModes);
257
258 //--------------------------------------------------------------------------
268 Eigen::Tensor<double, 3 > turbulencePPETensor1_cache(label NUmodes,
269 label NSUPmodes, label NPmodes, label nNutModes);
270
271 //--------------------------------------------------------------------------
281 Eigen::Tensor<double, 3 > turbulencePPETensor2(label NUmodes,
282 label NSUPmodes, label NPmodes, label nNutModes);
283
284 //--------------------------------------------------------------------------
294 Eigen::Tensor<double, 3 > turbulencePPETensor2_cache(label NUmodes,
295 label NSUPmodes, label NPmodes, label nNutModes);
296};
297
298#endif
299
300
301
302
303
304
305
306
307
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Eigen::MatrixXd coeffL2
The matrix of L2 projection coefficients for the eddy viscosity.
Eigen::MatrixXd btMatrix
Turbulent viscosity matrix.
Eigen::Tensor< double, 3 > cTotalPPETensor
Turbulent total viscosity tensor in the PPE equation.
label nNutModes
Number of viscoisty modes used for the projection.
Eigen::MatrixXd bTotalMatrix
Total B Matrix.
Eigen::Tensor< double, 3 > turbulenceTensor2_cache(label NUmodes, label NSUPmodes, label nNutModes)
Method to compute one of the turbulence eddy viscosity tensors using the cached procedure.
Eigen::MatrixXd btTurbulence(label NUmodes, label NSUPmodes)
bt added matrix for the turbulence treatement
autoPtr< volScalarField > _nut
Eddy viscosity field.
List< Eigen::MatrixXd > turbulenceTerm1(label NUmodes, label NSUPmodes, label nNutModes)
ct1 added matrix for the turbulence treatement
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
Eigen::Tensor< double, 3 > turbulenceTensor1_cache(label NUmodes, label NSUPmodes, label nNutModes)
ct1 tensor for the turbulence treatement using the cached procedure
List< Eigen::MatrixXd > ct2Matrix
Turbulent viscosity tensor.
List< Eigen::MatrixXd > cTotalMatrix
Total Turbulent tensor.
Eigen::Tensor< double, 3 > ct1PPETensor
Turbulent viscosity tensor in the PPE equation.
Eigen::Tensor< double, 3 > turbulencePPETensor2_cache(label NUmodes, label NSUPmodes, label NPmodes, label nNutModes)
ct2PPE added tensor for the turbulence treatement in the PPE method using the cached procedure
Eigen::Tensor< double, 3 > turbulencePPETensor1_cache(label NUmodes, label NSUPmodes, label NPmodes, label nNutModes)
ct1PPE added tensor for the turbulence treatement in the PPE method using the cached procedure
List< Eigen::MatrixXd > ct1Matrix
Turbulent viscosity tensor.
List< Eigen::MatrixXd > turbulenceTerm2(label NUmodes, label NSUPmodes, label nNutModes)
Method to compute one of the turbulence eddy viscosity tensors.
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes)
Project using a supremizer approach.
Eigen::Tensor< double, 3 > ct2PPETensor
Turbulent viscosity tensor in the PPE equation.
Eigen::Tensor< double, 3 > turbulencePPETensor1(label NUmodes, label NSUPmodes, label NPmodes, label nNutModes)
ct1PPE added tensor for the turbulence treatement in the PPE method
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::Tensor< double, 3 > turbulenceTensor1(label NUmodes, label NSUPmodes, label nNutModes)
ct1 tensor for the turbulence treatement
Eigen::Tensor< double, 3 > turbulenceTensor2(label NUmodes, label NSUPmodes, label nNutModes)
Method to compute one of the turbulence eddy viscosity tensors.
std::vector< std::shared_ptr< ithacaInterpolator > > rbfSplines
Interpolators for nut coefficient interpolation.
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes, label nNutModes)
Project using a PPE approach.
Eigen::VectorXd nutCoeff
The vector of L2 projection coefficients for the eddy viscosity snapshot.
volScalarModes nutModes
List of POD modes for eddy viscosity.
dictionary viscDict
The way to compute the eddy viscosity snapshots.
void truthSolve()
Perform a TruthSolve.
label NPmodes
Number of pressure modes used for the projection.
Definition steadyNS.H:146
steadyNS()
Null constructor.
Definition steadyNS.C:40
label NUmodes
Number of velocity modes used for the projection.
Definition steadyNS.H:143
label NSUPmodes
Number of supremizer modes used for the projection.
Definition steadyNS.H:149
Header file of the reductionProblem class.