FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

gaussLaplacianSchemes.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023     
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "gaussLaplacianScheme.H"
00027 #include <finiteVolume/fvMesh.H>
00028 
00029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace fv
00034 {
00035     makeFvLaplacianScheme(gaussLaplacianScheme)
00036 }
00037 }
00038 
00039 #define declareFvmLaplacianScalarGamma(Type)                                 \
00040                                                                              \
00041 template<>                                                                   \
00042 Foam::tmp<Foam::fvMatrix<Foam::Type> >                                       \
00043 Foam::fv::gaussLaplacianScheme<Foam::Type, Foam::scalar>::fvmLaplacian       \
00044 (                                                                            \
00045     const GeometricField<scalar, fvsPatchField, surfaceMesh>& gamma,         \
00046     GeometricField<Type, fvPatchField, volMesh>& vf                          \
00047 )                                                                            \
00048 {                                                                            \
00049     const fvMesh& mesh = this->mesh();                                       \
00050                                                                              \
00051     GeometricField<scalar, fvsPatchField, surfaceMesh> gammaMagSf            \
00052     (                                                                        \
00053         gamma*mesh.magSf()                                                   \
00054     );                                                                       \
00055                                                                              \
00056     tmp<fvMatrix<Type> > tfvm = fvmLaplacianUncorrected(gammaMagSf, vf);     \
00057     fvMatrix<Type>& fvm = tfvm();                                            \
00058                                                                              \
00059     if (this->tsnGradScheme_().corrected())                                  \
00060     {                                                                        \
00061         if (mesh.fluxRequired(vf.name()))                                    \
00062         {                                                                    \
00063             fvm.faceFluxCorrectionPtr() = new                                \
00064             GeometricField<Type, fvsPatchField, surfaceMesh>                 \
00065             (                                                                \
00066                 gammaMagSf*this->tsnGradScheme_().correction(vf)             \
00067             );                                                               \
00068                                                                              \
00069             fvm.source() -=                                                  \
00070                 mesh.V()*                                                    \
00071                 fvc::div                                                     \
00072                 (                                                            \
00073                     *fvm.faceFluxCorrectionPtr()                             \
00074                 )().internalField();                                         \
00075         }                                                                    \
00076         else                                                                 \
00077         {                                                                    \
00078             fvm.source() -=                                                  \
00079                 mesh.V()*                                                    \
00080                 fvc::div                                                     \
00081                 (                                                            \
00082                     gammaMagSf*this->tsnGradScheme_().correction(vf)         \
00083                 )().internalField();                                         \
00084         }                                                                    \
00085     }                                                                        \
00086                                                                              \
00087     return tfvm;                                                             \
00088 }                                                                            \
00089                                                                              \
00090                                                                              \
00091 template<>                                                                   \
00092 Foam::tmp<Foam::GeometricField<Foam::Type, Foam::fvPatchField, Foam::volMesh> >\
00093 Foam::fv::gaussLaplacianScheme<Foam::Type, Foam::scalar>::fvcLaplacian       \
00094 (                                                                            \
00095     const GeometricField<scalar, fvsPatchField, surfaceMesh>& gamma,         \
00096     const GeometricField<Type, fvPatchField, volMesh>& vf                    \
00097 )                                                                            \
00098 {                                                                            \
00099     const fvMesh& mesh = this->mesh();                                       \
00100                                                                              \
00101     tmp<GeometricField<Type, fvPatchField, volMesh> > tLaplacian             \
00102     (                                                                        \
00103         fvc::div(gamma*this->tsnGradScheme_().snGrad(vf)*mesh.magSf())       \
00104     );                                                                       \
00105                                                                              \
00106     tLaplacian().rename("laplacian(" + gamma.name() + ',' + vf.name() + ')');\
00107                                                                              \
00108     return tLaplacian;                                                       \
00109 }
00110 
00111 
00112 declareFvmLaplacianScalarGamma(scalar);
00113 declareFvmLaplacianScalarGamma(vector);
00114 declareFvmLaplacianScalarGamma(sphericalTensor);
00115 declareFvmLaplacianScalarGamma(symmTensor);
00116 declareFvmLaplacianScalarGamma(tensor);
00117 
00118 
00119 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines