vanDriestDelta.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) 2011-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 
26 #include "vanDriestDelta.H"
27 #include "wallFvPatch.H"
28 #include "fvPatchDistWave.H"
29 #include "WallLocationYPlus.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
40 }
41 }
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::LESModels::vanDriestDelta::calcDelta()
46 {
47  const fvMesh& mesh = momentumTransportModel_.mesh();
48 
50  const tmp<volScalarField> tnu = momentumTransportModel_.nu();
51  const volScalarField& nu = tnu();
52  tmp<volScalarField> nuSgs = momentumTransportModel_.nut();
53 
54  volScalarField yStar
55  (
56  IOobject
57  (
58  "yStar",
59  mesh.time().constant(),
60  mesh
61  ),
62  mesh,
64  );
65 
66  const fvPatchList& patches = mesh.boundary();
67  volScalarField::Boundary& yStarBf = yStar.boundaryFieldRef();
68 
70  {
71  if (isA<wallFvPatch>(patches[patchi]))
72  {
73  const fvPatchVectorField& Uw = U.boundaryField()[patchi];
74  const scalarField& nuw = nu.boundaryField()[patchi];
75  const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
76 
77  yStarBf[patchi] =
78  nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + vSmall);
79  }
80  }
81 
83  (
85  );
86 
87  WallLocationYPlus<wallPoint>::trackData td;
88  td.yPlusCutOff = yPlusCutOff_;
89 
90  fvPatchDistWave::calculateAndCorrect<WallLocationYPlus>
91  (
92  mesh,
93  mesh.boundaryMesh().findPatchIDs<wallPolyPatch>(),
94  minWallFaceFraction_,
95  2, // <-- roughly equivalent to old point-cell corrections
96  y,
97  yStar,
98  td
99  );
100 
101  delta_ = min
102  (
103  static_cast<const volScalarField&>(geometricDelta_()),
104  (kappa_/Cdelta_)*((scalar(1) + small) - exp(-y/yStar/Aplus_))*y
105  );
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110 
112 (
113  const word& name,
115  const dictionary& dict
116 )
117 :
119  geometricDelta_
120  (
121  LESdelta::New
122  (
123  IOobject::groupName("geometricDelta", turbulence.U().group()),
124  turbulence,
125  dict.optionalSubDict(type() + "Coeffs")
126  )
127  ),
128  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
129  Aplus_
130  (
131  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
132  (
133  "Aplus",
134  26.0
135  )
136  ),
137  Cdelta_
138  (
139  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
140  (
141  "Cdelta",
142  0.158
143  )
144  ),
145  calcInterval_
146  (
147  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<label>
148  (
149  "calcInterval",
150  1
151  )
152  ),
153  yPlusCutOff_
154  (
155  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
156  (
157  "yPlusCutOff",
158  500
159  )
160  ),
161  minWallFaceFraction_
162  (
163  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
164  (
165  "minWallFaceFraction",
166  0.1
167  )
168  )
169 {
170  delta_ = geometricDelta_();
171 }
172 
173 
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175 
177 {
178  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
179 
180  geometricDelta_().read(coeffsDict);
181  dict.readIfPresent<scalar>("kappa", kappa_);
182  coeffsDict.readIfPresent<scalar>("Aplus", Aplus_);
183  coeffsDict.readIfPresent<scalar>("Cdelta", Cdelta_);
184  coeffsDict.readIfPresent<label>("calcInterval", calcInterval_);
185 
186  calcDelta();
187 }
188 
189 
191 {
192  if (momentumTransportModel_.mesh().time().timeIndex() % calcInterval_ == 0)
193  {
194  geometricDelta_().correct();
195  calcDelta();
196  }
197 }
198 
199 
200 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Apply van Driest damping function to the specified geometric delta to improve near-wall behavior or L...
vanDriestDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, momentumTransportModel and dictionary.
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Abstract base class for LES deltas.
Definition: LESdelta.H:51
volScalarField delta_
Definition: LESdelta.H:58
const momentumTransportModel & momentumTransportModel_
Definition: LESdelta.H:56
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
A class for handling words, derived from string.
Definition: word.H:62
label patchi
const fvPatchList & patches
U
Definition: pEqn.H:72
defineTypeNameAndDebug(cubeRootVolDelta, 0)
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
const char *const group
Group name for atomic constants.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar exp(const dimensionedScalar &ds)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
fvPatchField< vector > fvPatchVectorField
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))