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-2024 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 interpolation 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  :
130  Co1_(readScalar(is)),
131  tScheme1_
132  (
134  ),
135  Co2_(readScalar(is)),
136  tScheme2_
137  (
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  :
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  (
194  const VolField<Type>& vf
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() != dimVolumetricFlux)
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().name(),
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().primitiveField())
237  *(0.5*mesh.time().deltaTValue());
238 
240 
242  (
243  vf.name() + "BlendingFactor",
244  scalar(1)
245  - max
246  (
247  min
248  (
249  (fvc::interpolate(Co) - Co1_)/(Co2_ - Co1_),
250  scalar(1)
251  ),
252  scalar(0)
253  )
254  );
255  }
256 
257 
258  //- Return the interpolation weighting factors
260  weights
261  (
262  const VolField<Type>& vf
263  ) const
264  {
266 
267  return
268  bf*tScheme1_().weights(vf)
269  + (scalar(1) - bf)*tScheme2_().weights(vf);
270  }
271 
272 
273  //- Return the face-interpolate of the given cell field
274  // with explicit correction
277  (
278  const VolField<Type>& vf
279  ) const
280  {
282 
283  return
284  bf*tScheme1_().interpolate(vf)
285  + (scalar(1) - bf)*tScheme2_().interpolate(vf);
286  }
287 
288 
289  //- Return true if this scheme uses an explicit correction
290  virtual bool corrected() const
291  {
292  return tScheme1_().corrected() || tScheme2_().corrected();
293  }
294 
295 
296  //- Return the explicit correction to the face-interpolate
297  // for the given field
299  correction
300  (
301  const VolField<Type>& vf
302  ) const
303  {
305 
306  if (tScheme1_().corrected())
307  {
308  if (tScheme2_().corrected())
309  {
310  return
311  (
312  bf
313  * tScheme1_().correction(vf)
314  + (scalar(1) - bf)
315  * tScheme2_().correction(vf)
316  );
317  }
318  else
319  {
320  return
321  (
322  bf
323  * tScheme1_().correction(vf)
324  );
325  }
326  }
327  else if (tScheme2_().corrected())
328  {
329  return
330  (
331  (scalar(1) - bf)
332  * tScheme2_().correction(vf)
333  );
334  }
335  else
336  {
337  return tmp<SurfaceField<Type>>
338  (
339  nullptr
340  );
341  }
342  }
343 
344 
345  // Member Operators
346 
347  //- Disallow default bitwise assignment
348  void operator=(const cellCoBlended&) = delete;
349 };
350 
351 
352 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
353 
354 } // End namespace Foam
355 
356 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357 
358 #endif
359 
360 // ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
void correctBoundaryConditions()
Correct boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
Base class for blended schemes to provide access to the blending factor surface field.
Two-scheme cell-based Courant number based blending interpolation scheme.
Definition: cellCoBlended.H:92
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
tmp< SurfaceField< Type > > interpolate(const VolField< Type > &vf) const
Return the face-interpolate of the given cell field.
tmp< surfaceScalarField > weights(const VolField< Type > &vf) const
Return the interpolation weighting factors.
virtual tmp< surfaceScalarField > blendingFactor(const VolField< Type > &vf) const
Return the face-based blending factor.
cellCoBlended(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
void operator=(const cellCoBlended &)=delete
Disallow default bitwise assignment.
virtual tmp< SurfaceField< Type > > correction(const VolField< Type > &vf) const
Return the explicit correction to the face-interpolate.
TypeName("cellCoBlended")
Runtime type information.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const dimensionSet dimMassFlux
const dimensionSet dimVolumetricFlux
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError