Go to the documentation of this file.
38 const DimensionedField<scalar, volMesh>& iF
41 mixedFvPatchScalarField(
p, iF),
44 psiName_(
"thermo:psi"),
46 accommodationCoeff_(1.0),
47 Twall_(
p.size(), 0.0),
52 valueFraction() = 0.0;
58 const smoluchowskiJumpTFvPatchScalarField& ptf,
60 const DimensionedField<scalar, volMesh>& iF,
61 const fvPatchFieldMapper& mapper
64 mixedFvPatchScalarField(ptf,
p, iF, mapper),
66 rhoName_(ptf.rhoName_),
67 psiName_(ptf.psiName_),
69 accommodationCoeff_(ptf.accommodationCoeff_),
78 const DimensionedField<scalar, volMesh>& iF,
79 const dictionary&
dict
82 mixedFvPatchScalarField(
p, iF),
83 UName_(
dict.lookupOrDefault<word>(
"U",
"U")),
84 rhoName_(
dict.lookupOrDefault<word>(
"rho",
"rho")),
85 psiName_(
dict.lookupOrDefault<word>(
"psi",
"thermo:psi")),
86 muName_(
dict.lookupOrDefault<word>(
"mu",
"thermo:mu")),
88 Twall_(
"Twall",
dict,
p.size()),
89 gamma_(
dict.lookupOrDefault<scalar>(
"gamma", 1.4))
93 mag(accommodationCoeff_) < SMALL
94 ||
mag(accommodationCoeff_) > 2.0
100 ) <<
"unphysical accommodationCoeff specified"
101 <<
"(0 < accommodationCoeff <= 1)" <<
endl
105 if (
dict.found(
"value"))
119 valueFraction() = 0.0;
125 const smoluchowskiJumpTFvPatchScalarField& ptpsf,
126 const DimensionedField<scalar, volMesh>& iF
129 mixedFvPatchScalarField(ptpsf, iF),
130 accommodationCoeff_(ptpsf.accommodationCoeff_),
131 Twall_(ptpsf.Twall_),
141 const fvPatchFieldMapper& m
144 mixedFvPatchScalarField::autoMap(m);
151 const fvPatchField<scalar>& ptf,
171 const fvPatchField<scalar>& ppsi =
177 const dictionary& thermophysicalProperties =
184 thermophysicalProperties.subDict(
"mixture").subDict(
"transport")
192 *2.0*gamma_/
Pr.
value()/(gamma_ + 1.0)
193 *(2.0 - accommodationCoeff_)/accommodationCoeff_
196 Field<scalar> aCoeff(prho.snGrad() - prho/C2);
197 Field<scalar> KEbyRho(0.5*
magSqr(pU));
199 valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C2));
203 mixedFvPatchScalarField::updateCoeffs();
212 writeEntryIfDifferent<word>(os,
"U",
"U", UName_);
213 writeEntryIfDifferent<word>(os,
"rho",
"rho", rhoName_);
214 writeEntryIfDifferent<word>(os,
"psi",
"thermo:psi", psiName_);
215 writeEntryIfDifferent<word>(os,
"mu",
"thermo:mu", muName_);
217 os.writeKeyword(
"accommodationCoeff")
219 Twall_.writeEntry(
"Twall", os);
220 os.writeKeyword(
"gamma")
222 writeEntry(
"value", os);
233 smoluchowskiJumpTFvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
virtual void write(Ostream &) const
Write.
virtual void operator=(const UList< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
List< label > labelList
A List of labels.
static const word dictName
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
dimensioned< scalar > mag(const dimensioned< Type > &)
const scalar piByTwo(0.5 *pi)
smoluchowskiJumpTFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar Pr("Pr", dimless, laminarTransport)
virtual void write(Ostream &) const
Write.
Macros for easy insertion into run-time selection tables.
Vector< scalar > vector
A scalar version of the templated Vector.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
errorManipArg< error, int > exit(error &err, const int errNo=1)
GeometricField< vector, fvPatchField, volMesh > volVectorField
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
fvPatchField< vector > fvPatchVectorField
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
dimensionedScalar sqrt(const dimensionedScalar &ds)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
stressControl lookup("compactNormalStress") >> compactNormalStress