fWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "turbulenceModel.H"
29 #include "v2f.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const fvPatch& p,
45 )
46 :
48 {}
49 
50 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
60 {}
61 
62 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
71 {}
72 
73 
75 (
76  const fWallFunctionFvPatchScalarField& v2wfpsf
77 )
78 :
80 {}
81 
82 
84 (
85  const fWallFunctionFvPatchScalarField& v2wfpsf,
87 )
88 :
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  if (updated())
98  {
99  return;
100  }
101 
102  const label patchi = patch().index();
103 
104  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
105  (
107  (
109  internalField().group()
110  )
111  );
112  const v2fBase& v2fModel = refCast<const v2fBase>(turbModel);
113 
115  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
116 
117  const scalarField& y = turbModel.y()[patchi];
118 
119  const tmp<volScalarField> tk = turbModel.k();
120  const volScalarField& k = tk();
121 
122  const tmp<volScalarField> tepsilon = turbModel.epsilon();
123  const volScalarField& epsilon = tepsilon();
124 
125  const tmp<volScalarField> tv2 = v2fModel.v2();
126  const volScalarField& v2 = tv2();
127 
128  const tmp<scalarField> tnuw = turbModel.nu(patchi);
129  const scalarField& nuw = tnuw();
130 
131  const scalar Cmu25 = pow025(nutw.Cmu());
132 
133  scalarField& f = *this;
134 
135  // Set f wall values
136  forAll(f, facei)
137  {
138  label celli = patch().faceCells()[facei];
139 
140  scalar uTau = Cmu25*sqrt(k[celli]);
141 
142  scalar yPlus = uTau*y[facei]/nuw[facei];
143 
144  if (yPlus > nutw.yPlusLam())
145  {
146  scalar N = 6.0;
147  scalar v2c = v2[celli];
148  scalar epsc = epsilon[celli];
149  scalar kc = k[celli];
150 
151  f[facei] = N*v2c*epsc/(sqr(kc) + rootVSmall);
152  f[facei] /= sqr(uTau) + rootVSmall;
153  }
154  else
155  {
156  f[facei] = 0.0;
157  }
158  }
159 
161 
162  // TODO: perform averaging for cells sharing more than one boundary face
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
169 (
172 );
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 } // End namespace RASModels
177 } // End namespace Foam
178 
179 // ************************************************************************* //
const char *const group
Group name for atomic constants.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:173
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:356
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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:61
dimensionedScalar pow025(const dimensionedScalar &ds)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
label k
Boltzmann constant.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
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.
scalar uTau
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields. ...
Definition: v2fBase.H:54
scalar y
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
const nearWallDist & y() const
Return the near wall distances.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:325
makePatchTypeField(fvPatchScalarField, fWallFunctionFvPatchScalarField)
labelList f(nPoints)
fWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
This boundary condition provides a turbulence damping function, f, wall function condition for low- a...
scalar yPlus
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:229
scalar epsilon
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:167
label N
Definition: createFields.H:22
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Namespace for OpenFOAM.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:332