43template<
class BasicTurbulenceModel>
54 const scalar Cpm = pm.
dict().
get<scalar>(
C.name());
56 for (
const label zonei : cellZoneIDs)
60 for (
const label celli :
cells)
69template<
class BasicTurbulenceModel>
73 const porosityModels::powerLawLopesdaCosta& pm
76 if (pm.dict().found(
C.
name()))
78 const labelList& cellZoneIDs = pm.cellZoneIDs();
81 const scalar Cpm = pm.dict().
get<scalar>(
C.
name());
83 for (
const label zonei : cellZoneIDs)
89 const label celli =
cells[i];
90 C[celli] = Cpm*Sigma[i];
97template<
class BasicTurbulenceModel>
106 const fv::explicitPorositySource& eps =
111 const porosityModels::powerLawLopesdaCosta& pm =
117 setPorosityCoefficient(Cmu_, pm);
118 Cmu_.correctBoundaryConditions();
119 setPorosityCoefficient(C1_, pm);
120 setPorosityCoefficient(C2_, pm);
121 setPorosityCoefficient(sigmak_, pm);
122 setPorosityCoefficient(sigmaEps_, pm);
123 sigmak_.correctBoundaryConditions();
124 sigmaEps_.correctBoundaryConditions();
126 setCdSigma(CdSigma_, pm);
127 setPorosityCoefficient(betap_, pm);
128 setPorosityCoefficient(betad_, pm);
129 setPorosityCoefficient(C4_, pm);
130 setPorosityCoefficient(C5_, pm);
137template<
class BasicTurbulenceModel>
140 this->nut_ = Cmu_*
sqr(k_)/epsilon_;
141 this->nut_.correctBoundaryConditions();
144 BasicTurbulenceModel::correctNut();
148template<
class BasicTurbulenceModel>
159template<
class BasicTurbulenceModel>
170 *(C4_*betap_*epsilon_()/k_()*magU3 - C5_*betad_*magU*epsilon_()),
178template<
class BasicTurbulenceModel>
179kEpsilonLopesdaCosta<BasicTurbulenceModel>::kEpsilonLopesdaCosta
181 const alphaField&
alpha,
187 const word& propertiesName,
344 IOobject::groupName(
"k", alphaRhoPhi.group()),
356 IOobject::groupName(
"epsilon", alphaRhoPhi.group()),
370 this->printCoeffs(
type);
379template<
class BasicTurbulenceModel>
391template<
class BasicTurbulenceModel>
394 if (!this->turbulence_)
400 const alphaField&
alpha = this->alpha_;
401 const rhoField&
rho = this->rho_;
407 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
423 epsilon_.boundaryFieldRef().updateCoeffs();
425 epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
431 tmp<fvScalarMatrix> epsEqn
440 + epsilonSource(magU, magU3)
444 epsEqn.ref().relax();
446 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
449 bound(epsilon_, this->epsilonMin_);
452 tmp<fvScalarMatrix> kEqn
461 + kSource(magU, magU3)
469 bound(k_, this->kMin_);
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
Bound the given scalar field if it has gone unbounded.
Graphite solid properties.
DimensionedField< scalar, volMesh > Internal
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Templated abstract base class for RAS turbulence models.
BasicTurbulenceModel::alphaField alphaField
volScalarField::Internal C1_
volScalarField::Internal C5_
BasicTurbulenceModel::rhoField rhoField
virtual tmp< fvScalarMatrix > kSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
void setCdSigma(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
volScalarField::Internal betad_
volScalarField::Internal betap_
volScalarField::Internal C4_
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
void setPorosityCoefficients()
volScalarField::Internal CdSigma_
volScalarField::Internal C2_
virtual void correctNut()
BasicTurbulenceModel::transportModel transportModel
virtual tmp< fvScalarMatrix > epsilonSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
void setPorosityCoefficient(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
bool get(const label i) const
Return bool value at specified position, always false for out-of-range access.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Generic dimensioned Type class.
const word & name() const noexcept
Return const reference to name.
Eddy viscosity turbulence model base class.
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
virtual tmp< volScalarField > nut() const
Applies an explicit porosity source to the momentum equation within a specified region.
const porosityModel & model() const noexcept
Access to the porosityModel.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
optionList(const optionList &)=delete
No copy construct.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
const labelList & cellZoneIDs() const
Return const access to the cell zone IDs.
const dictionary & dict() const
Return dictionary used for model construction.
const scalarField & Sigma() const
Return the porosity surface area per unit volume zone field.
Variant of the power law porosity model with spatially varying drag coefficient.
A class for managing temporary objects.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
zeroField Su(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.