cellCoBlended.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) 2015-2021 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 public:
113 
114  //- Runtime type information
115  TypeName("cellCoBlended");
116 
117 
118  // Constructors
119 
120  //- Construct from mesh and Istream.
121  // The name of the flux field is read from the Istream and looked-up
122  // from the mesh objectRegistry
124  (
125  const fvMesh& mesh,
126  Istream& is
127  )
128  :
129  surfaceInterpolationScheme<Type>(mesh),
130  Co1_(readScalar(is)),
131  tScheme1_
132  (
133  surfaceInterpolationScheme<Type>::New(mesh, is)
134  ),
135  Co2_(readScalar(is)),
136  tScheme2_
137  (
138  surfaceInterpolationScheme<Type>::New(mesh, is)
139  ),
140  faceFlux_
141  (
142  mesh.lookupObject<surfaceScalarField>(word(is))
143  )
144  {
145  if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
146  {
148  << "coefficients = " << Co1_ << " and " << Co2_
149  << " should be > 0 and Co2 > Co1"
150  << exit(FatalIOError);
151  }
152  }
153 
154 
155  //- Construct from mesh, faceFlux and Istream
157  (
158  const fvMesh& mesh,
159  const surfaceScalarField& faceFlux,
160  Istream& is
161  )
162  :
163  surfaceInterpolationScheme<Type>(mesh),
164  Co1_(readScalar(is)),
165  tScheme1_
166  (
167  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
168  ),
169  Co2_(readScalar(is)),
170  tScheme2_
171  (
172  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
173  ),
174  faceFlux_(faceFlux)
175  {
176  if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
177  {
179  << "coefficients = " << Co1_ << " and " << Co2_
180  << " should be > 0 and Co2 > Co1"
181  << exit(FatalIOError);
182  }
183  }
184 
185  //- Disallow default bitwise copy construction
186  cellCoBlended(const cellCoBlended&) = delete;
187 
188 
189  // Member Functions
190 
191  //- Return the face-based blending factor
193  (
195  ) const
196  {
197  const fvMesh& mesh = this->mesh();
198  tmp<surfaceScalarField> tUflux = faceFlux_;
199 
200  if (faceFlux_.dimensions() == dimMassFlux)
201  {
202  // Currently assume that the density field
203  // corresponding to the mass-flux is named "rho"
204  const volScalarField& rho =
205  mesh.objectRegistry::template lookupObject<volScalarField>
206  ("rho");
207 
208  tUflux = faceFlux_/fvc::interpolate(rho);
209  }
210  else if (faceFlux_.dimensions() != dimFlux)
211  {
213  << "dimensions of faceFlux are not correct"
214  << exit(FatalError);
215  }
216 
217  volScalarField Co
218  (
219  IOobject
220  (
221  "Co",
222  mesh.time().timeName(),
223  mesh
224  ),
225  mesh,
227  extrapolatedCalculatedFvPatchScalarField::typeName
228  );
229 
230  scalarField sumPhi
231  (
232  fvc::surfaceSum(mag(tUflux))().primitiveField()
233  );
234 
235  Co.primitiveFieldRef() =
236  (sumPhi/mesh.V().field())*(0.5*mesh.time().deltaTValue());
237  Co.correctBoundaryConditions();
238 
240  (
241  vf.name() + "BlendingFactor",
242  scalar(1)
243  - max
244  (
245  min
246  (
247  (fvc::interpolate(Co) - Co1_)/(Co2_ - Co1_),
248  scalar(1)
249  ),
250  scalar(0)
251  )
252  );
253  }
254 
255 
256  //- Return the interpolation weighting factors
259  (
261  ) const
262  {
264 
265  return
266  bf*tScheme1_().weights(vf)
267  + (scalar(1) - bf)*tScheme2_().weights(vf);
268  }
269 
270 
271  //- Return the face-interpolate of the given cell field
272  // with explicit correction
275  (
277  ) const
278  {
280 
281  return
282  bf*tScheme1_().interpolate(vf)
283  + (scalar(1) - bf)*tScheme2_().interpolate(vf);
284  }
285 
286 
287  //- Return true if this scheme uses an explicit correction
288  virtual bool corrected() const
289  {
290  return tScheme1_().corrected() || tScheme2_().corrected();
291  }
292 
293 
294  //- Return the explicit correction to the face-interpolate
295  // for the given field
298  (
300  ) const
301  {
303 
304  if (tScheme1_().corrected())
305  {
306  if (tScheme2_().corrected())
307  {
308  return
309  (
310  bf
311  * tScheme1_().correction(vf)
312  + (scalar(1) - bf)
313  * tScheme2_().correction(vf)
314  );
315  }
316  else
317  {
318  return
319  (
320  bf
321  * tScheme1_().correction(vf)
322  );
323  }
324  }
325  else if (tScheme2_().corrected())
326  {
327  return
328  (
329  (scalar(1) - bf)
330  * tScheme2_().correction(vf)
331  );
332  }
333  else
334  {
336  (
337  nullptr
338  );
339  }
340  }
341 
342 
343  // Member Operators
344 
345  //- Disallow default bitwise assignment
346  void operator=(const cellCoBlended&) = delete;
347 };
348 
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 } // End namespace Foam
353 
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355 
356 #endif
357 
358 // ************************************************************************* //
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
const word & name() const
Return name.
Definition: IOobject.H:303
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:323
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 tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< scalar >> &)
Return a temporary field constructed from name,.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const dimensionSet dimless
const fvMesh & mesh() const
Return mesh reference.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
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.
const dimensionSet dimFlux
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 successful.
Definition: doubleScalar.H:75
void operator=(const cellCoBlended &)=delete
Disallow default bitwise assignment.
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.
cellCoBlended(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
const dimensionSet dimMassFlux
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
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
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
Namespace for OpenFOAM.
IOerror FatalIOError