53 for (
int iter = 0; iter < 10; ++iter)
56 ypl = (
log(
max(E*ypl, scalar(1)))/kappa + Pc)*Sct/Sc;
67 const label patchi = patch().index();
71 9.24*(
pow(Sc_/Sct_, 0.75) - 1)*(1 + 0.28*
exp(-0.007*Sc_/Sct_));
72 const scalar Cmu25 =
pow025(wallCoeffs_.
Cmu());
73 const scalar kappa = wallCoeffs_.
kappa();
74 const scalar E = wallCoeffs_.
E();
90 const auto yStar = [&](
const label facei) -> scalar
94 Cmu25*
sqrt(kwc[facei])*ywc[facei]/nuwc[facei]
99 const auto yPlusVis = [&](
const label facei) -> scalar
101 return (Sc_*yStar(facei));
105 const auto yPlusLog = [&](
const label facei) -> scalar
109 Sct_*(
log(
max(E*yStar(facei), 1 + 1
e-4))/kappa + Pc)
114 auto&
yPlus = tyPlus.ref();
118 case blenderType::EXPONENTIAL:
123 const scalar yStarFace = yStar(facei);
125 0.01*
pow4(Sc_*yStarFace)/(1 + 5*
pow3(Sc_)*yStarFace);
126 const scalar invGamma = scalar(1)/
max(Gamma, ROOTVSMALL);
131 yPlusVis(facei)*
exp(-Gamma)
132 + yPlusLog(facei)*
exp(-invGamma)
138 case blenderType::STEPWISE:
140 static const scalar yStarLam =
146 if (yStar(facei) < yStarLam)
148 yPlus[facei] = yPlusVis(facei);
152 yPlus[facei] = yPlusLog(facei);
158 case blenderType::BINOMIAL:
165 pow(yPlusVis(facei),
n_) +
pow(yPlusLog(facei),
n_),
172 case blenderType::MAX:
176 yPlus[facei] =
max(yPlusVis(facei), yPlusLog(facei));
181 case blenderType::TANH:
185 const scalar yPlusVisFace = yPlusVis(facei);
186 const scalar yPlusLogFace = yPlusLog(facei);
187 const scalar b1 = yPlusVisFace + yPlusLogFace;
191 pow(yPlusVisFace, 1.2) +
pow(yPlusLogFace, 1.2),
194 const scalar phiTanh =
tanh(
pow4(0.1*yStar(facei)));
196 yPlus[facei] = phiTanh*b1 + (1 - phiTanh)*b2;
209 const label patchi = patch().index();
226 const auto& Csurf = *
this;
229 const scalar t = db().time().timeOutputValue();
234 const word& fldName = internalField().name();
240 return (tCstar - tCwc)/ta;
244void sorptionWallFunctionFvPatchScalarField::writeLocalEntries
250 os.writeEntryIfDifferent<
bool>(
"laminar",
false, laminar_);
251 os.writeEntry(
"Sc", Sc_);
252 os.writeEntry(
"Sct", Sct_);
253 os.writeEntryIfDifferent<scalar>(
"D", -1, D_);
255 os.writeEntryIfDifferent<
word>(
"k",
"k", kName_);
256 os.writeEntryIfDifferent<
word>(
"nu",
"nu", nuName_);
260 kAbsPtr_->writeData(
os);
273 fixedGradientFvPatchScalarField(
p, iF),
294 fixedGradientFvPatchScalarField(ptf,
p, iF, mapper),
296 laminar_(ptf.laminar_),
297 kAbsPtr_(ptf.kAbsPtr_.clone(patch().patch())),
302 nuName_(ptf.nuName_),
303 wallCoeffs_(ptf.wallCoeffs_)
314 fixedGradientFvPatchScalarField(
p, iF),
316 laminar_(
dict.getOrDefault<bool>(
"laminar", false)),
320 D_(
dict.getOrDefault<scalar>(
"D", -1)),
321 kName_(
dict.getOrDefault<
word>(
"k",
"k")),
322 nuName_(
dict.getOrDefault<
word>(
"nu",
"nu")),
330 <<
"Molecular diffusion coefficient cannot be non-positive. "
339 <<
"Adsorption or absorption coefficient is not set."
343 if (!this->readGradientEntry(
dict) || !this->readValueEntry(
dict))
345 extrapolateInternal();
353 const sorptionWallFunctionFvPatchScalarField& swfpsf
356 fixedGradientFvPatchScalarField(swfpsf),
357 wallFunctionBlenders(swfpsf),
358 laminar_(swfpsf.laminar_),
359 kAbsPtr_(swfpsf.kAbsPtr_.clone(patch().patch())),
363 kName_(swfpsf.kName_),
364 nuName_(swfpsf.nuName_),
365 wallCoeffs_(swfpsf.wallCoeffs_)
375 fixedGradientFvPatchScalarField(swfpsf, iF),
377 laminar_(swfpsf.laminar_),
378 kAbsPtr_(swfpsf.kAbsPtr_.clone(patch().patch())),
382 kName_(swfpsf.kName_),
383 nuName_(swfpsf.nuName_),
384 wallCoeffs_(swfpsf.wallCoeffs_)
395 fixedGradientFvPatchScalarField::autoMap(mapper);
399 kAbsPtr_->autoMap(mapper);
410 fixedGradientFvPatchScalarField::rmap(ptf, addr);
417 kAbsPtr_->rmap(swfptf.kAbsPtr_(), addr);
429 gradient() = flux()/
patch().deltaCoeffs();
431 fixedGradientFvPatchScalarField::updateCoeffs();
439 writeLocalEntries(
os);
450 sorptionWallFunctionFvPatchScalarField
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual void write(Ostream &) const
Write.
A FieldMapper for finite-volume patch fields.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
The sorptionWallFunction is a wall boundary condition to specify scalar/concentration gradient for tu...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
sorptionWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for managing temporary objects.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field.
The class wallFunctionBlenders is a base class that hosts common entries for various derived wall-fun...
void writeEntries(Ostream &) const
Write wall-function blending data as dictionary entries.
wallFunctionBlenders()
Default construct with default coefficients.
blenderType
Options for the blending treatment of viscous and inertial sublayers.
@ STEPWISE
"Stepwise switch (discontinuous)"
enum blenderType blender_
Blending treatment.
scalar n_
Binomial blending exponent being used when blenderType is blenderType::BINOMIAL.
scalar kappa() const noexcept
Return the object: kappa.
scalar E() const noexcept
Return the object: E.
void writeEntries(Ostream &) const
Write wall-function coefficients as dictionary entries.
scalar Cmu() const noexcept
Return the object: Cmu.
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,...
const std::string patch
OpenFOAM patch number as a std::string.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar exp(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
static scalar calcYStarLam(const scalar kappa, const scalar E, const scalar Sc, const scalar Sct, const scalar Pc)
Estimate the y* at the intersection of the two sublayers.
errorManipArg< error, int > exit(error &err, const int errNo=1)
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar pow025(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.