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-2023 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  const dictionary& dict
45 )
46 :
47  fixedValueFvPatchField<scalar>(p, iF, dict)
48 {}
49 
50 
52 (
54  const fvPatch& p,
56  const fieldMapper& mapper
57 )
58 :
59  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
60 {}
61 
62 
64 (
65  const v2WallFunctionFvPatchScalarField& v2wfpsf,
67 )
68 :
69  fixedValueFvPatchField<scalar>(v2wfpsf, iF)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
77  if (updated())
78  {
79  return;
80  }
81 
82  const label patchi = patch().index();
83 
84 
85  const momentumTransportModel& turbModel =
87 
90 
91  const scalarField& y = turbModel.yb()[patchi];
92 
93  const tmp<volScalarField> tk = turbModel.k();
94  const volScalarField& k = tk();
95 
96  const tmp<scalarField> tnuw = turbModel.nu(patchi);
97  const scalarField& nuw = tnuw();
98 
99  const scalar Cmu25 = pow025(nutw.Cmu());
100 
101  scalarField& v2 = *this;
102 
103  // Set v2 wall values
104  forAll(v2, facei)
105  {
106  label celli = patch().faceCells()[facei];
107 
108  scalar uTau = Cmu25*sqrt(k[celli]);
109 
110  scalar yPlus = uTau*y[facei]/nuw[facei];
111 
112  if (yPlus > nutw.yPlusLam())
113  {
114  scalar Cv2 = 0.193;
115  scalar Bv2 = -0.94;
116  v2[facei] = Cv2/nutw.kappa()*log(yPlus) + Bv2;
117  }
118  else
119  {
120  scalar Cv2 = 0.193;
121  v2[facei] = Cv2*pow4(yPlus);
122  }
123 
124  v2[facei] *= sqr(uTau);
125  }
126 
128 
129  // TODO: perform averaging for cells sharing more than one boundary face
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 
136 (
139 );
140 
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 
143 } // End namespace RASModels
144 } // End namespace Foam
145 
146 // ************************************************************************* //
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
This boundary condition provides a turbulence stress normal to streamlines wall function condition fo...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
v2WallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:362
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:143
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:356
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:374
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
Abstract base class for turbulence models (RAS, LES and laminar).
const volScalarField::Boundary & yb() const
Return the near wall distance.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
const Type & lookupType(const word &group=word::null) const
Lookup and return the object of the given Type.
A class for managing temporary objects.
Definition: tmp.H:55
const scalar yPlus
const scalar uTau
label patchi
makePatchTypeField(fvPatchScalarField, fWallFunctionFvPatchScalarField)
Namespace for OpenFOAM.
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & p