Loading...
Searching...
No Matches
LESdelta.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-2015 OpenFOAM Foundation
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "LESdelta.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
38}
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
43(
44 const word& name,
46)
47:
48 turbulenceModel_(turbulence),
49 delta_
50 (
52 (
53 name,
54 turbulence.mesh().time().timeName(),
56 IOobject::NO_READ,
57 IOobject::NO_WRITE
58 ),
61 )
62{}
63
64
65// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
66
68(
69 const word& name,
71 const dictionary& dict,
72 const word& lookupName
73)
74{
75 const word deltaType(dict.get<word>(lookupName));
76
77 Info<< "Selecting LES " << lookupName << " type " << deltaType << endl;
78
79 auto* ctorPtr = dictionaryConstructorTable(deltaType);
80
81 if (!ctorPtr)
82 {
84 (
85 dict,
86 "LESdelta",
87 deltaType,
88 *dictionaryConstructorTablePtr_
90 }
91
92 return autoPtr<LESdelta>(ctorPtr(name, turbulence, dict));
93}
94
95
97(
98 const word& name,
99 const turbulenceModel& turbulence,
100 const dictionary& dict,
101 const dictionaryConstructorTableType& additionalConstructors,
102 const word& lookupName
103)
104{
105 const word deltaType(dict.get<word>(lookupName));
106
107 Info<< "Selecting LES " << lookupName << " type " << deltaType << endl;
108
109 // First any additional ones
110 {
111 auto ctorIter = additionalConstructors.cfind(deltaType);
112
113 if (ctorIter.good())
114 {
115 return autoPtr<LESdelta>(ctorIter.val()(name, turbulence, dict));
116 }
117 }
118
119 auto* ctorPtr = dictionaryConstructorTable(deltaType);
120
121 if (!ctorPtr)
122 {
124 (
125 dict,
126 "LESdelta",
127 deltaType,
128 *dictionaryConstructorTablePtr_
129 );
130
131 if (!additionalConstructors.empty())
132 {
134 << " and " << additionalConstructors.sortedToc() << nl;
135 }
136
138 << exit(FatalIOError);
139 }
140
141 return autoPtr<LESdelta>(ctorPtr(name, turbulence, dict));
142}
143
144
145// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Abstract base class for LES deltas.
Definition LESdelta.H:50
LESdelta(const LESdelta &)=delete
No copy construct.
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:61
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition LESdelta.H:146
volScalarField delta_
Definition LESdelta.H:57
const turbulenceModel & turbulenceModel_
Definition LESdelta.H:55
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
compressible::turbulenceModel & turbulence
auto & name
word timeName
Definition getTimeIndex.H:3
Namespace for OpenFOAM.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict