v2WallFunctionFvPatchScalarField.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"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 (
42  const fvPatch& p,
44 )
45 :
47 {}
48 
49 
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
58 {}
59 
60 
62 (
64  const fvPatch& p,
66  const fvPatchFieldMapper& mapper
67 )
68 :
69  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
70 {}
71 
72 
74 (
75  const v2WallFunctionFvPatchScalarField& v2wfpsf,
77 )
78 :
80 {}
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
86 {
87  if (updated())
88  {
89  return;
90  }
91 
92  const label patchi = patch().index();
93 
94  const momentumTransportModel& turbModel =
96  (
98  (
99  momentumTransportModel::typeName,
100  internalField().group()
101  )
102  );
103 
105  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
106 
107  const scalarField& y = turbModel.y()[patchi];
108 
109  const tmp<volScalarField> tk = turbModel.k();
110  const volScalarField& k = tk();
111 
112  const tmp<scalarField> tnuw = turbModel.nu(patchi);
113  const scalarField& nuw = tnuw();
114 
115  const scalar Cmu25 = pow025(nutw.Cmu());
116 
117  scalarField& v2 = *this;
118 
119  // Set v2 wall values
120  forAll(v2, facei)
121  {
122  label celli = patch().faceCells()[facei];
123 
124  scalar uTau = Cmu25*sqrt(k[celli]);
125 
126  scalar yPlus = uTau*y[facei]/nuw[facei];
127 
128  if (yPlus > nutw.yPlusLam())
129  {
130  scalar Cv2 = 0.193;
131  scalar Bv2 = -0.94;
132  v2[facei] = Cv2/nutw.kappa()*log(yPlus) + Bv2;
133  }
134  else
135  {
136  scalar Cv2 = 0.193;
137  v2[facei] = Cv2*pow4(yPlus);
138  }
139 
140  v2[facei] *= sqr(uTau);
141  }
142 
144 
145  // TODO: perform averaging for cells sharing more than one boundary face
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
152 (
155 );
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace RASModels
160 } // End namespace Foam
161 
162 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dictionary dict
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
dimensionedScalar log(const dimensionedScalar &ds)
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
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:369
This boundary condition provides a turbulence stress normal to streamlines wall function condition fo...
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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.
const volScalarField::Boundary & y() const
Return the near wall distance.
Macros for easy insertion into run-time selection tables.
scalar uTau
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.
v2WallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
makePatchTypeField(fvPatchScalarField, fWallFunctionFvPatchScalarField)
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
label patchi
scalar yPlus
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow4(const dimensionedScalar &ds)
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:147
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:357