maxDeltaxyzCubeRootLESDelta.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) 2016-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace LESModels
36{
39 (
43 );
44}
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
50void Foam::LESModels::maxDeltaxyzCubeRootLESDelta::calcDelta()
51{
52 maxDeltaxyz_.calcDelta();
53 cubeRootVolDelta_.calcDelta();
54
55 delta_ =
56 max
57 (
58 static_cast<const volScalarField&>(maxDeltaxyz_),
59 static_cast<const volScalarField&>(cubeRootVolDelta_)
60 );
61
62 // Handle coupled boundaries
64}
65
66
67// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68
70(
71 const word& name,
72 const turbulenceModel& turbulence,
73 const dictionary& dict
74)
75:
77 maxDeltaxyz_
78 (
79 name + "maxDeltaxyz",
81 dict.subDict(typeName + "Coeffs")
82 ),
83 cubeRootVolDelta_
84 (
85 name + "cubeRootVolDelta",
87 dict.subDict(typeName + "Coeffs")
88 )
89{
90 calcDelta();
91}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
97{
98 maxDeltaxyz_.read(dict.subDict(typeName + "Coeffs"));
99 cubeRootVolDelta_.read(dict.subDict(typeName + "Coeffs"));
100
101 calcDelta();
102}
103
104
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.
void correctBoundaryConditions()
Correct boundary field.
void calcDelta()
Calculate the delta values.
void calcDelta()
Calculate the delta values.
Definition: maxDeltaxyz.C:46
Abstract base class for LES deltas.
Definition: LESdelta.H:54
volScalarField delta_
Definition: LESdelta.H:61
virtual bool read()
Re-read model coefficients if they have changed.
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
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
compressible::turbulenceModel & turbulence
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict