cellCoBlended< Type > Class Template Reference

Two-scheme cell-based Courant number based blending differencing scheme. More...

Inheritance diagram for cellCoBlended< Type >:
[legend]
Collaboration diagram for cellCoBlended< Type >:
[legend]

Public Member Functions

 TypeName ("cellCoBlended")
 Runtime type information. More...
 
 cellCoBlended (const fvMesh &mesh, Istream &is)
 Construct from mesh and Istream. More...
 
 cellCoBlended (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
 Construct from mesh, faceFlux and Istream. More...
 
virtual tmp< surfaceScalarFieldblendingFactor (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the face-based blending factor. More...
 
tmp< surfaceScalarFieldweights (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the interpolation weighting factors. More...
 
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the face-interpolate of the given cell field. More...
 
virtual bool corrected () const
 Return true if this scheme uses an explicit correction. More...
 
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the explicit correction to the face-interpolate. More...
 
- Public Member Functions inherited from surfaceInterpolationScheme< Type >
 TypeName ("surfaceInterpolationScheme")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (tmp, surfaceInterpolationScheme, Mesh,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
 
 declareRunTimeSelectionTable (tmp, surfaceInterpolationScheme, MeshFlux,(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData),(mesh, faceFlux, schemeData))
 
 surfaceInterpolationScheme (const fvMesh &mesh)
 Construct from mesh. More...
 
virtual ~surfaceInterpolationScheme ()=default
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
virtual tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Return the face-interpolate of the given cell field. More...
 
tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const surfaceVectorField &Sf, const tmp< GeometricField< Type, fvPatchField, volMesh >> &) const
 Return the face-interpolate of the given tmp cell field. More...
 
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const tmp< GeometricField< Type, fvPatchField, volMesh >> &) const
 Return the face-interpolate of the given tmp cell field. More...
 
- Public Member Functions inherited from refCount
constexpr refCount () noexcept
 Default construct, initializing count to 0. More...
 
int count () const noexcept
 Return the current reference count. More...
 
bool unique () const noexcept
 Return true if the reference count is zero. More...
 
void operator++ () noexcept
 Increment the reference count. More...
 
void operator++ (int) noexcept
 Increment the reference count. More...
 
void operator-- () noexcept
 Decrement the reference count. More...
 
void operator-- (int) noexcept
 Decrement the reference count. More...
 
- Public Member Functions inherited from blendedSchemeBase< Type >
 blendedSchemeBase ()
 Constructor. More...
 
virtual ~blendedSchemeBase ()=default
 Destructor. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from surfaceInterpolationScheme< Type >
static tmp< surfaceInterpolationScheme< Type > > New (const fvMesh &mesh, Istream &schemeData)
 Return new tmp interpolation scheme. More...
 
static tmp< surfaceInterpolationScheme< Type > > New (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
 Return new tmp interpolation scheme. More...
 
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
 Return the face-interpolate of the given cell field. More...
 
template<class SFType >
static tmp< GeometricField< typename innerProduct< typename SFType::value_type, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const SFType &Sf, const GeometricField< Type, fvPatchField, volMesh > &vf, const tmp< surfaceScalarField > &tlambdas)
 Return the face-interpolate of the given cell field. More...
 
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &)
 Return the face-interpolate of the given cell field. More...
 

Detailed Description

template<class Type>
class Foam::cellCoBlended< Type >

Two-scheme cell-based Courant number based blending differencing scheme.

This scheme is equivalent to the CoBlended scheme except that the Courant number is evaluated for cells using the same approach as use in the finite-volume solvers and then interpolated to the faces rather than being estimated directly at the faces based on the flux. This is a more consistent method for evaluating the Courant number but suffers from the need to interpolate which introduces a degree of freedom. However, the interpolation scheme for "Co" is run-time selected and may be specified in "interpolationSchemes" and "localMax" might be most appropriate.

Example of the cellCoBlended scheme specification using LUST for Courant numbers less than 1 and linearUpwind for Courant numbers greater than 10:

divSchemes
{
    .
    .
    div(phi,U)  Gauss cellCoBlended 1 LUST grad(U) 10 linearUpwind grad(U);
    .
    .
}

interpolationSchemes
{
    .
    .
    interpolate(Co) localMax;
    .
    .
}
See also
Foam::CoBlended Foam::localBlended
Source files

Definition at line 93 of file cellCoBlended.H.

Constructor & Destructor Documentation

◆ cellCoBlended() [1/2]

cellCoBlended ( const fvMesh mesh,
Istream is 
)
inline

Construct from mesh and Istream.

The name of the flux field is read from the Istream and looked-up from the mesh objectRegistry

Definition at line 137 of file cellCoBlended.H.

References Foam::exit(), Foam::FatalIOError, and FatalIOErrorInFunction.

Here is the call graph for this function:

◆ cellCoBlended() [2/2]

cellCoBlended ( const fvMesh mesh,
const surfaceScalarField faceFlux,
Istream is 
)
inline

Construct from mesh, faceFlux and Istream.

Definition at line 170 of file cellCoBlended.H.

References Foam::exit(), Foam::FatalIOError, and FatalIOErrorInFunction.

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

TypeName ( "cellCoBlended< Type >"  )

Runtime type information.

◆ blendingFactor()

virtual tmp<surfaceScalarField> blendingFactor ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const
inlinevirtual

Return the face-based blending factor.

Implements blendedSchemeBase< Type >.

Definition at line 203 of file cellCoBlended.H.

References TimeState::deltaTValue(), Foam::dimArea, Foam::dimDensity, Foam::dimless, Foam::dimVelocity, Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::fvc::interpolate(), Foam::mag(), Foam::max(), surfaceInterpolationScheme< Type >::mesh(), Foam::min(), rho, Foam::fvc::surfaceSum(), fvMesh::time(), Time::timeName(), fvMesh::V(), and Foam::Zero.

Referenced by cellCoBlended< Type >::correction(), cellCoBlended< Type >::interpolate(), and cellCoBlended< Type >::weights().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ weights()

tmp<surfaceScalarField> weights ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const
inlinevirtual

Return the interpolation weighting factors.

Implements surfaceInterpolationScheme< Type >.

Definition at line 272 of file cellCoBlended.H.

References cellCoBlended< Type >::blendingFactor().

Here is the call graph for this function:

◆ interpolate()

tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > interpolate ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const
inlinevirtual

Return the face-interpolate of the given cell field.

with explicit correction

Reimplemented from surfaceInterpolationScheme< Type >.

Definition at line 288 of file cellCoBlended.H.

References cellCoBlended< Type >::blendingFactor().

Here is the call graph for this function:

◆ corrected()

virtual bool corrected ( ) const
inlinevirtual

Return true if this scheme uses an explicit correction.

Reimplemented from surfaceInterpolationScheme< Type >.

Definition at line 301 of file cellCoBlended.H.

Referenced by cellCoBlended< Type >::correction().

Here is the caller graph for this function:

◆ correction()

virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > correction ( const GeometricField< Type, fvPatchField, volMesh > &  vf) const
inlinevirtual

Return the explicit correction to the face-interpolate.

for the given field

Reimplemented from surfaceInterpolationScheme< Type >.

Definition at line 311 of file cellCoBlended.H.

References cellCoBlended< Type >::blendingFactor(), and cellCoBlended< Type >::corrected().

Here is the call graph for this function:

The documentation for this class was generated from the following file: