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 -------------------------------------------------------------------------------
10 License
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 
26 Class
27  Foam::sorptionWallFunctionFvPatchScalarField
28 
29 Group
30  grpWallFunctions
31 
32 Description
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 
100 Usage
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 
147 SourceFiles
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 
163 namespace Foam
164 {
165 
166 /*---------------------------------------------------------------------------*\
167  Class sorptionWallFunctionFvPatchScalarField Declaration
168 \*---------------------------------------------------------------------------*/
169 
170 class sorptionWallFunctionFvPatchScalarField
171 :
172  public fixedGradientFvPatchScalarField,
173  private wallFunctionBlenders
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 
214 public:
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 and return a clone
255  virtual tmp<fvPatchScalarField> clone() const
256  {
257  return tmp<fvPatchScalarField>
258  (
260  );
261  }
262 
263  //- Construct as copy setting internal field reference
265  (
267  const DimensionedField<scalar, volMesh>&
268  );
269 
270  //- Construct and return a clone setting internal field reference
272  (
274  ) const
275  {
277  (
279  (
280  *this,
281  iF
282  )
283  );
284  }
285 
286 
287  // Member Functions
288 
289  // Mapping
290 
291  //- Map (and resize as needed) from self given a mapping object
292  virtual void autoMap(const fvPatchFieldMapper&);
293 
294  //- Reverse map the given fvPatchField onto this fvPatchField
295  virtual void rmap
296  (
297  const fvPatchScalarField&,
298  const labelList&
299  );
300 
301 
302  // Evaluation
303 
304  //- Update the coefficients associated with the patch field
305  virtual void updateCoeffs();
306 
307 
308  // I-O
309 
310  //- Write
311  virtual void write(Ostream&) const;
312 };
313 
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 
317 } // End namespace Foam
318 
319 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320 
321 #endif
322 
323 // ************************************************************************* //
sorptionWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A FieldMapper for finite-volume patch fields.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
TypeName("sorptionWallFunction")
Runtime type information.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.