All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
LimitedScheme.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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::LimitedScheme
26 
27 Description
28  Class to create NVD/TVD limited weighting-factors.
29 
30  The particular differencing scheme class is supplied as a template
31  argument, the weight function of which is called by the weight function
32  of this class for the internal faces as well as faces of coupled
33  patches (e.g. processor-processor patches). The weight function is
34  supplied the central-differencing weighting factor, the face-flux, the
35  cell and face gradients (from which the normalised variable
36  distribution may be created) and the cell centre distance.
37 
38  This code organisation is both neat and efficient, allowing for
39  convenient implementation of new schemes to run on parallelised cases.
40 
41 SourceFiles
42  LimitedScheme.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef LimitedScheme_H
47 #define LimitedScheme_H
48 
50 #include "LimitFuncs.H"
51 #include "NVDTVD.H"
52 #include "NVDVTVDV.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 /*---------------------------------------------------------------------------*\
60  Class LimitedScheme Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 template<class Type, class Limiter, template<class> class LimitFunc>
64 class LimitedScheme
65 :
67  public Limiter
68 {
69  // Private Member Functions
70 
71  //- Calculate the limiter
72  void calcLimiter
73  (
75  surfaceScalarField& limiterField
76  ) const;
77 
78 
79 public:
80 
81  //- Runtime type information
82  TypeName("LimitedScheme");
83 
84  typedef Limiter LimiterType;
85 
86  // Constructors
87 
88  //- Construct from mesh and faceFlux and limiter scheme
90  (
91  const fvMesh& mesh,
92  const surfaceScalarField& faceFlux,
93  const Limiter& weight
94  )
95  :
96  limitedSurfaceInterpolationScheme<Type>(mesh, faceFlux),
97  Limiter(weight)
98  {}
99 
100  //- Construct from mesh and Istream.
101  // The name of the flux field is read from the Istream and looked-up
102  // from the mesh objectRegistry
104  (
105  const fvMesh& mesh,
106  Istream& is
107  )
108  :
109  limitedSurfaceInterpolationScheme<Type>(mesh, is),
110  Limiter(is)
111  {}
112 
113  //- Construct from mesh, faceFlux and Istream
115  (
116  const fvMesh& mesh,
117  const surfaceScalarField& faceFlux,
118  Istream& is
119  )
120  :
121  limitedSurfaceInterpolationScheme<Type>(mesh, faceFlux),
122  Limiter(is)
123  {}
124 
125  //- Disallow default bitwise copy construction
126  LimitedScheme(const LimitedScheme&) = delete;
127 
128 
129  // Member Functions
130 
131  //- Return the interpolation weighting factors
133  (
135  ) const;
136 
137 
138  // Member Operators
139 
140  //- Disallow default bitwise assignment
141  void operator=(const LimitedScheme&) = delete;
142 };
143 
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 } // End namespace Foam
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 // Add the patch constructor functions to the hash tables
153 #define makeLimitedSurfaceInterpolationTypeScheme\
154 ( \
155  SS, \
156  LIMITER, \
157  NVDTVD, \
158  LIMFUNC, \
159  TYPE \
160 ) \
161  \
162 typedef LimitedScheme<TYPE, LIMITER<NVDTVD>, limitFuncs::LIMFUNC> \
163  LimitedScheme##TYPE##LIMITER##NVDTVD##LIMFUNC##_; \
164 defineTemplateTypeNameAndDebugWithName \
165  (LimitedScheme##TYPE##LIMITER##NVDTVD##LIMFUNC##_, #SS, 0); \
166  \
167 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
168 <LimitedScheme<TYPE, LIMITER<NVDTVD>, limitFuncs::LIMFUNC>> \
169  add##SS##LIMFUNC##TYPE##MeshConstructorToTable_; \
170  \
171 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
172 <LimitedScheme<TYPE, LIMITER<NVDTVD>, limitFuncs::LIMFUNC>> \
173  add##SS##LIMFUNC##TYPE##MeshFluxConstructorToTable_; \
174  \
175 limitedSurfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
176 <LimitedScheme<TYPE, LIMITER<NVDTVD>, limitFuncs::LIMFUNC>> \
177  add##SS##LIMFUNC##TYPE##MeshConstructorToLimitedTable_; \
178  \
179 limitedSurfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
180 <LimitedScheme<TYPE, LIMITER<NVDTVD>, limitFuncs::LIMFUNC>> \
181  add##SS##LIMFUNC##TYPE##MeshFluxConstructorToLimitedTable_;
182 
184 #define makeLimitedSurfaceInterpolationScheme(SS, LIMITER) \
185  \
186 makeLimitedSurfaceInterpolationTypeScheme(SS,LIMITER,NVDTVD,magSqr,scalar) \
187 makeLimitedSurfaceInterpolationTypeScheme(SS,LIMITER,NVDTVD,magSqr,vector) \
188 makeLimitedSurfaceInterpolationTypeScheme \
189 ( \
190  SS, \
191  LIMITER, \
192  NVDTVD, \
193  magSqr, \
194  sphericalTensor \
195 ) \
196 makeLimitedSurfaceInterpolationTypeScheme(SS,LIMITER,NVDTVD,magSqr,symmTensor)\
197 makeLimitedSurfaceInterpolationTypeScheme(SS,LIMITER,NVDTVD,magSqr,tensor)
198 
200 #define makeLimitedVSurfaceInterpolationScheme(SS, LIMITER) \
201 makeLimitedSurfaceInterpolationTypeScheme(SS,LIMITER,NVDVTVDV,null,vector)
202 
204 #define makeLLimitedSurfaceInterpolationTypeScheme\
205 ( \
206  SS, \
207  LLIMITER, \
208  LIMITER, \
209  NVDTVD, \
210  LIMFUNC, \
211  TYPE \
212 ) \
213  \
214 typedef LimitedScheme<TYPE, LLIMITER<LIMITER<NVDTVD>>, limitFuncs::LIMFUNC> \
215  LimitedScheme##TYPE##LLIMITER##LIMITER##NVDTVD##LIMFUNC##_; \
216 defineTemplateTypeNameAndDebugWithName \
217  (LimitedScheme##TYPE##LLIMITER##LIMITER##NVDTVD##LIMFUNC##_, #SS, 0); \
218  \
219 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
220 <LimitedScheme<TYPE, LLIMITER<LIMITER<NVDTVD>>, limitFuncs::LIMFUNC>> \
221  add##SS##LIMFUNC##TYPE##MeshConstructorToTable_; \
222  \
223 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
224 <LimitedScheme<TYPE, LLIMITER<LIMITER<NVDTVD>>, limitFuncs::LIMFUNC>> \
225  add##SS##LIMFUNC##TYPE##MeshFluxConstructorToTable_; \
226  \
227 limitedSurfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
228 <LimitedScheme<TYPE, LLIMITER<LIMITER<NVDTVD>>, limitFuncs::LIMFUNC>> \
229  add##SS##LIMFUNC##TYPE##MeshConstructorToLimitedTable_; \
230  \
231 limitedSurfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
232 <LimitedScheme<TYPE, LLIMITER<LIMITER<NVDTVD>>, limitFuncs::LIMFUNC>> \
233  add##SS##LIMFUNC##TYPE##MeshFluxConstructorToLimitedTable_;
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #ifdef NoRepository
239  #include "LimitedScheme.C"
240 #endif
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 #endif
245 
246 // ************************************************************************* //
LimitedScheme(const fvMesh &mesh, const surfaceScalarField &faceFlux, const Limiter &weight)
Construct from mesh and faceFlux and limiter scheme.
Definition: LimitedScheme.H:89
Class to create NVD/TVD limited weighting-factors.
Definition: LimitedScheme.H:63
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Generic GeometricField class.
void operator=(const LimitedScheme &)=delete
Disallow default bitwise assignment.
TypeName("LimitedScheme")
Runtime type information.
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Abstract base class for limited surface interpolation schemes.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.