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 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"
31 #include "wallDistData.H"
32 #include "wallPointYPlus.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace LESModels
40 {
42  addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
43 }
44 }
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::LESModels::vanDriestDelta::calcDelta()
49 {
50  const fvMesh& mesh = turbulenceModel_.mesh();
51 
53  const tmp<volScalarField> tnu = turbulenceModel_.nu();
54  const volScalarField& nu = tnu();
55  tmp<volScalarField> nuSgs = turbulenceModel_.nut();
56 
57  volScalarField ystar
58  (
59  IOobject
60  (
61  "ystar",
62  mesh.time().constant(),
63  mesh
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  scalar cutOff = wallPointYPlus::yPlusCutOff;
87  wallDistData<wallPointYPlus> y(mesh, ystar);
89 
90  delta_ = min
91  (
92  static_cast<const volScalarField&>(geometricDelta_()),
93  (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
94  );
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99 
100 Foam::LESModels::vanDriestDelta::vanDriestDelta
101 (
102  const word& name,
104  const dictionary& dict
105 )
106 :
108  geometricDelta_
109  (
111  (
112  IOobject::groupName("geometricDelta", turbulence.U().group()),
113  turbulence,
114  // Note: cannot use optionalSubDict - if no *Coeffs dict the
115  // code will get stuck in a loop attempting to read the delta entry
116  // - consider looking up "geometricDelta" instead of "delta"?
117  dict.subDict(type() + "Coeffs")
118  )
119  ),
120  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
121  Aplus_
122  (
123  dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
124  (
125  "Aplus",
126  26.0
127  )
128  ),
129  Cdelta_
130  (
131  dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
132  (
133  "Cdelta",
134  0.158
135  )
136  ),
137  calcInterval_
138  (
139  dict.optionalSubDict(type() + "Coeffs").getOrDefault<label>
140  (
141  "calcInterval",
142  1
143  )
144  )
145 {
146  delta_ = geometricDelta_();
147 }
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 {
154  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
155 
156  geometricDelta_().read(coeffsDict);
157  dict.readIfPresent<scalar>("kappa", kappa_);
158  coeffsDict.readIfPresent<scalar>("Aplus", Aplus_);
159  coeffsDict.readIfPresent<scalar>("Cdelta", Cdelta_);
160  coeffsDict.readIfPresent<label>("calcInterval", calcInterval_);
161 
162  calcDelta();
163 }
164 
165 
167 {
168  if (turbulenceModel_.mesh().time().timeIndex() % calcInterval_ == 0)
169  {
170  geometricDelta_().correct();
171  calcDelta();
172  }
173 }
174 
175 
176 // ************************************************************************* //
Foam::LESdelta::New
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict, const word &lookupName="delta")
Return a reference to the selected LES delta.
Definition: LESdelta.C:69
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::wallPointYPlus::yPlusCutOff
static scalar yPlusCutOff
The cut-off value for y+.
Definition: wallPointYPlus.H:84
wallDistData.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:61
wallFvPatch.H
Foam::LESModels::vanDriestDelta::correct
virtual void correct()
Definition: vanDriestDelta.C:166
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fvPatchList
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:47
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::LESdelta::turbulenceModel_
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:59
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
vanDriestDelta.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::LESModels::vanDriestDelta::read
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: vanDriestDelta.C:152
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
vanDriestDelta
Simple cube-root of cell volume delta used in incompressible LES models.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dictionary::read
bool read(Istream &is)
Read dictionary from Istream. Discards the header.
Definition: dictionaryIO.C:141
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::fvPatchVectorField
fvPatchField< vector > fvPatchVectorField
Definition: fvPatchFieldsFwd.H:43
U
U
Definition: pEqn.H:72
Foam::LESdelta
Abstract base class for LES deltas.
Definition: LESdelta.H:53
Foam::LESModels::defineTypeNameAndDebug
defineTypeNameAndDebug(cubeRootVolDelta, 0)
wallPointYPlus.H
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::turbulenceModel::mesh
const fvMesh & mesh() const
Definition: turbulenceModel.H:129
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::LESModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:144
y
scalar y
Definition: LISASMDCalcMethod1.H:14