PrandtlDelta.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "PrandtlDelta.H"
27 #include "wallDist.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36  defineTypeNameAndDebug(PrandtlDelta, 0);
37  addToRunTimeSelectionTable(LESdelta, PrandtlDelta, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::LESModels::PrandtlDelta::calcDelta()
45 {
46  delta_ = min
47  (
48  static_cast<const volScalarField&>(geometricDelta_()),
49  (kappa_/Cdelta_)*wallDist::New(turbulenceModel_.mesh()).y()
50  );
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 Foam::LESModels::PrandtlDelta::PrandtlDelta
57 (
58  const word& name,
59  const turbulenceModel& turbulence,
60  const dictionary& dict
61 )
62 :
63  LESdelta(name, turbulence),
64  geometricDelta_
65  (
67  (
68  name,
69  turbulence,
70  dict.optionalSubDict(type() + "Coeffs")
71  )
72  ),
73  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
74  Cdelta_
75  (
76  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
77  (
78  "Cdelta",
79  0.158
80  )
81  )
82 {
83  calcDelta();
84 }
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
90 {
91  const dictionary& coeffDict(dict.optionalSubDict(type() + "Coeffs"));
92 
93  geometricDelta_().read(coeffDict);
94  dict.readIfPresent<scalar>("kappa", kappa_);
95  coeffDict.readIfPresent<scalar>("Cdelta", Cdelta_);
96  calcDelta();
97 }
98 
99 
101 {
102  geometricDelta_().correct();
103 
104  if (turbulenceModel_.mesh().changing())
105  {
106  calcDelta();
107  }
108 }
109 
110 
111 // ************************************************************************* //
bool read(Istream &)
Read dictionary from Istream.
Definition: dictionaryIO.C:125
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Abstract base class for LES deltas.
Definition: LESdelta.H:50
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
scalar y
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
A class for handling words, derived from string.
Definition: word.H:59
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
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
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
defineTypeNameAndDebug(cubeRootVolDelta, 0)
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: PrandtlDelta.C:89
Namespace for OpenFOAM.