37 template<
class BasicMomentumTransportModel>
38 tmp<volScalarField::Internal>
39 SpalartAllmarasIDDES<BasicMomentumTransportModel>::IDDESalpha()
const 44 max(0.25 - this->y_()/IDDESDelta_.hmax(), scalar(-5))
49 template<
class BasicMomentumTransportModel>
50 tmp<volScalarField::Internal>
51 SpalartAllmarasIDDES<BasicMomentumTransportModel>::ft
64 template<
class BasicMomentumTransportModel>
65 tmp<volScalarField::Internal>
66 SpalartAllmarasIDDES<BasicMomentumTransportModel>::fl
74 tanh(
pow(
sqr(cl_)*rd(this->nu(), magGradU), 10))
79 template<
class BasicMomentumTransportModel>
80 tmp<volScalarField::Internal>
81 SpalartAllmarasIDDES<BasicMomentumTransportModel>::rd
98 )*
sqr(this->kappa_*this->y_())
106 template<
class BasicMomentumTransportModel>
107 tmp<volScalarField::Internal>
108 SpalartAllmarasIDDES<BasicMomentumTransportModel>::fd
116 1 -
tanh(
pow3(8*rd(this->nuEff(), magGradU)))
123 template<
class BasicMomentumTransportModel>
124 tmp<volScalarField::Internal>
136 modelName(
"expTerm"),
147 2*(
pos0(alpha)*
pow(expTerm, -11.09) +
neg(alpha)*
pow(expTerm, -9.0))
156 min(2*
pow(expTerm, -9.0), scalar(1))
163 max(1 - fd(magGradU), fStep)
171 1 -
max(ft(magGradU), fl(magGradU))
179 modelName(
"fRestore"),
180 max(fHill - 1, scalar(0))*fAmp
195 - this->Cb1_*this->fv2(chi, fv1)
196 /(this->Cw1_*
sqr(this->kappa_)*fwStar_)
208 fHyb*(1 + fRestore*Psi)*this->y_
209 + (1 - fHyb)*this->CDES_*Psi*this->
delta()
217 template<
class BasicMomentumTransportModel>
225 const viscosity& viscosity,
265 IDDESDelta_(refCast<IDDESDelta>(this->delta_()))
271 template<
class BasicMomentumTransportModel>
276 fwStar_.readIfPresent(this->coeffDict());
277 cl_.readIfPresent(this->coeffDict());
278 ct_.readIfPresent(this->coeffDict());
dimensionedScalar tanh(const dimensionedScalar &ds)
BasicMomentumTransportModel::rhoField rhoField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Generic dimensioned Type class.
SpalartAllmarasIDDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
virtual bool read()
Read model coefficients if they have changed.
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimLength
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from string.
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
BasicMomentumTransportModel::alphaField alphaField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual tmp< volScalarField::Internal > dTilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volTensorField::Internal &gradU) const
Length scale.
A class for managing temporary objects.