PrandtlDelta.H
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-2015 OpenFOAM Foundation
9 Copyright (C) 2019 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
27Class
28 Foam::PrandtlDelta
29
30Description
31 Apply Prandtl mixing-length based damping function to the specified
32 geometric delta to improve near-wall behavior or LES models.
33
34 \verbatim
35 delta = min(geometricDelta, (kappa/Cdelta)*y)
36 \endverbatim
37
38 Example specification in the turbulenceProperties dictionary:
39 \verbatim
40 delta Prandtl;
41
42 PrandtlCoeffs
43 {
44 delta cubeRootVol;
45
46 cubeRootVolCoeffs
47 {
48 deltaCoeff 1;
49 }
50
51 // Default coefficients
52 kappa 0.41;
53 Cdelta 0.158;
54 }
55 \endverbatim
56
57SourceFiles
58 PrandtlDelta.C
59
60\*---------------------------------------------------------------------------*/
61
62#ifndef PrandtlDelta_H
63#define PrandtlDelta_H
64
65#include "LESdelta.H"
66
67// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68
69namespace Foam
70{
71namespace LESModels
72{
73
74/*---------------------------------------------------------------------------*\
75 Class PrandtlDelta Declaration
76\*---------------------------------------------------------------------------*/
78class PrandtlDelta
79:
80 public LESdelta
81{
82 // Private data
83
84 autoPtr<LESdelta> geometricDelta_;
85 scalar kappa_;
86 scalar Cdelta_;
87
88
89 // Private Member Functions
90
91 //- No copy construct
92 PrandtlDelta(const PrandtlDelta&) = delete;
93
94 //- No copy assignment
95 void operator=(const PrandtlDelta&) = delete;
96
97 // Calculate the delta values
98 void calcDelta();
99
100
101public:
102
103 //- Runtime type information
104 TypeName("Prandtl");
105
106
107 // Constructors
108
109 //- Construct from name, turbulenceModel and dictionary
111 (
112 const word& name,
114 const dictionary&
115 );
116
117
118 //- Destructor
119 virtual ~PrandtlDelta() = default;
120
121
122 // Member Functions
123
124 //- Read the LESdelta dictionary
125 virtual void read(const dictionary&);
126
127 // Correct values
128 virtual void correct();
129};
130
131
132// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133
134} // End namespace LESModels
135} // End namespace Foam
136
137// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138
139#endif
140
141// ************************************************************************* //
TypeName("Prandtl")
Runtime type information.
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: PrandtlDelta.C:92
virtual ~PrandtlDelta()=default
Destructor.
Abstract base class for LES deltas.
Definition: LESdelta.H:54
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition: LESdelta.H:134
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73