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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd.
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 {
37  addToRunTimeSelectionTable(LESdelta, IDDESDelta, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
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().internalField();
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 
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  // Handle coupled boundaries
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137 
139 (
140  const word& name,
142  const dictionary& dict
143 )
144 :
146  hmax_
147  (
148  IOobject::groupName("hmax", turbulence.U().group()),
149  turbulence,
150  dict
151  ),
152  Cw_
153  (
154  dict.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cw", 0.15)
155  )
156 {
157  calcDelta();
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
164 {
165  const dictionary& coeffsDict(dict.subDict(type() + "Coeffs"));
166 
167  coeffsDict.readIfPresent<scalar>("Cw", Cw_);
168 
169  calcDelta();
170 }
171 
172 
174 {
175  if (turbulenceModel_.mesh().changing())
176  {
177  hmax_.correct();
178  calcDelta();
179  }
180 }
181 
182 
183 // ************************************************************************* //
wallDist.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::LESModels::IDDESDelta::IDDESDelta
IDDESDelta(const IDDESDelta &)
Disallow default bitwise copy construct and assignment.
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:59
Foam::LESModels::IDDESDelta::hmax_
maxDeltaxyz hmax_
Definition: IDDESDelta.H:58
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::LESModels::IDDESDelta::read
void read(const dictionary &)
Read the LESdelta dictionary.
Definition: IDDESDelta.C:163
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Foam::LESModels::IDDESDelta::Cw_
scalar Cw_
Definition: IDDESDelta.H:60
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
IDDESDelta
IDDESDelta used by the IDDES (improved low Re Spalart-Allmaras DES model) The min and max delta are c...
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::LESModels::IDDESDelta::correct
void correct()
Definition: IDDESDelta.C:173
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::LESdelta::turbulenceModel_
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:57
IDDESDelta.H
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::wallDist::n
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:168
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::LESModels::IDDESDelta::calcDelta
void calcDelta()
Calculate the delta values.
Definition: IDDESDelta.C:44
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::LESdelta
Abstract base class for LES deltas.
Definition: LESdelta.H:50
Foam::LESModels::defineTypeNameAndDebug
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::turbulenceModel::mesh
const fvMesh & mesh() const
Definition: turbulenceModel.H:127
Foam::LESModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::LESModels::IDDESDelta::hmax
const volScalarField & hmax() const
Return the hmax delta field.
Definition: IDDESDelta.H:102
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
y
scalar y
Definition: LISASMDCalcMethod1.H:14