LESdelta.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 Class
25  Foam::LESdelta
26 
27 Description
28  Abstract base class for LES deltas
29 
30 SourceFiles
31  LESdelta.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef LESdelta_H
36 #define LESdelta_H
37 
38 #include "turbulenceModel.H"
39 #include "volFields.H"
40 #include "runTimeSelectionTables.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class LESdelta Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 class LESdelta
52 {
53 
54 protected:
55 
56  // Protected data
57 
59 
61 
62 
63 private:
64 
65  // Private Member Functions
66 
67  // Disallow default bitwise copy construct and assignment
68  LESdelta(const LESdelta&);
69  void operator=(const LESdelta&);
70 
71 
72 public:
73 
74  //- Runtime type information
75  TypeName("LESdelta");
76 
77 
78  // Declare run-time constructor selection table
79 
81  (
82  autoPtr,
83  LESdelta,
84  dictionary,
85  (
86  const word& name,
88  const dictionary& dict
89  ),
90  (name, turbulence, dict)
91  );
92 
93 
94  // Constructors
95 
96  //- Construct from name, turbulenceModel and dictionary
97  LESdelta
98  (
99  const word& name,
100  const turbulenceModel& turbulence
101  );
102 
103 
104  // Selectors
105 
106  //- Return a reference to the selected LES delta
107  static autoPtr<LESdelta> New
108  (
109  const word& name,
110  const turbulenceModel& turbulence,
111  const dictionary& dict
112  );
113 
114  //- Return a reference to the selected LES delta
115  static autoPtr<LESdelta> New
116  (
117  const word& name,
118  const turbulenceModel& turbulence,
119  const dictionary& dict,
120  const dictionaryConstructorTable&
121  );
122 
123 
124  //- Destructor
125  virtual ~LESdelta()
126  {}
127 
128 
129  // Member Functions
130 
131  //- Return turbulenceModel reference
132  const turbulenceModel& turbulence() const
133  {
134  return turbulenceModel_;
135  }
136 
137  //- Read the LESdelta dictionary
138  virtual void read(const dictionary&) = 0;
139 
140  // Correct values
141  virtual void correct() = 0;
142 
143 
144  // Member Operators
146  virtual operator const volScalarField&() const
147  {
148  return delta_;
149  }
150 };
151 
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 } // End namespace Foam
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 #endif
160 
161 // ************************************************************************* //
dictionary dict
virtual ~LESdelta()
Destructor.
Definition: LESdelta.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Abstract base class for LES deltas.
Definition: LESdelta.H:50
TypeName("LESdelta")
Runtime type information.
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from string.
Definition: word.H:59
declareRunTimeSelectionTable(autoPtr, LESdelta, dictionary,(const word &name, const turbulenceModel &turbulence, const dictionary &dict),(name, turbulence, dict))
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict)
Return a reference to the selected LES delta.
Definition: LESdelta.C:66
virtual void read(const dictionary &)=0
Read the LESdelta dictionary.
volScalarField delta_
Definition: LESdelta.H:59
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition: LESdelta.H:131
Macros to ease declaration of run-time selection tables.
const turbulenceModel & turbulenceModel_
Definition: LESdelta.H:57
virtual void correct()=0
Namespace for OpenFOAM.