Loading...
Searching...
No Matches
SteadyNSTurbIntrusive.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 SteadyNSTurbIntrusive
26Description
27 Reduction class for a Stationary turbulent Navier-Stokes problem using fully intrusive approach.
28SourceFiles
29 SteadyNSTurbIntrusive.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef SteadyNSTurbIntrusive_H
38#define SteadyNSTurbIntrusive_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 <datatable.h>
49#include <bspline.h>
50#include <bsplinebuilder.h>
51#include <rbfspline.h>
52#include <spline.h>
53#include <Eigen/Dense>
54#include <unsupported/Eigen/NonLinearOptimization>
55#include <unsupported/Eigen/NumericalDiff>
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58/*---------------------------------------------------------------------------*\
59 Class SteadyNSturb Declaration
60\*---------------------------------------------------------------------------*/
61
63
66{
67
68
69 public:
70 // Constructors
73 SteadyNSTurbIntrusive(int argc, char* argv[]);
74
75 // Member Functions
77 PtrList<volScalarField> nutFields;
78
80 PtrList<volScalarField> nutModes;
81
83 Eigen::MatrixXd bMatrix;
84
86 Eigen::MatrixXd kMatrix;
87
89 Eigen::MatrixXd btMatrix;
90
92 Eigen::Tensor<double, 3 > convTensor;
93
95 Eigen::Tensor<double, 3 > ct2Tensor;
96
98 Eigen::Tensor<double, 3 > ct1Tensor;
99
101 Eigen::Tensor<double, 3 > cTotalTensor;
102
104 Eigen::MatrixXd bTotalMatrix;
105
108
110 autoPtr<volScalarField> _nut;
111 autoPtr<volScalarField> _nuTilda;
112
113 //--------------------------------------------------------------------------
120 void truthSolve(List<scalar> mu_now);
121
122 //--------------------------------------------------------------------------
128 void project(fileName folder, label nModes);
129
130
131 //--------------------------------------------------------------------------
138 Eigen::MatrixXd btTurbulence(label nModes);
139
140
141 //--------------------------------------------------------------------------
148 Eigen::Tensor<double, 3 > turbulenceTensor1(label nModes);
149
150
151 //--------------------------------------------------------------------------
158 Eigen::Tensor<double, 3 > turbulenceTensor2(label nModes);
159
160 //--------------------------------------------------------------------------
167 Eigen::MatrixXd diffusiveTerm(label nModes);
168
169 //--------------------------------------------------------------------------
176 Eigen::Tensor<double, 3> convectiveTerm(label nModes);
177
178 //--------------------------------------------------------------------------
187 Eigen::MatrixXd pressureGradientTerm(label nModes);
188
189 //--------------------------------------------------------------------------
197 List< Eigen::MatrixXd > bcVelocityVec(label nModes);
198
199 //--------------------------------------------------------------------------
206 List< Eigen::MatrixXd > bcVelocityMat(label nModes);
207
208
209};
210
211#endif
212
213
214
215
216
217
218
219
220
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Implementation of a parametrized full order steady turbulent Navier Stokes problem and preparation ...
Eigen::MatrixXd diffusiveTerm(label nModes)
Diffusive Term.
Eigen::MatrixXd bTotalMatrix
Total B Matrix.
PtrList< volScalarField > nutModes
List of POD modes for eddy viscosity.
Eigen::MatrixXd btTurbulence(label nModes)
bt added matrix for the turbulence treatement
Eigen::Tensor< double, 3 > turbulenceTensor2(label nModes)
Method to compute one of the turbulence eddy viscosity tensors.
Eigen::Tensor< double, 3 > convTensor
Convective tensor.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
List< Eigen::MatrixXd > bcVelocityVec(label nModes)
Boundary integral modes on boundary used by the penaly method.
Eigen::Tensor< double, 3 > cTotalTensor
Total Turbulent tensor.
Eigen::Tensor< double, 3 > convectiveTerm(label nModes)
The method for computing the convective term tensor.
autoPtr< volScalarField > _nut
Eddy viscosity field.
Eigen::MatrixXd kMatrix
Pressure Gradient matrix.
Eigen::MatrixXd btMatrix
Turbulent viscosity matrix.
List< Eigen::MatrixXd > bcVelocityMat(label nModes)
Boundary integral modes on boundary used by the penaly method.
Eigen::Tensor< double, 3 > turbulenceTensor1(label nModes)
Method to compute one of the turbulence eddy viscosity tensors.
Eigen::MatrixXd pressureGradientTerm(label nModes)
The method for computing the pressure gradient term with number of modes of pressure being equal to t...
label nModesOnline
Number of modes used in the online stage.
Eigen::MatrixXd bMatrix
Diffusive matrix.
Eigen::Tensor< double, 3 > ct2Tensor
Turbulent viscosity tensor.
Eigen::Tensor< double, 3 > ct1Tensor
Turbulent viscosity tensor.
autoPtr< volScalarField > _nuTilda
void project()
General projection operation.
void truthSolve()
Perform a TruthSolve.
Implementation of a parametrized full order steady NS problem and preparation of the the reduced ma...
Definition steadyNS.H:70
Header file of the reductionProblem class.