initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / ddtSchemes / ddtScheme / ddtScheme.H
blob31f8fa090d75eafd520384fc3777b2f6948893c4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::fv::ddtScheme
28 Description
29     Abstract base class for ddt schemes.
31 SourceFiles
32     ddtScheme.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef ddtScheme_H
37 #define ddtScheme_H
39 #include "tmp.H"
40 #include "dimensionedType.H"
41 #include "volFieldsFwd.H"
42 #include "surfaceFieldsFwd.H"
43 #include "typeInfo.H"
44 #include "runTimeSelectionTables.H"
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 namespace Foam
51 template<class Type>
52 class fvMatrix;
54 class fvMesh;
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 namespace fv
61 /*---------------------------------------------------------------------------*\
62                            Class ddtScheme Declaration
63 \*---------------------------------------------------------------------------*/
65 template<class Type>
66 class ddtScheme
68     public refCount
71 protected:
73     // Protected data
75         const fvMesh& mesh_;
78     // Private Member Functions
80         //- Disallow copy construct
81         ddtScheme(const ddtScheme&);
83         //- Disallow default bitwise assignment
84         void operator=(const ddtScheme&);
87 public:
89     //- Runtime type information
90     virtual const word& type() const = 0;
93     // Declare run-time constructor selection tables
95         declareRunTimeSelectionTable
96         (
97             tmp,
98             ddtScheme,
99             Istream,
100             (const fvMesh& mesh, Istream& schemeData),
101             (mesh, schemeData)
102         );
105     // Constructors
107         //- Construct from mesh
108         ddtScheme(const fvMesh& mesh)
109         :
110             mesh_(mesh)
111         {}
113         //- Construct from mesh and Istream
114         ddtScheme(const fvMesh& mesh, Istream&)
115         :
116             mesh_(mesh)
117         {}
120     // Selectors
122         //- Return a pointer to a new ddtScheme created on freestore
123         static tmp<ddtScheme<Type> > New
124         (
125             const fvMesh& mesh,
126             Istream& schemeData
127         );
130     // Destructor
132         virtual ~ddtScheme();
135     // Member Functions
137         //- Return mesh reference
138         const fvMesh& mesh() const
139         {
140             return mesh_;
141         }
143         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
144         (
145             const dimensioned<Type>&
146         ) = 0;
148         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
149         (
150             const GeometricField<Type, fvPatchField, volMesh>&
151         ) = 0;
153         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
154         (
155             const dimensionedScalar&,
156             const GeometricField<Type, fvPatchField, volMesh>&
157         ) = 0;
159         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
160         (
161             const volScalarField&,
162             const GeometricField<Type, fvPatchField, volMesh>&
163         ) = 0;
165         virtual tmp<fvMatrix<Type> > fvmDdt
166         (
167             GeometricField<Type, fvPatchField, volMesh>&
168         ) = 0;
170         virtual tmp<fvMatrix<Type> > fvmDdt
171         (
172             const dimensionedScalar&,
173             GeometricField<Type, fvPatchField, volMesh>&
174         ) = 0;
176         virtual tmp<fvMatrix<Type> > fvmDdt
177         (
178             const volScalarField&,
179             GeometricField<Type, fvPatchField, volMesh>&
180         ) = 0;
183         typedef GeometricField
184         <
185             typename flux<Type>::type,
186             fvsPatchField,
187             surfaceMesh
188         > fluxFieldType;
190         tmp<surfaceScalarField> fvcDdtPhiCoeff
191         (
192             const GeometricField<Type, fvPatchField, volMesh>& U,
193             const fluxFieldType& phi,
194             const fluxFieldType& phiCorr
195         );
197         tmp<surfaceScalarField> fvcDdtPhiCoeff
198         (
199             const GeometricField<Type, fvPatchField, volMesh>& U,
200             const fluxFieldType& phi
201         );
203         virtual tmp<fluxFieldType> fvcDdtPhiCorr
204         (
205             const volScalarField& rA,
206             const GeometricField<Type, fvPatchField, volMesh>& U,
207             const fluxFieldType& phi
208         ) = 0;
210         tmp<surfaceScalarField> fvcDdtPhiCoeff
211         (
212             const volScalarField& rho,
213             const GeometricField<Type, fvPatchField, volMesh>& rhoU,
214             const fluxFieldType& phi
215         );
217         virtual tmp<fluxFieldType> fvcDdtPhiCorr
218         (
219             const volScalarField& rA,
220             const volScalarField& rho,
221             const GeometricField<Type, fvPatchField, volMesh>& U,
222             const fluxFieldType& phi
223         ) = 0;
226         virtual tmp<surfaceScalarField> meshPhi
227         (
228             const GeometricField<Type, fvPatchField, volMesh>&
229         ) = 0;
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 } // End namespace fv
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 } // End namespace Foam
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 // Add the patch constructor functions to the hash tables
245 #define makeFvDdtTypeScheme(SS, Type)                                          \
246                                                                                \
247 defineNamedTemplateTypeNameAndDebug(SS<Type>, 0);                              \
248                                                                                \
249 ddtScheme<Type>::addIstreamConstructorToTable<SS<Type> >                       \
250     add##SS##Type##IstreamConstructorToTable_;
253 #define makeFvDdtScheme(SS)                                                    \
254                                                                                \
255 makeFvDdtTypeScheme(SS, scalar)                                                \
256 makeFvDdtTypeScheme(SS, vector)                                                \
257 makeFvDdtTypeScheme(SS, sphericalTensor)                                       \
258 makeFvDdtTypeScheme(SS, symmTensor)                                            \
259 makeFvDdtTypeScheme(SS, tensor)                                                \
260                                                                                \
261 template<>                                                                     \
262 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr                              \
263 (                                                                              \
264     const volScalarField& rA,                                                  \
265     const volScalarField& U,                                                   \
266     const surfaceScalarField& phi                                              \
267 )                                                                              \
268 {                                                                              \
269     notImplemented(#SS"<scalar>::fvcDdtPhiCorr");                              \
270     return surfaceScalarField::null();                                         \
271 }                                                                              \
272                                                                                \
273 template<>                                                                     \
274 tmp<surfaceScalarField> SS<scalar>::fvcDdtPhiCorr                              \
275 (                                                                              \
276     const volScalarField& rA,                                                  \
277     const volScalarField& rho,                                                 \
278     const volScalarField& U,                                                   \
279     const surfaceScalarField& phi                                              \
280 )                                                                              \
281 {                                                                              \
282     notImplemented(#SS"<scalar>::fvcDdtPhiCorr");                              \
283     return surfaceScalarField::null();                                         \
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 #ifdef NoRepository
290 #   include "ddtScheme.C"
291 #endif
293 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 #endif
297 // ************************************************************************* //