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-2024 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 {
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 
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,
132  const dictionary& dict
133 )
134 :
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& coeffDict(dict.optionalSubDict(type() + "Coeffs"));
160 
161  coeffDict.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 // ************************************************************************* //
scalar y
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
static wallDist & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
IDDESDelta used by the IDDES (improved low Re Spalart-Allmaras DES model) The min and max delta are c...
Definition: IDDESDelta.H:55
IDDESDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, mesh and IOdictionary.
Definition: IDDESDelta.C:129
void read(const dictionary &)
Read the LESdelta dictionary.
Definition: IDDESDelta.C:157
const volScalarField & hmax() const
Return the hmax delta field.
Definition: IDDESDelta.H:101
Abstract base class for LES deltas.
Definition: LESdelta.H:51
volScalarField delta_
Definition: LESdelta.H:58
const momentumTransportModel & momentumTransportModel_
Definition: LESdelta.H:56
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Abstract base class for turbulence models (RAS, LES and laminar).
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1023
const vectorField & faceCentres() const
const cellList & cells() const
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:173
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const cellShapeList & cells
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
defineTypeNameAndDebug(cubeRootVolDelta, 0)
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
const char *const group
Group name for atomic constants.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
List< cell > cellList
list of cells
Definition: cellList.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
static const char nl
Definition: Ostream.H:267
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
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))