37 namespace laminarModels
42 template<
class BasicTurbulenceModel>
51 const word& propertiesName,
69 lambda_(
"lambda",
dimTime, this->coeffDict_),
76 this->runTime_.timeName(),
86 this->printCoeffs(
type);
93 template<
class BasicTurbulenceModel>
98 nuM_.readIfPresent(this->coeffDict());
99 lambda_.readIfPresent(this->coeffDict());
108 template<
class BasicTurbulenceModel>
116 template<
class BasicTurbulenceModel>
120 return devRhoReff(this->U_);
124 template<
class BasicTurbulenceModel>
135 this->runTime_.timeName(),
140 this->alpha_*this->rho_*sigma_
141 - (this->alpha_*this->rho_*this->
nu())
148 template<
class BasicTurbulenceModel>
159 this->alpha_*this->rho_*this->nuM_*
fvc::grad(
U)
161 +
fvc::div(this->alpha_*this->rho_*sigma_)
168 template<
class BasicTurbulenceModel>
189 template<
class BasicTurbulenceModel>
193 const alphaField&
alpha = this->alpha_;
194 const rhoField&
rho = this->rho_;
210 this->runTime_.constant(),
234 sigmaEqn.
ref().relax();
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Templated abstract base class for laminar transport models.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const dimensionSet dimViscosity
Ignore writing from objectRegistry::writeObject()
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
BasicTurbulenceModel::transportModel transportModel
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
BasicTurbulenceModel::rhoField rhoField
BasicTurbulenceModel::alphaField alphaField
A class for handling words, derived from Foam::string.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Maxwell(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual bool read()
Read model coefficients if they have changed.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and related properties. ...
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
virtual void correct()
Correct the laminar transport.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
A class for managing temporary objects.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.