26 #include "JohnsonJacksonParticleThetaFvPatchScalarField.H" 28 #include "twoPhaseSystem.H" 37 JohnsonJacksonParticleThetaFvPatchScalarField
47 const DimensionedField<scalar, volMesh>& iF
50 mixedFvPatchScalarField(p, iF),
51 restitutionCoefficient_(
"restitutionCoefficient",
dimless, 0),
52 specularityCoefficient_(
"specularityCoefficient",
dimless, 0)
59 const JohnsonJacksonParticleThetaFvPatchScalarField& ptf,
61 const DimensionedField<scalar, volMesh>& iF,
62 const fvPatchFieldMapper& mapper
65 mixedFvPatchScalarField(ptf, p, iF, mapper),
66 restitutionCoefficient_(ptf.restitutionCoefficient_),
67 specularityCoefficient_(ptf.specularityCoefficient_)
76 const DimensionedField<scalar, volMesh>& iF,
77 const dictionary& dict
80 mixedFvPatchScalarField(p, iF),
81 restitutionCoefficient_
83 "restitutionCoefficient",
85 dict.
lookup(
"restitutionCoefficient")
87 specularityCoefficient_
89 "specularityCoefficient",
91 dict.
lookup(
"specularityCoefficient")
96 (restitutionCoefficient_.value() < 0)
97 || (restitutionCoefficient_.value() > 1)
101 <<
"The restitution coefficient has to be between 0 and 1" 107 (specularityCoefficient_.
value() < 0)
108 || (specularityCoefficient_.
value() > 1)
112 <<
"The specularity coefficient has to be between 0 and 1" 126 const JohnsonJacksonParticleThetaFvPatchScalarField& ptf
129 mixedFvPatchScalarField(ptf),
130 restitutionCoefficient_(ptf.restitutionCoefficient_),
131 specularityCoefficient_(ptf.specularityCoefficient_)
138 const JohnsonJacksonParticleThetaFvPatchScalarField& ptf,
139 const DimensionedField<scalar, volMesh>& iF
142 mixedFvPatchScalarField(ptf, iF),
143 restitutionCoefficient_(ptf.restitutionCoefficient_),
144 specularityCoefficient_(ptf.specularityCoefficient_)
152 const fvPatchFieldMapper& m
155 mixedFvPatchScalarField::autoMap(m);
165 mixedFvPatchScalarField::rmap(ptf, addr);
177 const twoPhaseSystem& fluid = db().lookupObject<twoPhaseSystem>
182 const phaseModel& phased
184 fluid.phase1().name() == internalField().group()
192 patch().lookupPatchField<volScalarField, scalar>
194 phased.volScalarField::name()
200 patch().lookupPatchField<volVectorField, vector>
208 patch().lookupPatchField<volScalarField, scalar>
216 patch().lookupPatchField<volScalarField, scalar>
230 .lookupObject<IOdictionary>
235 .subDict(
"kineticTheoryCoeffs")
240 if (restitutionCoefficient_.value() != 1.0)
244 *specularityCoefficient_.
value()
246 /(scalar(1) -
sqr(restitutionCoefficient_.value()));
248 this->refGrad() = 0.0;
255 *(scalar(1) -
sqr(restitutionCoefficient_.value()))
260 this->valueFraction() =
c/(
c + patch().deltaCoeffs());
267 this->refValue() = 0.0;
272 *specularityCoefficient_.
value()
279 this->valueFraction() = 0.0;
282 mixedFvPatchScalarField::updateCoeffs();
292 os.writeKeyword(
"restitutionCoefficient")
294 os.writeKeyword(
"specularityCoefficient")
296 writeEntry(
"value", os);
fvPatchField< vector > fvPatchVectorField
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void write(Ostream &) const
Write.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Macros for easy insertion into run-time selection tables.
JohnsonJacksonParticleThetaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual void write(Ostream &) const
Write.
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
const Type & value() const
Return const reference to value.
List< label > labelList
A List of labels.
errorManip< error > abort(error &err)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionedScalar c
Speed of light in a vacuum.
virtual void updateCoeffs()
Update the coefficients.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)