gaussLaplacianScheme.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  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 Class
25  Foam::fv::gaussLaplacianScheme
26 
27 Description
28  Basic second-order laplacian using face-gradients and Gauss' theorem.
29 
30 SourceFiles
31  gaussLaplacianScheme.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef gaussLaplacianScheme_H
36 #define gaussLaplacianScheme_H
37 
38 #include "laplacianScheme.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace fv
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class gaussLaplacianScheme Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 template<class Type, class GType>
56 :
57  public fv::laplacianScheme<Type, GType>
58 {
59  // Private Member Functions
60 
62  (
63  const surfaceVectorField& SfGammaCorr,
65  );
66 
67  //- Disallow default bitwise copy construct
69 
70  //- Disallow default bitwise assignment
71  void operator=(const gaussLaplacianScheme&);
72 
73 
74 public:
75 
76  //- Runtime type information
77  TypeName("Gauss");
78 
79 
80  // Constructors
81 
82  //- Construct null
84  :
85  laplacianScheme<Type, GType>(mesh)
86  {}
87 
88  //- Construct from Istream
90  :
91  laplacianScheme<Type, GType>(mesh, is)
92  {}
93 
94  //- Construct from mesh, interpolation and snGradScheme schemes
96  (
97  const fvMesh& mesh,
99  const tmp<snGradScheme<Type>>& sngs
100  )
101  :
102  laplacianScheme<Type, GType>(mesh, igs, sngs)
103  {}
104 
105 
106  //- Destructor
107  virtual ~gaussLaplacianScheme()
108  {}
109 
110 
111  // Member Functions
112 
114  (
115  const surfaceScalarField& gammaMagSf,
116  const surfaceScalarField& deltaCoeffs,
118  );
119 
121  (
123  );
124 
126  (
129  );
130 
132  (
135  );
136 };
137 
138 
139 // Use macros to emulate partial-specialisation of the the Laplacian functions
140 // for scalar diffusivity gamma
142 #define defineFvmLaplacianScalarGamma(Type) \
143  \
144 template<> \
145 tmp<fvMatrix<Type>> gaussLaplacianScheme<Type, scalar>::fvmLaplacian \
146 ( \
147  const GeometricField<scalar, fvsPatchField, surfaceMesh>&, \
148  const GeometricField<Type, fvPatchField, volMesh>& \
149 ); \
150  \
151 template<> \
152 tmp<GeometricField<Type, fvPatchField, volMesh>> \
153 gaussLaplacianScheme<Type, scalar>::fvcLaplacian \
154 ( \
155  const GeometricField<scalar, fvsPatchField, surfaceMesh>&, \
156  const GeometricField<Type, fvPatchField, volMesh>& \
157 );
158 
159 
165 
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace fv
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 } // End namespace Foam
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #ifdef NoRepository
178  #include "gaussLaplacianScheme.C"
179 #endif
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 #endif
184 
185 // ************************************************************************* //
static tmp< fvMatrix< Type > > fvmLaplacianUncorrected(const surfaceScalarField &gammaMagSf, const surfaceScalarField &deltaCoeffs, const GeometricField< Type, fvPatchField, volMesh > &)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Abstract base class for laplacian schemes.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual ~gaussLaplacianScheme()
Destructor.
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
labelList fv(nPoints)
Basic second-order laplacian using face-gradients and Gauss&#39; theorem.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcLaplacian(const GeometricField< Type, fvPatchField, volMesh > &)
Abstract base class for snGrad schemes.
Definition: snGradScheme.H:60
TypeName("Gauss")
Runtime type information.
#define defineFvmLaplacianScalarGamma(Type)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const fvMesh & mesh() const
Return mesh reference.
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars.
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< fvMatrix< Type > > fvmLaplacian(const GeometricField< GType, fvsPatchField, surfaceMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Namespace for OpenFOAM.