limiterBlended.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-2020 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 public:
73 
74  //- Runtime type information
75  TypeName("limiterBlended");
76 
77 
78  // Constructors
79 
80  //- Construct from mesh and Istream.
81  // The name of the flux field is read from the Istream and looked-up
82  // from the mesh objectRegistry
84  (
85  const fvMesh& mesh,
86  Istream& is
87  )
88  :
89  surfaceInterpolationScheme<Type>(mesh),
90  tLimitedScheme_
91  (
92  limitedSurfaceInterpolationScheme<Type>::New(mesh, is)
93  ),
94  tScheme1_
95  (
96  surfaceInterpolationScheme<Type>::New(mesh, is)
97  ),
98  tScheme2_
99  (
100  surfaceInterpolationScheme<Type>::New(mesh, is)
101  )
102  {}
103 
104  //- Construct from mesh, faceFlux and Istream
106  (
107  const fvMesh& mesh,
108  const surfaceScalarField& faceFlux,
109  Istream& is
110  )
111  :
112  surfaceInterpolationScheme<Type>(mesh),
113  tLimitedScheme_
114  (
115  limitedSurfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
116  ),
117  tScheme1_
118  (
119  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
120  ),
121  tScheme2_
122  (
123  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
124  )
125  {}
126 
127  //- Disallow default bitwise copy construction
128  limiterBlended(const limiterBlended&) = delete;
129 
130 
131  // Member Functions
132 
133  //- Return the interpolation weighting factors
135  (
137  ) const
138  {
139  surfaceScalarField blendingFactor
140  (
141  tLimitedScheme_().limiter(vf)
142  );
143 
144  return
145  blendingFactor*tScheme1_().weights(vf)
146  + (scalar(1) - blendingFactor)*tScheme2_().weights(vf);
147  }
148 
149  //- Return the face-interpolate of the given cell field
150  // with explicit correction
153  {
154  surfaceScalarField blendingFactor
155  (
156  tLimitedScheme_().limiter(vf)
157  );
158 
159  return
160  blendingFactor*tScheme1_().interpolate(vf)
161  + (scalar(1) - blendingFactor)*tScheme2_().interpolate(vf);
162  }
163 
164 
165  //- Return true if this scheme uses an explicit correction
166  virtual bool corrected() const
167  {
168  return tScheme1_().corrected() || tScheme2_().corrected();
169  }
170 
171 
172  //- Return the explicit correction to the face-interpolate
173  // for the given field
176  (
178  ) const
179  {
180  surfaceScalarField blendingFactor
181  (
182  tLimitedScheme_().limiter(vf)
183  );
184 
185  if (tScheme1_().corrected())
186  {
187  if (tScheme2_().corrected())
188  {
189  return
190  (
191  blendingFactor
192  * tScheme1_().correction(vf)
193  + (scalar(1) - blendingFactor)
194  * tScheme2_().correction(vf)
195  );
196  }
197  else
198  {
199  return
200  (
201  blendingFactor
202  * tScheme1_().correction(vf)
203  );
204  }
205  }
206  else if (tScheme2_().corrected())
207  {
208  return
209  (
210  (scalar(1) - blendingFactor)
211  * tScheme2_().correction(vf)
212  );
213  }
214  else
215  {
217  (
218  nullptr
219  );
220  }
221  }
222 
223 
224  // Member Operators
225 
226  //- Disallow default bitwise assignment
227  void operator=(const limiterBlended&) = delete;
228 };
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace Foam
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 #endif
238 
239 // ************************************************************************* //
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
void operator=(const limiterBlended &)=delete
Disallow default bitwise assignment.
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.
limiterBlended(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
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.
void limiter(surfaceScalarField &lambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
Blends two specified schemes using the limiter function provided by a limitedSurfaceInterpolationSche...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
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.