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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 
26 Class
27  Foam::limiterBlended
28 
29 Group
30  grpFvSurfaceInterpolationSchemes
31 
32 Description
33  Blends two specified schemes using the limiter function provided by a
34  limitedSurfaceInterpolationScheme.
35 
36  The limited scheme is specified first followed by the scheme to be scaled
37  by the limiter and then the scheme scaled by 1 - limiter e.g.
38 
39  div(phi,U) Gauss limiterBlended vanLeer linear linearUpwind grad(U);
40 
41 SourceFiles
42  limiterBlended.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef limiterBlended_H
47 #define limiterBlended_H
48 
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class limiterBlended Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class Type>
61 class limiterBlended
62 :
63  public surfaceInterpolationScheme<Type>
64 {
65  // Private Member Functions
66 
67  //- Limited scheme providing the limiter
69 
70  //- Scheme 1
72 
73  //- Scheme 2
75 
76 
77  //- No copy construct
78  limiterBlended(const limiterBlended&) = delete;
79 
80  //- No copy assignment
81  void operator=(const limiterBlended&) = delete;
82 
83 
84 public:
85 
86  //- Runtime type information
87  TypeName("limiterBlended");
88 
89 
90  // Constructors
91 
92  //- Construct from mesh and Istream.
93  // The name of the flux field is read from the Istream and looked-up
94  // from the mesh objectRegistry
96  (
97  const fvMesh& mesh,
98  Istream& is
99  )
100  :
102  tLimitedScheme_
103  (
105  ),
106  tScheme1_
107  (
109  ),
110  tScheme2_
111  (
113  )
114  {}
115 
116  //- Construct from mesh, faceFlux and Istream
118  (
119  const fvMesh& mesh,
120  const surfaceScalarField& faceFlux,
121  Istream& is
122  )
123  :
125  tLimitedScheme_
126  (
127  limitedSurfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
128  ),
129  tScheme1_
130  (
131  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
132  ),
133  tScheme2_
134  (
135  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
136  )
137  {}
138 
139 
140  // Member Functions
141 
142  //- Return the interpolation weighting factors
144  (
146  ) const
147  {
148  surfaceScalarField blendingFactor
149  (
150  tLimitedScheme_().limiter(vf)
151  );
152 
153  return
154  blendingFactor*tScheme1_().weights(vf)
155  + (scalar(1) - blendingFactor)*tScheme2_().weights(vf);
156  }
157 
158  //- Return the face-interpolate of the given cell field
159  // with explicit correction
162  {
163  surfaceScalarField blendingFactor
164  (
165  tLimitedScheme_().limiter(vf)
166  );
167 
168  return
169  blendingFactor*tScheme1_().interpolate(vf)
170  + (scalar(1) - blendingFactor)*tScheme2_().interpolate(vf);
171  }
172 
173 
174  //- Return true if this scheme uses an explicit correction
175  virtual bool corrected() const
176  {
177  return tScheme1_().corrected() || tScheme2_().corrected();
178  }
179 
180 
181  //- Return the explicit correction to the face-interpolate
182  // for the given field
185  (
187  ) const
188  {
189  surfaceScalarField blendingFactor
190  (
191  tLimitedScheme_().limiter(vf)
192  );
193 
194  if (tScheme1_().corrected())
195  {
196  if (tScheme2_().corrected())
197  {
198  return
199  (
200  blendingFactor
201  * tScheme1_().correction(vf)
202  + (scalar(1) - blendingFactor)
203  * tScheme2_().correction(vf)
204  );
205  }
206  else
207  {
208  return
209  (
210  blendingFactor
211  * tScheme1_().correction(vf)
212  );
213  }
214  }
215  else if (tScheme2_().corrected())
216  {
217  return
218  (
219  (scalar(1) - blendingFactor)
220  * tScheme2_().correction(vf)
221  );
222  }
223  else
224  {
226  (
227  nullptr
228  );
229  }
230  }
231 };
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 } // End namespace Foam
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #endif
241 
242 // ************************************************************************* //
Foam::surfaceInterpolationScheme::New
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Definition: surfaceInterpolationScheme.C:40
Foam::limiterBlended
Blends two specified schemes using the limiter function provided by a limitedSurfaceInterpolationSche...
Definition: limiterBlended.H:60
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::limiterBlended::corrected
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: limiterBlended.H:174
Foam::limiterBlended::weights
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Definition: limiterBlended.H:143
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::limiterBlended::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Definition: limiterBlended.H:160
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::limiterBlended::TypeName
TypeName("limiterBlended")
Runtime type information.
Foam::limitedSurfaceInterpolationScheme
Abstract base class for limited surface interpolation schemes.
Definition: limitedSurfaceInterpolationScheme.H:54
Foam::limiterBlended::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: limiterBlended.H:184
limitedSurfaceInterpolationScheme.H
Foam::surfaceInterpolationScheme
Abstract base class for surface interpolation schemes.
Definition: surfaceInterpolationScheme.H:57
Foam::surfaceInterpolationScheme::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: surfaceInterpolationScheme.H:144
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::limiter
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:37