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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "vanDriestDelta.H"
29 #include "wallFvPatch.H"
30 #include "wallDistData.H"
31 #include "wallPointYPlus.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace LESModels
39 {
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
63  ),
64  mesh,
65  dimensionedScalar("ystar", dimLength, GREAT)
66  );
67 
68  const fvPatchList& patches = mesh.boundary();
69  volScalarField::Boundary& ystarBf = ystar.boundaryFieldRef();
70 
71  forAll(patches, patchi)
72  {
73  if (isA<wallFvPatch>(patches[patchi]))
74  {
75  const fvPatchVectorField& Uw = U.boundaryField()[patchi];
76  const scalarField& nuw = nu.boundaryField()[patchi];
77  const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
78 
79  ystarBf[patchi] =
80  nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + VSMALL);
81  }
82  }
83 
84  scalar cutOff = wallPointYPlus::yPlusCutOff;
86  wallDistData<wallPointYPlus> y(mesh, ystar);
88 
89  delta_ = min
90  (
91  static_cast<const volScalarField&>(geometricDelta_()),
92  (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
93  );
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
98 
99 Foam::LESModels::vanDriestDelta::vanDriestDelta
100 (
101  const word& name,
103  const dictionary& dict
104 )
105 :
107  geometricDelta_
108  (
110  (
111  IOobject::groupName("geometricDelta", turbulence.U().group()),
112  turbulence,
113  // Note: cannot use optionalSubDict - if no *Coeffs dict the
114  // code will get stuck in a loop attempting to read the delta entry
115  // - consider looking up "geometricDelta" instead of "delta"?
116  dict.subDict(type() + "Coeffs")
117  )
118  ),
119  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
120  Aplus_
121  (
122  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
123  (
124  "Aplus",
125  26.0
126  )
127  ),
128  Cdelta_
129  (
130  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
131  (
132  "Cdelta",
133  0.158
134  )
135  ),
136  calcInterval_
137  (
138  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<label>
139  (
140  "calcInterval",
141  1
142  )
143  )
144 {
145  delta_ = geometricDelta_();
146 }
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 {
153  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
154 
155  geometricDelta_().read(coeffsDict);
156  dict.readIfPresent<scalar>("kappa", kappa_);
157  coeffsDict.readIfPresent<scalar>("Aplus", Aplus_);
158  coeffsDict.readIfPresent<scalar>("Cdelta", Cdelta_);
159  coeffsDict.readIfPresent<label>("calcInterval", calcInterval_);
160 
161  calcDelta();
162 }
163 
164 
166 {
167  if (turbulenceModel_.mesh().time().timeIndex() % calcInterval_ == 0)
168  {
169  geometricDelta_().correct();
170  calcDelta();
171  }
172 }
173 
174 
175 // ************************************************************************* //
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:62
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:53
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:62
wallFvPatch.H
Foam::LESModels::vanDriestDelta::correct
virtual void correct()
Definition: vanDriestDelta.C:165
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:290
Foam::fvPatchList
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:47
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::label
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:62
Foam::LESdelta::turbulenceModel_
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:60
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
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:151
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
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.
Definition: dictionaryIO.C:141
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:60
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::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
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::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:144
y
scalar y
Definition: LISASMDCalcMethod1.H:14