1 #ifndef ABLATELIBRARY_LES_H
2 #define ABLATELIBRARY_LES_H
4 #include "eos/transport/transportModel.hpp"
5 #include "finiteVolume/fluxCalculator/fluxCalculator.hpp"
6 #include "flowProcess.hpp"
7 #include "navierStokesTransport.hpp"
9 namespace ablate::finiteVolume::processes {
16 const std::string tkeField;
19 inline const static PetscReal c_k = 0.094;
20 inline const static PetscReal c_e = 1.048;
21 inline const static PetscReal c_p = 1004.0;
22 inline const static PetscReal scT = 1.00;
23 inline const static PetscReal prT = 1.00;
26 std::vector<PetscInt> numberComponents;
33 explicit LES(std::string tkeField = {});
50 static PetscErrorCode
LesMomentumFlux(PetscInt dim,
const PetscFVFaceGeom* fg,
const PetscInt uOff[],
const PetscInt uOff_x[],
const PetscScalar field[],
const PetscScalar grad[],
51 const PetscInt aOff[],
const PetscInt aOff_x[],
const PetscScalar aux[],
const PetscScalar gradAux[], PetscScalar flux[],
void* ctx);
60 static PetscErrorCode
LesEnergyFlux(PetscInt dim,
const PetscFVFaceGeom* fg,
const PetscInt uOff[],
const PetscInt uOff_x[],
const PetscScalar field[],
const PetscScalar grad[],
61 const PetscInt aOff[],
const PetscInt aOff_x[],
const PetscScalar aux[],
const PetscScalar gradAux[], PetscScalar flux[],
void* ctx);
71 static PetscErrorCode
LesTkeFlux(PetscInt dim,
const PetscFVFaceGeom* fg,
const PetscInt uOff[],
const PetscInt uOff_x[],
const PetscScalar field[],
const PetscScalar grad[],
72 const PetscInt aOff[],
const PetscInt aOff_x[],
const PetscScalar aux[],
const PetscScalar gradAux[], PetscScalar flux[],
void* ctx);
82 static PetscErrorCode
LesEvFlux(PetscInt dim,
const PetscFVFaceGeom* fg,
const PetscInt uOff[],
const PetscInt uOff_x[],
const PetscScalar field[],
const PetscScalar grad[],
const PetscInt aOff[],
83 const PetscInt aOff_x[],
const PetscScalar aux[],
const PetscScalar gradAux[], PetscScalar flux[],
void* ctx);
94 static PetscErrorCode
LesViscosity(PetscInt dim,
const PetscFVFaceGeom* fg,
const PetscScalar* densityField,
const PetscReal turbulence, PetscReal& mut);
Definition: finiteVolumeSolver.hpp:28
Definition: flowProcess.hpp:8
static PetscErrorCode LesViscosity(PetscInt dim, const PetscFVFaceGeom *fg, const PetscScalar *densityField, const PetscReal turbulence, PetscReal &mut)
Definition: les.cpp:210
static PetscErrorCode LesEnergyFlux(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void *ctx)
Definition: les.cpp:93
static PetscErrorCode LesMomentumFlux(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void *ctx)
Definition: les.cpp:55
LES(std::string tkeField={})
Definition: les.cpp:7
static PetscErrorCode LesTkeFlux(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void *ctx)
Definition: les.cpp:133
void Setup(ablate::finiteVolume::FiniteVolumeSolver &flow) override
Definition: les.cpp:9
static PetscErrorCode LesEvFlux(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void *ctx)
Definition: les.cpp:177