weightedFlux< Type > Class Template Reference

Weighted flux interpolation scheme class. More...

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

Public Member Functions

 TypeName ("weightedFlux")
 Runtime type information. More...
 
 weightedFlux (const fvMesh &mesh, Istream &is)
 Construct from Istream. More...
 
 weightedFlux (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
 Construct from faceFlux and Istream. More...
 
virtual ~weightedFlux ()
 Destructor. More...
 
tmp< surfaceScalarFieldweights (const GeometricField< Type, fvPatchField, volMesh > &) const
 Return the interpolation weighting factors. More...
 
const surfaceScalarFieldoDelta () const
 Return the distance between face and owner cell. More...
 
const surfaceScalarFieldnDelta () const
 Return the distance between face and neighbour cell. More...
 
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< Type, fvPatchField, volMesh > &vf) const
 Interpolate the cell values to faces. 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 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 > &) const
 Return the explicit correction to the face-interpolate. 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...
 

Protected Member Functions

void clearOut ()
 Clear all fields. 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::weightedFlux< Type >

Weighted flux interpolation scheme class.

This scheme is used to compute fluxes with variable diffusivity or conductivity, as e.g.

  • a thermal flux: lambda*grad(T)
  • a mass flux: D*grad(u)
  • an electric current: -sigma*grad(potential)

When using the Gauss theorem to compute a gradient, cell centred values need to be interpolated to the faces. Using this scheme, temperature (T) is weighted by thermal conductivity when being interpolated. Similarly, velocity is weighted by diffusivity (D) and the electric potential by the electric conductivity (sigma). Lambda, D or sigma are read from the object registry - the names need to be specified in fvSchemes as e.g.

    gradSchemes
    {
        grad(T)             Gauss weightedFlux lambda;
        grad(u)             Gauss weightedFlux D;
        grad(potential)     Gauss weightedFlux sigma;
    }

For more details, see equation 16 and 17 in

    Weber, N., Beckstein, P., Galindo, V., Starace, M. & Weier, T. (2018).
    Electro-vortex flow simulation using coupled meshes.
    Computers and Fluids 168, 101-109.
    doi:10.1016/j.compfluid.2018.03.047
    https://arxiv.org/pdf/1707.06546.pdf
Note
For support, contact Norbe.nosp@m.rt.W.nosp@m.eber@.nosp@m.hzdr.nosp@m..de
Source files

Definition at line 86 of file weightedFlux.H.

Constructor & Destructor Documentation

◆ weightedFlux() [1/2]

weightedFlux ( const fvMesh mesh,
Istream is 
)
inline

Construct from Istream.

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

Definition at line 135 of file weightedFlux.H.

◆ weightedFlux() [2/2]

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

Construct from faceFlux and Istream.

Definition at line 149 of file weightedFlux.H.

◆ ~weightedFlux()

~weightedFlux ( )
virtual

Destructor.

Definition at line 44 of file weightedFlux.C.

Member Function Documentation

◆ clearOut()

void clearOut ( )
protected

Clear all fields.

Definition at line 34 of file weightedFlux.C.

References Foam::deleteDemandDrivenData().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "weightedFlux< Type >"  )

Runtime type information.

◆ weights()

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

Return the interpolation weighting factors.

Implements surfaceInterpolationScheme< Type >.

Definition at line 171 of file weightedFlux.H.

References surfaceInterpolationScheme< Type >::mesh().

Here is the call graph for this function:

◆ oDelta()

const surfaceScalarField& oDelta ( ) const
inline

Return the distance between face and owner cell.

Definition at line 179 of file weightedFlux.H.

◆ nDelta()

const surfaceScalarField& nDelta ( ) const
inline

Return the distance between face and neighbour cell.

Definition at line 190 of file weightedFlux.H.

◆ interpolate()

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

Interpolate the cell values to faces.

Reimplemented from surfaceInterpolationScheme< Type >.

Definition at line 154 of file weightedFlux.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), fvsPatchField< Type >::coupled(), forAll, mesh, and Foam::New().

Here is the call graph for this function:

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