ddtScheme.H
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-2018 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::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 "tmp.H"
39 #include "dimensionedType.H"
40 #include "volFieldsFwd.H"
41 #include "surfaceFieldsFwd.H"
42 #include "typeInfo.H"
43 #include "runTimeSelectionTables.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 tmp<ddtScheme<Type>>::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  virtual ~ddtScheme();
131 
132 
133  // Member Functions
134 
135  //- Return mesh reference
136  const fvMesh& mesh() const
137  {
138  return mesh_;
139  }
140 
142  (
143  const dimensioned<Type>&
144  ) = 0;
145 
147  (
149  ) = 0;
150 
152  (
153  const dimensionedScalar&,
155  ) = 0;
156 
158  (
159  const volScalarField&,
161  ) = 0;
162 
164  (
165  const volScalarField& alpha,
166  const volScalarField& rho,
168  ) = 0;
169 
171  (
173  );
174 
175  virtual tmp<fvMatrix<Type>> fvmDdt
176  (
178  ) = 0;
179 
180  virtual tmp<fvMatrix<Type>> fvmDdt
181  (
182  const dimensionedScalar&,
184  ) = 0;
185 
186  virtual tmp<fvMatrix<Type>> fvmDdt
187  (
188  const volScalarField&,
190  ) = 0;
191 
192  virtual tmp<fvMatrix<Type>> fvmDdt
193  (
194  const volScalarField& alpha,
195  const volScalarField& rho,
197  ) = 0;
198 
199  typedef GeometricField
200  <
201  typename flux<Type>::type,
204  > fluxFieldType;
205 
207  (
209  const fluxFieldType& phi,
210  const fluxFieldType& phiCorr
211  );
212 
214  (
216  const fluxFieldType& phi,
217  const fluxFieldType& phiCorr,
218  const volScalarField& rho
219  );
220 
222  (
224  const fluxFieldType& phi
225  );
226 
228  (
230  const fluxFieldType& phi,
231  const volScalarField& rho
232  );
233 
235  (
238  ) = 0;
239 
241  (
243  const fluxFieldType& phi
244  ) = 0;
245 
247  (
248  const volScalarField& rho,
251  ) = 0;
252 
254  (
255  const volScalarField& rho,
257  const fluxFieldType& phi
258  ) = 0;
259 
261  (
263  ) = 0;
264 };
265 
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 } // End namespace fv
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 } // End namespace Foam
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 // Add the patch constructor functions to the hash tables
279 #define makeFvDdtTypeScheme(SS, Type) \
280  defineNamedTemplateTypeNameAndDebug(Foam::fv::SS<Foam::Type>, 0); \
281  \
282  namespace Foam \
283  { \
284  namespace fv \
285  { \
286  ddtScheme<Type>::addIstreamConstructorToTable<SS<Type>> \
287  add##SS##Type##IstreamConstructorToTable_; \
288  } \
289  }
291 #define makeFvDdtScheme(SS) \
292  \
293 makeFvDdtTypeScheme(SS, scalar) \
294 makeFvDdtTypeScheme(SS, vector) \
295 makeFvDdtTypeScheme(SS, sphericalTensor) \
296 makeFvDdtTypeScheme(SS, symmTensor) \
297 makeFvDdtTypeScheme(SS, tensor) \
298  \
299 namespace Foam \
300 { \
301 namespace fv \
302 { \
303  \
304 template<> \
305 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
306 ( \
307  const volScalarField& U, \
308  const surfaceScalarField& Uf \
309 ) \
310 { \
311  NotImplemented; \
312  return surfaceScalarField::null(); \
313 } \
314  \
315 template<> \
316 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
317 ( \
318  const volScalarField& U, \
319  const surfaceScalarField& phi \
320 ) \
321 { \
322  NotImplemented; \
323  return surfaceScalarField::null(); \
324 } \
325  \
326 template<> \
327 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
328 ( \
329  const volScalarField& rho, \
330  const volScalarField& U, \
331  const surfaceScalarField& Uf \
332 ) \
333 { \
334  NotImplemented; \
335  return surfaceScalarField::null(); \
336 } \
337  \
338 template<> \
339 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
340 ( \
341  const volScalarField& rho, \
342  const volScalarField& U, \
343  const surfaceScalarField& phi \
344 ) \
345 { \
346  NotImplemented; \
347  return surfaceScalarField::null(); \
348 } \
349  \
350 } \
351 }
352 
353 
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355 
356 #ifdef NoRepository
357  #include "ddtScheme.C"
358 #endif
359 
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 
362 #endif
363 
364 // ************************************************************************* //
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)=0
surfaceScalarField & phi
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:47
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)=0
Reference counter for various OpenFOAM components.
Definition: refCount.H:49
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)=0
declareRunTimeSelectionTable(tmp, ddtScheme, Istream,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > fluxFieldType
Definition: ddtScheme.H:203
Generic GeometricField class.
Generic dimensioned Type class.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:64
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)=0
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:135
virtual ~ddtScheme()
Destructor.
Definition: ddtScheme.C:90
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
autoPtr< surfaceVectorField > Uf
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)=0
virtual const word & type() const =0
Runtime type information.
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:153
U
Definition: pEqn.H:72
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:46
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Macros to ease declaration of run-time selection tables.
void operator=(const ddtScheme &)
Disallow default bitwise assignment.
A class for managing temporary objects.
Definition: PtrList.H:53
const fvMesh & mesh_
Definition: ddtScheme.H:73
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
ddtScheme(const ddtScheme &)
Disallow copy construct.
typeOfRank< typename pTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank) - 2 >::type type
Definition: products.H:115
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.