The sorptionWallFunction is a wall boundary condition to specify scalar/concentration gradient for turbulent and laminar flows.
More...
#include <sorptionWallFunctionFvPatchScalarField.H>


Public Member Functions | |
| TypeName ("sorptionWallFunction") | |
| Runtime type information. | |
| sorptionWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
| Construct from patch and internal field. | |
| sorptionWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
| Construct from patch, internal field and dictionary. | |
| sorptionWallFunctionFvPatchScalarField (const sorptionWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
| Construct by mapping given sorptionWallFunctionFvPatchScalarField onto a new patch. | |
| sorptionWallFunctionFvPatchScalarField (const sorptionWallFunctionFvPatchScalarField &) | |
| Construct as copy. | |
| sorptionWallFunctionFvPatchScalarField (const sorptionWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
| Construct as copy setting internal field reference. | |
| virtual tmp< fvPatchField< scalar > > | clone () const |
| Return a clone. | |
| virtual tmp< fvPatchField< scalar > > | clone (const DimensionedField< scalar, volMesh > &iF) const |
| Clone with an internal field reference. | |
| virtual void | autoMap (const fvPatchFieldMapper &) |
| Map (and resize as needed) from self given a mapping object. | |
| 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. | |
| virtual void | write (Ostream &) const |
| Write. | |
The sorptionWallFunction is a wall boundary condition to specify scalar/concentration gradient for turbulent and laminar flows.
The governing equation of the boundary condition is:
![\[ \nabla C = \frac{C^* - C_p}{\Delta_y} = \frac{F}{a \Delta_y}
\]](form_1004.png)
with
![\[ C^* = \frac{C_{surf}}{K}
\]](form_1005.png)
and with the mass-transfer coefficient is calculated for turbulent flows
![\[ a = \frac{C_\mu^{0.25} k_p^{0.5}}{y^+_{blended}}
\]](form_1006.png)
or for laminar-flow and molecular-diffusion-only states
![\[ a = \frac{D_m}{y_1}
\]](form_1007.png)
where
![]() | = | Gradient of concentration |
![]() | = | Wall-adjacent concentration |
![]() | = | Near-wall cell concentration |
![]() | = | First-cell centre wall distance |
![]() | = | Flux of concentration |
![]() | = | Mass-transfer coefficient |
![]() | = | Wall-surface concentration |
![]() | = | Adsorption or absorption/permeation coefficient |
![]() | = | Empirical model coefficient |
![]() | = | Turbulent kinetic energy in near-wall cell |
![]() | = | Non-dimensional blended near-wall cell height |
![]() | = | Molecular-diffusion coefficient |
![]() | = | First-cell centre wall distance |
Required fields:
x | Arbitrary scalar field, e.g. species, passive scalars etc.
Reference:
Standard model for exponential blending (tag:FDC):
Foat, T., Drodge, J., Charleson, A., Whatmore, B.,
Pownall, S., Glover, P., ... & Marr, A. (2022).
Predicting vapour transport from semi-volatile organic
compounds concealed within permeable packaging.
International Journal of Heat and Mass Transfer, 183, 122012.
DOI:10.1016/j.ijheatmasstransfer.2021.122012
Standard model for stepwise blending (tag:F):
Foat, T. (2021).
Modelling vapour transport in indoor environments for
improved detection of explosives using dogs.
Doctoral dissertation. University of Southampton.
URI:http://eprints.soton.ac.uk/id/eprint/456709
<patchName>
{
// Mandatory entries
type sorptionWallFunction;
Sc <scalar>;
Sct <scalar>;
kAbs <PatchFunction1<scalar>>;
// Optional entries
laminar <bool>;
D <scalar>;
kName <word>;
nuName <word>;
// Inherited entries
Cmu <scalar>;
kappa <scalar>;
E <scalar>;
blending <word>;
...
}
where the entries mean:
| Property | Description | Type | Reqd | Deflt |
|---|---|---|---|---|
type | Type name: sorptionWallFunction | word | yes | - |
Sc | Schmidt number | scalar | yes | - |
Sct | Turbulent Schmidt number | scalar | yes | - |
kAbs | Adsorption or absorption/permeation coefficient | PatchFunction1<scalar> | yes | - |
laminar | Flag to calculate mass-transfer coefficient under the laminar-flow or molecular-diffusion-only states | bool | no | false |
kName | Name of operand turbulent kinetic energy field | word | no | k |
nuName | Name of operand kinematic viscosity field | word | no | nu |
The inherited entries are elaborated in:
Definition at line 271 of file sorptionWallFunctionFvPatchScalarField.H.
| sorptionWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF ) |
Construct from patch and internal field.
Definition at line 262 of file sorptionWallFunctionFvPatchScalarField.C.
References p, and wallFunctionBlenders::wallFunctionBlenders().
Referenced by sorptionWallFunctionFvPatchScalarField(), sorptionWallFunctionFvPatchScalarField(), sorptionWallFunctionFvPatchScalarField(), and TypeName().


| sorptionWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const dictionary & | dict ) |
Construct from patch, internal field and dictionary.
Definition at line 302 of file sorptionWallFunctionFvPatchScalarField.C.
References dict, Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, Foam::New(), p, wallFunctionBlenders::STEPWISE, wallFunctionBlenders::wallFunctionBlenders(), and Foam::Zero.

| sorptionWallFunctionFvPatchScalarField | ( | const sorptionWallFunctionFvPatchScalarField & | ptf, |
| const fvPatch & | p, | ||
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const fvPatchFieldMapper & | mapper ) |
Construct by mapping given sorptionWallFunctionFvPatchScalarField onto a new patch.
Definition at line 281 of file sorptionWallFunctionFvPatchScalarField.C.
References clone(), p, sorptionWallFunctionFvPatchScalarField(), and wallFunctionBlenders::wallFunctionBlenders().

| sorptionWallFunctionFvPatchScalarField | ( | const sorptionWallFunctionFvPatchScalarField & | swfpsf | ) |
Construct as copy.
Definition at line 346 of file sorptionWallFunctionFvPatchScalarField.C.
References clone(), sorptionWallFunctionFvPatchScalarField(), and wallFunctionBlenders::wallFunctionBlenders().

| sorptionWallFunctionFvPatchScalarField | ( | const sorptionWallFunctionFvPatchScalarField & | swfpsf, |
| const DimensionedField< scalar, volMesh > & | iF ) |
Construct as copy setting internal field reference.
Definition at line 364 of file sorptionWallFunctionFvPatchScalarField.C.
References clone(), sorptionWallFunctionFvPatchScalarField(), and wallFunctionBlenders::wallFunctionBlenders().

| TypeName | ( | "sorptionWallFunction" | ) |
Runtime type information.
References sorptionWallFunctionFvPatchScalarField().

|
inlinevirtual |
Return a clone.
Definition at line 399 of file sorptionWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().
Referenced by sorptionWallFunctionFvPatchScalarField(), sorptionWallFunctionFvPatchScalarField(), and sorptionWallFunctionFvPatchScalarField().


|
inlinevirtual |
Clone with an internal field reference.
Definition at line 407 of file sorptionWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().

|
virtual |
Map (and resize as needed) from self given a mapping object.
Definition at line 385 of file sorptionWallFunctionFvPatchScalarField.C.
|
virtual |
Reverse map the given fvPatchField onto this fvPatchField.
Definition at line 399 of file sorptionWallFunctionFvPatchScalarField.C.
References Foam::refCast().

|
virtual |
Update the coefficients associated with the patch field.
Definition at line 417 of file sorptionWallFunctionFvPatchScalarField.C.
|
virtual |
Write.
Definition at line 430 of file sorptionWallFunctionFvPatchScalarField.C.
References os(), fixedGradientFvPatchField< Type >::write(), and fvPatchField< Type >::writeValueEntry().
