PrandtlDelta.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 "PrandtlDelta.H"
27 #include "wallDist.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
37  addToRunTimeSelectionTable(LESdelta, PrandtlDelta, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 {
46  delta_ = min
47  (
48  static_cast<const volScalarField&>(geometricDelta_()),
50  );
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const word& name,
60  const dictionary& dict
61 )
62 :
64  geometricDelta_
65  (
67  (
68  name,
69  turbulence,
70  dict.subDict(type() + "Coeffs")
71  )
72  ),
73  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
74  Cdelta_
75  (
76  dict.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
77  )
78 {
79  calcDelta();
80 }
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
86 {
87  const dictionary& coeffDict(dict.subDict(type() + "Coeffs"));
88 
89  geometricDelta_().read(coeffDict);
90  dict.readIfPresent<scalar>("kappa", kappa_);
91  coeffDict.readIfPresent<scalar>("Cdelta", Cdelta_);
92  calcDelta();
93 }
94 
95 
97 {
98  geometricDelta_().correct();
99 
100  if (turbulenceModel_.mesh().changing())
101  {
102  calcDelta();
103  }
104 }
105 
106 
107 // ************************************************************************* //
wallDist.H
Foam::LESModels::PrandtlDelta::Cdelta_
scalar Cdelta_
Definition: PrandtlDelta.H:82
PrandtlDelta.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::LESModels::PrandtlDelta::kappa_
scalar kappa_
Definition: PrandtlDelta.H:81
Foam::wallDist::y
const volScalarField & y() const
Return reference to cached distance-to-wall field.
Definition: wallDist.H:142
Foam::LESModels::PrandtlDelta::calcDelta
void calcDelta()
Definition: PrandtlDelta.C:44
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::LESModels::PrandtlDelta::geometricDelta_
autoPtr< LESdelta > geometricDelta_
Definition: PrandtlDelta.H:80
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:59
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
PrandtlDelta
Apply Prandtl mixing-length based damping function to the specified geometric delta to improve near-w...
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::LESdelta::turbulenceModel_
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:57
Foam::LESModels::PrandtlDelta::read
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: PrandtlDelta.C:85
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::PrandtlDelta::PrandtlDelta
PrandtlDelta(const PrandtlDelta &)
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::LESModels::PrandtlDelta::correct
virtual void correct()
Definition: PrandtlDelta.C:96
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::turbulenceModel::mesh
const fvMesh & mesh() const
Definition: turbulenceModel.H:127
Foam::LESModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
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