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-2023 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 public:
78 
79  //- Runtime type information
80  virtual const word& type() const = 0;
81 
82 
83  // Declare run-time constructor selection tables
84 
86  (
87  tmp,
88  ddtScheme,
89  Istream,
90  (const fvMesh& mesh, Istream& schemeData),
91  (mesh, schemeData)
92  );
93 
94 
95  // Constructors
96 
97  //- Construct from mesh
98  ddtScheme(const fvMesh& mesh)
99  :
100  mesh_(mesh)
101  {}
102 
103  //- Construct from mesh and Istream
104  ddtScheme(const fvMesh& mesh, Istream&)
105  :
106  mesh_(mesh)
107  {}
108 
109  //- Disallow default bitwise copy construction
110  ddtScheme(const ddtScheme&);
111 
112 
113  // Selectors
114 
115  //- Return a pointer to a new ddtScheme created on freestore
116  static tmp<ddtScheme<Type>> New
117  (
118  const fvMesh& mesh,
119  Istream& schemeData
120  );
121 
122 
123  //- Destructor
124  virtual ~ddtScheme();
125 
126 
127  // Member Functions
128 
129  //- Return mesh reference
130  const fvMesh& mesh() const
131  {
132  return mesh_;
133  }
134 
135  virtual tmp<VolField<Type>> fvcDdt
136  (
137  const dimensioned<Type>&
138  ) = 0;
139 
140  virtual tmp<VolField<Type>> fvcDdt
141  (
142  const VolField<Type>&
143  ) = 0;
144 
145  virtual tmp<VolField<Type>> fvcDdt
146  (
147  const dimensionedScalar&,
148  const VolField<Type>&
149  ) = 0;
150 
151  virtual tmp<VolField<Type>> fvcDdt
152  (
153  const volScalarField&,
154  const VolField<Type>&
155  ) = 0;
156 
157  virtual tmp<VolField<Type>> fvcDdt
158  (
159  const volScalarField& alpha,
160  const volScalarField& rho,
161  const VolField<Type>&
162  ) = 0;
163 
165  (
166  const SurfaceField<Type>&
167  );
168 
169  virtual tmp<fvMatrix<Type>> fvmDdt
170  (
171  const VolField<Type>&
172  ) = 0;
173 
174  virtual tmp<fvMatrix<Type>> fvmDdt
175  (
176  const dimensionedScalar&,
177  const VolField<Type>&
178  ) = 0;
179 
180  virtual tmp<fvMatrix<Type>> fvmDdt
181  (
182  const volScalarField&,
183  const VolField<Type>&
184  ) = 0;
185 
186  virtual tmp<fvMatrix<Type>> fvmDdt
187  (
188  const volScalarField& alpha,
189  const volScalarField& rho,
190  const VolField<Type>& vf
191  ) = 0;
192 
194 
196  (
197  const VolField<Type>& U,
198  const fluxFieldType& phi,
199  const fluxFieldType& phiCorr
200  );
201 
203  (
204  const VolField<Type>& U,
205  const fluxFieldType& phi,
206  const fluxFieldType& phiCorr,
207  const volScalarField& rho
208  );
209 
211  (
212  const VolField<Type>& U,
213  const fluxFieldType& phi
214  );
215 
217  (
218  const VolField<Type>& U,
219  const fluxFieldType& phi,
220  const volScalarField& rho
221  );
222 
224  (
225  const VolField<Type>& U,
226  const SurfaceField<Type>& Uf
227  ) = 0;
228 
230  (
231  const VolField<Type>& U,
232  const fluxFieldType& phi
233  ) = 0;
234 
236  (
237  const volScalarField& rho,
238  const VolField<Type>& U,
239  const SurfaceField<Type>& rhoUf
240  ) = 0;
241 
243  (
244  const volScalarField& rho,
245  const VolField<Type>& U,
246  const fluxFieldType& phi
247  ) = 0;
248 
250  (
251  const volScalarField& alpha,
252  const volScalarField& rho,
253  const VolField<Type>& U,
254  const SurfaceField<Type>& Uf
255  ) = 0;
256 
258  (
259  const volScalarField& alpha,
260  const volScalarField& rho,
261  const VolField<Type>& U,
262  const fluxFieldType& phi
263  ) = 0;
264 
266  (
267  const VolField<Type>&
268  ) = 0;
269 
270  virtual tmp<scalarField> meshPhi
271  (
272  const VolField<Type>&,
273  const label patchi
274  ) = 0;
275 
276 
277  // Member Operators
278 
279  //- Disallow default bitwise assignment
280  void operator=(const ddtScheme&) = delete;
281 };
282 
283 
284 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 
286 } // End namespace fv
287 
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 
290 } // End namespace Foam
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 // Add the patch constructor functions to the hash tables
295 
296 #define makeFvDdtTypeScheme(SS, Type) \
297  defineNamedTemplateTypeNameAndDebug(Foam::fv::SS<Foam::Type>, 0); \
298  \
299  namespace Foam \
300  { \
301  namespace fv \
302  { \
303  ddtScheme<Type>::addIstreamConstructorToTable<SS<Type>> \
304  add##SS##Type##IstreamConstructorToTable_; \
305  } \
306  }
307 
308 #define makeFvDdtScheme(SS) \
309  \
310 makeFvDdtTypeScheme(SS, scalar) \
311 makeFvDdtTypeScheme(SS, vector) \
312 makeFvDdtTypeScheme(SS, sphericalTensor) \
313 makeFvDdtTypeScheme(SS, symmTensor) \
314 makeFvDdtTypeScheme(SS, tensor) \
315  \
316 namespace Foam \
317 { \
318 namespace fv \
319 { \
320  \
321 template<> \
322 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
323 ( \
324  const volScalarField& U, \
325  const surfaceScalarField& Uf \
326 ) \
327 { \
328  NotImplemented; \
329  return surfaceScalarField::null(); \
330 } \
331  \
332 template<> \
333 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
334 ( \
335  const volScalarField& U, \
336  const surfaceScalarField& phi \
337 ) \
338 { \
339  NotImplemented; \
340  return surfaceScalarField::null(); \
341 } \
342  \
343 template<> \
344 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
345 ( \
346  const volScalarField& rho, \
347  const volScalarField& U, \
348  const surfaceScalarField& rhoUf \
349 ) \
350 { \
351  NotImplemented; \
352  return surfaceScalarField::null(); \
353 } \
354  \
355 template<> \
356 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
357 ( \
358  const volScalarField& rho, \
359  const volScalarField& U, \
360  const surfaceScalarField& phi \
361 ) \
362 { \
363  NotImplemented; \
364  return surfaceScalarField::null(); \
365 } \
366  \
367 template<> \
368 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
369 ( \
370  const volScalarField& alpha, \
371  const volScalarField& rho, \
372  const volScalarField& U, \
373  const surfaceScalarField& Uf \
374 ) \
375 { \
376  NotImplemented; \
377  return surfaceScalarField::null(); \
378 } \
379  \
380 template<> \
381 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
382 ( \
383  const volScalarField& alpha, \
384  const volScalarField& rho, \
385  const volScalarField& U, \
386  const surfaceScalarField& phi \
387 ) \
388 { \
389  NotImplemented; \
390  return surfaceScalarField::null(); \
391 } \
392  \
393 } \
394 }
395 
396 
397 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
398 
399 #ifdef NoRepository
400  #include "ddtScheme.C"
401 #endif
402 
403 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
404 
405 #endif
406 
407 // ************************************************************************* //
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Abstract base class for ddt schemes.
Definition: ddtScheme.H:67
SurfaceField< typename flux< Type >::type > fluxFieldType
Definition: ddtScheme.H:192
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)=0
virtual const word & type() const =0
Runtime type information.
virtual ~ddtScheme()
Destructor.
Definition: ddtScheme.C:89
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)=0
const fvMesh & mesh_
Definition: ddtScheme.H:73
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:129
virtual tmp< surfaceScalarField > fvcDdtPhiCoeff(const VolField< Type > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:152
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:45
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)=0
declareRunTimeSelectionTable(tmp, ddtScheme, Istream,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)=0
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const VolField< Type > &U, const fluxFieldType &phi)=0
void operator=(const ddtScheme &)=delete
Disallow default bitwise assignment.
ddtScheme(const fvMesh &mesh)
Construct from mesh.
Definition: ddtScheme.H:97
Reference counter for various OpenFOAM components.
Definition: refCount.H:50
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
labelList fv(nPoints)
Macros to ease declaration of run-time selection tables.