CoBlended.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) 2012-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::CoBlended
28 
29 Group
30  grpFvSurfaceInterpolationSchemes
31 
32 Description
33  Two-scheme Courant number based blending differencing scheme.
34 
35  Similar to localBlended but uses a blending factor computed from the
36  face-based Courant number and the lower and upper Courant number limits
37  supplied:
38  \f[
39  weight = 1 - clamp((Co - Co1)/(Co2 - Co1), zero_one{})
40  \f]
41  where
42  \vartable
43  Co1 | Courant number below which scheme1 is used
44  Co2 | Courant number above which scheme2 is used
45  \endvartable
46 
47  The weight applies to the first scheme and 1-weight to the second scheme.
48 
49 Usage
50  Example of the CoBlended scheme specification using LUST for Courant numbers
51  less than 1 and linearUpwind for Courant numbers greater than 10:
52  \verbatim
53  divSchemes
54  {
55  .
56  .
57  div(phi,U) Gauss CoBlended 1 LUST grad(U) 10 linearUpwind grad(U);
58  .
59  .
60  }
61  \endverbatim
62 
63 SourceFiles
64  CoBlended.C
65 
66 \*---------------------------------------------------------------------------*/
67 
68 #ifndef CoBlended_H
69 #define CoBlended_H
70 
72 #include "blendedSchemeBase.H"
73 #include "surfaceInterpolate.H"
74 
75 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 
77 namespace Foam
78 {
79 
80 /*---------------------------------------------------------------------------*\
81  Class CoBlended Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 template<class Type>
85 class CoBlended
86 :
87  public surfaceInterpolationScheme<Type>,
88  public blendedSchemeBase<Type>
89 {
90  // Private data
91 
92  //- Courant number below which scheme1 is used
93  const scalar Co1_;
94 
95  //- Scheme 1
97 
98  //- Courant number above which scheme2 is used
99  const scalar Co2_;
100 
101  //- Scheme 2
103 
104  //- The face-flux used to compute the face Courant number
105  const surfaceScalarField& faceFlux_;
106 
107 
108  // Private Member Functions
109 
110  //- No copy construct
111  CoBlended(const CoBlended&) = delete;
112 
113  //- No copy assignment
114  void operator=(const CoBlended&) = delete;
115 
116 
117 public:
118 
119  //- Runtime type information
120  TypeName("CoBlended");
121 
122 
123  // Constructors
124 
125  //- Construct from mesh and Istream.
126  // The name of the flux field is read from the Istream and looked-up
127  // from the mesh objectRegistry
128  CoBlended
129  (
130  const fvMesh& mesh,
131  Istream& is
132  )
133  :
135  Co1_(readScalar(is)),
136  tScheme1_
137  (
139  ),
140  Co2_(readScalar(is)),
141  tScheme2_
142  (
144  ),
145  faceFlux_
146  (
148  )
149  {
150  if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
151  {
153  << "coefficients = " << Co1_ << " and " << Co2_
154  << " should be > 0 and Co2 > Co1"
155  << exit(FatalIOError);
156  }
157  }
158 
159 
160  //- Construct from mesh, faceFlux and Istream
161  CoBlended
162  (
163  const fvMesh& mesh,
164  const surfaceScalarField& faceFlux,
165  Istream& is
166  )
167  :
168  surfaceInterpolationScheme<Type>(mesh),
169  Co1_(readScalar(is)),
170  tScheme1_
171  (
172  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
173  ),
174  Co2_(readScalar(is)),
175  tScheme2_
176  (
177  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
178  ),
179  faceFlux_(faceFlux)
180  {
181  if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
182  {
184  << "coefficients = " << Co1_ << " and " << Co2_
185  << " should be > 0 and Co2 > Co1"
187  }
188  }
189 
190 
191  // Member Functions
192 
193  //- Return the face-based blending factor
194  virtual tmp<surfaceScalarField> blendingFactor
195  (
196  const GeometricField<Type, fvPatchField, volMesh>& vf
197  ) const
198  {
199  const fvMesh& mesh = this->mesh();
200  tmp<surfaceScalarField> tUflux = faceFlux_;
201 
202  if (faceFlux_.dimensions() == dimMass/dimTime)
203  {
204  // Currently assume that the density field
205  // corresponding to the mass-flux is named "rho"
206  const volScalarField& rho =
207  mesh.objectRegistry::template lookupObject<volScalarField>
208  ("rho");
209 
210  tUflux = faceFlux_/fvc::interpolate(rho);
211  }
212  else if (faceFlux_.dimensions() != dimVolume/dimTime)
213  {
215  << "dimensions of faceFlux are not correct"
216  << exit(FatalError);
217  }
218 
219  return tmp<surfaceScalarField>
220  (
222  (
223  vf.name() + "BlendingFactor",
224  scalar(1)
225  - clamp
226  (
227  (
229  *mag(tUflux)/mesh.magSf()
230  - Co1_
231  )/(Co2_ - Co1_),
232  zero_one{}
233  )
234  )
235  );
236  }
237 
238 
239  //- Return the interpolation weighting factors
240  tmp<surfaceScalarField>
241  weights
242  (
243  const GeometricField<Type, fvPatchField, volMesh>& vf
244  ) const
245  {
247 
248  return
249  bf*tScheme1_().weights(vf)
250  + (scalar(1) - bf)*tScheme2_().weights(vf);
251  }
252 
253 
254  //- Return the face-interpolate of the given cell field
255  // with explicit correction
256  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
258  (
259  const GeometricField<Type, fvPatchField, volMesh>& vf
260  ) const
261  {
263 
264  return
265  bf*tScheme1_().interpolate(vf)
266  + (scalar(1) - bf)*tScheme2_().interpolate(vf);
267  }
268 
269 
270  //- Return true if this scheme uses an explicit correction
271  virtual bool corrected() const
272  {
273  return tScheme1_().corrected() || tScheme2_().corrected();
274  }
275 
276 
277  //- Return the explicit correction to the face-interpolate
278  // for the given field
280  correction
281  (
283  ) const
284  {
286 
287  if (tScheme1_().corrected())
288  {
289  if (tScheme2_().corrected())
290  {
291  return
292  (
293  bf
294  * tScheme1_().correction(vf)
295  + (scalar(1) - bf)
296  * tScheme2_().correction(vf)
297  );
298  }
299  else
300  {
301  return
302  (
303  bf
304  * tScheme1_().correction(vf)
305  );
306  }
307  }
308  else if (tScheme2_().corrected())
309  {
310  return
311  (
312  (scalar(1) - bf)
313  * tScheme2_().correction(vf)
314  );
315  }
316  else
317  {
318  return nullptr;
319  }
320  }
321 };
322 
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
326 } // End namespace Foam
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 #endif
331 
332 // ************************************************************************* //
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: CoBlended.H:317
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const fvMesh & mesh() const
Return mesh reference.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
A class for handling words, derived from Foam::string.
Definition: word.H:63
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Base class for blended schemes to provide access to the blending factor surface field.
Two-scheme Courant number based blending differencing scheme.
Definition: CoBlended.H:88
TypeName("CoBlended")
Runtime type information.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:65
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
Definition: CoBlended.H:221
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: CoBlended.H:304
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Definition: CoBlended.H:289
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Abstract base class for surface interpolation schemes.
virtual const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:61
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Definition: CoBlended.H:270
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
const dimensionSet & dimensions() const noexcept
Return dimensions.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...