Loading...
Searching...
No Matches
sorptionWallFunctionFvPatchScalarField.C
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-2023 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
26\*---------------------------------------------------------------------------*/
27
30#include "fvPatchFieldMapper.H"
31#include "volFields.H"
32#include "wallDist.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
38
39// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
40
41//- Estimate the y* at the intersection of the two sublayers
42static scalar calcYStarLam
43(
44 const scalar kappa,
45 const scalar E,
46 const scalar Sc,
47 const scalar Sct,
48 const scalar Pc
49)
50{
51 scalar ypl = 11;
52
53 for (int iter = 0; iter < 10; ++iter)
54 {
55 // (F:Eq. 5.5)
56 ypl = (log(max(E*ypl, scalar(1)))/kappa + Pc)*Sct/Sc;
57 }
58
59 return ypl;
60}
61
62
63// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
64
65tmp<scalarField> sorptionWallFunctionFvPatchScalarField::yPlus() const
66{
67 const label patchi = patch().index();
68
69 // Calculate the empirical constant given by (Jayatilleke, 1966) (FDC:Eq. 6)
70 const scalar Pc =
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();
75
76 // Calculate fields of interest
77 const auto& k = db().lookupObject<volScalarField>(kName_);
78 tmp<scalarField> tkwc = k.boundaryField()[patchi].patchInternalField();
79 const scalarField& kwc = tkwc.cref();
80
81 const auto& nu = db().lookupObject<volScalarField>(nuName_);
82 tmp<scalarField> tnuwc = nu.boundaryField()[patchi].patchInternalField();
83 const scalarField& nuwc = tnuwc.cref();
84
85 const volScalarField& y = wallDist::New(internalField().mesh()).y();
86 tmp<scalarField> tywc = y.boundaryField()[patchi].patchInternalField();
87 const scalarField& ywc = tywc.cref();
88
89 // (FDC:Eq. 3)
90 const auto yStar = [&](const label facei) -> scalar
91 {
92 return
93 (
94 Cmu25*sqrt(kwc[facei])*ywc[facei]/nuwc[facei]
95 );
96 };
97
98 // (FDC:Eq. 4)
99 const auto yPlusVis = [&](const label facei) -> scalar
100 {
101 return (Sc_*yStar(facei));
102 };
103
104 // (FDC:Eq. 5)
105 const auto yPlusLog = [&](const label facei) -> scalar
106 {
107 return
108 (
109 Sct_*(log(max(E*yStar(facei), 1 + 1e-4))/kappa + Pc)
110 );
111 };
112
113 auto tyPlus = tmp<scalarField>::New(patch().size());
114 auto& yPlus = tyPlus.ref();
115
116 switch (blender_)
117 {
118 case blenderType::EXPONENTIAL:
119 {
120 forAll(yPlus, facei)
121 {
122 // (FDC:Eq. 2)
123 const scalar yStarFace = yStar(facei);
124 const scalar Gamma =
125 0.01*pow4(Sc_*yStarFace)/(1 + 5*pow3(Sc_)*yStarFace);
126 const scalar invGamma = scalar(1)/max(Gamma, ROOTVSMALL);
127
128 // (FDC:Eq. 1)
129 yPlus[facei] =
130 (
131 yPlusVis(facei)*exp(-Gamma)
132 + yPlusLog(facei)*exp(-invGamma)
133 );
134 }
135 break;
136 }
137
138 case blenderType::STEPWISE:
139 {
140 static const scalar yStarLam =
141 calcYStarLam(kappa, E, Sc_, Sct_, Pc);
142
143 forAll(yPlus, facei)
144 {
145 // (F:Eq. 5.3)
146 if (yStar(facei) < yStarLam)
147 {
148 yPlus[facei] = yPlusVis(facei);
149 }
150 else
151 {
152 yPlus[facei] = yPlusLog(facei);
153 }
154 }
155 break;
156 }
157
158 case blenderType::BINOMIAL:
159 {
160 forAll(yPlus, facei)
161 {
162 yPlus[facei] =
163 pow
164 (
165 pow(yPlusVis(facei), n_) + pow(yPlusLog(facei), n_),
166 scalar(1)/n_
167 );
168 }
169 break;
170 }
171
172 case blenderType::MAX:
173 {
174 forAll(yPlus, facei)
175 {
176 yPlus[facei] = max(yPlusVis(facei), yPlusLog(facei));
177 }
178 break;
179 }
180
181 case blenderType::TANH:
182 {
183 forAll(yPlus, facei)
184 {
185 const scalar yPlusVisFace = yPlusVis(facei);
186 const scalar yPlusLogFace = yPlusLog(facei);
187 const scalar b1 = yPlusVisFace + yPlusLogFace;
188 const scalar b2 =
189 pow
190 (
191 pow(yPlusVisFace, 1.2) + pow(yPlusLogFace, 1.2),
192 1.0/1.2
193 );
194 const scalar phiTanh = tanh(pow4(0.1*yStar(facei)));
195
196 yPlus[facei] = phiTanh*b1 + (1 - phiTanh)*b2;
197 }
198 break;
199 }
200 }
201
202 return tyPlus;
203}
204
205
206tmp<scalarField> sorptionWallFunctionFvPatchScalarField::flux() const
207{
208 // Calculate fields of interest
209 const label patchi = patch().index();
210
211 const auto& k = db().lookupObject<volScalarField>(kName_);
212 tmp<scalarField> tkwc = k.boundaryField()[patchi].patchInternalField();
213
214 const volScalarField& y = wallDist::New(internalField().mesh()).y();
215 tmp<scalarField> tywc = y.boundaryField()[patchi].patchInternalField();
216
217
218 // Calculate mass-transfer coefficient field (FDC:Eqs. 8-9)
220 laminar_
221 ? D_/tywc
222 : pow025(wallCoeffs_.Cmu())*sqrt(tkwc)/yPlus();
223
224
225 // Calculate wall-surface concentration
226 const auto& Csurf = *this;
227
228 // Calculate wall-adjacent concentration
229 const scalar t = db().time().timeOutputValue();
230 tmp<scalarField> tkAbs = kAbsPtr_->value(t);
231 tmp<scalarField> tCstar = Csurf/tkAbs;
232
233 // Calculate near-wall cell concentration
234 const word& fldName = internalField().name();
235 const auto& C = db().lookupObject<volScalarField>(fldName);
236 tmp<scalarField> tCwc = C.boundaryField()[patchi].patchInternalField();
237
238 // Return adsorption or absorption/permeation flux
239 // normalised by the mass-transfer coefficient (FDC:Fig. 1)
240 return (tCstar - tCwc)/ta;
241}
242
243
244void sorptionWallFunctionFvPatchScalarField::writeLocalEntries
245(
246 Ostream& os
247) const
248{
250 os.writeEntryIfDifferent<bool>("laminar", false, laminar_);
251 os.writeEntry("Sc", Sc_);
252 os.writeEntry("Sct", Sct_);
253 os.writeEntryIfDifferent<scalar>("D", -1, D_);
254 wallCoeffs_.writeEntries(os);
255 os.writeEntryIfDifferent<word>("k", "k", kName_);
256 os.writeEntryIfDifferent<word>("nu", "nu", nuName_);
257
258 if (kAbsPtr_)
259 {
260 kAbsPtr_->writeData(os);
261 }
263
264
265// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
266
268(
269 const fvPatch& p,
271)
272:
273 fixedGradientFvPatchScalarField(p, iF),
275 laminar_(false),
276 kAbsPtr_(nullptr),
277 Sc_(1),
278 Sct_(1),
279 D_(-1),
280 kName_("k"),
281 nuName_("nu"),
282 wallCoeffs_()
283{}
284
285
287(
289 const fvPatch& p,
291 const fvPatchFieldMapper& mapper
292)
293:
294 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
296 laminar_(ptf.laminar_),
297 kAbsPtr_(ptf.kAbsPtr_.clone(patch().patch())),
298 Sc_(ptf.Sc_),
299 Sct_(ptf.Sct_),
300 D_(ptf.D_),
301 kName_(ptf.kName_),
302 nuName_(ptf.nuName_),
303 wallCoeffs_(ptf.wallCoeffs_)
304{}
305
306
308(
309 const fvPatch& p,
311 const dictionary& dict
312)
313:
314 fixedGradientFvPatchScalarField(p, iF), // Bypass dictionary constructor
315 wallFunctionBlenders(dict, blenderType::STEPWISE, scalar(2)),
316 laminar_(dict.getOrDefault<bool>("laminar", false)),
317 kAbsPtr_(PatchFunction1<scalar>::New(p.patch(), "kAbs", dict)),
318 Sc_(dict.getCheck<scalar>("Sc", scalarMinMax::ge(0))),
319 Sct_(dict.getCheck<scalar>("Sct", scalarMinMax::ge(0))),
320 D_(dict.getOrDefault<scalar>("D", -1)),
321 kName_(dict.getOrDefault<word>("k", "k")),
322 nuName_(dict.getOrDefault<word>("nu", "nu")),
323 wallCoeffs_(dict)
324{
325 if (laminar_)
326 {
327 if (D_ < 0)
328 {
330 << "Molecular diffusion coefficient cannot be non-positive. "
331 << "D = " << D_
332 << exit(FatalIOError);
333 }
334 }
335
336 if (!kAbsPtr_)
337 {
339 << "Adsorption or absorption coefficient is not set."
340 << exit(FatalIOError);
341 }
342
343 if (!this->readGradientEntry(dict) || !this->readValueEntry(dict))
344 {
345 extrapolateInternal();
346 gradient() = Zero;
347 }
348}
349
350
352(
353 const sorptionWallFunctionFvPatchScalarField& swfpsf
354)
355:
356 fixedGradientFvPatchScalarField(swfpsf),
357 wallFunctionBlenders(swfpsf),
358 laminar_(swfpsf.laminar_),
359 kAbsPtr_(swfpsf.kAbsPtr_.clone(patch().patch())),
360 Sc_(swfpsf.Sc_),
361 Sct_(swfpsf.Sct_),
362 D_(swfpsf.D_),
363 kName_(swfpsf.kName_),
364 nuName_(swfpsf.nuName_),
365 wallCoeffs_(swfpsf.wallCoeffs_)
366{}
367
368
370(
373)
374:
375 fixedGradientFvPatchScalarField(swfpsf, iF),
376 wallFunctionBlenders(swfpsf),
377 laminar_(swfpsf.laminar_),
378 kAbsPtr_(swfpsf.kAbsPtr_.clone(patch().patch())),
379 Sc_(swfpsf.Sc_),
380 Sct_(swfpsf.Sct_),
381 D_(swfpsf.D_),
382 kName_(swfpsf.kName_),
383 nuName_(swfpsf.nuName_),
384 wallCoeffs_(swfpsf.wallCoeffs_)
386
387
388// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
389
391(
392 const fvPatchFieldMapper& mapper
393)
394{
395 fixedGradientFvPatchScalarField::autoMap(mapper);
396
397 if (kAbsPtr_)
398 {
399 kAbsPtr_->autoMap(mapper);
400 }
401}
402
403
405(
406 const fvPatchScalarField& ptf,
407 const labelList& addr
408)
409{
410 fixedGradientFvPatchScalarField::rmap(ptf, addr);
411
412 const auto& swfptf =
414
415 if (kAbsPtr_)
416 {
417 kAbsPtr_->rmap(swfptf.kAbsPtr_(), addr);
418 }
419}
420
421
423{
424 if (updated())
425 {
426 return;
427 }
428
429 gradient() = flux()/patch().deltaCoeffs();
431 fixedGradientFvPatchScalarField::updateCoeffs();
432}
433
434
436{
438
439 writeLocalEntries(os);
440
443
444
445// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446
448(
450 sorptionWallFunctionFvPatchScalarField
451);
452
453
454// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
455
456} // End namespace Foam
457
458// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
scalar y
label k
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Definition C.H:49
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,...
Definition Ostream.H:59
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,...
Definition dictionary.H:133
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.
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
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.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field.
Definition wallDist.H:220
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.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
scalar yPlus
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.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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.
Definition List.H:62
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
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).
Definition zero.H:127
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)
Definition errorManip.H:125
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar pow025(const dimensionedScalar &ds)
volScalarField & nu
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299