Loading...
Searching...
No Matches
sorptionWallFunctionFvPatchScalarField.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::sorptionWallFunctionFvPatchScalarField
28
29Group
30 grpWallFunctions
31
32Description
33 The \c sorptionWallFunction is a wall boundary condition to
34 specify scalar/concentration gradient for turbulent and laminar flows.
35
36 The governing equation of the boundary condition is:
37
38 \f[
39 \nabla C = \frac{C^* - C_p}{\Delta_y} = \frac{F}{a \Delta_y}
40 \f]
41
42 with
43
44 \f[
45 C^* = \frac{C_{surf}}{K}
46 \f]
47
48 and with the mass-transfer coefficient is calculated for turbulent flows
49
50 \f[
51 a = \frac{C_\mu^{0.25} k_p^{0.5}}{y^+_{blended}}
52 \f]
53
54 or for laminar-flow and molecular-diffusion-only states
55
56 \f[
57 a = \frac{D_m}{y_1}
58 \f]
59
60 where
61 \vartable
62 \nabla C | Gradient of concentration
63 C^* | Wall-adjacent concentration
64 C_p | Near-wall cell concentration
65 \Delta_y | First-cell centre wall distance
66 F | Flux of concentration
67 a | Mass-transfer coefficient
68 C_{surf} | Wall-surface concentration
69 K | Adsorption or absorption/permeation coefficient
70 C_\mu | Empirical model coefficient
71 k_p | Turbulent kinetic energy in near-wall cell
72 y^+_{blended} | Non-dimensional blended near-wall cell height
73 D_m | Molecular-diffusion coefficient
74 y_1 | First-cell centre wall distance
75 \endvartable
76
77 Required fields:
78 \verbatim
79 x | Arbitrary scalar field, e.g. species, passive scalars etc.
80 \endverbatim
81
82 Reference:
83 \verbatim
84 Standard model for exponential blending (tag:FDC):
85 Foat, T., Drodge, J., Charleson, A., Whatmore, B.,
86 Pownall, S., Glover, P., ... & Marr, A. (2022).
87 Predicting vapour transport from semi-volatile organic
88 compounds concealed within permeable packaging.
89 International Journal of Heat and Mass Transfer, 183, 122012.
90 DOI:10.1016/j.ijheatmasstransfer.2021.122012
91
92 Standard model for stepwise blending (tag:F):
93 Foat, T. (2021).
94 Modelling vapour transport in indoor environments for
95 improved detection of explosives using dogs.
96 Doctoral dissertation. University of Southampton.
97 URI:http://eprints.soton.ac.uk/id/eprint/456709
98 \endverbatim
99
100Usage
101 Example of the boundary condition specification:
102 \verbatim
103 <patchName>
104 {
105 // Mandatory entries
106 type sorptionWallFunction;
107 Sc <scalar>;
108 Sct <scalar>;
109 kAbs <PatchFunction1<scalar>>;
110
111 // Optional entries
112 laminar <bool>;
113 D <scalar>;
114 kName <word>;
115 nuName <word>;
116
117 // Inherited entries
118 Cmu <scalar>;
119 kappa <scalar>;
120 E <scalar>;
121 blending <word>;
122 ...
123 }
124 \endverbatim
125
126 where the entries mean:
127 \table
128 Property | Description | Type | Reqd | Deflt
129 type | Type name: sorptionWallFunction | word | yes | -
130 Sc | Schmidt number | scalar | yes | -
131 Sct | Turbulent Schmidt number | scalar | yes | -
132 kAbs | Adsorption or absorption/permeation coefficient <!--
133 --> | PatchFunction1<scalar> | yes | -
134 laminar | Flag to calculate mass-transfer coefficient under the <!--
135 --> laminar-flow or molecular-diffusion-only states <!--
136 --> | bool | no | false
137 kName | Name of operand turbulent kinetic energy field | word | no | k
138 nuName | Name of operand kinematic viscosity field | word | no | nu
139 \endtable
140
141 The inherited entries are elaborated in:
142 - \link fixedGradientFvPatchField.H \endlink
143 - \link wallFunctionCoefficients.H \endlink
144 - \link wallFunctionBlenders.H \endlink
145 - \link PatchFunction1.H \endlink
146
147SourceFiles
148 sorptionWallFunctionFvPatchScalarField.C
149
150\*---------------------------------------------------------------------------*/
151
152#ifndef Foam_sorptionWallFunctionFvPatchScalarFields_H
153#define Foam_sorptionWallFunctionFvPatchScalarFields_H
154
155#include "fvPatchFields.H"
158#include "wallFunctionBlenders.H"
159#include "PatchFunction1.H"
160
161// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162
163namespace Foam
164{
165
166/*---------------------------------------------------------------------------*\
167 Class sorptionWallFunctionFvPatchScalarField Declaration
168\*---------------------------------------------------------------------------*/
169
171:
172 public fixedGradientFvPatchScalarField,
174{
175 //- Flag to calculate mass-transfer coefficient
176 //- under the laminar-flow or molecular-diffusion-only states
177 bool laminar_;
178
179 //- Adsorption or absorption/permeation coefficient
180 autoPtr<PatchFunction1<scalar>> kAbsPtr_;
181
182 //- Schmidt number
183 scalar Sc_;
184
185 //- Turbulent Schmidt number
186 scalar Sct_;
187
188 //- Molecular diffusion coefficient
189 scalar D_;
190
191 //- Name of operand turbulent kinetic energy field
192 word kName_;
193
194 //- Name of operand kinematic viscosity field
195 word nuName_;
196
197 //- Standard wall-function coefficients
198 wallFunctionCoefficients wallCoeffs_;
199
200
201 // Private Member Functions
202
203 //- Return the non-dimensional near-wall cell height field
204 //- blended between the viscous and inertial sublayers
205 tmp<scalarField> yPlus() const;
206
207 //- Return the flux normalised by the mass-transfer coefficient
208 tmp<scalarField> flux() const;
209
210 //- Write local wall-function variables
211 void writeLocalEntries(Ostream&) const;
212
213
214public:
215
216 //- Runtime type information
217 TypeName("sorptionWallFunction");
218
219
220 // Constructors
221
222 //- Construct from patch and internal field
224 (
225 const fvPatch&,
226 const DimensionedField<scalar, volMesh>&
227 );
228
229 //- Construct from patch, internal field and dictionary
231 (
232 const fvPatch&,
233 const DimensionedField<scalar, volMesh>&,
234 const dictionary&
235 );
236
237 //- Construct by mapping given
238 //- sorptionWallFunctionFvPatchScalarField onto
239 //- a new patch
241 (
243 const fvPatch&,
244 const DimensionedField<scalar, volMesh>&,
245 const fvPatchFieldMapper&
246 );
247
248 //- Construct as copy
250 (
252 );
253
254 //- Construct as copy setting internal field reference
256 (
258 const DimensionedField<scalar, volMesh>&
259 );
260
261 //- Return a clone
262 virtual tmp<fvPatchField<scalar>> clone() const
263 {
264 return fvPatchField<scalar>::Clone(*this);
265 }
266
267 //- Clone with an internal field reference
268 virtual tmp<fvPatchField<scalar>> clone
269 (
270 const DimensionedField<scalar, volMesh>& iF
271 ) const
272 {
273 return fvPatchField<scalar>::Clone(*this, iF);
274 }
275
276
277 // Member Functions
278
279 // Mapping
280
281 //- Map (and resize as needed) from self given a mapping object
282 virtual void autoMap(const fvPatchFieldMapper&);
283
284 //- Reverse map the given fvPatchField onto this fvPatchField
285 virtual void rmap
286 (
287 const fvPatchScalarField&,
288 const labelList&
289 );
290
291
292 // Evaluation
293
294 //- Update the coefficients associated with the patch field
295 virtual void updateCoeffs();
296
297
298 // I-O
299
300 //- Write
301 virtual void write(Ostream&) const;
302};
303
304
305// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306
307} // End namespace Foam
308
309// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310
311#endif
312
313// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
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 tmp< fvPatchField< scalar > > clone(const DimensionedField< scalar, volMesh > &iF) const
Clone with an internal field reference.
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.
TypeName("sorptionWallFunction")
Runtime type information.
A class for managing temporary objects.
Definition tmp.H:75
The class wallFunctionBlenders is a base class that hosts common entries for various derived wall-fun...
wallFunctionBlenders()
Default construct with default coefficients.
Class to host the wall-function coefficients being used in the wall function boundary conditions.
A class for handling words, derived from Foam::string.
Definition word.H:66
scalar yPlus
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
fvPatchField< scalar > fvPatchScalarField
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68