Loading...
Searching...
No Matches
UnsteadyNSTurbIntrusive.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 UnsteadyNSTurbIntrusive
26Description
27 Reduction class for a NON-Stationary turbulent NS problem using fully intrusive approach.
28SourceFiles
29 UnsteadyNSTurbIntrusive.C
30
31\*---------------------------------------------------------------------------*/
32
37
38#ifndef UnsteadyNSTurbIntrusive_H
39#define UnsteadyNSTurbIntrusive_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 UnsteadyNSTurbIntrusive(int argc, char* argv[]);
76
77 // Specific variable for the unstationary case
78 // Member Functions
80 PtrList<volScalarField> nutFields;
81
84
86 Eigen::MatrixXd bMatrix;
87
89 Eigen::MatrixXd kMatrix;
90
92 Eigen::MatrixXd btMatrix;
93
95 Eigen::Tensor<double, 3 > convTensor;
96
98 Eigen::Tensor<double, 3 > ct2Tensor;
99
101 Eigen::Tensor<double, 3 > ct1Tensor;
102
104 Eigen::Tensor<double, 3 > ct1PPETensor;
105
107 Eigen::Tensor<double, 3 > ct2PPETensor;
108
110 Eigen::Tensor<double, 3 > cTotalTensor;
111
113 Eigen::Tensor<double, 3 > cTotalPPETensor;
114
116 Eigen::MatrixXd bTotalMatrix;
117
120
122 autoPtr<volScalarField> _nut;
123
124 //--------------------------------------------------------------------------
130 void truthSolve(List<scalar> mu_now);
131
132 //--------------------------------------------------------------------------
138 void project(fileName folder, label nModes);
139
140 //--------------------------------------------------------------------------
148 void projectPPE(fileName folder, label nUModes, label nPModes);
149
150 //--------------------------------------------------------------------------
157 Eigen::MatrixXd btTurbulence(label nModes);
158
159
160 //--------------------------------------------------------------------------
167 Eigen::Tensor<double, 3 > turbulenceTensor1(label nModes);
168
169 //--------------------------------------------------------------------------
176 Eigen::Tensor<double, 3 > turbulenceTensor2(label nModes);
177
178 //--------------------------------------------------------------------------
186 Eigen::Tensor<double, 3 > turbulencePPETensor1(label nUModes, label nPModes);
187
188 //--------------------------------------------------------------------------
196 Eigen::Tensor<double, 3 > turbulencePPETensor2(label nUModes, label nPModes);
197
198 //--------------------------------------------------------------------------
205 Eigen::MatrixXd diffusiveTerm(label nModes);
206
207 //--------------------------------------------------------------------------
214 Eigen::Tensor<double, 3> convectiveTerm(label nModes);
215
216 //--------------------------------------------------------------------------
225 Eigen::MatrixXd pressureGradientTerm(label nModes);
226
227 //--------------------------------------------------------------------------
237 Eigen::MatrixXd pressureGradientTerm(label nUModes, label nPModes);
238
239 //--------------------------------------------------------------------------
246 Eigen::MatrixXd laplacianPressure(label nPModes);
247
248 //--------------------------------------------------------------------------
256 Eigen::Tensor<double, 3> divMomentum(label nUModes, label nPModes);
257
258 //--------------------------------------------------------------------------
266 Eigen::MatrixXd pressureBC1(label nUModes, label nPModes);
267
268 //--------------------------------------------------------------------------
276 Eigen::Tensor<double, 3 > pressureBC2(label nUModes, label nPModes);
277
278 //--------------------------------------------------------------------------
286 Eigen::MatrixXd pressureBC3(label nUModes, label nPModes);
287
288 //--------------------------------------------------------------------------
296 Eigen::MatrixXd pressureBC4(label nUModes, label nPModes);
297
298 //--------------------------------------------------------------------------
306 List< Eigen::MatrixXd > bcVelocityVec(label nModes);
307
308 //--------------------------------------------------------------------------
315 List< Eigen::MatrixXd > bcVelocityMat(label nModes);
316
317};
318
319// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320
321#endif
322
323
324
325
326
327
328
329
330
331
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
Eigen::Tensor< double, 3 > convectiveTerm(label nModes)
The method for computing the convective term tensor.
Eigen::MatrixXd pressureBC4(label nUModes, label nPModes)
Term N° 4 given by the additional boundary condition using a PPE approach.
Eigen::MatrixXd kMatrix
Pressure Gradient matrix.
Eigen::Tensor< double, 3 > cTotalTensor
Total Turbulent tensor.
Eigen::Tensor< double, 3 > turbulencePPETensor1(label nUModes, label nPModes)
Method to compute one of the turbulence eddy viscosity tensors in the PPE approach.
List< Eigen::MatrixXd > bcVelocityMat(label nModes)
Boundary integral modes on boundary used by the penaly method.
Eigen::MatrixXd pressureGradientTerm(label nModes)
The method for computing the pressure gradient term with number of modes of pressure being equal to t...
Eigen::Tensor< double, 3 > ct1PPETensor
Turbulent viscosity tensor PPE approach.
Eigen::Tensor< double, 3 > turbulenceTensor1(label nModes)
Method to compute one of the turbulence tensors.
Eigen::MatrixXd btTurbulence(label nModes)
bt added matrix for the turbulence treatement
Eigen::Tensor< double, 3 > divMomentum(label nUModes, label nPModes)
Divergence of convective term (PPE approach)
Eigen::Tensor< double, 3 > turbulenceTensor2(label nModes)
Method to compute one of the turbulence tensors.
Eigen::Tensor< double, 3 > ct2PPETensor
Turbulent viscosity tensor PPE approach.
Eigen::Tensor< double, 3 > pressureBC2(label nUModes, label nPModes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Eigen::Tensor< double, 3 > turbulencePPETensor2(label nUModes, label nPModes)
Method to compute one of the turbulence eddy viscosity tensors in the PPE approach.
Eigen::MatrixXd pressureBC3(label nUModes, label nPModes)
Term N° 3 given by the additional boundary condition using a PPE approach.
Eigen::Tensor< double, 3 > convTensor
Convective tensor.
void projectPPE(fileName folder, label nUModes, label nPModes)
The projection method for computing the reduced order matrices with the usage of Poisson equation for...
Eigen::Tensor< double, 3 > ct2Tensor
Turbulent viscosity tensor.
Eigen::Tensor< double, 3 > cTotalPPETensor
Total Turbulent tensor PPE approach.
Eigen::MatrixXd btMatrix
Turbulent viscosity matrix.
Eigen::MatrixXd bTotalMatrix
Total B Matrix.
Eigen::MatrixXd laplacianPressure(label nPModes)
Laplacian of pressure term (PPE approach)
Eigen::Tensor< double, 3 > ct1Tensor
Turbulent viscosity tensor.
UnsteadyNSTurbIntrusive()
Construct Null.
autoPtr< volScalarField > _nut
Eddy viscosity field.
Eigen::MatrixXd bMatrix
Diffusive matrix.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
label nModesOnline
Number of modes used in the online stage.
volScalarModes nutModes
List of POD modes for eddy viscosity.
Eigen::MatrixXd diffusiveTerm(label nModes)
Diffusive Term.
Eigen::MatrixXd pressureBC1(label nUModes, label nPModes)
Term N° 1 given by the additional boundary condition using a PPE approach.
List< Eigen::MatrixXd > bcVelocityVec(label nModes)
Boundary integral modes on boundary used by the penaly method.
void project()
General projection operation.
void truthSolve()
Perform a TruthSolve.
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
Definition unsteadyNS.H:62
Header file of the steadyNS class.
Header file of the unsteadyNS class.