Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2021 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 "PrandtlDelta.H"
30#include "wallDist.H"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace LESModels
38{
41}
42}
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
47void Foam::LESModels::PrandtlDelta::calcDelta()
48{
49 delta_ = min
50 (
51 static_cast<const volScalarField&>(geometricDelta_()),
52 (kappa_/Cdelta_)*wallDist::New(turbulenceModel_.mesh()).y()
53 );
54}
55
56
57// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58
59Foam::LESModels::PrandtlDelta::PrandtlDelta
60(
61 const word& name,
63 const dictionary& dict
64)
65:
67 geometricDelta_
68 (
70 (
71 IOobject::groupName("geometricDelta", turbulence.U().group()),
73 dict.optionalSubDict(type() + "Coeffs")
74 )
75 ),
76 kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
77 Cdelta_
78 (
79 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
80 (
81 "Cdelta",
82 0.158
83 )
84 )
86 calcDelta();
87}
88
89
90// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91
93{
94 const dictionary& coeffDict(dict.optionalSubDict(type() + "Coeffs"));
95
96 geometricDelta_().read(coeffDict);
97 dict.readIfPresent<scalar>("kappa", kappa_);
98 coeffDict.readIfPresent<scalar>("Cdelta", Cdelta_);
99 calcDelta();
100}
101
102
104{
105 geometricDelta_().correct();
106
107 if (turbulenceModel_.mesh().changing())
108 {
109 calcDelta();
110 }
111}
112
113
114// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Abstract base class for LES deltas.
Definition LESdelta.H:50
LESdelta(const LESdelta &)=delete
No copy construct.
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict, const word &lookupName="delta")
Return a reference to the selected LES delta.
Definition LESdelta.C:61
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition LESdelta.H:146
volScalarField delta_
Definition LESdelta.H:57
const turbulenceModel & turbulenceModel_
Definition LESdelta.H:55
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Abstract base class for turbulence models (RAS, LES and laminar).
const fvMesh & mesh() const
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field.
Definition wallDist.H:220
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
Namespace for LES SGS models.
Namespace for OpenFOAM.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dictionary dict