initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / laplacianSchemes / laplacianScheme / laplacianScheme.H
blobb36d87df9e1a35595976daf4d820c2e94fcd9a4f
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::laplacianScheme
28 Description
29     Abstract base class for laplacian schemes.
31 SourceFiles
32     laplacianScheme.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef laplacianScheme_H
37 #define laplacianScheme_H
39 #include "tmp.H"
40 #include "volFieldsFwd.H"
41 #include "surfaceFieldsFwd.H"
42 #include "linear.H"
43 #include "correctedSnGrad.H"
44 #include "typeInfo.H"
45 #include "runTimeSelectionTables.H"
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 namespace Foam
52 template<class Type>
53 class fvMatrix;
55 class fvMesh;
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 namespace fv
62 /*---------------------------------------------------------------------------*\
63                            Class laplacianScheme Declaration
64 \*---------------------------------------------------------------------------*/
66 template<class Type, class GType>
67 class laplacianScheme
69     public refCount
72 protected:
74     // Protected data
76         const fvMesh& mesh_;
77         tmp<surfaceInterpolationScheme<GType> > tinterpGammaScheme_;
78         tmp<snGradScheme<Type> > tsnGradScheme_;
81     // Private Member Functions
83         //- Disallow copy construct
84         laplacianScheme(const laplacianScheme&);
86         //- Disallow default bitwise assignment
87         void operator=(const laplacianScheme&);
90 public:
92     //- Runtime type information
93     virtual const word& type() const = 0;
96     // Declare run-time constructor selection tables
98         declareRunTimeSelectionTable
99         (
100             tmp,
101             laplacianScheme,
102             Istream,
103             (const fvMesh& mesh, Istream& schemeData),
104             (mesh, schemeData)
105         );
108     // Constructors
110         //- Construct from mesh
111         laplacianScheme(const fvMesh& mesh)
112         :
113             mesh_(mesh),
114             tinterpGammaScheme_(new linear<GType>(mesh)),
115             tsnGradScheme_(new correctedSnGrad<Type>(mesh))
116         {}
118         //- Construct from mesh and Istream
119         laplacianScheme(const fvMesh& mesh, Istream& is)
120         :
121             mesh_(mesh),
122             tinterpGammaScheme_(NULL),
123             tsnGradScheme_(NULL)
124         {
125             tinterpGammaScheme_ = tmp<surfaceInterpolationScheme<GType> >
126             (
127                 surfaceInterpolationScheme<GType>::New(mesh, is)
128             );
130             tsnGradScheme_ = tmp<snGradScheme<Type> >
131             (
132                 snGradScheme<Type>::New(mesh, is)
133             );
134         }
136         //- Construct from mesh, interpolation and snGradScheme schemes
137         laplacianScheme
138         (
139             const fvMesh& mesh,
140             const tmp<surfaceInterpolationScheme<GType> >& igs,
141             const tmp<snGradScheme<Type> >& sngs
142         )
143         :
144             mesh_(mesh),
145             tinterpGammaScheme_(igs),
146             tsnGradScheme_(sngs)
147         {}
150     // Selectors
152         //- Return a pointer to a new laplacianScheme created on freestore
153         static tmp<laplacianScheme<Type, GType> > New
154         (
155             const fvMesh& mesh,
156             Istream& schemeData
157         );
160     // Destructor
162         virtual ~laplacianScheme();
165     // Member Functions
167         //- Return mesh reference
168         const fvMesh& mesh() const
169         {
170             return mesh_;
171         }
173         virtual tmp<fvMatrix<Type> > fvmLaplacian
174         (
175             const GeometricField<GType, fvsPatchField, surfaceMesh>&,
176             GeometricField<Type, fvPatchField, volMesh>&
177         ) = 0;
179         virtual tmp<fvMatrix<Type> > fvmLaplacian
180         (
181             const GeometricField<GType, fvPatchField, volMesh>&,
182             GeometricField<Type, fvPatchField, volMesh>&
183         );
185         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcLaplacian
186         (
187             const GeometricField<Type, fvPatchField, volMesh>&
188         ) = 0;
190         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcLaplacian
191         (
192             const GeometricField<GType, fvsPatchField, surfaceMesh>&,
193             const GeometricField<Type, fvPatchField, volMesh>&
194         ) = 0;
196         virtual tmp<GeometricField<Type, fvPatchField, volMesh> > fvcLaplacian
197         (
198             const GeometricField<GType, fvPatchField, volMesh>&,
199             const GeometricField<Type, fvPatchField, volMesh>&
200         );
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 } // End namespace fv
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 } // End namespace Foam
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 // Add the patch constructor functions to the hash tables
216 #define makeFvLaplacianTypeScheme(SS, Type, GType)                             \
217                                                                                \
218     typedef SS<Type, GType> SS##Type##GType;                                   \
219     defineNamedTemplateTypeNameAndDebug(SS##Type##GType, 0);                   \
220                                                                                \
221     laplacianScheme<Type, GType>::                                             \
222         addIstreamConstructorToTable<SS<Type, GType> >                         \
223     add##SS##Type##GType##IstreamConstructorToTable_;
226 #define makeFvLaplacianScheme(SS)                                              \
227                                                                                \
228 makeFvLaplacianTypeScheme(SS, scalar, scalar)                                  \
229 makeFvLaplacianTypeScheme(SS, scalar, symmTensor)                              \
230 makeFvLaplacianTypeScheme(SS, scalar, tensor)                                  \
231 makeFvLaplacianTypeScheme(SS, vector, scalar)                                  \
232 makeFvLaplacianTypeScheme(SS, sphericalTensor, scalar)                         \
233 makeFvLaplacianTypeScheme(SS, symmTensor, scalar)                              \
234 makeFvLaplacianTypeScheme(SS, tensor, scalar)                                  \
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 #ifdef NoRepository
240 #   include "laplacianScheme.C"
241 #endif
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 #endif
247 // ************************************************************************* //