cubeRootVolDelta.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 "cubeRootVolDelta.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
36  addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
44 {
45  const fvMesh& mesh = turbulenceModel_.mesh();
46 
47  label nD = mesh.nGeometricD();
48 
49  if (nD == 3)
50  {
51  delta_.internalField() = 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_.internalField() = 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  // Handle coupled boundaries
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
88 (
89  const word& name,
91  const dictionary& dict
92 )
93 :
95  deltaCoeff_
96  (
97  dict.subDict(type() + "Coeffs").lookupOrDefault<scalar>("deltaCoeff", 1)
98  )
99 {
100  calcDelta();
101 }
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
107 {
108  dict.subDict(type() + "Coeffs").readIfPresent<scalar>
109  (
110  "deltaCoeff",
111  deltaCoeff_
112  );
113 
114  calcDelta();
115 }
116 
117 
119 {
120  if (turbulenceModel_.mesh().changing())
121  {
122  calcDelta();
123  }
124 }
125 
126 
127 // ************************************************************************* //
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::LESModels::cubeRootVolDelta::correct
virtual void correct()
Definition: cubeRootVolDelta.C:118
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:59
cubeRootVolDelta.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::directions
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:63
Foam::LESModels::cubeRootVolDelta::cubeRootVolDelta
cubeRootVolDelta(const cubeRootVolDelta &)
Disallow default bitwise copy construct and assignment.
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
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
cubeRootVolDelta
Simple cube-root of cell volume delta used in LES models.
Foam::LESdelta::turbulenceModel_
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:57
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
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::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::cubeRootVolDelta::read
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: cubeRootVolDelta.C:106
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
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::LESModels::cubeRootVolDelta::deltaCoeff_
scalar deltaCoeff_
Definition: cubeRootVolDelta.H:56
Foam::Vector< label >
Foam::LESModels::cubeRootVolDelta::calcDelta
void calcDelta()
Definition: cubeRootVolDelta.C:43
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::turbulenceModel::mesh
const fvMesh & mesh() const
Definition: turbulenceModel.H:127
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::LESModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
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