laplacianScheme.H 8.45 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
Henry Weller's avatar
Henry Weller committed
5
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
6 7 8 9 10
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

11 12 13 14
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
15 16 17 18 19 20 21

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
22
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67

Class
    Foam::fv::laplacianScheme

Description
    Abstract base class for laplacian schemes.

SourceFiles
    laplacianScheme.C

\*---------------------------------------------------------------------------*/

#ifndef laplacianScheme_H
#define laplacianScheme_H

#include "tmp.H"
#include "volFieldsFwd.H"
#include "surfaceFieldsFwd.H"
#include "linear.H"
#include "correctedSnGrad.H"
#include "typeInfo.H"
#include "runTimeSelectionTables.H"

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace Foam
{

template<class Type>
class fvMatrix;

class fvMesh;

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace fv
{

/*---------------------------------------------------------------------------*\
                           Class laplacianScheme Declaration
\*---------------------------------------------------------------------------*/

template<class Type, class GType>
class laplacianScheme
:
68
    public tmp<laplacianScheme<Type, GType>>::refCount
69 70 71 72 73 74 75
{

protected:

    // Protected data

        const fvMesh& mesh_;
76 77
        tmp<surfaceInterpolationScheme<GType>> tinterpGammaScheme_;
        tmp<snGradScheme<Type>> tsnGradScheme_;
78 79


80 81
private:

82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125
    // Private Member Functions

        //- Disallow copy construct
        laplacianScheme(const laplacianScheme&);

        //- Disallow default bitwise assignment
        void operator=(const laplacianScheme&);


public:

    //- Runtime type information
    virtual const word& type() const = 0;


    // Declare run-time constructor selection tables

        declareRunTimeSelectionTable
        (
            tmp,
            laplacianScheme,
            Istream,
            (const fvMesh& mesh, Istream& schemeData),
            (mesh, schemeData)
        );


    // Constructors

        //- Construct from mesh
        laplacianScheme(const fvMesh& mesh)
        :
            mesh_(mesh),
            tinterpGammaScheme_(new linear<GType>(mesh)),
            tsnGradScheme_(new correctedSnGrad<Type>(mesh))
        {}

        //- Construct from mesh and Istream
        laplacianScheme(const fvMesh& mesh, Istream& is)
        :
            mesh_(mesh),
            tinterpGammaScheme_(NULL),
            tsnGradScheme_(NULL)
        {
126
            tinterpGammaScheme_ = tmp<surfaceInterpolationScheme<GType>>
127 128 129 130
            (
                surfaceInterpolationScheme<GType>::New(mesh, is)
            );

131
            tsnGradScheme_ = tmp<snGradScheme<Type>>
132 133 134 135 136 137 138 139 140
            (
                snGradScheme<Type>::New(mesh, is)
            );
        }

        //- Construct from mesh, interpolation and snGradScheme schemes
        laplacianScheme
        (
            const fvMesh& mesh,
141 142
            const tmp<surfaceInterpolationScheme<GType>>& igs,
            const tmp<snGradScheme<Type>>& sngs
143 144 145 146 147 148 149 150 151 152 153
        )
        :
            mesh_(mesh),
            tinterpGammaScheme_(igs),
            tsnGradScheme_(sngs)
        {}


    // Selectors

        //- Return a pointer to a new laplacianScheme created on freestore
154
        static tmp<laplacianScheme<Type, GType>> New
155 156 157 158 159 160
        (
            const fvMesh& mesh,
            Istream& schemeData
        );


161 162
    //- Destructor
    virtual ~laplacianScheme();
163 164 165 166 167 168 169 170 171 172


    // Member Functions

        //- Return mesh reference
        const fvMesh& mesh() const
        {
            return mesh_;
        }

173
        virtual tmp<fvMatrix<Type>> fvmLaplacian
174 175
        (
            const GeometricField<GType, fvsPatchField, surfaceMesh>&,
176
            const GeometricField<Type, fvPatchField, volMesh>&
177 178
        ) = 0;

179
        virtual tmp<fvMatrix<Type>> fvmLaplacian
180 181
        (
            const GeometricField<GType, fvPatchField, volMesh>&,
182
            const GeometricField<Type, fvPatchField, volMesh>&
183 184
        );

185
        virtual tmp<GeometricField<Type, fvPatchField, volMesh>> fvcLaplacian
186 187 188 189
        (
            const GeometricField<Type, fvPatchField, volMesh>&
        ) = 0;

190
        virtual tmp<GeometricField<Type, fvPatchField, volMesh>> fvcLaplacian
191 192 193 194 195
        (
            const GeometricField<GType, fvsPatchField, surfaceMesh>&,
            const GeometricField<Type, fvPatchField, volMesh>&
        ) = 0;

196
        virtual tmp<GeometricField<Type, fvPatchField, volMesh>> fvcLaplacian
197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215
        (
            const GeometricField<GType, fvPatchField, volMesh>&,
            const GeometricField<Type, fvPatchField, volMesh>&
        );
};


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

} // End namespace fv

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

} // End namespace Foam

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

// Add the patch constructor functions to the hash tables

216
#define makeFvLaplacianTypeScheme(SS, GType, Type)                             \
217
    typedef Foam::fv::SS<Foam::Type, Foam::GType> SS##Type##GType;             \
218 219
    defineNamedTemplateTypeNameAndDebug(SS##Type##GType, 0);                   \
                                                                               \
220 221 222 223 224 225 226
    namespace Foam                                                             \
    {                                                                          \
        namespace fv                                                           \
        {                                                                      \
            typedef SS<Type, GType> SS##Type##GType;                           \
                                                                               \
            laplacianScheme<Type, GType>::                                     \
227
                addIstreamConstructorToTable<SS<Type, GType>>                  \
228 229 230
                add##SS##Type##GType##IstreamConstructorToTable_;              \
        }                                                                      \
    }
231 232 233 234 235 236 237


#define makeFvLaplacianScheme(SS)                                              \
                                                                               \
makeFvLaplacianTypeScheme(SS, scalar, scalar)                                  \
makeFvLaplacianTypeScheme(SS, symmTensor, scalar)                              \
makeFvLaplacianTypeScheme(SS, tensor, scalar)                                  \
238 239 240 241 242 243 244 245 246 247 248 249
makeFvLaplacianTypeScheme(SS, scalar, vector)                                  \
makeFvLaplacianTypeScheme(SS, symmTensor, vector)                              \
makeFvLaplacianTypeScheme(SS, tensor, vector)                                  \
makeFvLaplacianTypeScheme(SS, scalar, sphericalTensor)                         \
makeFvLaplacianTypeScheme(SS, symmTensor, sphericalTensor)                     \
makeFvLaplacianTypeScheme(SS, tensor, sphericalTensor)                         \
makeFvLaplacianTypeScheme(SS, scalar, symmTensor)                              \
makeFvLaplacianTypeScheme(SS, symmTensor, symmTensor)                          \
makeFvLaplacianTypeScheme(SS, tensor, symmTensor)                              \
makeFvLaplacianTypeScheme(SS, scalar, tensor)                                  \
makeFvLaplacianTypeScheme(SS, symmTensor, tensor)                              \
makeFvLaplacianTypeScheme(SS, tensor, tensor)
250 251 252 253 254 255 256 257 258 259 260 261

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#ifdef NoRepository
#   include "laplacianScheme.C"
#endif

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#endif

// ************************************************************************* //