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 -------------------------------------------------------------------------------
11 License
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"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace LESModels
39 {
40  defineTypeNameAndDebug(IDDESDelta, 0);
41  addToRunTimeSelectionTable(LESdelta, IDDESDelta, dictionary);
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::LESModels::IDDESDelta::calcDelta()
49 {
50  const volScalarField& hmax = hmaxPtr_();
51  const fvMesh& mesh = turbulenceModel_.mesh();
52 
53  // Wall-normal vectors
54  const volVectorField& n = wallDist::New(mesh).n();
55 
56  tmp<volScalarField> tfaceToFacenMax
57  (
58  new volScalarField
59  (
60  IOobject
61  (
62  "faceToFaceMax",
63  mesh.time().timeName(),
64  mesh,
67  ),
68  mesh,
70  )
71  );
72 
73  scalarField& faceToFacenMax = tfaceToFacenMax.ref().primitiveFieldRef();
74 
75  const cellList& cells = mesh.cells();
76  const vectorField& faceCentres = mesh.faceCentres();
77 
78  forAll(cells, celli)
79  {
80  scalar maxDelta = 0.0;
81  const labelList& cFaces = cells[celli];
82  const vector nci = n[celli];
83 
84  forAll(cFaces, cFacei)
85  {
86  label facei = cFaces[cFacei];
87  const point& fci = faceCentres[facei];
88 
89  forAll(cFaces, cFacej)
90  {
91  label facej = cFaces[cFacej];
92  const point& fcj = faceCentres[facej];
93  scalar ndfc = nci & (fcj - fci);
94 
95  if (ndfc > maxDelta)
96  {
97  maxDelta = ndfc;
98  }
99  }
100  }
101 
102  faceToFacenMax[celli] = maxDelta;
103  }
104 
105 
106  label nD = mesh.nGeometricD();
107 
108  if (nD == 2)
109  {
111  << "Case is 2D, LES is not strictly applicable" << nl
112  << endl;
113  }
114  else if (nD != 3)
115  {
117  << "Case must be either 2D or 3D" << exit(FatalError);
118  }
119 
121  min
122  (
123  max
124  (
125  max
126  (
127  Cw_*wallDist::New(mesh).y(),
128  Cw_*hmax
129  ),
130  tfaceToFacenMax
131  ),
132  hmax
133  );
134 
135  // Handle coupled boundaries
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 
142 Foam::LESModels::IDDESDelta::IDDESDelta
143 (
144  const word& name,
146  const dictionary& dict
147 )
148 :
150  hmaxPtr_(nullptr),
151  Cw_
152  (
153  dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
154  (
155  "Cw",
156  0.15
157  )
158  )
159 {
160  if (dict.optionalSubDict(type() + "Coeffs").found("hmax"))
161  {
162  // User-defined hmax
163  hmaxPtr_ =
165  (
166  IOobject::groupName("hmax", turbulence.U().group()),
167  turbulence,
168  dict.optionalSubDict(type() + "Coeffs"),
169  "hmax"
170  );
171  }
172  else
173  {
174  Info<< "Employing " << maxDeltaxyz::typeName << " for hmax" << endl;
175 
176  hmaxPtr_.reset
177  (
178  new maxDeltaxyz
179  (
180  IOobject::groupName("hmax", turbulence.U().group()),
181  turbulence,
182  dict.optionalSubDict(type() + "Coeffs")
183  )
184  );
185  }
187  calcDelta();
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
194 {
195  const dictionary& coeffsDict(dict.optionalSubDict(type() + "Coeffs"));
197  coeffsDict.readIfPresent<scalar>("Cw", Cw_);
198 
199  calcDelta();
200 }
201 
202 
204 {
205  if (turbulenceModel_.mesh().changing())
206  {
207  hmaxPtr_->correct();
208  calcDelta();
209  }
210 }
211 
212 
213 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:300
const fvMesh & mesh() const
dictionary dict
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:152
const volVectorField & U() const
Access function to velocity field.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition: LESdelta.H:146
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Abstract base class for LES deltas.
Definition: LESdelta.H:49
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
void read(const dictionary &)
Read the LESdelta dictionary.
Definition: IDDESDelta.C:186
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
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
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Vector< scalar > vector
Definition: vector.H:57
volScalarField delta_
Definition: LESdelta.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
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...
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const volScalarField & hmax() const
Return the hmax delta field.
Definition: IDDESDelta.H:122
vector point
Point is a vector.
Definition: point.H:37
#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.
Nothing to be read.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
Delta calculated by taking the maximum distance between the cell centre and any face centre...
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:55
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127