Go to the documentation of this file.
30 #include "twoPhaseSystem.H"
39 JohnsonJacksonParticleThetaFvPatchScalarField
52 mixedFvPatchScalarField(
p, iF),
53 restitutionCoefficient_(
"restitutionCoefficient",
dimless,
Zero),
54 specularityCoefficient_(
"specularityCoefficient",
dimless,
Zero)
67 mixedFvPatchScalarField(ptf,
p, iF, mapper),
68 restitutionCoefficient_(ptf.restitutionCoefficient_),
69 specularityCoefficient_(ptf.specularityCoefficient_)
82 mixedFvPatchScalarField(
p, iF),
83 restitutionCoefficient_(
"restitutionCoefficient",
dimless,
dict),
84 specularityCoefficient_(
"specularityCoefficient",
dimless,
dict)
88 (restitutionCoefficient_.
value() < 0)
89 || (restitutionCoefficient_.
value() > 1)
93 <<
"The restitution coefficient has to be between 0 and 1"
99 (specularityCoefficient_.
value() < 0)
100 || (specularityCoefficient_.
value() > 1)
104 <<
"The specularity coefficient has to be between 0 and 1"
121 mixedFvPatchScalarField(ptf),
122 restitutionCoefficient_(ptf.restitutionCoefficient_),
123 specularityCoefficient_(ptf.specularityCoefficient_)
134 mixedFvPatchScalarField(ptf, iF),
135 restitutionCoefficient_(ptf.restitutionCoefficient_),
136 specularityCoefficient_(ptf.specularityCoefficient_)
147 mixedFvPatchScalarField::autoMap(m);
157 mixedFvPatchScalarField::rmap(ptf, addr);
169 const twoPhaseSystem&
fluid = db().lookupObject<twoPhaseSystem>
174 const phaseModel& phased
184 patch().lookupPatchField<volScalarField, scalar>
186 phased.volScalarField::name()
192 patch().lookupPatchField<volVectorField, vector>
200 patch().lookupPatchField<volScalarField, scalar>
208 patch().lookupPatchField<volScalarField, scalar>
222 .lookupObject<IOdictionary>
227 .subDict(
"kineticTheoryCoeffs")
231 if (restitutionCoefficient_.value() != 1.0)
235 *specularityCoefficient_.value()
237 /(scalar(1) -
sqr(restitutionCoefficient_.value()));
239 this->refGrad() = 0.0;
246 *(scalar(1) -
sqr(restitutionCoefficient_.value()))
251 this->valueFraction() =
c/(
c +
patch().deltaCoeffs());
258 this->refValue() = 0.0;
263 *specularityCoefficient_.value()
270 this->valueFraction() = 0;
273 mixedFvPatchScalarField::updateCoeffs();
283 os.
writeEntry(
"restitutionCoefficient", restitutionCoefficient_);
284 os.
writeEntry(
"specularityCoefficient", specularityCoefficient_);
285 writeEntry(
"value",
os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
virtual void write(Ostream &) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Class which solves the volume fraction equations for two phases.
constexpr const char *const group
static constexpr const zero Zero
JohnsonJacksonParticleThetaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
const dimensionedScalar alpha
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
const Type & value() const
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionedScalar kappa
Robin condition for the particulate granular temperature.
const phaseModel & phase2() const
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Generic templated field type.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const phaseModel & phase1() const
virtual void write(Ostream &) const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
errorManip< error > abort(error &err)
virtual void autoMap(const fvPatchFieldMapper &)
const word & name() const
fvPatchField< vector > fvPatchVectorField
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr scalar pi(M_PI)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & writeEntry(const keyType &key, const T &value)
const dimensionedScalar c
Foam::fvPatchFieldMapper.
static word groupName(StringType base, const word &group)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void updateCoeffs()
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimless