cubeRootVolDelta.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 "cubeRootVolDelta.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35  defineTypeNameAndDebug(cubeRootVolDelta, 0);
36  addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 void Foam::LESModels::cubeRootVolDelta::calcDelta()
44 {
45  const fvMesh& mesh = momentumTransportModel_.mesh();
46 
47  label nD = mesh.nGeometricD();
48 
49  if (nD == 3)
50  {
51  delta_.primitiveFieldRef() = deltaCoeff_*pow(mesh.V(), 1.0/3.0);
52  }
53  else if (nD == 2)
54  {
56  << "Case is 2D, LES is not strictly applicable\n"
57  << endl;
58 
59  const Vector<label>& directions = mesh.geometricD();
60 
61  scalar thickness = 0.0;
62  for (direction dir=0; dir<directions.nComponents; dir++)
63  {
64  if (directions[dir] == -1)
65  {
66  thickness = mesh.bounds().span()[dir];
67  break;
68  }
69  }
70 
71  delta_.primitiveFieldRef() = deltaCoeff_*sqrt(mesh.V()/thickness);
72  }
73  else
74  {
76  << "Case is not 3D or 2D, LES is not applicable"
77  << exit(FatalError);
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
85 (
86  const word& name,
87  const momentumTransportModel& turbulence,
88  const dictionary& dict
89 )
90 :
91  LESdelta(name, turbulence),
92  deltaCoeff_
93  (
94  dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
95  (
96  "deltaCoeff",
97  1
98  )
99  )
100 {
101  calcDelta();
102 }
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {
109  dict.optionalSubDict(type() + "Coeffs").readIfPresent<scalar>
110  (
111  "deltaCoeff",
112  deltaCoeff_
113  );
114 
115  calcDelta();
116 }
117 
118 
120 {
121  if (momentumTransportModel_.mesh().changing())
122  {
123  calcDelta();
124  }
125 }
126 
127 
128 // ************************************************************************* //
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
virtual void read(const dictionary &)
Read the LESdelta dictionary.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
uint8_t direction
Definition: direction.H:45
Abstract base class for LES deltas.
Definition: LESdelta.H:50
addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Macros for easy insertion into run-time selection tables.
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
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Abstract base class for turbulence models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
cubeRootVolDelta(const word &name, const momentumTransportModel &turbulence, const dictionary &)
Construct from name, momentumTransportModel and dictionary.
Namespace for OpenFOAM.