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

13 14 15 16
    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.
17 18 19 20 21 22 23

    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
24
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
25 26 27 28

Class
    Foam::fv::laplacianScheme

29 30 31
Group
    grpFvLaplacianSchemes

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 71 72
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
:
73
    public refCount
74 75 76 77 78 79 80
{

protected:

    // Protected data

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


85 86
private:

87 88
    // Private Member Functions

89 90
        //- No copy construct
        laplacianScheme(const laplacianScheme&) = delete;
91

92 93
        //- No copy assignment
        void operator=(const laplacianScheme&) = delete;
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 126 127


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

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

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


    // Selectors

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


166
    //- Destructor
167
    virtual ~laplacianScheme() = default;
168 169 170 171 172 173 174 175 176 177


    // Member Functions

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

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

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

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

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

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

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


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

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

#ifdef NoRepository
259
    #include "laplacianScheme.C"
260 261 262 263 264 265 266
#endif

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

#endif

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