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-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"
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  const dictionary& dict
46 )
47 :
48  fixedValueFvPatchField<scalar>(p, iF, dict)
49 {}
50 
51 
53 (
55  const fvPatch& p,
57  const fieldMapper& mapper
58 )
59 :
60  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
61 {}
62 
63 
65 (
66  const fWallFunctionFvPatchScalarField& v2wfpsf,
68 )
69 :
70  fixedValueFvPatchField<scalar>(v2wfpsf, iF)
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
78  if (updated())
79  {
80  return;
81  }
82 
83  const label patchi = patch().index();
84 
85  const momentumTransportModel& turbModel =
87 
88  const v2fBase& v2fModel = refCast<const v2fBase>(turbModel);
89 
92 
93  const scalarField& y = turbModel.yb()[patchi];
94 
95  const tmp<volScalarField> tk = turbModel.k();
96  const volScalarField& k = tk();
97 
98  const tmp<volScalarField> tepsilon = turbModel.epsilon();
99  const volScalarField& epsilon = tepsilon();
100 
101  const tmp<volScalarField> tv2 = v2fModel.v2();
102  const volScalarField& v2 = tv2();
103 
104  const tmp<scalarField> tnuw = turbModel.nu(patchi);
105  const scalarField& nuw = tnuw();
106 
107  const scalar Cmu25 = pow025(nutw.Cmu());
108 
109  scalarField& f = *this;
110 
111  // Set f wall values
112  forAll(f, facei)
113  {
114  label celli = patch().faceCells()[facei];
115 
116  scalar uTau = Cmu25*sqrt(k[celli]);
117 
118  scalar yPlus = uTau*y[facei]/nuw[facei];
119 
120  if (yPlus > nutw.yPlusLam())
121  {
122  scalar N = 6.0;
123  scalar v2c = v2[celli];
124  scalar epsc = epsilon[celli];
125  scalar kc = k[celli];
126 
127  f[facei] = N*v2c*epsc/(sqr(kc) + rootVSmall);
128  f[facei] /= sqr(uTau) + rootVSmall;
129  }
130  else
131  {
132  f[facei] = 0.0;
133  }
134  }
135 
137 
138  // TODO: perform averaging for cells sharing more than one boundary face
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 
145 (
148 );
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 } // End namespace RASModels
153 } // End namespace Foam
154 
155 // ************************************************************************* //
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 damping function, f, wall function condition for low- a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields.
Definition: v2fBase.H:55
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.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
static scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
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
const scalar epsilon
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 sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
labelList f(nPoints)
dictionary dict
volScalarField & p