cellCoBlended.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) 2015-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::cellCoBlended
28 
29 Group
30  grpFvSurfaceInterpolationSchemes
31 
32 Description
33  Two-scheme cell-based Courant number based blending differencing scheme.
34 
35  This scheme is equivalent to the CoBlended scheme except that the Courant
36  number is evaluated for cells using the same approach as use in the
37  finite-volume solvers and then interpolated to the faces rather than being
38  estimated directly at the faces based on the flux. This is a more
39  consistent method for evaluating the Courant number but suffers from the
40  need to interpolate which introduces a degree of freedom. However, the
41  interpolation scheme for "Co" is run-time selected and may be specified in
42  "interpolationSchemes" and "localMax" might be most appropriate.
43 
44  Example of the cellCoBlended scheme specification using LUST for Courant
45  numbers less than 1 and linearUpwind for Courant numbers greater than 10:
46  \verbatim
47  divSchemes
48  {
49  .
50  .
51  div(phi,U) Gauss cellCoBlended 1 LUST grad(U) 10 linearUpwind grad(U);
52  .
53  .
54  }
55 
56  interpolationSchemes
57  {
58  .
59  .
60  interpolate(Co) localMax;
61  .
62  .
63  }
64  \endverbatim
65 
66 See also
67  Foam::CoBlended
68  Foam::localBlended
69 
70 SourceFiles
71  cellCoBlended.C
72 
73 \*---------------------------------------------------------------------------*/
74 
75 #ifndef cellCoBlended_H
76 #define cellCoBlended_H
77 
79 #include "blendedSchemeBase.H"
80 #include "surfaceInterpolate.H"
82 #include "fvcSurfaceIntegrate.H"
83 
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86 namespace Foam
87 {
88 
89 /*---------------------------------------------------------------------------*\
90  Class cellCoBlended Declaration
91 \*---------------------------------------------------------------------------*/
92 
93 template<class Type>
94 class cellCoBlended
95 :
96  public surfaceInterpolationScheme<Type>,
97  public blendedSchemeBase<Type>
98 {
99  // Private data
100 
101  //- Courant number below which scheme1 is used
102  const scalar Co1_;
103 
104  //- Scheme 1
106 
107  //- Courant number above which scheme2 is used
108  const scalar Co2_;
109 
110  //- Scheme 2
112 
113  //- The face-flux used to compute the face Courant number
114  const surfaceScalarField& faceFlux_;
115 
116 
117  // Private Member Functions
118 
119  //- No copy construct
120  cellCoBlended(const cellCoBlended&) = delete;
121 
122  //- No copy assignment
123  void operator=(const cellCoBlended&) = delete;
124 
125 
126 public:
127 
128  //- Runtime type information
129  TypeName("cellCoBlended");
130 
131 
132  // Constructors
133 
134  //- Construct from mesh and Istream.
135  // The name of the flux field is read from the Istream and looked-up
136  // from the mesh objectRegistry
138  (
139  const fvMesh& mesh,
140  Istream& is
141  )
142  :
144  Co1_(readScalar(is)),
145  tScheme1_
146  (
148  ),
149  Co2_(readScalar(is)),
150  tScheme2_
151  (
153  ),
154  faceFlux_
155  (
157  )
158  {
159  if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
160  {
162  << "coefficients = " << Co1_ << " and " << Co2_
163  << " should be > 0 and Co2 > Co1"
164  << exit(FatalIOError);
165  }
166  }
167 
168 
169  //- Construct from mesh, faceFlux and Istream
171  (
172  const fvMesh& mesh,
173  const surfaceScalarField& faceFlux,
174  Istream& is
175  )
176  :
178  Co1_(readScalar(is)),
179  tScheme1_
180  (
181  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
182  ),
183  Co2_(readScalar(is)),
184  tScheme2_
185  (
186  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
187  ),
188  faceFlux_(faceFlux)
189  {
190  if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
191  {
193  << "coefficients = " << Co1_ << " and " << Co2_
194  << " should be > 0 and Co2 > Co1"
195  << exit(FatalIOError);
196  }
197  }
198 
199 
200  // Member Functions
201 
202  //- Return the face-based blending factor
204  (
206  ) const
207  {
208  const fvMesh& mesh = this->mesh();
209  tmp<surfaceScalarField> tUflux = faceFlux_;
210 
211  if (faceFlux_.dimensions() == dimMass/dimTime)
212  {
213  // Currently assume that the density field
214  // corresponding to the mass-flux is named "rho"
215  const volScalarField& rho =
216  mesh.objectRegistry::template lookupObject<volScalarField>
217  ("rho");
218 
219  tUflux = faceFlux_/fvc::interpolate(rho);
220  }
221  else if (faceFlux_.dimensions() != dimVolume/dimTime)
222  {
224  << "dimensions of faceFlux are not correct"
225  << exit(FatalError);
226  }
227 
228  volScalarField Co
229  (
230  IOobject
231  (
232  "Co",
233  mesh.time().timeName(),
234  mesh
235  ),
236  mesh,
239  );
240 
241  scalarField sumPhi
242  (
243  fvc::surfaceSum(mag(tUflux))().primitiveField()
244  );
245 
246  Co.primitiveFieldRef() =
247  (sumPhi/mesh.V().field())*(0.5*mesh.time().deltaTValue());
248  Co.correctBoundaryConditions();
249 
251  (
253  (
254  vf.name() + "BlendingFactor",
255  scalar(1)
256  - clamp
257  (
258  (fvc::interpolate(Co) - Co1_)/(Co2_ - Co1_),
259  zero_one{}
260  )
261  )
262  );
263  }
264 
265 
266  //- Return the interpolation weighting factors
267  tmp<surfaceScalarField>
268  weights
269  (
270  const GeometricField<Type, fvPatchField, volMesh>& vf
271  ) const
272  {
274 
275  return
276  bf*tScheme1_().weights(vf)
277  + (scalar(1) - bf)*tScheme2_().weights(vf);
278  }
279 
280 
281  //- Return the face-interpolate of the given cell field
282  // with explicit correction
283  tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
285  (
286  const GeometricField<Type, fvPatchField, volMesh>& vf
287  ) const
288  {
290 
291  return
292  bf*tScheme1_().interpolate(vf)
293  + (scalar(1) - bf)*tScheme2_().interpolate(vf);
294  }
295 
296 
297  //- Return true if this scheme uses an explicit correction
298  virtual bool corrected() const
299  {
300  return tScheme1_().corrected() || tScheme2_().corrected();
301  }
302 
303 
304  //- Return the explicit correction to the face-interpolate
305  // for the given field
306  virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
307  correction
308  (
310  ) const
311  {
313 
314  if (tScheme1_().corrected())
315  {
316  if (tScheme2_().corrected())
317  {
318  return
319  (
320  bf
321  * tScheme1_().correction(vf)
322  + (scalar(1) - bf)
323  * tScheme2_().correction(vf)
324  );
325  }
326  else
327  {
328  return
329  (
330  bf
331  * tScheme1_().correction(vf)
332  );
333  }
334  }
335  else if (tScheme2_().corrected())
336  {
337  return
338  (
339  (scalar(1) - bf)
340  * tScheme2_().correction(vf)
341  );
342  }
343  else
344  {
345  return nullptr;
346  }
347  }
348 };
349 
350 
351 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
352 
353 } // End namespace Foam
354 
355 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
356 
357 #endif
358 
359 // ************************************************************************* //
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
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...
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
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< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
TypeName("cellCoBlended")
Runtime type information.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
const dimensionSet dimless
Dimensionless.
const fvMesh & mesh() const
Return mesh reference.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
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.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
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
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Abstract base class for surface interpolation schemes.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Two-scheme cell-based Courant number based blending differencing scheme.
Definition: cellCoBlended.H:89
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 ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127