PLIC Class Reference

Piecewise-Linear Interface Calculation (PLIC) corrected scheme is a surface interpolation scheme for flux calculation in advection of a bounded variable, e.g. phase fraction and for interface capturing in the volume of fluid (VoF) method. More...

Inheritance diagram for PLIC:
Collaboration diagram for PLIC:

Public Member Functions

 TypeName ("PLIC")
 Runtime type information. More...
 
 PLIC (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
 Construct from faceFlux and Istream. More...
 
virtual tmp< surfaceScalarFieldinterpolate (const volScalarField &vf) const
 Return the face-interpolate of the given cell field. More...
 
void operator= (const PLIC &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from MPLIC
 TypeName ("MPLIC")
 Runtime type information. More...
 
 MPLIC (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
 Construct from faceFlux and Istream. More...
 
virtual tmp< surfaceScalarFieldweights (const GeometricField< scalar, fvPatchField, volMesh > &) const
 Return the interpolation weighting factors. More...
 
virtual tmp< surfaceScalarFieldinterpolate (const GeometricField< scalar, fvPatchField, volMesh > &vf) const
 Return the face-interpolate of the given cell field. More...
 
void operator= (const MPLIC &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from surfaceInterpolationScheme< scalar >
 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...
 
 surfaceInterpolationScheme (const surfaceInterpolationScheme &)
 Disallow default bitwise copy construction. More...
 
virtual ~surfaceInterpolationScheme ()
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > interpolate (const tmp< GeometricField< scalar, fvPatchField, volMesh >> &) const
 Return the face-interpolate of the given tmp cell field. More...
 
virtual tmp< GeometricField< typename innerProduct< vector, scalar >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const surfaceVectorField &Sf, const GeometricField< scalar, fvPatchField, volMesh > &vf) const
 Return the face-interpolate of the given cell field. More...
 
tmp< GeometricField< typename innerProduct< vector, scalar >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const surfaceVectorField &Sf, const tmp< GeometricField< scalar, fvPatchField, volMesh >> &) const
 Return the face-interpolate of the given tmp cell field. More...
 
virtual bool corrected () const
 Return true if this scheme uses an explicit correction. More...
 
virtual tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > correction (const GeometricField< scalar, fvPatchField, volMesh > &) const
 Return the explicit correction to the face-interpolate. More...
 
void operator= (const surfaceInterpolationScheme &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from refCount
int count () const
 Return the current reference count. More...
 
bool unique () const
 Return true if the reference count is zero. More...
 
void operator++ ()
 Increment the reference count. More...
 
void operator++ (int)
 Increment the reference count. More...
 
void operator-- ()
 Decrement the reference count. More...
 
void operator-- (int)
 Decrement the reference count. More...
 

Protected Attributes

const surfaceScalarFieldphi_
 
tmp< surfaceInterpolationScheme< scalar > > tScheme_
 Base scheme to which the compression is applied. More...
 
- Protected Attributes inherited from MPLIC
const surfaceScalarFieldphi_
 

Additional Inherited Members

- Static Public Member Functions inherited from surfaceInterpolationScheme< scalar >
static tmp< surfaceInterpolationScheme< scalar > > New (const fvMesh &mesh, Istream &schemeData)
 Return new tmp interpolation scheme. More...
 
static tmp< surfaceInterpolationScheme< scalar > > New (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
 Return new tmp interpolation scheme. More...
 
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< scalar, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
 Return the face-interpolate of the given cell field. More...
 
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > interpolate (const GeometricField< scalar, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &)
 Return the face-interpolate of the given cell field. More...
 
static tmp< GeometricField< typename innerProduct< typename SFType::value_type, scalar >::type, fvsPatchField, surfaceMesh > > dotInterpolate (const SFType &Sf, const GeometricField< scalar, fvPatchField, volMesh > &vf, const tmp< surfaceScalarField > &tlambdas)
 Return the face-interpolate of the given cell field. More...
 
- Protected Member Functions inherited from MPLIC
void setCellAlphaf (const label celli, const scalarField &phi, scalarField &alphaf, boolList &correctedFaces, const DynamicList< scalar > &cellAlphaf, const fvMesh &mesh) const
 Set alphaPhi for the faces of the given cell. More...
 
tmp< surfaceScalarFieldsurfaceAlpha (const volScalarField &alpha, const surfaceScalarField &phi, scalarField &spicedTvff, const bool unweighted, const scalar tol, const bool isMPLIC=true) const
 Return alpha interpolation. More...
 
- Protected Member Functions inherited from refCount
 refCount ()
 Construct null initializing count to 0. More...
 

Detailed Description

Piecewise-Linear Interface Calculation (PLIC) corrected scheme is a surface interpolation scheme for flux calculation in advection of a bounded variable, e.g. phase fraction and for interface capturing in the volume of fluid (VoF) method.

The interface is represented by single cuts which split each cell to match the volume fraction of the phase in the cell. The cut planes are oriented according to the point field of the local phase fraction. The phase fraction at each cell face - the interpolated value - is then calculated from the face area on either side of the cut. For cases where the single-cut does not accurately represent the cell volume fraction the specified default scheme is used, e.g. interfaceCompression.

Example:

divSchemes
{
    .
    .
    div(phi,alpha)     Gauss PLIC interfaceCompression vanLeer 1;
    .
    .
}
See also
Foam::PLICU Foam::MPLIC Foam::MPLICU Foam::interfaceCompression
Source files

Definition at line 77 of file PLIC.H.

Constructor & Destructor Documentation

◆ PLIC()

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

Construct from faceFlux and Istream.

Definition at line 102 of file PLIC.H.

References PLIC::interpolate(), and PLIC::operator=().

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

TypeName ( "PLIC"  )

Runtime type information.

◆ interpolate()

Foam::tmp< Foam::surfaceScalarField > interpolate ( const volScalarField vf) const
virtual

Return the face-interpolate of the given cell field.

Reimplemented in PLICU.

Definition at line 43 of file PLIC.C.

References Foam::e, Foam::interpolate(), DimensionedField< Type, GeoMesh >::mesh(), and tmp< T >::ref().

Referenced by PLIC::PLIC().

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

◆ operator=()

void operator= ( const PLIC )
delete

Disallow default bitwise assignment.

Referenced by PLIC::PLIC().

Here is the caller graph for this function:

Member Data Documentation

◆ phi_

const surfaceScalarField& phi_
protected

Definition at line 86 of file PLIC.H.

◆ tScheme_

tmp<surfaceInterpolationScheme<scalar> > tScheme_
protected

Base scheme to which the compression is applied.

Definition at line 89 of file PLIC.H.


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