PrandtlDelta.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) 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 "PrandtlDelta.H"
30 #include "wallDist.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace LESModels
38 {
40  addToRunTimeSelectionTable(LESdelta, PrandtlDelta, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::LESModels::PrandtlDelta::calcDelta()
48 {
49  delta_ = min
50  (
51  static_cast<const volScalarField&>(geometricDelta_()),
52  (kappa_/Cdelta_)*wallDist::New(turbulenceModel_.mesh()).y()
53  );
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
59 Foam::LESModels::PrandtlDelta::PrandtlDelta
60 (
61  const word& name,
63  const dictionary& dict
64 )
65 :
67  geometricDelta_
68  (
70  (
71  name,
72  turbulence,
73  dict.optionalSubDict(type() + "Coeffs")
74  )
75  ),
76  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
77  Cdelta_
78  (
79  dict.optionalSubDict(type() + "Coeffs").getOrDefault<scalar>
80  (
81  "Cdelta",
82  0.158
83  )
84  )
85 {
86  calcDelta();
87 }
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  const dictionary& coeffDict(dict.optionalSubDict(type() + "Coeffs"));
95 
96  geometricDelta_().read(coeffDict);
97  dict.readIfPresent<scalar>("kappa", kappa_);
98  coeffDict.readIfPresent<scalar>("Cdelta", Cdelta_);
99  calcDelta();
100 }
101 
102 
104 {
105  geometricDelta_().correct();
106 
107  if (turbulenceModel_.mesh().changing())
108  {
109  calcDelta();
110  }
111 }
112 
113 
114 // ************************************************************************* //
wallDist.H
Foam::LESdelta::New
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:69
PrandtlDelta.H
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::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::LESdelta::delta_
volScalarField delta_
Definition: LESdelta.H:62
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
PrandtlDelta
Apply Prandtl mixing-length based damping function to the specified geometric delta to improve near-w...
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
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::LESModels::PrandtlDelta::read
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Definition: PrandtlDelta.C:92
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dictionary::read
bool read(Istream &is)
Read dictionary from Istream.
Definition: dictionaryIO.C:141
Foam::LESdelta
Abstract base class for LES deltas.
Definition: LESdelta.H:53
Foam::LESModels::defineTypeNameAndDebug
defineTypeNameAndDebug(cubeRootVolDelta, 0)
Foam::LESModels::PrandtlDelta::correct
virtual void correct()
Definition: PrandtlDelta.C:103
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::LESModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(LESfluidThermoCompressibleTurbulenceModel, SmagorinskyLESfluidThermoCompressibleTurbulenceModel, dictionary)
y
scalar y
Definition: LISASMDCalcMethod1.H:14