51 mixedFvPatchScalarField(
p, iF,
dict, false),
52 restitutionCoefficient_
56 specularityCoefficient_
61 if (restitutionCoefficient_ < 0 || restitutionCoefficient_ > 1)
64 <<
"The restitution coefficient has to be between 0 and 1"
68 if (specularityCoefficient_ < 0 || specularityCoefficient_ > 1)
71 <<
"The specularity coefficient has to be between 0 and 1"
91 mixedFvPatchScalarField(ptf,
p, iF, mapper),
92 restitutionCoefficient_(ptf.restitutionCoefficient_),
93 specularityCoefficient_(ptf.specularityCoefficient_)
105 mixedFvPatchScalarField(ptf, iF),
106 restitutionCoefficient_(ptf.restitutionCoefficient_),
107 specularityCoefficient_(ptf.specularityCoefficient_)
132 patch().lookupPatchField<volScalarField, scalar>
134 phase.volScalarField::name()
140 patch().lookupPatchField<volVectorField, vector>
148 patch().lookupPatchField<volScalarField, scalar>
152 Foam::typedName<RASModels::kineticTheoryModel>(
"gs0"),
160 patch().lookupPatchField<volScalarField, scalar>
164 Foam::typedName<RASModels::kineticTheoryModel>(
"kappa"),
173 if (restitutionCoefficient_ != 1.0)
177 *specularityCoefficient_
179 /(scalar(1) -
sqr(restitutionCoefficient_));
181 this->refGrad() = 0.0;
188 *(scalar(1) -
sqr(restitutionCoefficient_))
190 /
max(4*
kappa*phase.alphaMax(), small)
193 this->valueFraction() =
c/(
c + patch().deltaCoeffs());
200 this->refValue() = 0.0;
205 *specularityCoefficient_
210 /
max(6*
kappa*phase.alphaMax(), small);
212 this->valueFraction() = 0;
215 mixedFvPatchScalarField::updateCoeffs();
225 writeEntry(os,
"restitutionCoefficient", restitutionCoefficient_);
226 writeEntry(os,
"specularityCoefficient", specularityCoefficient_);
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
static word groupName(Name name, const word &group)
Robin condition for the particulate granular temperature.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients.
JohnsonJacksonParticleThetaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base class for field mapping.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Class to represent a system of phases and model interfacial transfers between them.
static const word propertiesName
Default name of the phase properties dictionary.
const phaseModelList & phases() const
Return the phase models.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const char *const group
Group name for atomic constants.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
errorManip< error > abort(error &err)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const unitConversion unitFraction