All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
LESdelta.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 "LESdelta.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(LESdelta, 0);
34  defineRunTimeSelectionTable(LESdelta, dictionary);
35 }
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const word& name,
42  const momentumTransportModel& turbulence
43 )
44 :
45  momentumTransportModel_(turbulence),
46  delta_
47  (
48  IOobject
49  (
50  name,
51  turbulence.mesh().time().timeName(),
52  turbulence.mesh(),
55  ),
56  turbulence.mesh(),
57  dimensionedScalar(name, dimLength, small),
58  calculatedFvPatchScalarField::typeName
59  )
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
64 
66 (
67  const word& name,
68  const momentumTransportModel& turbulence,
69  const dictionary& dict
70 )
71 {
72  const word deltaType(dict.lookup("delta"));
73 
74  Info<< "Selecting LES delta type " << deltaType << endl;
75 
76  dictionaryConstructorTable::iterator cstrIter =
77  dictionaryConstructorTablePtr_->find(deltaType);
78 
79  if (cstrIter == dictionaryConstructorTablePtr_->end())
80  {
82  << "Unknown LESdelta type "
83  << deltaType << nl << nl
84  << "Valid LESdelta types are :" << endl
85  << dictionaryConstructorTablePtr_->sortedToc()
86  << exit(FatalError);
87  }
88 
89  return autoPtr<LESdelta>(cstrIter()(name, turbulence, dict));
90 }
91 
92 
94 (
95  const word& name,
96  const momentumTransportModel& turbulence,
97  const dictionary& dict,
98  const dictionaryConstructorTable& additionalConstructors
99 )
100 {
101  const word deltaType(dict.lookup("delta"));
102 
103  Info<< "Selecting LES delta type " << deltaType << endl;
104 
105  // First on additional ones
106  dictionaryConstructorTable::const_iterator cstrIter =
107  additionalConstructors.find(deltaType);
108 
109  if (cstrIter != additionalConstructors.end())
110  {
111  return autoPtr<LESdelta>(cstrIter()(name, turbulence, dict));
112  }
113  else
114  {
115  dictionaryConstructorTable::const_iterator cstrIter =
116  dictionaryConstructorTablePtr_->find(deltaType);
117 
118  if (cstrIter == dictionaryConstructorTablePtr_->end())
119  {
121  << "Unknown LESdelta type "
122  << deltaType << nl << nl
123  << "Valid LESdelta types are :" << endl
124  << additionalConstructors.sortedToc()
125  << " and "
126  << dictionaryConstructorTablePtr_->sortedToc()
127  << exit(FatalError);
128  return autoPtr<LESdelta>();
129  }
130  else
131  {
132  return autoPtr<LESdelta>(cstrIter()(name, turbulence, dict));
133  }
134  }
135 }
136 
137 
138 // ************************************************************************* //
dictionary dict
static autoPtr< LESdelta > New(const word &name, const momentumTransportModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
Definition: LESdelta.C:66
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
LESdelta(const word &name, const momentumTransportModel &turbulence)
Construct from name, momentumTransportModel and dictionary.
Definition: LESdelta.C:40
A class for handling words, derived from string.
Definition: word.H:59
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Abstract base class for turbulence models (RAS, LES and laminar).
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812