FreeFOAM The Cross-Platform CFD Toolkit
ddtScheme.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) 1991-2010 OpenCFD Ltd.
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::fv::ddtScheme
26 
27 Description
28  Abstract base class for ddt schemes.
29 
30 SourceFiles
31  ddtScheme.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef ddtScheme_H
36 #define ddtScheme_H
37 
38 #include <OpenFOAM/tmp.H>
42 #include <OpenFOAM/typeInfo.H>
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 template<class Type>
51 class fvMatrix;
52 
53 class fvMesh;
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace fv
58 {
59 
60 /*---------------------------------------------------------------------------*\
61  Class ddtScheme Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 template<class Type>
65 class ddtScheme
66 :
67  public refCount
68 {
69 
70 protected:
71 
72  // Protected data
73 
74  const fvMesh& mesh_;
75 
76 
77  // Private Member Functions
78 
79  //- Disallow copy construct
80  ddtScheme(const ddtScheme&);
81 
82  //- Disallow default bitwise assignment
83  void operator=(const ddtScheme&);
84 
85 
86 public:
87 
88  //- Runtime type information
89  virtual const word& type() const = 0;
90 
91 
92  // Declare run-time constructor selection tables
93 
95  (
96  tmp,
97  ddtScheme,
98  Istream,
99  (const fvMesh& mesh, Istream& schemeData),
100  (mesh, schemeData)
101  );
102 
103 
104  // Constructors
105 
106  //- Construct from mesh
107  ddtScheme(const fvMesh& mesh)
108  :
109  mesh_(mesh)
110  {}
111 
112  //- Construct from mesh and Istream
113  ddtScheme(const fvMesh& mesh, Istream&)
114  :
115  mesh_(mesh)
116  {}
117 
118 
119  // Selectors
120 
121  //- Return a pointer to a new ddtScheme created on freestore
122  static tmp<ddtScheme<Type> > New
123  (
124  const fvMesh& mesh,
125  Istream& schemeData
126  );
127 
128 
129  // Destructor
130 
131  virtual ~ddtScheme();
132 
133 
134  // Member Functions
135 
136  //- Return mesh reference
137  const fvMesh& mesh() const
138  {
139  return mesh_;
140  }
141 
143  (
144  const dimensioned<Type>&
145  ) = 0;
146 
148  (
150  ) = 0;
151 
153  (
154  const dimensionedScalar&,
156  ) = 0;
157 
159  (
160  const volScalarField&,
162  ) = 0;
163 
164  virtual tmp<fvMatrix<Type> > fvmDdt
165  (
167  ) = 0;
168 
169  virtual tmp<fvMatrix<Type> > fvmDdt
170  (
171  const dimensionedScalar&,
173  ) = 0;
174 
175  virtual tmp<fvMatrix<Type> > fvmDdt
176  (
177  const volScalarField&,
179  ) = 0;
180 
181 
182  typedef GeometricField
183  <
184  typename flux<Type>::type,
188 
190  (
192  const fluxFieldType& phi,
193  const fluxFieldType& phiCorr
194  );
195 
197  (
199  const fluxFieldType& phi
200  );
201 
203  (
204  const volScalarField& rA,
206  const fluxFieldType& phi
207  ) = 0;
208 
210  (
211  const volScalarField& rho,
213  const fluxFieldType& phi
214  );
215 
217  (
218  const volScalarField& rA,
219  const volScalarField& rho,
221  const fluxFieldType& phi
222  ) = 0;
223 
224 
226  (
228  ) = 0;
229 };
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace fv
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 // Add the patch constructor functions to the hash tables
243 
244 #define makeFvDdtTypeScheme(SS, Type) \
245  \
246 defineNamedTemplateTypeNameAndDebug(SS<Type>, 0); \
247  \
248 ddtScheme<Type>::addIstreamConstructorToTable<SS<Type> > \
249  add##SS##Type##IstreamConstructorToTable_;
250 
251 
252 #define makeFvDdtScheme(SS) \
253  \
254 makeFvDdtTypeScheme(SS, scalar) \
255 makeFvDdtTypeScheme(SS, vector) \
256 makeFvDdtTypeScheme(SS, sphericalTensor) \
257 makeFvDdtTypeScheme(SS, symmTensor) \
258 makeFvDdtTypeScheme(SS, tensor) \
259  \
260 template<> \
261 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
262 ( \
263  const volScalarField& rA, \
264  const volScalarField& U, \
265  const surfaceScalarField& phi \
266 ) \
267 { \
268  notImplemented(#SS"<scalar>::fvcDdtPhiCorr"); \
269  return surfaceScalarField::null(); \
270 } \
271  \
272 template<> \
273 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
274 ( \
275  const volScalarField& rA, \
276  const volScalarField& rho, \
277  const volScalarField& U, \
278  const surfaceScalarField& phi \
279 ) \
280 { \
281  notImplemented(#SS"<scalar>::fvcDdtPhiCorr"); \
282  return surfaceScalarField::null(); \
283 }
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 #ifdef NoRepository
289 # include <finiteVolume/ddtScheme.C>
290 #endif
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 #endif
295 
296 // ************************ vim: set sw=4 sts=4 et: ************************ //