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 | 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 "cubeRootVolDelta.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
39  addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
47 {
48  const fvMesh& mesh = turbulenceModel_.mesh();
49 
50  label nD = mesh.nGeometricD();
51 
52  if (nD == 3)
53  {
54  delta_.primitiveFieldRef() = deltaCoeff_ * cbrt(mesh.V());
55  }
56  else if (nD == 2)
57  {
59  << "Case is 2D, LES is not strictly applicable\n"
60  << endl;
61 
62  const Vector<label>& directions = mesh.geometricD();
63 
64  scalar thickness = 0.0;
65  for (direction dir=0; dir<directions.nComponents; dir++)
66  {
67  if (directions[dir] == -1)
68  {
69  thickness = mesh.bounds().span()[dir];
70  break;
71  }
72  }
73 
74  delta_.primitiveFieldRef() = deltaCoeff_*sqrt(mesh.V()/thickness);
75  }
76  else
77  {
78  if (debug)
79  {
81  << "Case is not 3D or 2D, LES is not applicable"
82  << exit(FatalError);
83  }
84  }
85 
86  // Handle coupled boundaries
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
93 Foam::LESModels::cubeRootVolDelta::cubeRootVolDelta
94 (
95  const word& name,
97  const dictionary& dict
98 )
99 :
101  deltaCoeff_
102  (
103  dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
104  (
105  "deltaCoeff",
106  1
107  )
108  )
109 {
110  calcDelta();
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  dict.optionalSubDict(type() + "Coeffs").readIfPresent<scalar>
119  (
120  "deltaCoeff",
121  deltaCoeff_
122  );
123 
124  calcDelta();
125 }
126 
127 
129 {
130  if (turbulenceModel_.mesh().changing())
131  {
132  calcDelta();
133  }
134 }
135 
136 
137 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::LESModels::cubeRootVolDelta::correct
virtual void correct()
Definition: cubeRootVolDelta.C:128
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:62
cubeRootVolDelta.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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:67
cubeRootVolDelta
Simple cube-root of cell volume delta used in LES models.
Foam::LESdelta::turbulenceModel_
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:60
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
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:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModels::cubeRootVolDelta::read
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: cubeRootVolDelta.C:116
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::LESdelta
Abstract base class for LES deltas.
Definition: LESdelta.H:53
Foam::LESModels::defineTypeNameAndDebug
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Foam::Vector< label >
Foam::LESModels::cubeRootVolDelta::calcDelta
void calcDelta()
Calculate the delta values.
Definition: cubeRootVolDelta.C:46
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::turbulenceModel::mesh
const fvMesh & mesh() const
Definition: turbulenceModel.H:129
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::LESModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303