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-2021 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 "momentumTransportModel.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 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
59 {}
60 
61 
63 (
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper
68 )
69 :
70  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
71 {}
72 
73 
75 (
76  const fWallFunctionFvPatchScalarField& v2wfpsf,
78 )
79 :
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
87 {
88  if (updated())
89  {
90  return;
91  }
92 
93  const label patchi = patch().index();
94 
95  const momentumTransportModel& turbModel =
97  (
99  (
100  momentumTransportModel::typeName,
101  internalField().group()
102  )
103  );
104  const v2fBase& v2fModel = refCast<const v2fBase>(turbModel);
105 
107  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
108 
109  const scalarField& y = turbModel.y()[patchi];
110 
111  const tmp<volScalarField> tk = turbModel.k();
112  const volScalarField& k = tk();
113 
114  const tmp<volScalarField> tepsilon = turbModel.epsilon();
115  const volScalarField& epsilon = tepsilon();
116 
117  const tmp<volScalarField> tv2 = v2fModel.v2();
118  const volScalarField& v2 = tv2();
119 
120  const tmp<scalarField> tnuw = turbModel.nu(patchi);
121  const scalarField& nuw = tnuw();
122 
123  const scalar Cmu25 = pow025(nutw.Cmu());
124 
125  scalarField& f = *this;
126 
127  // Set f wall values
128  forAll(f, facei)
129  {
130  label celli = patch().faceCells()[facei];
131 
132  scalar uTau = Cmu25*sqrt(k[celli]);
133 
134  scalar yPlus = uTau*y[facei]/nuw[facei];
135 
136  if (yPlus > nutw.yPlusLam())
137  {
138  scalar N = 6.0;
139  scalar v2c = v2[celli];
140  scalar epsc = epsilon[celli];
141  scalar kc = k[celli];
142 
143  f[facei] = N*v2c*epsc/(sqr(kc) + rootVSmall);
144  f[facei] /= sqr(uTau) + rootVSmall;
145  }
146  else
147  {
148  f[facei] = 0.0;
149  }
150  }
151 
153 
154  // TODO: perform averaging for cells sharing more than one boundary face
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
161 (
164 );
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 } // End namespace RASModels
169 } // End namespace Foam
170 
171 // ************************************************************************* //
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
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:174
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:355
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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
const nearWallDist & y() const
Return the near wall distances.
label k
Boltzmann constant.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
scalar uTau
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields. ...
Definition: v2fBase.H:54
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:337
makePatchTypeField(fvPatchScalarField, fWallFunctionFvPatchScalarField)
labelList f(nPoints)
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
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 void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:209
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:147
label N
Definition: createFields.H:9
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
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:343