Loading...
Searching...
No Matches
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"
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_();
52
53 // Wall-normal vectors
55
56 auto tfaceToFacenMax = volScalarField::New
57 (
58 "faceToFaceMax",
60 mesh,
62 );
63
64 scalarField& faceToFacenMax = tfaceToFacenMax.ref().primitiveFieldRef();
65
66 const cellList& cells = mesh.cells();
67 const vectorField& faceCentres = mesh.faceCentres();
68
69 forAll(cells, celli)
70 {
71 scalar maxDelta = 0.0;
72 const labelList& cFaces = cells[celli];
73 const vector nci = n[celli];
74
75 forAll(cFaces, cFacei)
76 {
77 label facei = cFaces[cFacei];
78 const point& fci = faceCentres[facei];
79
80 forAll(cFaces, cFacej)
81 {
82 label facej = cFaces[cFacej];
83 const point& fcj = faceCentres[facej];
84 scalar ndfc = nci & (fcj - fci);
85
86 if (ndfc > maxDelta)
87 {
88 maxDelta = ndfc;
89 }
90 }
91 }
92
93 faceToFacenMax[celli] = maxDelta;
94 }
95
96
97 label nD = mesh.nGeometricD();
98
99 if (nD == 2)
100 {
102 << "Case is 2D, LES is not strictly applicable" << nl
103 << endl;
104 }
105 else if (nD != 3)
106 {
108 << "Case must be either 2D or 3D" << exit(FatalError);
109 }
110
111 delta_.primitiveFieldRef() =
112 min
113 (
114 max
115 (
116 max
117 (
118 Cw_*wallDist::New(mesh).y(),
119 Cw_*hmax
120 ),
121 tfaceToFacenMax
122 ),
123 hmax
124 );
125
126 // Handle coupled boundaries
127 delta_.correctBoundaryConditions();
128}
129
130
131// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
132
133Foam::LESModels::IDDESDelta::IDDESDelta
134(
135 const word& name,
137 const dictionary& dict
138)
139:
141 hmaxPtr_(nullptr),
142 Cw_
143 (
144 dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
145 (
146 "Cw",
147 0.15
148 )
149 )
150{
151 if (dict.optionalSubDict(type() + "Coeffs").found("hmax"))
152 {
153 // User-defined hmax
154 hmaxPtr_ =
156 (
157 IOobject::groupName("hmax", turbulence.U().group()),
159 dict.optionalSubDict(type() + "Coeffs"),
160 "hmax"
161 );
162 }
163 else
164 {
165 Info<< "Employing " << maxDeltaxyz::typeName << " for hmax" << endl;
166
167 hmaxPtr_.reset
168 (
169 new maxDeltaxyz
170 (
171 IOobject::groupName("hmax", turbulence.U().group()),
173 dict.optionalSubDict(type() + "Coeffs")
174 )
175 );
176 }
178 calcDelta();
179}
180
181
182// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
183
185{
186 const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
188 coeffsDict.readIfPresent<scalar>("Cw", Cw_);
189
190 calcDelta();
191}
192
193
195{
196 if (turbulenceModel_.mesh().changing())
197 {
198 hmaxPtr_->correct();
199 calcDelta();
200 }
201}
202
203
204// ************************************************************************* //
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.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
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:122
void read(const dictionary &)
Read the LESdelta dictionary.
Definition IDDESDelta.C:177
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...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Delta calculated by taking the maximum distance between the cell centre and any face centre....
Abstract base class for turbulence models (RAS, LES and laminar).
const fvMesh & mesh() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition wallDist.C:168
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for LES SGS models.
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:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299