vanDriestDelta.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "vanDriestDelta.H"
30 #include "wallFvPatch.H"
32 #include "wallDistAddressing.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace LESModels
39 {
40  defineTypeNameAndDebug(vanDriestDelta, 0);
41  addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
42 }
43 }
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::LESModels::vanDriestDelta::calcDelta()
48 {
49  const fvMesh& mesh = turbulenceModel_.mesh();
50 
52  const tmp<volScalarField> tnu = turbulenceModel_.nu();
53  const volScalarField& nu = tnu();
54  tmp<volScalarField> nuSgs = turbulenceModel_.nut();
55 
56  volScalarField ystar
57  (
58  IOobject
59  (
60  "ystar",
61  mesh.time().constant(),
62  mesh.thisDb(),
64  ),
65  mesh,
66  dimensionedScalar("ystar", dimLength, GREAT)
67  );
68 
69  const fvPatchList& patches = mesh.boundary();
70  volScalarField::Boundary& ystarBf = ystar.boundaryFieldRef();
71 
72  forAll(patches, patchi)
73  {
74  if (isA<wallFvPatch>(patches[patchi]))
75  {
76  const fvPatchVectorField& Uw = U.boundaryField()[patchi];
77  const scalarField& nuw = nu.boundaryField()[patchi];
78  const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
79 
80  ystarBf[patchi] =
81  nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + VSMALL);
82  }
83  }
84 
85  // Construct wall transporter
86  const auto& wDist = wallDistAddressing::New(mesh);
87 
88  // Get distance to nearest wall
89  const auto& y = wDist.y();
90 
91  // Get ystar from nearest wall
92  wDist.map(ystar, mapDistribute::transform());
93 
94  // Calculate y/ystar (stored in ystar!) and do the clipping
95  constexpr scalar yPlusCutOff = 500;
96  // Allow for some precision loss from transformation/interpolation of GREAT
97  // (= unvisited value)(though ystar is scalar so should not be transformed)
98  constexpr scalar fuzzyGREAT = 0.5*GREAT;
99 
100  ystar.dimensions().reset(y.dimensions()/ystar.dimensions());
101  forAll(y, celli)
102  {
103  const scalar yPlus = y[celli]/ystar[celli];
104  if (y[celli] > fuzzyGREAT || (yPlus > yPlusCutOff))
105  {
106  // unvisited : y is GREAT, ystar is 1.0
107  ystar[celli] = GREAT;
108  }
109  else
110  {
111  ystar[celli] = yPlus;
112  }
113  }
114 
115  forAll(y.boundaryField(), patchi)
116  {
117  const auto& yp = y.boundaryField()[patchi];
118  auto& ystarp = ystar.boundaryFieldRef()[patchi];
119 
120  forAll(yp, i)
121  {
122  const scalar yPlus = yp[i]/ystarp[i];
123  if (yp[i] > fuzzyGREAT || (yPlus > yPlusCutOff))
124  {
125  ystarp[i] = GREAT;
126  }
127  else
128  {
129  ystarp[i] = yPlus;
130  }
131  }
132  }
133  ystar.correctBoundaryConditions();
134 
135  // Note: y/ystar stored in ystar
136  delta_ = min
137  (
138  static_cast<const volScalarField&>(geometricDelta_()),
139  //(kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
140  (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-ystar/Aplus_))*y
141  );
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
146 
147 Foam::LESModels::vanDriestDelta::vanDriestDelta
148 (
149  const word& name,
151  const dictionary& dict
152 )
153 :
155  geometricDelta_
156  (
157  LESdelta::New
158  (
159  IOobject::groupName("geometricDelta", turbulence.U().group()),
160  turbulence,
161  // Note: cannot use optionalSubDict - if no *Coeffs dict the
162  // code will get stuck in a loop attempting to read the delta entry
163  // - consider looking up "geometricDelta" instead of "delta"?
164  dict.subDict(type() + "Coeffs")
165  )
166  ),
167  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
168  Aplus_
169  (
170  dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
171  (
172  "Aplus",
173  26.0
174  )
175  ),
176  Cdelta_
177  (
178  dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
179  (
180  "Cdelta",
181  0.158
182  )
183  )
184 {
185  calcInterval_ = 1;
186  const dictionary& coeffsDict = dict.optionalSubDict(type() + "Coeffs");
187  if (!coeffsDict.readIfPresent("calcInterval", calcInterval_))
188  {
189  coeffsDict.readIfPresent("updateInterval", calcInterval_);
190  }
192  delta_ = geometricDelta_();
193 }
194 
195 
196 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197 
199 {
200  const dictionary& coeffsDict = dict.optionalSubDict(type() + "Coeffs");
201 
202  geometricDelta_().read(coeffsDict);
203  dict.readIfPresent<scalar>("kappa", kappa_);
204  coeffsDict.readIfPresent<scalar>("Aplus", Aplus_);
205  coeffsDict.readIfPresent<scalar>("Cdelta", Cdelta_);
206  calcInterval_ = 1;
207  if (!coeffsDict.readIfPresent<label>("calcInterval", calcInterval_))
208  {
209  coeffsDict.readIfPresent("updateInterval", calcInterval_);
210  }
211 
212  calcDelta();
213 }
214 
215 
217 {
218  if
219  (
220  (calcInterval_ > 0)
221  && (turbulenceModel_.mesh().time().timeIndex() % calcInterval_) == 0
222  )
223  {
224  geometricDelta_().correct();
225  calcDelta();
226  }
227 }
228 
229 
230 // ************************************************************************* //
const fvMesh & mesh() const
dictionary dict
fvPatchField< vector > fvPatchVectorField
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const volVectorField & U() const
Access function to velocity field.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static const wallDistAddressing & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
Abstract base class for LES deltas.
Definition: LESdelta.H:49
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar y
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
volScalarField delta_
Definition: LESdelta.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(cubeRootVolDelta, 0)
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar yPlus
const polyBoundaryMesh & patches
virtual void read(const dictionary &)
Read the LESdelta dictionary.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:55
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
volScalarField & nu
Do not request registration (bool: false)
Namespace for OpenFOAM.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:59