41 const DimensionedField<vector, volMesh>& iF
44 partialSlipFvPatchVectorField(
p, iF),
47 psiName_(
"thermo:psi"),
50 accommodationCoeff_(1.0),
51 Uwall_(
p.size(), Zero),
59 const maxwellSlipUFvPatchVectorField& mspvf,
61 const DimensionedField<vector, volMesh>& iF,
62 const fvPatchFieldMapper& mapper
65 partialSlipFvPatchVectorField(mspvf,
p, iF, mapper),
67 rhoName_(mspvf.rhoName_),
68 psiName_(mspvf.psiName_),
69 muName_(mspvf.muName_),
70 tauMCName_(mspvf.tauMCName_),
71 accommodationCoeff_(mspvf.accommodationCoeff_),
73 thermalCreep_(mspvf.thermalCreep_),
74 curvature_(mspvf.curvature_)
81 const DimensionedField<vector, volMesh>& iF,
85 partialSlipFvPatchVectorField(
p, iF),
86 TName_(
dict.getOrDefault<word>(
"T",
"T")),
87 rhoName_(
dict.getOrDefault<word>(
"rho",
"rho")),
88 psiName_(
dict.getOrDefault<word>(
"psi",
"thermo:psi")),
89 muName_(
dict.getOrDefault<word>(
"mu",
"thermo:mu")),
90 tauMCName_(
dict.getOrDefault<word>(
"tauMC",
"tauMC")),
91 accommodationCoeff_(
dict.
get<scalar>(
"accommodationCoeff")),
92 Uwall_(
"Uwall",
dict,
p.size()),
93 thermalCreep_(
dict.getOrDefault(
"thermalCreep", true)),
94 curvature_(
dict.getOrDefault(
"curvature", true))
98 mag(accommodationCoeff_) < SMALL
99 ||
mag(accommodationCoeff_) > 2.0
103 <<
"unphysical accommodationCoeff_ specified"
104 <<
"(0 < accommodationCoeff_ <= 1)" <<
endl
105 <<
exit(FatalIOError);
108 if (this->readValueEntry(
dict))
110 const auto* hasRefValue =
dict.findEntry(
"refValue", keyType::LITERAL);
111 const auto* hasFrac =
dict.findEntry(
"valueFraction", keyType::LITERAL);
113 if (hasRefValue && hasFrac)
115 this->refValue().assign(*hasRefValue,
p.size());
116 this->valueFraction().assign(*hasFrac,
p.size());
120 this->refValue() = *
this;
121 this->valueFraction() = scalar(1);
129 const maxwellSlipUFvPatchVectorField& mspvf,
130 const DimensionedField<vector, volMesh>& iF
133 partialSlipFvPatchVectorField(mspvf, iF),
134 TName_(mspvf.TName_),
135 rhoName_(mspvf.rhoName_),
136 psiName_(mspvf.psiName_),
137 muName_(mspvf.muName_),
138 tauMCName_(mspvf.tauMCName_),
139 accommodationCoeff_(mspvf.accommodationCoeff_),
140 Uwall_(mspvf.Uwall_),
141 thermalCreep_(mspvf.thermalCreep_),
142 curvature_(mspvf.curvature_)
162 * (2.0 - accommodationCoeff_)/accommodationCoeff_
165 Field<scalar> pnu(pmu/prho);
166 valueFraction() = (1.0/(1.0 +
patch().deltaCoeffs()*C1*pnu));
174 label patchi = this->
patch().index();
176 Field<vector> gradpT(
fvc::grad(vsfT)().boundaryField()[patchi]);
179 refValue() -= 3.0*pnu/(4.0*pT)*
transform(
I -
n*
n, gradpT);
191 partialSlipFvPatchVectorField::updateCoeffs();
204 os.
writeEntry(
"accommodationCoeff", accommodationCoeff_);
Macros for easy insertion into run-time selection tables.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Maxwell slip boundary condition including thermal creep and surface curvature terms that can be optio...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
maxwellSlipUFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
A class for handling words, derived from Foam::string.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Calculate the gradient of the given field.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
constexpr scalar piByTwo(0.5 *M_PI)
const std::string patch
OpenFOAM patch number as a std::string.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
static const Identity< scalar > I
Ostream & endl(Ostream &os)
Add newline and flush stream.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField