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-2020 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 "wallDistData.H"
29 #include "wallPointYPlus.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38  defineTypeNameAndDebug(vanDriestDelta, 0);
39  addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
40 }
41 }
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::LESModels::vanDriestDelta::calcDelta()
46 {
47  const fvMesh& mesh = momentumTransportModel_.mesh();
48 
49  const volVectorField& U = momentumTransportModel_.U();
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 
69  forAll(patches, patchi)
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 
82  scalar cutOff = wallPointYPlus::yPlusCutOff;
84  wallDistData<wallPointYPlus> y(mesh, ystar);
86 
87  delta_ = min
88  (
89  static_cast<const volScalarField&>(geometricDelta_()),
90  (kappa_/Cdelta_)*((scalar(1) + small) - exp(-y/ystar/Aplus_))*y
91  );
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96 
98 (
99  const word& name,
100  const momentumTransportModel& turbulence,
101  const dictionary& dict
102 )
103 :
104  LESdelta(name, turbulence),
105  geometricDelta_
106  (
108  (
109  IOobject::groupName("geometricDelta", turbulence.U().group()),
110  turbulence,
111  dict.optionalSubDict(type() + "Coeffs")
112  )
113  ),
114  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
115  Aplus_
116  (
117  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
118  (
119  "Aplus",
120  26.0
121  )
122  ),
123  Cdelta_
124  (
125  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
126  (
127  "Cdelta",
128  0.158
129  )
130  ),
131  calcInterval_
132  (
133  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<label>
134  (
135  "calcInterval",
136  1
137  )
138  )
139 {
140  delta_ = geometricDelta_();
141 }
142 
143 
144 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145 
147 {
148  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
149 
150  geometricDelta_().read(coeffsDict);
151  dict.readIfPresent<scalar>("kappa", kappa_);
152  coeffsDict.readIfPresent<scalar>("Aplus", Aplus_);
153  coeffsDict.readIfPresent<scalar>("Cdelta", Cdelta_);
154  coeffsDict.readIfPresent<label>("calcInterval", calcInterval_);
155 
156  calcDelta();
157 }
158 
159 
161 {
162  if (momentumTransportModel_.mesh().time().timeIndex() % calcInterval_ == 0)
163  {
164  geometricDelta_().correct();
165  calcDelta();
166  }
167 }
168 
169 
170 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
#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
fvPatchField< vector > fvPatchVectorField
static autoPtr< LESdelta > New(const word &name, const momentumTransportModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
Definition: LESdelta.C:66
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Abstract base class for LES deltas.
Definition: LESdelta.H:50
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
dimensionedScalar sqrt(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
bool read(Istream &, const bool keepHeader=false)
Read dictionary from Istream, optionally keeping the header.
Definition: dictionaryIO.C:104
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
scalar y
vanDriestDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, momentumTransportModel and dictionary.
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1040
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
static scalar yPlusCutOff
cut-off value for y+
virtual void read(const dictionary &)
Read the LESdelta dictionary.
dimensioned< scalar > mag(const dimensioned< Type > &)
Namespace for OpenFOAM.