IDDESDelta.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "IDDESDelta.H"
28 #include "wallDist.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36  defineTypeNameAndDebug(IDDESDelta, 0);
37  addToRunTimeSelectionTable(LESdelta, IDDESDelta, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::LESModels::IDDESDelta::calcDelta()
45 {
46  const volScalarField& hmax = hmax_;
47  const fvMesh& mesh = turbulenceModel_.mesh();
48 
49  // Wall-normal vectors
50  const volVectorField& n = wallDist::New(mesh).n();
51 
52  tmp<volScalarField> tfaceToFacenMax
53  (
54  new volScalarField
55  (
56  IOobject
57  (
58  "faceToFaceMax",
59  mesh.time().timeName(),
60  mesh,
63  ),
64  mesh,
65  dimensionedScalar("zero", dimLength, 0.0)
66  )
67  );
68 
69  scalarField& faceToFacenMax = tfaceToFacenMax.ref().primitiveFieldRef();
70 
71  const cellList& cells = mesh.cells();
72  const vectorField& faceCentres = mesh.faceCentres();
73 
74  forAll(cells, celli)
75  {
76  scalar maxDelta = 0.0;
77  const labelList& cFaces = cells[celli];
78  const vector nci = n[celli];
79 
80  forAll(cFaces, cFacei)
81  {
82  label facei = cFaces[cFacei];
83  const point& fci = faceCentres[facei];
84 
85  forAll(cFaces, cFacej)
86  {
87  label facej = cFaces[cFacej];
88  const point& fcj = faceCentres[facej];
89  scalar ndfc = nci & (fcj - fci);
90 
91  if (ndfc > maxDelta)
92  {
93  maxDelta = ndfc;
94  }
95  }
96  }
97 
98  faceToFacenMax[celli] = maxDelta;
99  }
100 
101 
102  label nD = mesh.nGeometricD();
103 
104  if (nD == 2)
105  {
107  << "Case is 2D, LES is not strictly applicable" << nl
108  << endl;
109  }
110  else if (nD != 3)
111  {
113  << "Case must be either 2D or 3D" << exit(FatalError);
114  }
115 
116  delta_.primitiveFieldRef() =
117  min
118  (
119  max
120  (
121  max
122  (
123  Cw_*wallDist::New(mesh).y(),
124  Cw_*hmax
125  ),
126  tfaceToFacenMax
127  ),
128  hmax
129  );
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134 
135 Foam::LESModels::IDDESDelta::IDDESDelta
136 (
137  const word& name,
138  const turbulenceModel& turbulence,
139  const dictionary& dict
140 )
141 :
142  LESdelta(name, turbulence),
143  hmax_
144  (
145  IOobject::groupName("hmax", turbulence.U().group()),
146  turbulence,
147  dict
148  ),
149  Cw_
150  (
151  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
152  (
153  "Cw",
154  0.15
155  )
156  )
157 {
158  calcDelta();
159 }
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
165 {
166  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
167 
168  coeffsDict.readIfPresent<scalar>("Cw", Cw_);
169 
170  calcDelta();
171 }
172 
173 
175 {
176  if (turbulenceModel_.mesh().changing())
177  {
178  hmax_.correct();
179  calcDelta();
180  }
181 }
182 
183 
184 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:176
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:170
const volVectorField & U() const
Access function to velocity field.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Abstract base class for LES deltas.
Definition: LESdelta.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void read(const dictionary &)
Read the LESdelta dictionary.
Definition: IDDESDelta.C:164
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
scalar y
dynamicFvMesh & mesh
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const char nl
Definition: Ostream.H:265
const Mesh & mesh() const
Return mesh.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
defineTypeNameAndDebug(cubeRootVolDelta, 0)
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< cell > cellList
list of cells
Definition: cellList.H:42
Namespace for OpenFOAM.