IDDESDelta.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) 2016-2020 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 "IDDESDelta.H"
31#include "wallDist.H"
32#include "maxDeltaxyz.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace LESModels
39{
42}
43}
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48void Foam::LESModels::IDDESDelta::calcDelta()
49{
50 const volScalarField& hmax = hmaxPtr_();
51 const fvMesh& mesh = turbulenceModel_.mesh();
52
53 // Wall-normal vectors
54 const volVectorField& n = wallDist::New(mesh).n();
55
56 tmp<volScalarField> tfaceToFacenMax
57 (
59 (
60 IOobject
61 (
62 "faceToFaceMax",
63 mesh.time().timeName(),
64 mesh,
67 ),
68 mesh,
70 )
71 );
72
73 scalarField& faceToFacenMax = tfaceToFacenMax.ref().primitiveFieldRef();
74
75 const cellList& cells = mesh.cells();
76 const vectorField& faceCentres = mesh.faceCentres();
77
78 forAll(cells, celli)
79 {
80 scalar maxDelta = 0.0;
81 const labelList& cFaces = cells[celli];
82 const vector nci = n[celli];
83
84 forAll(cFaces, cFacei)
85 {
86 label facei = cFaces[cFacei];
87 const point& fci = faceCentres[facei];
88
89 forAll(cFaces, cFacej)
90 {
91 label facej = cFaces[cFacej];
92 const point& fcj = faceCentres[facej];
93 scalar ndfc = nci & (fcj - fci);
94
95 if (ndfc > maxDelta)
96 {
97 maxDelta = ndfc;
98 }
99 }
100 }
101
102 faceToFacenMax[celli] = maxDelta;
103 }
104
105
106 label nD = mesh.nGeometricD();
107
108 if (nD == 2)
109 {
111 << "Case is 2D, LES is not strictly applicable" << nl
112 << endl;
113 }
114 else if (nD != 3)
115 {
117 << "Case must be either 2D or 3D" << exit(FatalError);
118 }
119
121 min
122 (
123 max
124 (
125 max
126 (
127 Cw_*wallDist::New(mesh).y(),
128 Cw_*hmax
129 ),
130 tfaceToFacenMax
131 ),
132 hmax
133 );
134
135 // Handle coupled boundaries
137}
138
139
140// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141
143(
144 const word& name,
145 const turbulenceModel& turbulence,
146 const dictionary& dict
147)
148:
150 hmaxPtr_(nullptr),
151 Cw_
152 (
153 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
154 (
155 "Cw",
156 0.15
157 )
158 )
159{
160 if (dict.optionalSubDict(type() + "Coeffs").found("hmax"))
161 {
162 // User-defined hmax
163 hmaxPtr_ =
165 (
168 dict.optionalSubDict(type() + "Coeffs"),
169 "hmax"
170 );
171 }
172 else
173 {
174 Info<< "Employing " << maxDeltaxyz::typeName << " for hmax" << endl;
175
176 hmaxPtr_.reset
177 (
178 new maxDeltaxyz
179 (
182 dict.optionalSubDict(type() + "Coeffs")
183 )
184 );
185 }
186
187 calcDelta();
188}
189
190
191// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192
194{
195 const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
196
197 coeffsDict.readIfPresent<scalar>("Cw", Cw_);
198
199 calcDelta();
200}
201
202
204{
205 if (turbulenceModel_.mesh().changing())
206 {
207 hmaxPtr_->correct();
208 calcDelta();
209 }
210}
211
212
213// ************************************************************************* //
scalar y
bool found
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
void correctBoundaryConditions()
Correct boundary field.
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:282
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const volScalarField & hmax() const
Return the hmax delta field.
Definition: IDDESDelta.H:109
Abstract base class for LES deltas.
Definition: LESdelta.H:54
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition: LESdelta.H:134
volScalarField delta_
Definition: LESdelta.H:61
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:59
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:577
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:883
const vectorField & faceCentres() const
const cellList & cells() const
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
const fvMesh & mesh() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
compressible::turbulenceModel & turbulence
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
List< label > labelList
A List of labels.
Definition: List.H:66
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.