38#ifndef UnsteadyNSTurb_H
39#define UnsteadyNSTurb_H
41#include "singlePhaseTransportModel.H"
42#include "turbulentTransportModel.H"
43#include "pimpleControl.H"
45#include "IOporosityModelList.H"
46#include "IOMRFZoneList.H"
47#include "fixedFluxPressureFvPatchScalarField.H"
53#include <bsplinebuilder.h>
89 std::vector<SPLINTER::DataTable*>
samples;
161 autoPtr<surfaceVectorField>
_Uf;
285 void truthSolve(List<scalar> mu_now, std::string& offlinepath);
329 Eigen::VectorXd initSnapInd, Eigen::VectorXd timeSnap);
341 List < Eigen::MatrixXd >
velParCoeff(Eigen::MatrixXd
A, Eigen::MatrixXd G);
360 Eigen::VectorXd initSnapInd, Eigen::VectorXd timeSnap);
376 Eigen::VectorXd par,
double timeSnap);
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.
label NUmodes
Number of velocity modes used for the projection.
label NSUPmodes
Number of supremizer modes used for the projection.
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
Header file of the steadyNS class.
Header file of the unsteadyNS class.