35 incompressibleMomentumTransportModel
42 namespace incompressible
52 RASincompressibleMomentumTransportModel,
67 tmp<volScalarField> LamBremhorstKE::Rt()
const
73 tmp<volScalarField> LamBremhorstKE::fMu(
const volScalarField& Rt)
const
75 tmp<volScalarField>
Ry(
sqrt(
k_)*
y()/nu());
76 return sqr(scalar(1) -
exp(-0.0165*
Ry))*(scalar(1) + 20.5/(Rt + small));
80 tmp<volScalarField> LamBremhorstKE::f1(
const volScalarField& fMu)
const
82 return scalar(1) +
pow3(0.05/(fMu + small));
86 tmp<volScalarField> LamBremhorstKE::f2(
const volScalarField& Rt)
const
88 return scalar(1) -
exp(-
sqr(Rt));
111 const geometricOneField&
alpha,
112 const geometricOneField&
rho,
120 eddyViscosity<incompressible::RASModel>
172 this->groupName(
"k"),
185 this->groupName(
"epsilon"),
197 if (
type == typeName)
233 tmp<volTensorField> tgradU =
fvc::grad(U_);
244 tmp<fvScalarMatrix> epsEqn
254 epsEqn.ref().relax();
260 tmp<fvScalarMatrix> kEqn
makeMomentumTransportModelTypes(geometricOneField, geometricOneField, incompressibleMomentumTransportModel)
Bound the given scalar field where it is below the specified minimum.
void updateCoeffs()
Update the boundary condition coefficients.
void correctBoundaryConditions()
Correct boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static dimensioned< scalar > lookupOrAddToDict(const word &, dictionary &, const dimensionSet &dims=dimless, const scalar &defaultValue=pTraits< scalar >::zero)
Construct from dictionary, with default value.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
virtual void correctNut()
Correct the eddy-viscosity nut.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
dimensionedScalar sigmaEps_
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
LamBremhorstKE(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
virtual bool read()
Read RASProperties dictionary.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
dimensionedScalar exp(const dimensionedScalar &ds)
VolField< vector > volVectorField
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
SurfaceField< scalar > surfaceScalarField
tensor Ry(const scalar &omega)
Rotational transformation tensor about the y-axis by omega radians.
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.name(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating face flux\n"<< endl;surfaceScalarField phi(IOobject("phi", runTime.name(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar(mesh.Sf().dimensions() *U.dimensions(), 0));autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))