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-2020 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 = momentumTransportModel_.mesh();
48 
49  // Wall-normal vectors
50  const volVectorField& n = wallDist::New(mesh).n();
51 
52  tmp<volScalarField> tfaceToFacenMax
53  (
55  (
56  "faceToFaceMax",
57  mesh,
59  )
60  );
61 
62  scalarField& faceToFacenMax = tfaceToFacenMax.ref().primitiveFieldRef();
63 
64  const cellList& cells = mesh.cells();
65  const vectorField& faceCentres = mesh.faceCentres();
66 
67  forAll(cells, celli)
68  {
69  scalar maxDelta = 0.0;
70  const labelList& cFaces = cells[celli];
71  const vector nci = n[celli];
72 
73  forAll(cFaces, cFacei)
74  {
75  label facei = cFaces[cFacei];
76  const point& fci = faceCentres[facei];
77 
78  forAll(cFaces, cFacej)
79  {
80  label facej = cFaces[cFacej];
81  const point& fcj = faceCentres[facej];
82  scalar ndfc = nci & (fcj - fci);
83 
84  if (ndfc > maxDelta)
85  {
86  maxDelta = ndfc;
87  }
88  }
89  }
90 
91  faceToFacenMax[celli] = maxDelta;
92  }
93 
94 
95  label nD = mesh.nGeometricD();
96 
97  if (nD == 2)
98  {
100  << "Case is 2D, LES is not strictly applicable" << nl
101  << endl;
102  }
103  else if (nD != 3)
104  {
106  << "Case must be either 2D or 3D" << exit(FatalError);
107  }
108 
109  delta_.primitiveFieldRef() =
110  min
111  (
112  max
113  (
114  max
115  (
116  Cw_*wallDist::New(mesh).y(),
117  Cw_*hmax
118  ),
119  tfaceToFacenMax
120  ),
121  hmax
122  );
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
129 (
130  const word& name,
131  const momentumTransportModel& turbulence,
132  const dictionary& dict
133 )
134 :
135  LESdelta(name, turbulence),
136  hmax_
137  (
138  IOobject::groupName("hmax", turbulence.U().group()),
139  turbulence,
140  dict
141  ),
142  Cw_
143  (
144  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
145  (
146  "Cw",
147  0.15
148  )
149  )
150 {
151  calcDelta();
152 }
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
158 {
159  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
160 
161  coeffsDict.readIfPresent<scalar>("Cw", Cw_);
162 
163  calcDelta();
164 }
165 
166 
168 {
169  if (momentumTransportModel_.mesh().changing())
170  {
171  hmax_.correct();
172  calcDelta();
173  }
174 }
175 
176 
177 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:162
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
Abstract base class for LES deltas.
Definition: LESdelta.H:50
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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:157
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
scalar y
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1040
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)
IDDESDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, mesh and IOdictionary.
Definition: IDDESDelta.C:129
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:260
const Mesh & mesh() const
Return mesh.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< cell > cellList
list of cells
Definition: cellList.H:42
Namespace for OpenFOAM.