/
LESdelta.C
138 lines (114 loc) · 4.14 KB
/
LESdelta.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "LESdelta.H"
#include "calculatedFvPatchFields.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(LESdelta, 0);
defineRunTimeSelectionTable(LESdelta, dictionary);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::LESdelta::LESdelta
(
const word& name,
const turbulenceModel& turbulence
)
:
turbulenceModel_(turbulence),
delta_
(
IOobject
(
name,
turbulence.mesh().time().timeName(),
turbulence.mesh(),
IOobject::NO_READ,
IOobject::NO_WRITE
),
turbulence.mesh(),
dimensionedScalar(name, dimLength, small),
calculatedFvPatchScalarField::typeName
)
{}
// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
Foam::autoPtr<Foam::LESdelta> Foam::LESdelta::New
(
const word& name,
const turbulenceModel& turbulence,
const dictionary& dict
)
{
const word deltaType(dict.lookup("delta"));
Info<< "Selecting LES delta type " << deltaType << endl;
dictionaryConstructorTable::iterator cstrIter =
dictionaryConstructorTablePtr_->find(deltaType);
if (cstrIter == dictionaryConstructorTablePtr_->end())
{
FatalErrorInFunction
<< "Unknown LESdelta type "
<< deltaType << nl << nl
<< "Valid LESdelta types are :" << endl
<< dictionaryConstructorTablePtr_->sortedToc()
<< exit(FatalError);
}
return autoPtr<LESdelta>(cstrIter()(name, turbulence, dict));
}
Foam::autoPtr<Foam::LESdelta> Foam::LESdelta::New
(
const word& name,
const turbulenceModel& turbulence,
const dictionary& dict,
const dictionaryConstructorTable& additionalConstructors
)
{
const word deltaType(dict.lookup("delta"));
Info<< "Selecting LES delta type " << deltaType << endl;
// First on additional ones
dictionaryConstructorTable::const_iterator cstrIter =
additionalConstructors.find(deltaType);
if (cstrIter != additionalConstructors.end())
{
return autoPtr<LESdelta>(cstrIter()(name, turbulence, dict));
}
else
{
dictionaryConstructorTable::const_iterator cstrIter =
dictionaryConstructorTablePtr_->find(deltaType);
if (cstrIter == dictionaryConstructorTablePtr_->end())
{
FatalErrorInFunction
<< "Unknown LESdelta type "
<< deltaType << nl << nl
<< "Valid LESdelta types are :" << endl
<< additionalConstructors.sortedToc()
<< " and "
<< dictionaryConstructorTablePtr_->sortedToc()
<< exit(FatalError);
return autoPtr<LESdelta>();
}
else
{
return autoPtr<LESdelta>(cstrIter()(name, turbulence, dict));
}
}
}
// ************************************************************************* //