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-2024 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 "scalarField.H"
41 #include "volFieldsFwd.H"
42 #include "surfaceFieldsFwd.H"
43 #include "typeInfo.H"
44 #include "runTimeSelectionTables.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 template<class Type>
52 class fvMatrix;
53 
54 class fvMesh;
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace fv
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class ddtScheme Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 template<class Type>
66 class ddtScheme
67 :
68  public tmp<ddtScheme<Type>>::refCount
69 {
70 
71 protected:
72 
73  // Protected data
74 
75  const fvMesh& mesh_;
76 
77 
78 public:
79 
80  //- Runtime type information
81  virtual const word& type() const = 0;
82 
83 
84  // Declare run-time constructor selection tables
85 
87  (
88  tmp,
89  ddtScheme,
90  Istream,
91  (const fvMesh& mesh, Istream& schemeData),
92  (mesh, schemeData)
93  );
94 
95 
96  // Constructors
97 
98  //- Construct from mesh
99  ddtScheme(const fvMesh& mesh)
100  :
101  mesh_(mesh)
102  {}
103 
104  //- Construct from mesh and Istream
105  ddtScheme(const fvMesh& mesh, Istream&)
106  :
107  mesh_(mesh)
108  {}
109 
110  //- Disallow default bitwise copy construction
111  ddtScheme(const ddtScheme&) = delete;
112 
113 
114  // Selectors
115 
116  //- Return a pointer to a new ddtScheme created on freestore
117  static tmp<ddtScheme<Type>> New
118  (
119  const fvMesh& mesh,
120  Istream& schemeData
121  );
122 
123 
124  //- Destructor
125  virtual ~ddtScheme();
126 
127 
128  // Member Functions
129 
130  //- Return mesh reference
131  const fvMesh& mesh() const
132  {
133  return mesh_;
134  }
135 
136  virtual tmp<VolField<Type>> fvcDdt
137  (
138  const dimensioned<Type>&
139  ) = 0;
140 
141  virtual tmp<VolField<Type>> fvcDdt
142  (
143  const VolField<Type>&
144  ) = 0;
145 
146  virtual tmp<VolField<Type>> fvcDdt
147  (
148  const dimensionedScalar&,
149  const VolField<Type>&
150  ) = 0;
151 
152  virtual tmp<VolField<Type>> fvcDdt
153  (
154  const volScalarField&,
155  const VolField<Type>&
156  ) = 0;
157 
158  virtual tmp<VolField<Type>> fvcDdt
159  (
160  const volScalarField& alpha,
161  const volScalarField& rho,
162  const VolField<Type>&
163  ) = 0;
164 
166  (
167  const SurfaceField<Type>&
168  );
169 
170  virtual tmp<fvMatrix<Type>> fvmDdt
171  (
172  const VolField<Type>&
173  ) = 0;
174 
175  virtual tmp<fvMatrix<Type>> fvmDdt
176  (
177  const dimensionedScalar&,
178  const VolField<Type>&
179  ) = 0;
180 
181  virtual tmp<fvMatrix<Type>> fvmDdt
182  (
183  const volScalarField&,
184  const VolField<Type>&
185  ) = 0;
186 
187  virtual tmp<fvMatrix<Type>> fvmDdt
188  (
189  const volScalarField& alpha,
190  const volScalarField& rho,
191  const VolField<Type>& vf
192  ) = 0;
193 
195 
197  (
198  const VolField<Type>& U,
199  const fluxFieldType& phi,
200  const fluxFieldType& phiCorr
201  );
202 
204  (
205  const VolField<Type>& U,
206  const fluxFieldType& phi,
207  const fluxFieldType& phiCorr,
208  const volScalarField& rho
209  );
210 
212  (
213  const VolField<Type>& U,
214  const fluxFieldType& phi
215  );
216 
218  (
219  const VolField<Type>& U,
220  const fluxFieldType& phi,
221  const volScalarField& rho
222  );
223 
225  (
226  const VolField<Type>& U,
227  const SurfaceField<Type>& Uf
228  ) = 0;
229 
231  (
232  const VolField<Type>& U,
233  const fluxFieldType& phi
234  ) = 0;
235 
237  (
238  const volScalarField& rho,
239  const VolField<Type>& U,
240  const SurfaceField<Type>& rhoUf
241  ) = 0;
242 
244  (
245  const volScalarField& rho,
246  const VolField<Type>& U,
247  const fluxFieldType& phi
248  ) = 0;
249 
251  (
252  const volScalarField& alpha,
253  const volScalarField& rho,
254  const VolField<Type>& U,
255  const SurfaceField<Type>& Uf
256  ) = 0;
257 
259  (
260  const volScalarField& alpha,
261  const volScalarField& rho,
262  const VolField<Type>& U,
263  const fluxFieldType& phi
264  ) = 0;
265 
267  (
268  const VolField<Type>&
269  ) = 0;
270 
271  virtual tmp<scalarField> meshPhi
272  (
273  const VolField<Type>&,
274  const label patchi
275  ) = 0;
276 
277 
278  // Member Operators
279 
280  //- Disallow default bitwise assignment
281  void operator=(const ddtScheme&) = delete;
282 };
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 } // End namespace fv
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 } // End namespace Foam
292 
293 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
294 
295 // Add the patch constructor functions to the hash tables
296 
297 #define makeFvDdtTypeScheme(SS, Type) \
298  defineNamedTemplateTypeNameAndDebug(Foam::fv::SS<Foam::Type>, 0); \
299  \
300  namespace Foam \
301  { \
302  namespace fv \
303  { \
304  ddtScheme<Type>::addIstreamConstructorToTable<SS<Type>> \
305  add##SS##Type##IstreamConstructorToTable_; \
306  } \
307  }
308 
309 #define makeFvDdtScheme(SS) \
310  \
311 makeFvDdtTypeScheme(SS, scalar) \
312 makeFvDdtTypeScheme(SS, vector) \
313 makeFvDdtTypeScheme(SS, sphericalTensor) \
314 makeFvDdtTypeScheme(SS, symmTensor) \
315 makeFvDdtTypeScheme(SS, tensor) \
316  \
317 namespace Foam \
318 { \
319 namespace fv \
320 { \
321  \
322 template<> \
323 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
324 ( \
325  const volScalarField& U, \
326  const surfaceScalarField& Uf \
327 ) \
328 { \
329  NotImplemented; \
330  return surfaceScalarField::null(); \
331 } \
332  \
333 template<> \
334 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
335 ( \
336  const volScalarField& U, \
337  const surfaceScalarField& phi \
338 ) \
339 { \
340  NotImplemented; \
341  return surfaceScalarField::null(); \
342 } \
343  \
344 template<> \
345 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
346 ( \
347  const volScalarField& rho, \
348  const volScalarField& U, \
349  const surfaceScalarField& rhoUf \
350 ) \
351 { \
352  NotImplemented; \
353  return surfaceScalarField::null(); \
354 } \
355  \
356 template<> \
357 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
358 ( \
359  const volScalarField& rho, \
360  const volScalarField& U, \
361  const surfaceScalarField& phi \
362 ) \
363 { \
364  NotImplemented; \
365  return surfaceScalarField::null(); \
366 } \
367  \
368 template<> \
369 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
370 ( \
371  const volScalarField& alpha, \
372  const volScalarField& rho, \
373  const volScalarField& U, \
374  const surfaceScalarField& Uf \
375 ) \
376 { \
377  NotImplemented; \
378  return surfaceScalarField::null(); \
379 } \
380  \
381 template<> \
382 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
383 ( \
384  const volScalarField& alpha, \
385  const volScalarField& rho, \
386  const volScalarField& U, \
387  const surfaceScalarField& phi \
388 ) \
389 { \
390  NotImplemented; \
391  return surfaceScalarField::null(); \
392 } \
393  \
394 } \
395 }
396 
397 
398 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399 
400 #ifdef NoRepository
401  #include "ddtScheme.C"
402 #endif
403 
404 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
405 
406 #endif
407 
408 // ************************************************************************* //
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:99
Abstract base class for ddt schemes.
Definition: ddtScheme.H:68
SurfaceField< typename flux< Type >::type > fluxFieldType
Definition: ddtScheme.H:193
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:74
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:130
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:98
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.
Basic run-time type information using word as the type's name. Used to enhance the standard RTTI to c...