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