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-2022 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 
136  (
137  const dimensioned<Type>&
138  ) = 0;
139 
141  (
143  ) = 0;
144 
146  (
147  const dimensionedScalar&,
149  ) = 0;
150 
152  (
153  const volScalarField&,
155  ) = 0;
156 
158  (
159  const volScalarField& alpha,
160  const volScalarField& rho,
162  ) = 0;
163 
165  (
167  );
168 
169  virtual tmp<fvMatrix<Type>> fvmDdt
170  (
172  ) = 0;
173 
174  virtual tmp<fvMatrix<Type>> fvmDdt
175  (
176  const dimensionedScalar&,
178  ) = 0;
179 
180  virtual tmp<fvMatrix<Type>> fvmDdt
181  (
182  const volScalarField&,
184  ) = 0;
185 
186  virtual tmp<fvMatrix<Type>> fvmDdt
187  (
188  const volScalarField& alpha,
189  const volScalarField& rho,
191  ) = 0;
192 
193  typedef GeometricField
194  <
195  typename flux<Type>::type,
198  > fluxFieldType;
199 
201  (
203  const fluxFieldType& phi,
204  const fluxFieldType& phiCorr
205  );
206 
208  (
210  const fluxFieldType& phi,
211  const fluxFieldType& phiCorr,
212  const volScalarField& rho
213  );
214 
216  (
218  const fluxFieldType& phi
219  );
220 
222  (
224  const fluxFieldType& phi,
225  const volScalarField& rho
226  );
227 
229  (
232  ) = 0;
233 
235  (
237  const fluxFieldType& phi
238  ) = 0;
239 
241  (
242  const volScalarField& rho,
245  ) = 0;
246 
248  (
249  const volScalarField& rho,
251  const fluxFieldType& phi
252  ) = 0;
253 
255  (
257  ) = 0;
258 
259  virtual tmp<scalarField> meshPhi
260  (
262  const label patchi
263  ) = 0;
264 
265 
266  // Member Operators
267 
268  //- Disallow default bitwise assignment
269  void operator=(const ddtScheme&) = delete;
270 };
271 
272 
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 
275 } // End namespace fv
276 
277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 
279 } // End namespace Foam
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
283 // Add the patch constructor functions to the hash tables
285 #define makeFvDdtTypeScheme(SS, Type) \
286  defineNamedTemplateTypeNameAndDebug(Foam::fv::SS<Foam::Type>, 0); \
287  \
288  namespace Foam \
289  { \
290  namespace fv \
291  { \
292  ddtScheme<Type>::addIstreamConstructorToTable<SS<Type>> \
293  add##SS##Type##IstreamConstructorToTable_; \
294  } \
295  }
297 #define makeFvDdtScheme(SS) \
298  \
299 makeFvDdtTypeScheme(SS, scalar) \
300 makeFvDdtTypeScheme(SS, vector) \
301 makeFvDdtTypeScheme(SS, sphericalTensor) \
302 makeFvDdtTypeScheme(SS, symmTensor) \
303 makeFvDdtTypeScheme(SS, tensor) \
304  \
305 namespace Foam \
306 { \
307 namespace fv \
308 { \
309  \
310 template<> \
311 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
312 ( \
313  const volScalarField& U, \
314  const surfaceScalarField& Uf \
315 ) \
316 { \
317  NotImplemented; \
318  return surfaceScalarField::null(); \
319 } \
320  \
321 template<> \
322 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
323 ( \
324  const volScalarField& U, \
325  const surfaceScalarField& phi \
326 ) \
327 { \
328  NotImplemented; \
329  return surfaceScalarField::null(); \
330 } \
331  \
332 template<> \
333 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
334 ( \
335  const volScalarField& rho, \
336  const volScalarField& U, \
337  const surfaceScalarField& Uf \
338 ) \
339 { \
340  NotImplemented; \
341  return surfaceScalarField::null(); \
342 } \
343  \
344 template<> \
345 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
346 ( \
347  const volScalarField& rho, \
348  const volScalarField& U, \
349  const surfaceScalarField& phi \
350 ) \
351 { \
352  NotImplemented; \
353  return surfaceScalarField::null(); \
354 } \
355  \
356 } \
357 }
358 
359 
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 
362 #ifdef NoRepository
363  #include "ddtScheme.C"
364 #endif
365 
366 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
367 
368 #endif
369 
370 // ************************************************************************* //
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)=0
void operator=(const ddtScheme &)=delete
Disallow default bitwise assignment.
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:46
U
Definition: pEqn.H:72
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
ddtScheme(const fvMesh &mesh)
Construct from mesh.
Definition: ddtScheme.H:97
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:197
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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:129
virtual ~ddtScheme()
Destructor.
Definition: ddtScheme.C:90
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
phi
Definition: correctPhi.H:3
autoPtr< surfaceVectorField > Uf
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)=0
virtual const word & type() const =0
Runtime type information.
virtual tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:153
label patchi
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:95
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
const fvMesh & mesh_
Definition: ddtScheme.H:73
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.