Go to the documentation of this file.
30 #include "twoPhaseSystem.H"
39 JohnsonJacksonParticleSlipFvPatchVectorField
53 partialSlipFvPatchVectorField(
p, iF),
54 specularityCoefficient_(
"specularityCoefficient",
dimless,
Zero)
67 partialSlipFvPatchVectorField(ptf,
p, iF, mapper),
68 specularityCoefficient_(ptf.specularityCoefficient_)
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=
107 partialSlipFvPatchVectorField(ptf),
108 specularityCoefficient_(ptf.specularityCoefficient_)
119 partialSlipFvPatchVectorField(ptf, iF),
120 specularityCoefficient_(ptf.specularityCoefficient_)
131 partialSlipFvPatchVectorField::autoMap(m);
141 partialSlipFvPatchVectorField::rmap(ptf, addr);
160 fluid.phase1().name() == internalField().group()
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);
virtual void write(Ostream &) const
Write.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
JohnsonJacksonParticleSlipFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
A class for handling words, derived from Foam::string.
Class which solves the volume fraction equations for two phases.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
static constexpr const zero Zero
Global zero (0)
Partial slip boundary condition for the particulate velocity.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const Type & value() const
Return const reference to value.
Field< vector > vectorField
Specialisation of Field<T> for vector.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
errorManip< error > abort(error &err)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
fvPatchField< vector > fvPatchVectorField
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
constexpr scalar pi(M_PI)
const std::string patch
OpenFOAM patch number as a std::string.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
const dimensionedScalar c
Speed of light in a vacuum.
virtual void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void updateCoeffs()
Update the coefficients.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...