atmNutWallFunctionFvPatchScalarField.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 CENER
9  Copyright (C) 2020-2022 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 "bound.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 {
45  const label patchi = patch().index();
46 
47  const auto& turbModel = db().lookupObject<turbulenceModel>
48  (
50  (
52  internalField().group()
53  )
54  );
55 
56  const scalarField& y = turbModel.y()[patchi];
57 
58  const tmp<volScalarField> tk = turbModel.k();
59  const volScalarField& k = tk();
60 
61  const tmp<scalarField> tnuw = turbModel.nu(patchi);
62  const scalarField& nuw = tnuw();
63 
64  auto tnutw = tmp<scalarField>::New(*this);
65  auto& nutw = tnutw.ref();
66 
67  const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
68  const scalar kappa = wallCoeffs_.kappa();
69 
70  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
71  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
72 
73  const scalar t = db().time().timeOutputValue();
74  const scalarField z0(z0_->value(t));
75 
76  #ifdef FULLDEBUG
77  for (const scalar z : z0)
78  {
79  if (z < VSMALL)
80  {
82  << "z0 field can only contain positive values. "
83  << "Please check input field z0."
84  << exit(FatalError);
85  }
86  }
87  #endif
88 
89  const labelList& faceCells = patch().faceCells();
90 
91  forAll(nutw, facei)
92  {
93  const label celli = faceCells[facei];
94 
95  // (RH:Eq. 6)
96  const scalar Edash = (y[facei] + z0[facei])/(z0[facei] + z0Min_);
97 
98  // (RH:Eq. 6)
99  const scalar uStarU = magUp[facei]*kappa/log(max(Edash, 1 + SMALL));
100 
101  // (RH:Eq. 7)
102  const scalar uStarK = Cmu25*Foam::sqrt(k[celli]);
103 
104  // (SBJM:Eq. 7; SM:Eq. 25)
105  const scalar tauw = uStarU*uStarK;
106 
107  nutw[facei] =
108  max(tauw*y[facei]/(max(magUp[facei], SMALL)) - nuw[facei], 0.0);
109  }
110 
111  return tnutw;
112 }
113 
114 
116 (
117  Ostream& os
118 ) const
119 {
120  os.writeEntryIfDifferent<scalar>("z0Min", SMALL, z0Min_);
121 
122  if (z0_)
123  {
124  z0_->writeData(os);
125  }
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
130 
132 (
133  const fvPatch& p,
134  const DimensionedField<scalar, volMesh>& iF
135 )
136 :
138  z0Min_(SMALL),
139  z0_(nullptr)
140 {}
141 
142 
144 (
146  const fvPatch& p,
148  const fvPatchFieldMapper& mapper
149 )
150 :
152  z0Min_(ptf.z0Min_),
153  z0_(ptf.z0_.clone(p.patch()))
154 {}
155 
156 
158 (
159  const fvPatch& p,
161  const dictionary& dict
162 )
163 :
165  z0Min_
166  (
167  dict.getCheckOrDefault<scalar>
168  (
169  "z0Min",
170  SMALL,
171  scalarMinMax::ge(0)
172  )
173  ),
174  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
175 {}
176 
177 
179 (
181 )
182 :
184  z0Min_(rwfpsf.z0Min_),
185  z0_(rwfpsf.z0_.clone(this->patch().patch()))
186 {}
187 
188 
190 (
193 )
194 :
196  z0Min_(rwfpsf.z0Min_),
197  z0_(rwfpsf.z0_.clone(this->patch().patch()))
198 {}
199 
200 
201 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
202 
204 (
205  const fvPatchFieldMapper& m
206 )
207 {
208  nutkWallFunctionFvPatchScalarField::autoMap(m);
209 
210  if (z0_)
211  {
212  z0_->autoMap(m);
213  }
214 }
215 
216 
218 (
219  const fvPatchScalarField& ptf,
220  const labelList& addr
221 )
222 {
223  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
224 
225  const auto& nrwfpsf =
226  refCast<const atmNutWallFunctionFvPatchScalarField>(ptf);
227 
228  if (z0_)
229  {
230  z0_->rmap(nrwfpsf.z0_(), addr);
231  }
232 }
233 
234 
236 {
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
246 (
248  atmNutWallFunctionFvPatchScalarField
249 );
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 } // End namespace Foam
254 
255 // ************************************************************************* //
atmNutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dictionary dict
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
dimensionedScalar pow025(const dimensionedScalar &ds)
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.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
label k
Boltzmann constant.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
scalar tauw
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
scalar magUp
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on the turbulent kinetic energy, (i.e. k) for for low- and high-Reynolds number applications.
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
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.
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
This boundary condition provides a wall constraint on the turbulent viscosity (i.e. nut) based on the turbulent kinetic energy (i.e. k) and velocity (i.e. U) for atmospheric boundary layer modelling.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
Definition: fvPatchField.C:221
Bound the given scalar field if it has gone unbounded.
OBJstream os(runTime.globalPath()/outputName)
scalar Cmu() const noexcept
Return the object: Cmu.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
List< label > labelList
A List of labels.
Definition: List.H:62
scalar kappa() const noexcept
Return the object: kappa.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.