Go to the documentation of this file.
31 #include "twoPhaseSystem.H"
36 Foam::RASModels::kineticTheoryModel::kineticTheoryModel
43 const transportModel& phase,
44 const word& propertiesName,
67 kineticTheoryModels::viscosityModel::
New
74 kineticTheoryModels::conductivityModel::
New
81 kineticTheoryModels::radialModel::
New
86 granularPressureModel_
88 kineticTheoryModels::granularPressureModel::
New
93 frictionalStressModel_
95 kineticTheoryModels::frictionalStressModel::
New
101 equilibrium_(coeffDict_.
get<
bool>(
"equilibrium")),
103 alphaMax_(
"alphaMax",
dimless, coeffDict_),
104 alphaMinFriction_(
"alphaMinFriction",
dimless, coeffDict_),
105 residualAlpha_(
"residualAlpha",
dimless, coeffDict_),
112 IOobject::groupName(
"Theta", phase.
name()),
125 IOobject::groupName(
"lambda", phase.
name()),
139 IOobject::groupName(
"gs0", phase.
name()),
153 IOobject::groupName(
"kappa", phase.
name()),
167 IOobject::groupName(
"nuFric", phase.
name()),
177 if (
type == typeName)
198 RASModel<EddyDiffusivity<phaseCompressibleTurbulenceModel>>
202 coeffDict().readEntry(
"equilibrium", equilibrium_);
203 e_.readIfPresent(coeffDict());
204 alphaMax_.readIfPresent(coeffDict());
205 alphaMinFriction_.readIfPresent(coeffDict());
207 viscosityModel_->read();
208 conductivityModel_->read();
209 radialModel_->read();
210 granularPressureModel_->read();
211 frictionalStressModel_->read();
247 return tmp<volSymmTensorField>
271 tmp<volScalarField> tpPrime
274 *granularPressureModel_->granularPressureCoeffPrime
277 radialModel_->g0(alpha_, alphaMinFriction_, alphaMax_),
278 radialModel_->g0prime(alpha_, alphaMinFriction_, alphaMax_),
282 + frictionalStressModel_->frictionalPressurePrime
290 volScalarField::Boundary& bpPrime =
291 tpPrime.ref().boundaryFieldRef();
295 if (!bpPrime[patchi].
coupled())
297 bpPrime[patchi] == 0;
315 return devRhoReff(U_);
325 return tmp<volSymmTensorField>
367 const twoPhaseSystem&
fluid = refCast<const twoPhaseSystem>(phase_.fluid());
378 tmp<volScalarField> tda(phase_.d());
381 tmp<volTensorField> tgradU(
fvc::grad(U_));
386 gs0_ = radialModel_->g0(
alpha, alphaMinFriction_, alphaMax_);
391 nut_ = viscosityModel_->nu(
alpha, Theta_, gs0_,
rho, da, e_);
396 lambda_ = (4.0/3.0)*
sqr(
alpha)*da*gs0_*(1.0 + e_)*ThetaSqrt/sqrtPi;
401 rho*(2.0*nut_*
D + (lambda_ - (2.0/3.0)*nut_)*
tr(
D)*
I)
410 *
rho*gs0_*(1.0/da)*ThetaSqrt/sqrtPi
416 fluid.lookupSubModel<dragModel>
419 fluid.otherPhase(phase_)
431 *sqrtPi*(ThetaSqrt + ThetaSmallSqrt)
438 granularPressureModel_->granularPressureCoeff
448 kappa_ = conductivityModel_->kappa(
alpha, Theta_, gs0_,
rho, da, e_);
470 +
fvm::Sp(J2/(Theta_ + ThetaSmall), Theta_)
489 (sqrtPi/(3.0*(3.0 - e_)))
490 *(1.0 + 0.4*(1.0 + e_)*(3.0*e_ - 1.0)*
alpha*gs0_)
491 +1.6*
alpha*gs0_*(1.0 + e_)/sqrtPi
498 4.0*da*
rho*(1.0 + e_)*
alpha*gs0_/(3.0*sqrtPi) - 2.0*
K3/3.0
521 *(2.0*
K3*trD2 +
K2*tr2D)
530 kappa_ = conductivityModel_->kappa(
alpha, Theta_, gs0_,
rho, da, e_);
538 nut_ = viscosityModel_->nu(
alpha, Theta_, gs0_,
rho, da, e_);
543 lambda_ = (4.0/3.0)*
sqr(
alpha)*da*gs0_*(1.0 + e_)*ThetaSqrt/sqrtPi;
548 frictionalStressModel_->frictionalPressure
556 nuFric_ = frictionalStressModel_->nu
567 nuFric_ =
min(nuFric_, maxNut_ - nut_);
573 Info<< typeName <<
':' <<
nl
574 <<
" max(Theta) = " <<
max(Theta_).value() <<
nl
575 <<
" max(nut) = " <<
max(nut_).value() <<
endl;
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for managing temporary objects.
virtual tmp< volScalarField > k() const
static constexpr const zero Zero
static options & New(const fvMesh &mesh)
const dimensionedScalar alpha
bool read(const char *buf, int32_t &val)
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
virtual tmp< surfaceScalarField > pPrimef() const
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Ostream & endl(Ostream &os)
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
label min(const labelHashSet &set, label minValue=labelMax)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvMatrix< scalar > fvScalarMatrix
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
static const SymmTensor I
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label max(const labelHashSet &set, label maxValue=labelMin)
const dimensionSet dimViscosity
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
virtual tmp< volScalarField > epsilon() const
virtual tmp< volScalarField > pPrime() const
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< vector, fvPatchField, volMesh > volVectorField
virtual tmp< volSymmTensorField > devRhoReff() const
virtual ~kineticTheoryModel()
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
constexpr scalar pi(M_PI)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar sqrt(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
const dimensionedScalar & D
virtual void printCoeffs(const word &type)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
word name(const expressions::valueTypeCode typeCode)
virtual tmp< volScalarField > omega() const
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
static word groupName(StringType base, const word &group)
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
virtual tmp< volSymmTensorField > R() const
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
const dimensionSet dimless
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)