atmAlphatkWallFunctionFvPatchScalarField.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) 2020 ENERCON GmbH
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "turbulenceModel.H"
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
44 
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
52  if (!isA<wallFvPatch>(patch()))
53  {
55  << "Invalid wall function specification" << nl
56  << " Patch type for patch " << patch().name()
57  << " must be wall" << nl
58  << " Current patch type is " << patch().type() << nl << endl
59  << abort(FatalError);
60  }
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
68 (
69  const fvPatch& p,
71 )
72 :
73  fixedValueFvPatchScalarField(p, iF),
74  Cmu_(0.09),
75  kappa_(0.41),
76  Pr_(db().time(), "Pr"),
77  Prt_(nullptr),
78  z0_(nullptr)
79 {
80  checkType();
81 }
82 
83 
86 (
88  const fvPatch& p,
90  const fvPatchFieldMapper& mapper
91 )
92 :
93  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
94  Cmu_(ptf.Cmu_),
95  kappa_(ptf.kappa_),
96  Pr_(ptf.Pr_),
97  Prt_(ptf.Prt_.clone(p.patch())),
98  z0_(ptf.z0_.clone(p.patch()))
99 {
100  checkType();
101 }
102 
103 
106 (
107  const fvPatch& p,
109  const dictionary& dict
110 )
111 :
112  fixedValueFvPatchScalarField(p, iF, dict),
113  Cmu_
114  (
115  dict.getCheckOrDefault<scalar>
116  (
117  "Cmu",
118  0.09,
120  )
121  ),
122  kappa_
123  (
124  dict.getCheckOrDefault<scalar>
125  (
126  "kappa",
127  0.41,
129  )
130  ),
131  Pr_(TimeFunction1<scalar>(db().time(), "Pr", dict)),
132  Prt_(PatchFunction1<scalar>::New(p.patch(), "Prt", dict)),
133  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
134 {
135  checkType();
136 }
137 
138 
141 (
143 )
144 :
145  fixedValueFvPatchScalarField(wfpsf),
146  Cmu_(wfpsf.Cmu_),
147  kappa_(wfpsf.kappa_),
148  Pr_(wfpsf.Pr_),
149  Prt_(wfpsf.Prt_.clone(this->patch().patch())),
150  z0_(wfpsf.z0_.clone(this->patch().patch()))
151 {
152  checkType();
153 }
154 
155 
158 (
161 )
162 :
163  fixedValueFvPatchScalarField(wfpsf, iF),
164  Cmu_(wfpsf.Cmu_),
165  kappa_(wfpsf.kappa_),
166  Pr_(wfpsf.Pr_),
167  Prt_(wfpsf.Prt_.clone(this->patch().patch())),
168  z0_(wfpsf.z0_.clone(this->patch().patch()))
169 {
170  checkType();
171 }
172 
173 
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175 
177 {
178  if (updated())
179  {
180  return;
181  }
182 
183  const label patchi = patch().index();
184 
185  // Retrieve turbulence properties from model
186  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
187  (
189  (
191  internalField().group()
192  )
193  );
194 
195  const scalarField& y = turbModel.y()[patchi];
196 
197  const tmp<scalarField> tnuw = turbModel.nu(patchi);
198  const scalarField& nuw = tnuw();
199 
200  const tmp<volScalarField> tk = turbModel.k();
201  const volScalarField& k = tk();
202 
203  const scalar Cmu25 = pow025(Cmu_);
204 
205  const scalar t = db().time().timeOutputValue();
206  const scalar Pr = Pr_.value(t);
207 
208  #ifdef FULLDEBUG
209  if (Pr < VSMALL)
210  {
212  << "Pr cannot be negative or zero. "
213  << "Please check input Pr = " << Pr
214  << exit(FatalIOError);
215  }
216  #endif
217 
218  const scalarField Prt(Prt_->value(t));
219  const scalarField z0(z0_->value(t));
220 
221  #ifdef FULLDEBUG
222  forAll(Prt, i)
223  {
224  if (Prt[i] < VSMALL || z0[i] < VSMALL)
225  {
227  << "Elements of input surface fields can only be positive. "
228  << "Please check input fields z0 and Prt."
229  << exit(FatalIOError);
230  }
231  }
232  #endif
233 
234  const labelUList& faceCells = patch().faceCells();
235 
236  scalarField& alphatw = *this;
237 
238  forAll(alphatw, facei)
239  {
240  const label celli = faceCells[facei];
241 
242  const scalar uStar = Cmu25*Foam::sqrt(k[celli]);
243  const scalar Edash = (y[facei] + z0[facei])/(z0[facei] + 1e-4);
244 
245  // Update turbulent thermal conductivity
246  alphatw[facei] =
247  uStar*kappa_*y[facei]/(Prt[facei]*log(max(Edash, 1 + 1e-4)))
248  + nuw[facei]/Pr;
249  }
250 
251  // lower bound values to avoid unrealistic
252  // negative temperatures on the ground
253  alphatw = max(alphatw, scalar(0.01));
254 
256 }
257 
258 
260 (
261  const fvPatchFieldMapper& m
262 )
263 {
264  fixedValueFvPatchScalarField::autoMap(m);
265  Prt_->autoMap(m);
266  z0_->autoMap(m);
267 }
268 
269 
271 (
272  const fvPatchScalarField& ptf,
273  const labelList& addr
274 )
275 {
276  fixedValueFvPatchScalarField::rmap(ptf, addr);
277 
279  refCast<const atmAlphatkWallFunctionFvPatchScalarField>(ptf);
280 
281  z0_->rmap(nrwfpsf.z0_(), addr);
282  Prt_->rmap(nrwfpsf.Prt_(), addr);
283 }
284 
285 
287 {
289  os.writeEntry("Cmu", Cmu_);
290  os.writeEntry("kappa", kappa_);
291  Pr_.writeData(os);
292  Prt_->writeData(os);
293  z0_->writeData(os);
294  writeEntry("value", os);
295 }
296 
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
301 (
304 );
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 } // End namespace Foam
309 
310 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::atmAlphatkWallFunctionFvPatchScalarField::atmAlphatkWallFunctionFvPatchScalarField
atmAlphatkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmAlphatkWallFunctionFvPatchScalarField.C:68
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: atmAlphatkWallFunctionFvPatchScalarField.C:176
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
wallFvPatch.H
Foam::TimeFunction1< scalar >
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
fvPatchFieldMapper.H
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
Foam::atmAlphatkWallFunctionFvPatchScalarField::Cmu_
const scalar Cmu_
Empirical model constant.
Definition: atmAlphatkWallFunctionFvPatchScalarField.H:158
Foam::atmAlphatkWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmAlphatkWallFunctionFvPatchScalarField.C:286
Foam::atmAlphatkWallFunctionFvPatchScalarField::Prt_
autoPtr< PatchFunction1< scalar > > Prt_
Turbulent Prandtl number field.
Definition: atmAlphatkWallFunctionFvPatchScalarField.H:167
atmAlphatkWallFunctionFvPatchScalarField.H
Foam::TimeFunction1::value
virtual Type value(const scalar x) const
Return value as a function of (scalar) independent variable.
Definition: TimeFunction1.C:97
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::atmAlphatkWallFunctionFvPatchScalarField::checkType
virtual void checkType()
Check the type of the patch.
Definition: atmAlphatkWallFunctionFvPatchScalarField.C:50
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::atmAlphatkWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Solution parameters.
Definition: atmAlphatkWallFunctionFvPatchScalarField.H:175
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Prt
dimensionedScalar Prt("Prt", dimless, laminarTransport)
Foam::dictionary::getCheckOrDefault
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:203
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::atmAlphatkWallFunctionFvPatchScalarField::z0_
autoPtr< PatchFunction1< scalar > > z0_
Surface roughness length [m].
Definition: atmAlphatkWallFunctionFvPatchScalarField.H:170
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::PatchFunction1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: PatchFunction1.H:63
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::atmAlphatkWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmAlphatkWallFunctionFvPatchScalarField.C:271
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:313
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::MinMax::ge
static MinMax< T > ge(const T &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:31
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::UList< label >
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::atmAlphatkWallFunctionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmAlphatkWallFunctionFvPatchScalarField.C:260
Foam::atmAlphatkWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the kinematic turbulent thermal conductivity,...
Definition: atmAlphatkWallFunctionFvPatchScalarField.H:149
Foam::TimeFunction1::writeData
virtual void writeData(Ostream &os) const
Write in dictionary format.
Definition: TimeFunction1.C:132
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::atmAlphatkWallFunctionFvPatchScalarField::Pr_
TimeFunction1< scalar > Pr_
Molecular Prandtl number.
Definition: atmAlphatkWallFunctionFvPatchScalarField.H:164
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
turbulenceModel.H
Foam::atmAlphatkWallFunctionFvPatchScalarField::kappa_
const scalar kappa_
von Kármán constant
Definition: atmAlphatkWallFunctionFvPatchScalarField.H:161
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::atmAlphatkWallFunctionFvPatchScalarField::maxIters_
static label maxIters_
Definition: atmAlphatkWallFunctionFvPatchScalarField.H:176
y
scalar y
Definition: LISASMDCalcMethod1.H:14