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 | Copyright (C) 2011-2015 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 {
39  addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
40 }
41 }
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 {
47  const fvMesh& mesh = turbulenceModel_.mesh();
48 
51  const volScalarField& nu = tnu();
53 
54  volScalarField ystar
55  (
56  IOobject
57  (
58  "ystar",
59  mesh.time().constant(),
60  mesh
61  ),
62  mesh,
63  dimensionedScalar("ystar", dimLength, GREAT)
64  );
65 
66  const fvPatchList& patches = mesh.boundary();
68  {
69  if (isA<wallFvPatch>(patches[patchi]))
70  {
71  const fvPatchVectorField& Uw = U.boundaryField()[patchi];
72  const scalarField& nuw = nu.boundaryField()[patchi];
73  const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
74 
75  ystar.boundaryField()[patchi] =
76  nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + VSMALL);
77  }
78  }
79 
80  scalar cutOff = wallPointYPlus::yPlusCutOff;
84 
85  delta_ = min
86  (
87  static_cast<const volScalarField&>(geometricDelta_()),
88  (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
89  );
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
96 (
97  const word& name,
99  const dictionary& dict
100 )
101 :
103  geometricDelta_
104  (
106  (
107  IOobject::groupName("geometricDelta", turbulence.U().group()),
108  turbulence,
109  dict.subDict(type() + "Coeffs")
110  )
111  ),
112  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
113  Aplus_
114  (
115  dict.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
116  ),
117  Cdelta_
118  (
119  dict.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
120  ),
121  calcInterval_
122  (
123  dict.subDict(type() + "Coeffs").lookupOrDefault<label>
124  (
125  "calcInterval",
126  1
127  )
128  )
129 {
130  delta_ = geometricDelta_();
131 }
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
137 {
138  const dictionary& coeffsDict(dict.subDict(type() + "Coeffs"));
139 
140  geometricDelta_().read(coeffsDict);
141  dict.readIfPresent<scalar>("kappa", kappa_);
142  coeffsDict.readIfPresent<scalar>("Aplus", Aplus_);
143  coeffsDict.readIfPresent<scalar>("Cdelta", Cdelta_);
144  coeffsDict.readIfPresent<label>("calcInterval", calcInterval_);
145 
146  calcDelta();
147 }
148 
149 
151 {
152  if (turbulenceModel_.mesh().time().timeIndex() % calcInterval_ == 0)
153  {
154  geometricDelta_().correct();
155  calcDelta();
156  }
157 }
158 
159 
160 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:200
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::wallPointYPlus::yPlusCutOff
static scalar yPlusCutOff
cut-off value for y+
Definition: wallPointYPlus.H:82
wallDistData.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::LESModels::vanDriestDelta::kappa_
scalar kappa_
Definition: vanDriestDelta.H:57
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::LESModels::vanDriestDelta::Aplus_
scalar Aplus_
Definition: vanDriestDelta.H:58
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:59
wallFvPatch.H
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::LESModels::vanDriestDelta::correct
virtual void correct()
Definition: vanDriestDelta.C:150
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::LESModels::vanDriestDelta::calcDelta
void calcDelta()
Definition: vanDriestDelta.C:45
U
U
Definition: pEqn.H:46
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::wallDistData
Wall distance calculation. Like wallDist but also transports extra data (template argument).
Definition: wallDistData.H:57
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::LESdelta::turbulenceModel_
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:57
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
vanDriestDelta.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::LESModels::vanDriestDelta::read
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: vanDriestDelta.C:136
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
vanDriestDelta
Simple cube-root of cell volume delta used in incompressible LES models.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::vanDriestDelta::vanDriestDelta
vanDriestDelta(const vanDriestDelta &)
Disallow default bitwise copy construct and assignment.
Foam::LESdelta
Abstract base class for LES deltas.
Definition: LESdelta.H:50
Foam::LESModels::defineTypeNameAndDebug
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Foam::LESdelta::New
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
Definition: LESdelta.C:66
Foam::LESModels::vanDriestDelta::geometricDelta_
autoPtr< LESdelta > geometricDelta_
Definition: vanDriestDelta.H:56
wallPointYPlus.H
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::turbulenceModel::mesh
const fvMesh & mesh() const
Definition: turbulenceModel.H:127
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::LESModels::vanDriestDelta::Cdelta_
scalar Cdelta_
Definition: vanDriestDelta.H:59
Foam::LESModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:142
y
scalar y
Definition: LISASMDCalcMethod1.H:14