initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / laplacianSchemes / gaussLaplacianScheme / gaussLaplacianSchemes.C
blob5dab5ceef08bac4f124e260c26e50c55c9d6d8f9
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
24     
25 \*---------------------------------------------------------------------------*/
27 #include "gaussLaplacianScheme.H"
28 #include "fvMesh.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace fv
36     makeFvLaplacianScheme(gaussLaplacianScheme)
40 #define declareFvmLaplacianScalarGamma(Type)                                 \
41                                                                              \
42 template<>                                                                   \
43 Foam::tmp<Foam::fvMatrix<Foam::Type> >                                       \
44 Foam::fv::gaussLaplacianScheme<Foam::Type, Foam::scalar>::fvmLaplacian       \
45 (                                                                            \
46     const GeometricField<scalar, fvsPatchField, surfaceMesh>& gamma,         \
47     GeometricField<Type, fvPatchField, volMesh>& vf                          \
48 )                                                                            \
49 {                                                                            \
50     const fvMesh& mesh = this->mesh();                                       \
51                                                                              \
52     GeometricField<scalar, fvsPatchField, surfaceMesh> gammaMagSf            \
53     (                                                                        \
54         gamma*mesh.magSf()                                                   \
55     );                                                                       \
56                                                                              \
57     tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected(gammaMagSf, vf);     \
58     fvMatrix<Type>& fvm = tfvm();                                            \
59                                                                              \
60     if (this->tsnGradScheme_().corrected())                                  \
61     {                                                                        \
62         if (mesh.fluxRequired(vf.name()))                                    \
63         {                                                                    \
64             fvm.faceFluxCorrectionPtr() = new                                \
65             GeometricField<Type, fvsPatchField, surfaceMesh>                 \
66             (                                                                \
67                 gammaMagSf*this->tsnGradScheme_().correction(vf)             \
68             );                                                               \
69                                                                              \
70             fvm.source() -=                                                  \
71                 mesh.V()*                                                    \
72                 fvc::div                                                     \
73                 (                                                            \
74                     *fvm.faceFluxCorrectionPtr()                             \
75                 )().internalField();                                         \
76         }                                                                    \
77         else                                                                 \
78         {                                                                    \
79             fvm.source() -=                                                  \
80                 mesh.V()*                                                    \
81                 fvc::div                                                     \
82                 (                                                            \
83                     gammaMagSf*this->tsnGradScheme_().correction(vf)         \
84                 )().internalField();                                         \
85         }                                                                    \
86     }                                                                        \
87                                                                              \
88     return tfvm;                                                             \
89 }                                                                            \
90                                                                              \
91                                                                              \
92 template<>                                                                   \
93 Foam::tmp<Foam::GeometricField<Foam::Type, Foam::fvPatchField, Foam::volMesh> >\
94 Foam::fv::gaussLaplacianScheme<Foam::Type, Foam::scalar>::fvcLaplacian       \
95 (                                                                            \
96     const GeometricField<scalar, fvsPatchField, surfaceMesh>& gamma,         \
97     const GeometricField<Type, fvPatchField, volMesh>& vf                    \
98 )                                                                            \
99 {                                                                            \
100     const fvMesh& mesh = this->mesh();                                       \
101                                                                              \
102     tmp<GeometricField<Type, fvPatchField, volMesh> > tLaplacian             \
103     (                                                                        \
104         fvc::div(gamma*this->tsnGradScheme_().snGrad(vf)*mesh.magSf())       \
105     );                                                                       \
106                                                                              \
107     tLaplacian().rename("laplacian(" + gamma.name() + ',' + vf.name() + ')');\
108                                                                              \
109     return tLaplacian;                                                       \
113 declareFvmLaplacianScalarGamma(scalar);
114 declareFvmLaplacianScalarGamma(vector);
115 declareFvmLaplacianScalarGamma(sphericalTensor);
116 declareFvmLaplacianScalarGamma(symmTensor);
117 declareFvmLaplacianScalarGamma(tensor);
120 // ************************************************************************* //