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-2022 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 {
35 }
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const word& name,
43 )
44 :
45  momentumTransportModel_(turbulence),
46  delta_
47  (
48  IOobject
49  (
50  name,
51  turbulence.mesh().time().name(),
52  turbulence.mesh(),
53  IOobject::NO_READ,
54  IOobject::NO_WRITE
55  ),
56  turbulence.mesh(),
58  calculatedFvPatchScalarField::typeName
59  )
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
64 
66 (
67  const word& name,
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,
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 // ************************************************************************* //
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Abstract base class for LES deltas.
Definition: LESdelta.H:51
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
LESdelta(const word &name, const momentumTransportModel &turbulence)
Construct from name, momentumTransportModel and dictionary.
Definition: LESdelta.C:40
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
const dimensionSet dimLength
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
static const char nl
Definition: Ostream.H:266
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))