atmOmegaWallFunctionFvPatchScalarField.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,
43  scalarField& omega0
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 kappa = wallCoeffs_.kappa();
65 
66  const scalar t = db().time().timeOutputValue();
67  const scalarField z0(z0_->value(t));
68 
69  #ifdef FULLDEBUG
70  for (const scalar z : z0)
71  {
72  if (z < VSMALL)
73  {
75  << "z0 field can only contain positive values. "
76  << "Please check input field z0."
77  << exit(FatalError);
78  }
79  }
80  #endif
81 
83 
84  // Set omega and G
85  forAll(nutw, facei)
86  {
87  const label celli = faceCells[facei];
88 
89  const scalar w = cornerWeights[facei];
90 
91  omega0[celli] +=
92  w*sqrt(k[celli])/(Cmu25*kappa*(y[facei] + z0[facei]));
93 
94  G0[celli] +=
95  w
96  *(nutw[facei] + nuw[facei])
97  *magGradUw[facei]
98  *Cmu25*sqrt(k[celli])
99  /(kappa*(y[facei] + z0[facei]));
100  }
101 }
102 
103 
105 (
106  Ostream& os
107 ) const
108 {
109  if (z0_)
110  {
111  z0_->writeData(os);
112  }
113 
114  wallCoeffs_.writeEntries(os);
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119 
122 (
123  const fvPatch& p,
125 )
126 :
128  z0_(nullptr)
129 {}
130 
131 
134 (
136  const fvPatch& p,
138  const fvPatchFieldMapper& mapper
139 )
140 :
142  z0_(ptf.z0_.clone(p.patch()))
143 {}
144 
145 
148 (
149  const fvPatch& p,
151  const dictionary& dict
152 )
153 :
155  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
156 {}
157 
158 
161 (
163 )
164 :
166  z0_(owfpsf.z0_.clone(this->patch().patch()))
167 {}
168 
169 
172 (
175 )
176 :
178  z0_(owfpsf.z0_.clone(this->patch().patch()))
179 {}
180 
181 
182 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
183 
185 (
186  const fvPatchFieldMapper& m
187 )
188 {
190 
191  if (z0_)
192  {
193  z0_->autoMap(m);
194  }
195 }
196 
197 
199 (
200  const fvPatchScalarField& ptf,
201  const labelList& addr
202 )
203 {
205 
206  const auto& atmpsf =
207  refCast<const atmOmegaWallFunctionFvPatchScalarField>(ptf);
208 
209  if (z0_)
210  {
211  z0_->rmap(atmpsf.z0_(), addr);
212  }
213 }
214 
215 
217 (
218  Ostream& os
219 ) const
220 {
222  writeLocalEntries(os);
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 namespace Foam
230 {
232  (
234  atmOmegaWallFunctionFvPatchScalarField
235  );
236 }
237 
238 
239 // ************************************************************************* //
dictionary dict
const objectRegistry & db() const
The associated objectRegistry.
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
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
atmOmegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
This boundary condition provides a wall constraint on the specific dissipation rate (i...
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
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
void writeLocalEntries(Ostream &) const
Write local wall function variables.
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
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
fvPatchField< scalar > fvPatchScalarField
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 calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
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.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
OBJstream os(runTime.globalPath()/outputName)
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
scalar Cmu() const noexcept
Return the object: Cmu.
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion) ...
Definition: TimeStateI.H:24
This boundary condition provides a wall function for the specific dissipation rate (i...
autoPtr< PatchFunction1< scalar > > z0_
Surface roughness length field [m].
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.
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 Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.