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