32 template<
class BasicMomentumTransportModel>
38 scalar kMin = this->kMin_.value();
57 template<
class BasicMomentumTransportModel>
71 if (isA<wallFvPatch>(curPatch))
79 this->U_.boundaryField()[
patchi].snGrad()
85 this->mesh_.magSf().boundaryField()[
patchi];
91 = (Sf[facei]/magSf[facei])*snGradUw[facei];
103 template<
class BasicMomentumTransportModel>
120 template<
class BasicMomentumTransportModel>
123 const word& modelName,
132 BasicMomentumTransportModel
157 this->groupName(
"R"),
158 this->runTime_.
name(),
170 this->groupName(
"nut"),
171 this->runTime_.
name(),
179 if (couplingFactor_.value() < 0 || couplingFactor_.value() > 1)
182 <<
"couplingFactor = " << couplingFactor_
183 <<
" is not in range 0 - 1" <<
nl
191 template<
class BasicMomentumTransportModel>
198 template<
class BasicMomentumTransportModel>
206 template<
class BasicMomentumTransportModel>
211 tk.
ref().rename(
"k");
216 template<
class BasicMomentumTransportModel>
222 this->groupName(
"devTau"),
223 this->alpha_*this->rho_*R_
224 - (this->alpha_*this->rho_*this->nu())
230 template<
class BasicMomentumTransportModel>
231 template<
class RhoFieldType>
235 const RhoFieldType&
rho,
270 template<
class BasicMomentumTransportModel>
277 return DivDevRhoReff(this->rho_,
U);
281 template<
class BasicMomentumTransportModel>
289 return DivDevRhoReff(
rho,
U);
293 template<
class BasicMomentumTransportModel>
300 template<
class BasicMomentumTransportModel>
#define R(A, B, C, D, E, F, K, M)
#define forAll(list, i)
Loop across all elements in list.
Generic GeometricBoundaryField class.
Generic GeometricField class.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
BasicMomentumTransportModel::alphaField alphaField
void boundNormalStress(volSymmTensorField &R) const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual void validate()
Validate the turbulence fields after construction.
ReynoldsStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
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 tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volSymmTensorField > sigma() const
Return the Reynolds stress tensor [m^2/s^2].
virtual tmp< fvSymmTensorMatrix > RSource() const
Source term for the R equation.
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
void correctWallShearStress(volSymmTensorField &R) const
BasicMomentumTransportModel::rhoField rhoField
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
Abstract base class for all fluid physical properties.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the snGrad of the given volField.
const fvPatchList & patches
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
word name(const bool)
Return a word representation of a bool.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimTime
const dimensionSet dimVolume
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.