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 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 
31 #include "turbulenceModel.H"
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 
49  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
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(nutw.Cmu());
64  const scalar Cmu75 = pow(nutw.Cmu(), 0.75);
65 
66  const scalar t = db().time().timeOutputValue();
67  const scalarField z0(z0_->value(t));
68 
69  #ifdef FULLDEBUG
70  for (const auto& 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 
82  const labelUList& faceCells = patch.faceCells();
83 
84  // Set epsilon and G
85  forAll(nutw, facei)
86  {
87  const label celli = faceCells[facei];
88 
89  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
90 
91  const scalar w = cornerWeights[facei];
92 
93  // (PGVB:Eq. 7, RH:Eq. 8)
94  scalar epsilonc =
95  w*Cmu75*pow(k[celli], 1.5)/(nutw.kappa()*(y[facei] + z0[facei]));
96 
97  scalar Gc =
98  w
99  *(nutw[facei] + nuw[facei])
100  *magGradUw[facei]
101  *Cmu25*sqrt(k[celli])
102  /(nutw.kappa()*(y[facei] + z0[facei]));
103 
104  if (lowReCorrection_ && yPlus < nutw.yPlusLam())
105  {
106  epsilonc = w*2.0*k[celli]*nuw[facei]/sqr(y[facei] + z0[facei]);
107  Gc = 0;
108  }
109 
110  epsilon0[celli] += epsilonc;
111 
112  G0[celli] += Gc;
113  }
114 }
115 
116 
117 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
118 
121 (
122  const fvPatch& p,
124 )
125 :
127  z0_(nullptr)
128 {}
129 
130 
133 (
135  const fvPatch& p,
137  const fvPatchFieldMapper& mapper
138 )
139 :
140  epsilonWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
141  z0_(ptf.z0_.clone(p.patch()))
142 {}
143 
144 
147 (
148  const fvPatch& p,
150  const dictionary& dict
151 )
152 :
154  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
155 {}
156 
157 
160 (
162 )
163 :
165  z0_(ewfpsf.z0_.clone(this->patch().patch()))
166 {}
167 
168 
171 (
174 )
175 :
177  z0_(ewfpsf.z0_.clone(this->patch().patch()))
178 {}
179 
180 
181 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182 
184 (
185  const fvPatchFieldMapper& m
186 )
187 {
188  epsilonWallFunctionFvPatchScalarField::autoMap(m);
189  z0_->autoMap(m);
190 }
191 
192 
194 (
195  const fvPatchScalarField& ptf,
196  const labelList& addr
197 )
198 {
199  epsilonWallFunctionFvPatchScalarField::rmap(ptf, addr);
200 
202  refCast<const atmEpsilonWallFunctionFvPatchScalarField>(ptf);
203 
204  z0_->rmap(atmpsf.z0_(), addr);
205 }
206 
207 
209 (
210  Ostream& os
211 ) const
212 {
214  z0_->writeData(os);
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 namespace Foam
221 {
223  (
226  );
227 }
228 
229 
230 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::atmEpsilonWallFunctionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmEpsilonWallFunctionFvPatchScalarField.C:184
Foam::nutWallFunctionFvPatchScalarField::Cmu
scalar Cmu() const
Return Cmu.
Definition: nutWallFunctionFvPatchScalarField.H:369
Foam::atmEpsilonWallFunctionFvPatchScalarField::z0_
autoPtr< PatchFunction1< scalar > > z0_
Surface roughness length field [m].
Definition: atmEpsilonWallFunctionFvPatchScalarField.H:136
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
nutWallFunctionFvPatchScalarField.H
fvMatrix.H
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
Foam::constant::electromagnetic::epsilon0
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::Field< scalar >
Foam::constant::electromagnetic::G0
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
Foam::atmEpsilonWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent kinetic energy dissipation rate (...
Definition: atmEpsilonWallFunctionFvPatchScalarField.H:127
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::epsilonWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent kinetic energy dissipation rate,...
Definition: epsilonWallFunctionFvPatchScalarField.H:259
atmEpsilonWallFunctionFvPatchScalarField.H
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
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:123
os
OBJstream os(runTime.globalPath()/outputName)
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:60
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::nutWallFunctionFvPatchScalarField::yPlusLam
static scalar yPlusLam(const scalar kappa, const scalar E)
Estimate the y+ at the intersection of the two sublayers.
Definition: nutWallFunctionFvPatchScalarField.C:250
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:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::nutWallFunctionFvPatchScalarField::kappa
scalar kappa() const
Return kappa.
Definition: nutWallFunctionFvPatchScalarField.H:375
Foam::UList< label >
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::atmEpsilonWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmEpsilonWallFunctionFvPatchScalarField.C:194
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::GeometricField< scalar, fvPatchField, volMesh >
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
Foam::atmEpsilonWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmEpsilonWallFunctionFvPatchScalarField.C:209
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
turbulenceModel.H
Foam::atmEpsilonWallFunctionFvPatchScalarField::calculate
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
Definition: atmEpsilonWallFunctionFvPatchScalarField.C:38
Foam::atmEpsilonWallFunctionFvPatchScalarField::atmEpsilonWallFunctionFvPatchScalarField
atmEpsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmEpsilonWallFunctionFvPatchScalarField.C:121
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
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::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:144
y
scalar y
Definition: LISASMDCalcMethod1.H:14