30#include "twoPhaseSystem.H"
39 JohnsonJacksonParticleSlipFvPatchVectorField
46Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
47JohnsonJacksonParticleSlipFvPatchVectorField
50 const DimensionedField<vector, volMesh>& iF
53 partialSlipFvPatchVectorField(
p, iF),
54 specularityCoefficient_(
"specularityCoefficient",
dimless, 0)
58Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
59JohnsonJacksonParticleSlipFvPatchVectorField
61 const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
63 const DimensionedField<vector, volMesh>& iF,
64 const fvPatchFieldMapper& mapper
67 partialSlipFvPatchVectorField(ptf,
p, iF, mapper),
68 specularityCoefficient_(ptf.specularityCoefficient_)
72Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
73JohnsonJacksonParticleSlipFvPatchVectorField
76 const DimensionedField<vector, volMesh>& iF,
77 const dictionary&
dict
80 partialSlipFvPatchVectorField(
p, iF),
81 specularityCoefficient_(
"specularityCoefficient",
dimless,
dict)
85 (specularityCoefficient_.value() < 0)
86 || (specularityCoefficient_.value() > 1)
90 <<
"The specularity coefficient has to be between 0 and 1"
94 fvPatchVectorField::operator=
101Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
102JohnsonJacksonParticleSlipFvPatchVectorField
104 const JohnsonJacksonParticleSlipFvPatchVectorField& ptf
107 partialSlipFvPatchVectorField(ptf),
108 specularityCoefficient_(ptf.specularityCoefficient_)
112Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
113JohnsonJacksonParticleSlipFvPatchVectorField
115 const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
116 const DimensionedField<vector, volMesh>& iF
119 partialSlipFvPatchVectorField(ptf, iF),
120 specularityCoefficient_(ptf.specularityCoefficient_)
128 const fvPatchFieldMapper& m
131 partialSlipFvPatchVectorField::autoMap(m);
141 partialSlipFvPatchVectorField::rmap(ptf, addr);
153 const twoPhaseSystem&
fluid = db().lookupObject<twoPhaseSystem>
158 const phaseModel& phased
168 patch().lookupPatchField<volScalarField, scalar>
170 phased.volScalarField::name()
176 patch().lookupPatchField<volScalarField, scalar>
184 patch().lookupPatchField<volScalarField, scalar>
192 patch().lookupPatchField<volScalarField, scalar>
202 db().foundObject<volScalarField>(ThetaName)
203 ?
patch().lookupPatchField<volScalarField, scalar>(ThetaName)
213 .lookupObject<IOdictionary>
218 .subDict(
"kineticTheoryCoeffs")
227 *specularityCoefficient_.value()
232 this->valueFraction() =
c/(
c +
patch().deltaCoeffs());
234 partialSlipFvPatchVectorField::updateCoeffs();
244 os.
writeEntry(
"specularityCoefficient", specularityCoefficient_);
245 writeEntry(
"value",
os);
Macros for easy insertion into run-time selection tables.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
const Type & value() const
Return const reference to value.
virtual bool write()
Write the output fields.
const word & name() const
const phaseModel & phase1() const
Return phase model 1.
const phaseModel & phase2() const
Return phase model 2.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
constexpr const char *const group
Group name for atomic constants.
constexpr scalar pi(M_PI)
const dimensionedScalar c
Speed of light in a vacuum.
const std::string patch
OpenFOAM patch number as a std::string.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)