atmEpsilonWallFunctionFvPatchScalarField.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-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 
32 #include "fvMatrix.H"
34 
35 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
36 
38 (
39  const turbulenceModel& turbModel,
40  const List<scalar>& cornerWeights,
41  const fvPatch& patch,
42  scalarField& G0,
44 )
45 {
46  const label patchi = patch.index();
47 
48  const tmp<scalarField> tnutw = turbModel.nut(patchi);
49  const scalarField& nutw = tnutw();
50 
51  const scalarField& y = turbModel.y()[patchi];
52 
53  const tmp<scalarField> tnuw = turbModel.nu(patchi);
54  const scalarField& nuw = tnuw();
55 
56  const tmp<volScalarField> tk = turbModel.k();
57  const volScalarField& k = tk();
58 
59  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
60 
61  const scalarField magGradUw(mag(Uw.snGrad()));
62 
63  const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
64  const scalar Cmu75 = pow(wallCoeffs_.Cmu(), 0.75);
65  const scalar kappa = wallCoeffs_.kappa();
66  const scalar yPlusLam = wallCoeffs_.yPlusLam();
67 
68  const scalar t = db().time().timeOutputValue();
69  const scalarField z0(z0_->value(t));
70 
71  #ifdef FULLDEBUG
72  for (const auto& z : z0)
73  {
74  if (z < VSMALL)
75  {
77  << "z0 field can only contain positive values. "
78  << "Please check input field z0."
79  << exit(FatalError);
80  }
81  }
82  #endif
83 
85 
86  // Set epsilon and G
87  forAll(nutw, facei)
88  {
89  const label celli = faceCells[facei];
90 
91  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
92 
93  const scalar w = cornerWeights[facei];
94 
95  // (PGVB:Eq. 7, RH:Eq. 8)
96  scalar epsilonc =
97  w*Cmu75*pow(k[celli], 1.5)/(kappa*(y[facei] + z0[facei]));
98 
99  scalar Gc =
100  w
101  *(nutw[facei] + nuw[facei])
102  *magGradUw[facei]
103  *Cmu25*sqrt(k[celli])
104  /(kappa*(y[facei] + z0[facei]));
105 
106  if (lowReCorrection_ && yPlus < yPlusLam)
107  {
108  epsilonc = w*2.0*k[celli]*nuw[facei]/sqr(y[facei] + z0[facei]);
109  Gc = 0;
110  }
111 
112  epsilon0[celli] += epsilonc;
114  G0[celli] += Gc;
115  }
116 }
117 
118 
120 (
121  Ostream& os
122 ) const
123 {
124  os.writeEntryIfDifferent<bool>("lowReCorrection", false, lowReCorrection_);
125 
126  if (z0_)
127  {
128  z0_->writeData(os);
129  }
130 
131  wallCoeffs_.writeEntries(os);
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
136 
139 (
140  const fvPatch& p,
142 )
143 :
145  z0_(nullptr)
146 {}
147 
148 
151 (
153  const fvPatch& p,
155  const fvPatchFieldMapper& mapper
156 )
157 :
159  z0_(ptf.z0_.clone(p.patch()))
160 {}
161 
162 
165 (
166  const fvPatch& p,
168  const dictionary& dict
169 )
170 :
172  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
173 {}
174 
175 
178 (
180 )
181 :
183  z0_(ewfpsf.z0_.clone(this->patch().patch()))
184 {}
185 
186 
189 (
192 )
193 :
195  z0_(ewfpsf.z0_.clone(this->patch().patch()))
196 {}
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
202 (
203  const fvPatchFieldMapper& m
204 )
205 {
207 
208  if (z0_)
209  {
210  z0_->autoMap(m);
211  }
212 }
213 
214 
216 (
217  const fvPatchScalarField& ptf,
218  const labelList& addr
219 )
220 {
222 
223  const auto& atmpsf =
224  refCast<const atmEpsilonWallFunctionFvPatchScalarField>(ptf);
225  if (z0_)
226  {
227  z0_->rmap(atmpsf.z0_(), addr);
228  }
229 }
230 
231 
233 (
234  Ostream& os
235 ) const
236 {
238  writeLocalEntries(os);
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 namespace Foam
246 {
248  (
250  atmEpsilonWallFunctionFvPatchScalarField
251  );
252 }
253 
254 
255 // ************************************************************************* //
dictionary dict
const objectRegistry & db() const
The associated objectRegistry.
scalar yPlusLam() const noexcept
Return the object: yPlusLam.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:236
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:213
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const volVectorField & U() const
Access function to velocity field.
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
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedSymmTensor sqr(const dimensionedVector &dv)
atmEpsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
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
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
label k
Boltzmann constant.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
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.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
scalar y
fvPatchField< scalar > fvPatchScalarField
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
const Time & time() const noexcept
Return time registry.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:107
A FieldMapper for finite-volume patch fields.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:299
const nearWallDist & y() const
Return the near wall distances.
This boundary condition provides wall functions for the turbulent kinetic energy dissipation rate (i...
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...
void writeLocalEntries(Ostream &) const
Write local wall function variables.
OBJstream os(runTime.globalPath()/outputName)
scalar Cmu() const noexcept
Return the object: Cmu.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
This boundary condition provides a wall constraint on the turbulent kinetic energy dissipation rate (...
autoPtr< PatchFunction1< scalar > > z0_
Surface roughness length field [m].
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion) ...
Definition: TimeStateI.H:24
scalar yPlus
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.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
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
label index() const noexcept
The index of this patch in the boundary mesh.
Definition: fvPatch.H:218
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const bool lowReCorrection_
Apply low-Re correction term (default = no)
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.