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...
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< surfaceScalarField > | interpolate (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< surfaceScalarField > | weights (const VolField< scalar > &) const |
Return the interpolation weighting factors. More... | |
virtual tmp< surfaceScalarField > | interpolate (const VolField< scalar > &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 &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~surfaceInterpolationScheme () |
Destructor. More... | |
const fvMesh & | mesh () const |
Return mesh reference. More... | |
virtual tmp< SurfaceField< scalar > > | interpolate (const VolField< scalar > &) const |
Return the face-interpolate of the given cell field. More... | |
tmp< SurfaceField< scalar > > | interpolate (const tmp< VolField< scalar >> &) const |
Return the face-interpolate of the given tmp cell field. More... | |
virtual tmp< SurfaceField< typename innerProduct< vector, scalar >::type > > | dotInterpolate (const surfaceVectorField &Sf, const VolField< scalar > &vf) const |
Return the face-interpolate of the given cell field. More... | |
tmp< SurfaceField< typename innerProduct< vector, scalar >::type > > | dotInterpolate (const surfaceVectorField &Sf, const tmp< VolField< scalar >> &) const |
Return the face-interpolate of the given tmp cell field. More... | |
tmp< SurfaceField< typename innerProduct< vector, scalar >::type > > | dotInterpolate (const surfaceVectorField &Sf, const VolField< scalar > &) const |
virtual bool | corrected () const |
Return true if this scheme uses an explicit correction. More... | |
virtual tmp< SurfaceField< scalar > > | correction (const VolField< scalar > &) 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... | |
void | operator= (const refCount &)=delete |
Disallow bitwise assignment. More... | |
Protected Attributes | |
const surfaceScalarField & | phi_ |
tmp< surfaceInterpolationScheme< scalar > > | tScheme_ |
Base scheme to which the compression is applied. More... | |
Protected Attributes inherited from MPLIC | |
const surfaceScalarField & | phi_ |
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< SurfaceField< scalar > > | interpolate (const VolField< scalar > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &) |
Return the face-interpolate of the given cell field. More... | |
static tmp< SurfaceField< scalar > > | interpolate (const VolField< scalar > &, const tmp< surfaceScalarField > &) |
Return the face-interpolate of the given cell field. More... | |
static tmp< SurfaceField< typename innerProduct< typename SFType::value_type, scalar >::type > > | dotInterpolate (const SFType &Sf, const VolField< scalar > &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< surfaceScalarField > | surfaceAlpha (const volScalarField &alpha, const surfaceScalarField &phi, scalarField &splicedTvff, const bool unweighted, const scalar tol, const bool isMPLIC=true) const |
Return alpha interpolation. More... | |
Protected Member Functions inherited from refCount | |
refCount () | |
Construct null initialising count to 0. More... | |
refCount (const refCount &)=delete | |
Disallow copy. More... | |
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; . . }
|
inline |
TypeName | ( | "PLIC" | ) |
Runtime type information.
|
virtual |
Return the face-interpolate of the given cell field.
Reimplemented in PLICU.
Definition at line 42 of file PLIC.C.
References Foam::e, DimensionedField< Type, GeoMesh >::mesh(), PLIC::phi_, tmp< T >::ref(), SlicedGeometricField< Type, PatchField, SlicedPatchField, GeoMesh >::splice(), MPLIC::surfaceAlpha(), and PLIC::tScheme_.
|
delete |
Disallow default bitwise assignment.
|
protected |
Definition at line 86 of file PLIC.H.
Referenced by PLIC::interpolate(), and PLICU::interpolate().
|
protected |
Base scheme to which the compression is applied.
Definition at line 89 of file PLIC.H.
Referenced by PLIC::interpolate(), and PLICU::interpolate().