relaxedNonOrthoGaussLaplacianScheme.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 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
26Class
27 Foam::fv::relaxedNonOrthoGaussLaplacianScheme
28
29Description
30 Basic second-order laplacian using face-gradients and Gauss' theorem.
31
32Usage
33 Minimal example by using \c system/fvSchemes:
34 \verbatim
35 laplacianSchemes
36 {
37 laplacian(<term>) relaxedNonOrthoGauss <other options>;
38 }
39 \endverbatim
40
41 and by using \c system/fvSolution:
42 \verbatim
43 relaxationFactors
44 {
45 equations
46 {
47 <term> <relaxation factor>;
48 }
49 }
50 \endverbatim
51
52SourceFiles
53 relaxedNonOrthoGaussLaplacianScheme.C
54
55\*---------------------------------------------------------------------------*/
56
57#ifndef relaxedNonOrthoGaussLaplacianScheme_H
58#define relaxedNonOrthoGaussLaplacianScheme_H
59
60#include "laplacianScheme.H"
61
62// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63
64namespace Foam
65{
66
67// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68
69namespace fv
70{
71
72/*---------------------------------------------------------------------------*\
73 Class relaxedNonOrthoGaussLaplacianScheme Declaration
74\*---------------------------------------------------------------------------*/
75
76template<class Type, class GType>
78:
79 public fv::laplacianScheme<Type, GType>
80{
81 // Private Member Functions
82
84 (
85 const surfaceVectorField& SfGammaCorr,
87 );
88
89 //- No copy construct
91 (
93 ) = delete;
94
95 //- No copy assignment
96 void operator=(const relaxedNonOrthoGaussLaplacianScheme&) = delete;
97
98
99public:
100
101 //- Runtime type information
102 TypeName("relaxedNonOrthoGauss");
103
104
105 // Constructors
106
107 //- Construct null
109 :
110 laplacianScheme<Type, GType>(mesh)
111 {}
112
113 //- Construct from Istream
115 :
116 laplacianScheme<Type, GType>(mesh, is)
117 {}
118
119 //- Construct from mesh, interpolation and snGradScheme schemes
121 (
122 const fvMesh& mesh,
124 const tmp<snGradScheme<Type>>& sngs
125 )
126 :
127 laplacianScheme<Type, GType>(mesh, igs, sngs)
128 {}
129
130
131 //- Destructor
132 virtual ~relaxedNonOrthoGaussLaplacianScheme() = default;
133
134
135 // Member Functions
136
138 (
139 const surfaceScalarField& gammaMagSf,
140 const surfaceScalarField& deltaCoeffs,
142 );
143
145 (
147 );
148
150 (
153 );
154
156 (
159 );
160};
161
162
163// Use macros to emulate partial-specialisation of the Laplacian functions
164// for scalar diffusivity gamma
166#define defineFvmLaplacianScalarGamma(Type) \
167 \
168template<> \
169tmp<fvMatrix<Type>> \
170relaxedNonOrthoGaussLaplacianScheme<Type, scalar>::fvmLaplacian \
171( \
172 const GeometricField<scalar, fvsPatchField, surfaceMesh>&, \
173 const GeometricField<Type, fvPatchField, volMesh>& \
174); \
175 \
176template<> \
177tmp<GeometricField<Type, fvPatchField, volMesh>> \
178relaxedNonOrthoGaussLaplacianScheme<Type, scalar>::fvcLaplacian \
179( \
180 const GeometricField<scalar, fvsPatchField, surfaceMesh>&, \
181 const GeometricField<Type, fvPatchField, volMesh>& \
182);
183
190
191
192// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193
194} // End namespace fv
195
196// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197
198} // End namespace Foam
199
200// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201
202#ifdef NoRepository
204#endif
205
206// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207
208#endif
209
210// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for laplacian schemes.
const fvMesh & mesh() const
Return mesh reference.
Basic second-order laplacian using face-gradients and Gauss' theorem.
relaxedNonOrthoGaussLaplacianScheme(const fvMesh &mesh, const tmp< surfaceInterpolationScheme< GType > > &igs, const tmp< snGradScheme< Type > > &sngs)
Construct from mesh, interpolation and snGradScheme schemes.
TypeName("relaxedNonOrthoGauss")
Runtime type information.
relaxedNonOrthoGaussLaplacianScheme(const fvMesh &mesh, Istream &is)
Construct from Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcLaplacian(const GeometricField< Type, fvPatchField, volMesh > &)
relaxedNonOrthoGaussLaplacianScheme(const fvMesh &mesh)
Construct null.
tmp< fvMatrix< Type > > fvmLaplacian(const GeometricField< GType, fvsPatchField, surfaceMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
static tmp< fvMatrix< Type > > fvmLaplacianUncorrected(const surfaceScalarField &gammaMagSf, const surfaceScalarField &deltaCoeffs, const GeometricField< Type, fvPatchField, volMesh > &)
virtual ~relaxedNonOrthoGaussLaplacianScheme()=default
Destructor.
Abstract base class for runtime selected snGrad surface normal gradient schemes.
Definition: snGradScheme.H:79
Abstract base class for surface interpolation schemes.
A class for managing temporary objects.
Definition: tmp.H:65
#define defineFvmLaplacianScalarGamma(Type)
Namespace for OpenFOAM.
labelList fv(nPoints)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73