initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / fvm / fvmLaplacian.C
blob4a60eadd3d3211dbef5ab58cf55e79cdd831b263
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 \*---------------------------------------------------------------------------*/
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvMatrix.H"
30 #include "laplacianScheme.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 namespace fvm
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 template<class Type>
45 tmp<fvMatrix<Type> >
46 laplacian
48     GeometricField<Type, fvPatchField, volMesh>& vf,
49     const word& name
52     surfaceScalarField Gamma
53     (
54         IOobject
55         (
56             "1",
57             vf.time().constant(),
58             vf.mesh(),
59             IOobject::NO_READ
60         ),
61         vf.mesh(),
62         dimensionedScalar("1", dimless, 1.0)
63     );
65     return fvm::laplacian(Gamma, vf, name);
69 template<class Type>
70 tmp<fvMatrix<Type> >
71 laplacian
73     GeometricField<Type, fvPatchField, volMesh>& vf
76     surfaceScalarField Gamma
77     (
78         IOobject
79         (
80             "1",
81             vf.time().constant(),
82             vf.mesh(),
83             IOobject::NO_READ
84         ),
85         vf.mesh(),
86         dimensionedScalar("1", dimless, 1.0)
87     );
89     return fvm::laplacian
90     (
91         Gamma,
92         vf,
93         "laplacian(" + vf.name() + ')'
94     );
98 template<class Type>
99 tmp<fvMatrix<Type> >
100 laplacian
102     const zeroField&,
103     GeometricField<Type, fvPatchField, volMesh>& vf,
104     const word& name
107     return tmp<fvMatrix<Type> >
108     (
109         new fvMatrix<Type>(vf, dimensionSet(0, 0, -2, 0, 0))
110     );
114 template<class Type>
115 tmp<fvMatrix<Type> >
116 laplacian
118     const zeroField&,
119     GeometricField<Type, fvPatchField, volMesh>& vf
122     return tmp<fvMatrix<Type> >
123     (
124         new fvMatrix<Type>(vf, dimensionSet(0, 0, -2, 0, 0))
125     );
129 template<class Type>
130 tmp<fvMatrix<Type> >
131 laplacian
133     const oneField&,
134     GeometricField<Type, fvPatchField, volMesh>& vf,
135     const word& name
138     return fvm::laplacian(vf, name);
142 template<class Type>
143 tmp<fvMatrix<Type> >
144 laplacian
146     const oneField&,
147     GeometricField<Type, fvPatchField, volMesh>& vf
150     return fvm::laplacian(vf);
154 template<class Type, class GType>
155 tmp<fvMatrix<Type> >
156 laplacian
158     const dimensioned<GType>& gamma,
159     GeometricField<Type, fvPatchField, volMesh>& vf,
160     const word& name
163     GeometricField<GType, fvsPatchField, surfaceMesh> Gamma
164     (
165         IOobject
166         (
167             gamma.name(),
168             vf.instance(),
169             vf.mesh(),
170             IOobject::NO_READ
171         ),
172         vf.mesh(),
173         gamma
174     );
176     return fvm::laplacian(Gamma, vf, name);
180 template<class Type, class GType>
181 tmp<fvMatrix<Type> >
182 laplacian
184     const dimensioned<GType>& gamma,
185     GeometricField<Type, fvPatchField, volMesh>& vf
188     GeometricField<GType, fvsPatchField, surfaceMesh> Gamma
189     (
190         IOobject
191         (
192             gamma.name(),
193             vf.instance(),
194             vf.mesh(),
195             IOobject::NO_READ
196         ),
197         vf.mesh(),
198         gamma
199     );
201     return fvm::laplacian(Gamma, vf);
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 template<class Type, class GType>
208 tmp<fvMatrix<Type> >
209 laplacian
211     const GeometricField<GType, fvPatchField, volMesh>& gamma,
212     GeometricField<Type, fvPatchField, volMesh>& vf,
213     const word& name
216     return fv::laplacianScheme<Type, GType>::New
217     (
218         vf.mesh(),
219         vf.mesh().laplacianScheme(name)
220     )().fvmLaplacian(gamma, vf);
224 template<class Type, class GType>
225 tmp<fvMatrix<Type> >
226 laplacian
228     const tmp<GeometricField<GType, fvPatchField, volMesh> >& tgamma,
229     GeometricField<Type, fvPatchField, volMesh>& vf,
230     const word& name
233     tmp<fvMatrix<Type> > Laplacian(fvm::laplacian(tgamma(), vf, name));
234     tgamma.clear();
235     return Laplacian;
239 template<class Type, class GType>
240 tmp<fvMatrix<Type> >
241 laplacian
243     const GeometricField<GType, fvPatchField, volMesh>& gamma,
244     GeometricField<Type, fvPatchField, volMesh>& vf
247     return fvm::laplacian
248     (
249         gamma,
250         vf,
251         "laplacian(" + gamma.name() + ',' + vf.name() + ')'
252     );
256 template<class Type, class GType>
257 tmp<fvMatrix<Type> >
258 laplacian
260     const tmp<GeometricField<GType, fvPatchField, volMesh> >& tgamma,
261     GeometricField<Type, fvPatchField, volMesh>& vf
264     tmp<fvMatrix<Type> > Laplacian(fvm::laplacian(tgamma(), vf));
265     tgamma.clear();
266     return Laplacian;
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 template<class Type, class GType>
273 tmp<fvMatrix<Type> >
274 laplacian
276     const GeometricField<GType, fvsPatchField, surfaceMesh>& gamma,
277     GeometricField<Type, fvPatchField, volMesh>& vf,
278     const word& name
281     return fv::laplacianScheme<Type, GType>::New
282     (
283         vf.mesh(),
284         vf.mesh().laplacianScheme(name)
285     )().fvmLaplacian(gamma, vf);
289 template<class Type, class GType>
290 tmp<fvMatrix<Type> >
291 laplacian
293     const tmp<GeometricField<GType, fvsPatchField, surfaceMesh> >& tgamma,
294     GeometricField<Type, fvPatchField, volMesh>& vf,
295     const word& name
298     tmp<fvMatrix<Type> > tLaplacian = fvm::laplacian(tgamma(), vf, name);
299     tgamma.clear();
300     return tLaplacian;
304 template<class Type, class GType>
305 tmp<fvMatrix<Type> >
306 laplacian
308     const GeometricField<GType, fvsPatchField, surfaceMesh>& gamma,
309     GeometricField<Type, fvPatchField, volMesh>& vf
312     return fvm::laplacian
313     (
314         gamma,
315         vf,
316         "laplacian(" + gamma.name() + ',' + vf.name() + ')'
317     );
321 template<class Type, class GType>
322 tmp<fvMatrix<Type> >
323 laplacian
325     const tmp<GeometricField<GType, fvsPatchField, surfaceMesh> >& tGamma,
326     GeometricField<Type, fvPatchField, volMesh>& vf
329     tmp<fvMatrix<Type> > tfvm(fvm::laplacian(tGamma(), vf));
330     tGamma.clear();
331     return tfvm;
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 } // End namespace fvm
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341 } // End namespace Foam
343 // ************************************************************************* //