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-------------------------------------------------------------------------------
11License
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
37namespace Foam
38{
39namespace LESModels
40{
43}
44}
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48void 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
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
101(
102 const word& name,
103 const turbulenceModel& turbulence,
104 const dictionary& dict
105)
106:
108 geometricDelta_
109 (
111 (
112 IOobject::groupName("geometricDelta", turbulence.U().group()),
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// ************************************************************************* //
scalar y
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Abstract base class for LES deltas.
Definition: LESdelta.H:54
volScalarField delta_
Definition: LESdelta.H:61
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:59
virtual bool read()
Re-read model coefficients if they have changed.
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const fvMesh & mesh() const
static scalar yPlusCutOff
The cut-off value for y+.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
compressible::turbulenceModel & turbulence
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:64
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
fvPatchField< vector > fvPatchVectorField
volScalarField & nu
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333