limiterBlended.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-2017 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::limiterBlended
26 
27 Description
28  Blends two specified schemes using the limiter function provided by a
29  limitedSurfaceInterpolationScheme.
30 
31  The limited scheme is specified first followed by the scheme to be scaled
32  by the limiter and then the scheme scaled by 1 - limiter e.g.
33 
34  div(phi,U) Gauss limiterBlended vanLeer linear linearUpwind grad(U);
35 
36 SourceFiles
37  limiterBlended.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef limiterBlended_H
42 #define limiterBlended_H
43 
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class limiterBlended Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 template<class Type>
56 class limiterBlended
57 :
58  public surfaceInterpolationScheme<Type>
59 {
60  // Private Member Functions
61 
62  //- Limited scheme providing the limiter
64 
65  //- Scheme 1
67 
68  //- Scheme 2
70 
71 
72  //- Disallow default bitwise copy construct
74 
75  //- Disallow default bitwise assignment
76  void operator=(const limiterBlended&);
77 
78 
79 public:
80 
81  //- Runtime type information
82  TypeName("limiterBlended");
83 
84 
85  // Constructors
86 
87  //- Construct from mesh and Istream.
88  // The name of the flux field is read from the Istream and looked-up
89  // from the mesh objectRegistry
91  (
92  const fvMesh& mesh,
93  Istream& is
94  )
95  :
96  surfaceInterpolationScheme<Type>(mesh),
97  tLimitedScheme_
98  (
99  limitedSurfaceInterpolationScheme<Type>::New(mesh, is)
100  ),
101  tScheme1_
102  (
103  surfaceInterpolationScheme<Type>::New(mesh, is)
104  ),
105  tScheme2_
106  (
107  surfaceInterpolationScheme<Type>::New(mesh, is)
108  )
109  {}
110 
111  //- Construct from mesh, faceFlux and Istream
113  (
114  const fvMesh& mesh,
115  const surfaceScalarField& faceFlux,
116  Istream& is
117  )
118  :
119  surfaceInterpolationScheme<Type>(mesh),
120  tLimitedScheme_
121  (
122  limitedSurfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
123  ),
124  tScheme1_
125  (
126  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
127  ),
128  tScheme2_
129  (
130  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
131  )
132  {}
133 
134 
135  // Member Functions
136 
137  //- Return the interpolation weighting factors
139  (
141  ) const
142  {
143  surfaceScalarField blendingFactor
144  (
145  tLimitedScheme_().limiter(vf)
146  );
147 
148  return
149  blendingFactor*tScheme1_().weights(vf)
150  + (scalar(1) - blendingFactor)*tScheme2_().weights(vf);
151  }
152 
153  //- Return the face-interpolate of the given cell field
154  // with explicit correction
157  {
158  surfaceScalarField blendingFactor
159  (
160  tLimitedScheme_().limiter(vf)
161  );
162 
163  return
164  blendingFactor*tScheme1_().interpolate(vf)
165  + (scalar(1) - blendingFactor)*tScheme2_().interpolate(vf);
166  }
167 
168 
169  //- Return true if this scheme uses an explicit correction
170  virtual bool corrected() const
171  {
172  return tScheme1_().corrected() || tScheme2_().corrected();
173  }
174 
175 
176  //- Return the explicit correction to the face-interpolate
177  // for the given field
180  (
182  ) const
183  {
184  surfaceScalarField blendingFactor
185  (
186  tLimitedScheme_().limiter(vf)
187  );
188 
189  if (tScheme1_().corrected())
190  {
191  if (tScheme2_().corrected())
192  {
193  return
194  (
195  blendingFactor
196  * tScheme1_().correction(vf)
197  + (scalar(1) - blendingFactor)
198  * tScheme2_().correction(vf)
199  );
200  }
201  else
202  {
203  return
204  (
205  blendingFactor
206  * tScheme1_().correction(vf)
207  );
208  }
209  }
210  else if (tScheme2_().corrected())
211  {
212  return
213  (
214  (scalar(1) - blendingFactor)
215  * tScheme2_().correction(vf)
216  );
217  }
218  else
219  {
221  (
222  nullptr
223  );
224  }
225  }
226 };
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace Foam
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #endif
236 
237 // ************************************************************************* //
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
void limiter(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Blends two specified schemes using the limiter function provided by a limitedSurfaceInterpolationSche...
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
Abstract base class for surface interpolation schemes.
Namespace for OpenFOAM.
TypeName("limiterBlended")
Runtime type information.