laplacianScheme.H 8.49 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

Class
    Foam::fv::laplacianScheme

27 28 29
Group
    grpFvLaplacianSchemes

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 68 69 70
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
:
71
    public tmp<laplacianScheme<Type, GType>>::refCount
72 73 74 75 76 77 78
{

protected:

    // Protected data

        const fvMesh& mesh_;
79 80
        tmp<surfaceInterpolationScheme<GType>> tinterpGammaScheme_;
        tmp<snGradScheme<Type>> tsnGradScheme_;
81 82


83 84
private:

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),
126 127
            tinterpGammaScheme_(nullptr),
            tsnGradScheme_(nullptr)
128
        {
129
            tinterpGammaScheme_ = tmp<surfaceInterpolationScheme<GType>>
130 131 132 133
            (
                surfaceInterpolationScheme<GType>::New(mesh, is)
            );

134
            tsnGradScheme_ = tmp<snGradScheme<Type>>
135 136 137 138 139 140 141 142 143
            (
                snGradScheme<Type>::New(mesh, is)
            );
        }

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


    // Selectors

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


164 165
    //- Destructor
    virtual ~laplacianScheme();
166 167 168 169 170 171 172 173 174 175


    // Member Functions

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

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

182
        virtual tmp<fvMatrix<Type>> fvmLaplacian
183 184
        (
            const GeometricField<GType, fvPatchField, volMesh>&,
185
            const GeometricField<Type, fvPatchField, volMesh>&
186 187
        );

188
        virtual tmp<GeometricField<Type, fvPatchField, volMesh>> fvcLaplacian
189 190 191 192
        (
            const GeometricField<Type, fvPatchField, volMesh>&
        ) = 0;

193
        virtual tmp<GeometricField<Type, fvPatchField, volMesh>> fvcLaplacian
194 195 196 197 198
        (
            const GeometricField<GType, fvsPatchField, surfaceMesh>&,
            const GeometricField<Type, fvPatchField, volMesh>&
        ) = 0;

199
        virtual tmp<GeometricField<Type, fvPatchField, volMesh>> fvcLaplacian
200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218
        (
            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

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


#define makeFvLaplacianScheme(SS)                                              \
                                                                               \
makeFvLaplacianTypeScheme(SS, scalar, scalar)                                  \
makeFvLaplacianTypeScheme(SS, symmTensor, scalar)                              \
makeFvLaplacianTypeScheme(SS, tensor, scalar)                                  \
241 242 243 244 245 246 247 248 249 250 251 252
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)
253 254 255 256

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

#ifdef NoRepository
257
    #include "laplacianScheme.C"
258 259 260 261 262 263 264
#endif

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

#endif

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