Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::fv::ddtScheme
29
30Group
31 grpFvDdtSchemes
32
33Description
34 Abstract base class for ddt schemes.
35
36SourceFiles
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"
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace Foam
54{
55
56template<class Type>
57class fvMatrix;
58
59class fvMesh;
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63namespace fv
64{
66/*---------------------------------------------------------------------------*\
67 Class ddtScheme Declaration
68\*---------------------------------------------------------------------------*/
69
70class ddtSchemeBase
71{
72public:
74 //- Flag to use experimental ddtCorr from org version
75 //- Default is off for backwards compatibility
76 static bool experimentalDdtCorr;
77
79 {}
80};
82
83template<class Type>
84class ddtScheme
85:
86 public refCount,
87 public ddtSchemeBase
88{
89
90protected:
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
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
111public:
113 //- Runtime type information
114 virtual const word& type() const = 0;
115
116
117 // Declare run-time constructor selection tables
118
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)
141 mesh_(mesh),
142 ddtPhiCoeff_(-1)
143 {}
144
145
146 // Selectors
147
148 //- Return a pointer to a new ddtScheme created on freestore
150 (
151 const fvMesh& mesh,
152 Istream& schemeData
153 );
154
155
156 //- Destructor
157 virtual ~ddtScheme() = default;
158
159
160 // Member Functions
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
183
185 (
186 const volScalarField&,
188 ) = 0;
191 (
192 const volScalarField& alpha,
193 const volScalarField& rho,
195 ) = 0;
196
198 (
200 );
201
203 (
205 ) = 0;
208 (
209 const dimensionedScalar&,
211 ) = 0;
212
214 (
215 const volScalarField&,
217 ) = 0;
220 (
221 const volScalarField& alpha,
222 const volScalarField& rho,
224 ) = 0;
225
226 typedef GeometricField
227 <
228 typename flux<Type>::type,
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
258 const fluxFieldType& phi
259 );
260
262 (
264 const fluxFieldType& phi,
265 const volScalarField& rho
266 );
267
269 (
272 ) = 0;
273
275 (
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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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 \
327makeFvDdtTypeScheme(SS, scalar) \
328makeFvDdtTypeScheme(SS, vector) \
329makeFvDdtTypeScheme(SS, sphericalTensor) \
330makeFvDdtTypeScheme(SS, symmTensor) \
331makeFvDdtTypeScheme(SS, tensor) \
332 \
333namespace Foam \
334{ \
335namespace fv \
336{ \
337 \
338template<> \
339tmp<surfaceScalarField> SS<scalar>::fvcDdtUfCorr \
340( \
341 const volScalarField& U, \
342 const surfaceScalarField& Uf \
343) \
344{ \
345 NotImplemented; \
346 return surfaceScalarField::null(); \
347} \
348 \
349template<> \
350tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr \
351( \
352 const volScalarField& U, \
353 const surfaceScalarField& phi \
354) \
355{ \
356 NotImplemented; \
357 return surfaceScalarField::null(); \
358} \
359 \
360template<> \
361tmp<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 \
372template<> \
373tmp<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// ************************************************************************* //
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.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static bool experimentalDdtCorr
Flag to use experimental ddtCorr from org version Default is off for backwards compatibility.
Definition ddtScheme.H:73
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)=0
ddtScheme(const ddtScheme &)=delete
No copy construct.
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)=0
virtual ~ddtScheme()=default
Destructor.
const fvMesh & mesh_
Definition ddtScheme.H:91
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
Definition ddtScheme.C:303
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &)=0
Definition ddtScheme.C:85
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const GeometricField< Type, fvPatchField, volMesh > &)=0
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition ddtScheme.C:43
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)=0
virtual tmp< fvMatrix< Type > > fvmDdt(const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &)=0
ddtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
Definition ddtScheme.H:149
declareRunTimeSelectionTable(tmp, ddtScheme, Istream,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)=0
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)=0
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)=0
virtual tmp< fvMatrix< Type > > fvmDdt(const volScalarField &, const GeometricField< Type, fvPatchField, volMesh > &)=0
tmp< surfaceScalarField > fvcDdtPhiCoeffExperimental(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition ddtScheme.C:196
const fvMesh & mesh() const
Definition ddtScheme.H:179
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr, const volScalarField &rho)
Definition ddtScheme.C:277
void operator=(const ddtScheme &)=delete
No copy assignment.
virtual tmp< fvMatrix< Type > > fvmDdt(const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &vf)=0
Definition ddtScheme.C:98
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)=0
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)=0
virtual tmp< fluxFieldType > fvcDdtUfCorr(const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)=0
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fvcDdt(const GeometricField< Type, fvsPatchField, surfaceMesh > &)
Definition ddtScheme.C:111
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &rhoU, const fluxFieldType &phi, const volScalarField &rho)
Definition ddtScheme.C:333
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition ddtScheme.C:122
virtual const word & type() const =0
Runtime type information.
ddtScheme(const fvMesh &mesh)
Construct from mesh.
Definition ddtScheme.H:140
GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > fluxFieldType
Definition ddtScheme.H:247
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
typeOfRank< typenamepTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank) -2 >::type type
Definition products.H:155
constexpr refCount() noexcept
Default construct, initializing count to 0.
Definition refCount.H:63
Mesh data needed to do the Finite Volume discretisation.
Definition surfaceMesh.H:47
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
Namespace for finite-volume.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & alpha
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
Basic run-time type information using word as the type's name. Used to enhance the standard RTTI to c...
Forwards and collection of common volume field types.