37 JohnsonJacksonParticleSlipFvPatchVectorField
48 const DimensionedField<vector, volMesh>& iF
51 partialSlipFvPatchVectorField(p, iF),
52 specularityCoefficient_(
"specularityCoefficient",
dimless, 0)
60 const DimensionedField<vector, volMesh>& iF,
61 const dictionary& dict
64 partialSlipFvPatchVectorField(p, iF),
65 specularityCoefficient_
67 "specularityCoefficient",
69 dict.
lookup(
"specularityCoefficient")
74 (specularityCoefficient_.value() < 0)
75 || (specularityCoefficient_.value() > 1)
79 <<
"The specularity coefficient has to be between 0 and 1" 93 const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
95 const DimensionedField<vector, volMesh>& iF,
96 const fvPatchFieldMapper& mapper
99 partialSlipFvPatchVectorField(ptf, p, iF, mapper),
100 specularityCoefficient_(ptf.specularityCoefficient_)
107 const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
108 const DimensionedField<vector, volMesh>& iF
111 partialSlipFvPatchVectorField(ptf, iF),
112 specularityCoefficient_(ptf.specularityCoefficient_)
120 const fvPatchFieldMapper& m
123 partialSlipFvPatchVectorField::autoMap(m);
133 partialSlipFvPatchVectorField::rmap(ptf, addr);
145 const phaseSystem& fluid =
146 db().lookupObject<phaseSystem>(
"phaseProperties");
148 const phaseModel& phase
150 fluid.phases()[internalField().group()]
156 patch().lookupPatchField<volScalarField, scalar>
158 phase.volScalarField::name()
164 patch().lookupPatchField<volScalarField, scalar>
172 patch().lookupPatchField<volScalarField, scalar>
180 patch().lookupPatchField<volScalarField, scalar>
190 db().foundObject<volScalarField>(ThetaName)
191 ? patch().lookupPatchField<volScalarField, scalar>(ThetaName)
201 .lookupObject<IOdictionary>
206 .subDict(
"kineticTheoryCoeffs")
216 *specularityCoefficient_.value()
221 this->valueFraction() =
c/(
c + patch().deltaCoeffs());
223 partialSlipFvPatchVectorField::updateCoeffs();
233 writeEntry(os,
"specularityCoefficient", specularityCoefficient_);
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
fvPatchField< vector > fvPatchVectorField
JohnsonJacksonParticleSlipFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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 > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
virtual void write(Ostream &) const
Write.
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
Macros for easy insertion into run-time selection tables.
stressControl lookup("compactNormalStress") >> compactNormalStress
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
List< label > labelList
A List of labels.
errorManip< error > abort(error &err)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void updateCoeffs()
Update the coefficients.
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)