ENH: DEShybrid - code refactoring/simplification

This commit is contained in:
Andrew Heather 2022-08-25 09:59:01 +01:00
parent 32507b3251
commit 493bfdbdc4
2 changed files with 14 additions and 37 deletions

View File

@ -120,8 +120,7 @@ SourceFiles
#include "surfaceInterpolate.H"
#include "fvcGrad.H"
#include "blendedSchemeBase.H"
#include "turbulentTransportModel.H"
#include "turbulentFluidThermoModel.H"
#include "turbulenceModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -393,38 +392,26 @@ public:
{
const fvMesh& mesh = this->mesh();
typedef compressible::turbulenceModel cmpModel;
typedef incompressible::turbulenceModel icoModel;
// Retrieve LES delta from the mesh database
const auto& delta =
mesh.lookupObject<const volScalarField>(deltaName_);
// Lookup the LES delta from the mesh database
const volScalarField& delta = this->mesh().template
lookupObject<const volScalarField>(deltaName_);
// Retrieve turbulence model from the mesh database
const auto* modelPtr =
mesh.cfindObject<turbulenceModel>
(
turbulenceModel::propertiesName
);
// Could avoid the compressible/incompressible case by looking
// up all fields from the database - but retrieving from model
// ensures consistent fields are being employed e.g. for multiphase
// where group name is used
if (mesh.foundObject<icoModel>(icoModel::propertiesName))
if (modelPtr)
{
const icoModel& model =
mesh.lookupObject<icoModel>(icoModel::propertiesName);
const auto& model = *modelPtr;
return calcBlendingFactor
(
vf, model.nut(), model.nu(), model.U(), delta
);
}
else if (mesh.foundObject<cmpModel>(cmpModel::propertiesName))
{
const cmpModel& model =
mesh.lookupObject<cmpModel>(cmpModel::propertiesName);
return calcBlendingFactor
(
vf, model.nut(), model.mu()/model.rho(), model.U(), delta
);
}
FatalErrorInFunction
<< "Scheme requires a turbulence model to be present. "

View File

@ -1,19 +1,9 @@
EXE_INC = \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude \
-I$(LIB_SRC)/transportModels \
-I$(LIB_SRC)/transportModels/compressible/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/incompressible/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude
-I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude
LIB_LIBS = \
-lfiniteVolume \
-lmeshTools \
-lcompressibleTransportModels \
-lturbulenceModels \
-lincompressibleTurbulenceModels \
-lincompressibleTransportModels \
-lcompressibleTurbulenceModels \
-lfluidThermophysicalModels
-lturbulenceModels