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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9  Copyright (C) 2017 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::fv::ddtScheme
29 
30 Group
31  grpFvDdtSchemes
32 
33 Description
34  Abstract base class for ddt schemes.
35 
36 SourceFiles
37  ddtScheme.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef ddtScheme_H
42 #define ddtScheme_H
43 
44 #include "tmp.H"
45 #include "dimensionedType.H"
46 #include "volFieldsFwd.H"
47 #include "surfaceFieldsFwd.H"
48 #include "typeInfo.H"
49 #include "runTimeSelectionTables.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 template<class Type>
57 class fvMatrix;
58 
59 class fvMesh;
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace fv
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class ddtScheme Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class ddtSchemeBase
71 {
72 public:
73 
74  //- Flag to use experimental ddtCorr from org version
75  //- Default is off for backwards compatibility
76  static bool experimentalDdtCorr;
77 
79  {}
80 };
81 
82 
83 template<class Type>
84 class ddtScheme
85 :
86  public refCount,
87  public ddtSchemeBase
88 {
89 
90 protected:
91 
92  // Protected data
93 
94  const fvMesh& mesh_;
95 
96  //- Input for fvcDdtPhiCoeff
97  // If set to -1 (default) the code will calculate the coefficient
98  // If > 0 the coupling coeff is set to this value
99  scalar ddtPhiCoeff_;
100 
101 
102  // Private Member Functions
103 
104  //- No copy construct
105  ddtScheme(const ddtScheme&) = delete;
106 
107  //- No copy assignment
108  void operator=(const ddtScheme&) = delete;
109 
110 
111 public:
112 
113  //- Runtime type information
114  virtual const word& type() const = 0;
115 
116 
117  // Declare run-time constructor selection tables
118 
120  (
121  tmp,
122  ddtScheme,
123  Istream,
124  (const fvMesh& mesh, Istream& schemeData),
125  (mesh, schemeData)
126  );
127 
128 
129  // Constructors
130 
131  //- Construct from mesh
132  ddtScheme(const fvMesh& mesh)
133  :
134  mesh_(mesh),
135  ddtPhiCoeff_(-1)
136  {}
137 
138  //- Construct from mesh and Istream
139  ddtScheme(const fvMesh& mesh, Istream& is)
140  :
141  mesh_(mesh),
142  ddtPhiCoeff_(-1)
143  {}
144 
145 
146  // Selectors
147 
148  //- Return a pointer to a new ddtScheme created on freestore
149  static tmp<ddtScheme<Type>> New
150  (
151  const fvMesh& mesh,
152  Istream& schemeData
153  );
154 
155 
156  //- Destructor
157  virtual ~ddtScheme() = default;
158 
159 
160  // Member Functions
161 
162  //- Return mesh reference
163  const fvMesh& mesh() const
164  {
165  return mesh_;
166  }
167 
169  (
170  const dimensioned<Type>&
171  ) = 0;
172 
174  (
176  ) = 0;
177 
179  (
180  const dimensionedScalar&,
182  ) = 0;
183 
185  (
186  const volScalarField&,
188  ) = 0;
189 
191  (
192  const volScalarField& alpha,
193  const volScalarField& rho,
195  ) = 0;
196 
198  (
200  );
201 
202  virtual tmp<fvMatrix<Type>> fvmDdt
203  (
205  ) = 0;
206 
207  virtual tmp<fvMatrix<Type>> fvmDdt
208  (
209  const dimensionedScalar&,
211  ) = 0;
212 
213  virtual tmp<fvMatrix<Type>> fvmDdt
214  (
215  const volScalarField&,
217  ) = 0;
218 
219  virtual tmp<fvMatrix<Type>> fvmDdt
220  (
221  const volScalarField& alpha,
222  const volScalarField& rho,
224  ) = 0;
225 
226  typedef GeometricField
227  <
228  typename flux<Type>::type,
231  > fluxFieldType;
232 
234  (
236  const fluxFieldType& phi,
237  const fluxFieldType& phiCorr
238  );
239 
241  (
243  const fluxFieldType& phi,
244  const fluxFieldType& phiCorr
245  );
246 
248  (
250  const fluxFieldType& phi,
251  const fluxFieldType& phiCorr,
252  const volScalarField& rho
253  );
254 
256  (
258  const fluxFieldType& phi
259  );
260 
262  (
264  const fluxFieldType& phi,
265  const volScalarField& rho
266  );
267 
269  (
272  ) = 0;
273 
275  (
277  const fluxFieldType& phi
278  ) = 0;
279 
281  (
282  const volScalarField& rho,
285  ) = 0;
286 
288  (
289  const volScalarField& rho,
291  const fluxFieldType& phi
292  ) = 0;
293 
295  (
297  ) = 0;
298 };
299 
300 
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302 
303 } // End namespace fv
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 } // End namespace Foam
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 // Add the patch constructor functions to the hash tables
312 
313 #define makeFvDdtTypeScheme(SS, Type) \
314  defineNamedTemplateTypeNameAndDebug(Foam::fv::SS<Foam::Type>, 0); \
315  \
316  namespace Foam \
317  { \
318  namespace fv \
319  { \
320  ddtScheme<Type>::addIstreamConstructorToTable<SS<Type>> \
321  add##SS##Type##IstreamConstructorToTable_; \
322  } \
323  }
324 
325 #define makeFvDdtScheme(SS) \
326  \
327 makeFvDdtTypeScheme(SS, scalar) \
328 makeFvDdtTypeScheme(SS, vector) \
329 makeFvDdtTypeScheme(SS, sphericalTensor) \
330 makeFvDdtTypeScheme(SS, symmTensor) \
331 makeFvDdtTypeScheme(SS, tensor) \
332  \
333 namespace Foam \
334 { \
335 namespace fv \
336 { \
337  \
338 template<> \
339 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
340 ( \
341  const volScalarField& U, \
342  const surfaceScalarField& Uf \
343 ) \
344 { \
345  NotImplemented; \
346  return surfaceScalarField::null(); \
347 } \
348  \
349 template<> \
350 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
351 ( \
352  const volScalarField& U, \
353  const surfaceScalarField& phi \
354 ) \
355 { \
356  NotImplemented; \
357  return surfaceScalarField::null(); \
358 } \
359  \
360 template<> \
361 tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
362 ( \
363  const volScalarField& rho, \
364  const volScalarField& U, \
365  const surfaceScalarField& Uf \
366 ) \
367 { \
368  NotImplemented; \
369  return surfaceScalarField::null(); \
370 } \
371  \
372 template<> \
373 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
374 ( \
375  const volScalarField& rho, \
376  const volScalarField& U, \
377  const surfaceScalarField& phi \
378 ) \
379 { \
380  NotImplemented; \
381  return surfaceScalarField::null(); \
382 } \
383  \
384 } \
385 }
386 
387 
388 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
389 
390 #ifdef NoRepository
391  #include "ddtScheme.C"
392 #endif
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 #endif
397 
398 // ************************************************************************* //
volFieldsFwd.H
Foam::fv::ddtScheme::ddtScheme
ddtScheme(const ddtScheme &)=delete
No copy construct.
Foam::fv::ddtScheme::ddtScheme
ddtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
Definition: ddtScheme.H:138
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
typeInfo.H
Foam::fv::ddtScheme::fvcDdtPhiCoeffExperimental
tmp< surfaceScalarField > fvcDdtPhiCoeffExperimental(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:218
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::refCount
Reference counter for various OpenFOAM components.
Definition: refCount.H:50
Foam::fv::ddtScheme::mesh_
const fvMesh & mesh_
Definition: ddtScheme.H:93
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::surfaceMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:49
Foam::fv::ddtSchemeBase::experimentalDdtCorr
static bool experimentalDdtCorr
Definition: ddtScheme.H:75
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
rho
rho
Definition: readInitialConditions.H:88
ddtScheme.C
Foam::fv::ddtScheme::operator=
void operator=(const ddtScheme &)=delete
No copy assignment.
Foam::fv::ddtScheme::meshPhi
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)=0
Foam::fv::ddtScheme::ddtScheme
ddtScheme(const fvMesh &mesh)
Construct from mesh.
Definition: ddtScheme.H:131
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::fv::ddtScheme::fvcDdtPhiCoeff
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:144
Foam::fv::ddtScheme::type
virtual const word & type() const =0
Runtime type information.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fv::ddtScheme::fvcDdtPhiCorr
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)=0
Foam::fv::ddtScheme::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:162
fv
labelList fv(nPoints)
U
U
Definition: pEqn.H:72
tmp.H
Foam::fv::ddtSchemeBase
Non-templated base class for ddt schemes.
Definition: ddtScheme.H:69
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::fv::ddtScheme::fvcDdt
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)=0
Foam::fv::ddtScheme::~ddtScheme
virtual ~ddtScheme()=default
Destructor.
surfaceFieldsFwd.H
Foam::fv::ddtSchemeBase::ddtSchemeBase
ddtSchemeBase()
Definition: ddtScheme.H:77
Foam::fv::ddtScheme
Abstract base class for ddt schemes.
Definition: ddtScheme.H:83
Foam::fv::ddtScheme::ddtPhiCoeff_
scalar ddtPhiCoeff_
Input for fvcDdtPhiCoeff.
Definition: ddtScheme.H:98
Foam::fv::ddtScheme::fvcDdtUfCorr
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)=0
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::fv::ddtScheme::New
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:50
Foam::fv::ddtScheme::declareRunTimeSelectionTable
declareRunTimeSelectionTable(tmp, ddtScheme, Istream,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
Foam::fv::ddtScheme::fluxFieldType
GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > fluxFieldType
Definition: ddtScheme.H:230
Foam::innerProduct< vector, Type >::type
typeOfRank< typename pTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank) - 2 >::type type
Definition: products.H:149
dimensionedType.H
Foam::fv::ddtScheme::fvmDdt
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)=0