surfaceIteratorIso Class Reference

Finds the isovalue that matches the volume fraction. More...

Public Member Functions

 surfaceIteratorIso (const fvMesh &mesh, scalarField &pointVal, const scalar tol)
 Construct from fvMesh and a scalarField. More...
 
bool isASurfaceCell (const scalar alpha1) const
 Determine if a cell is a surface cell. More...
 
label vofCutCell (const label celli, const scalar alpha1, const scalar tol, const label maxIter)
 
const pointsubCellCentre () const
 The centre point of cutted volume. More...
 
scalar subCellVolume () const
 The volume of cutted volume. More...
 
const pointsurfaceCentre () const
 The centre of cutting isosurface. More...
 
const vectorsurfaceArea () const
 The area vector of cutting isosurface. More...
 
scalar VolumeOfFluid () const
 Volume of Fluid for cellI (subCellVolume_/mesh_.V()[cellI]) More...
 
scalar cutValue () const
 The cutValue. More...
 
label cellStatus ()
 The cellStatus. More...
 
const DynamicList< point > & facePoints ()
 The points of the cutting isosurface in sorted order. More...
 

Detailed Description

Finds the isovalue that matches the volume fraction.

Reference:

    Roenby, J., Bredmose, H. and Jasak, H. (2016).
    A computational method for sharp interface advection
    Royal Society Open Science, 3
    doi 10.1098/rsos.160405

Author Johan Roenby, DHI, all rights reserved.

Source files

Definition at line 64 of file surfaceIteratorIso.H.

Constructor & Destructor Documentation

◆ surfaceIteratorIso()

surfaceIteratorIso ( const fvMesh mesh,
scalarField pointVal,
const scalar  tol 
)

Construct from fvMesh and a scalarField.

The scalarField size should equal the number of mesh points

Definition at line 33 of file surfaceIteratorIso.C.

Member Function Documentation

◆ isASurfaceCell()

bool isASurfaceCell ( const scalar  alpha1) const
inline

Determine if a cell is a surface cell.

Definition at line 100 of file surfaceIteratorIso.H.

References alpha1.

◆ vofCutCell()

Foam::label vofCutCell ( const label  celli,
const scalar  alpha1,
const scalar  tol,
const label  maxIter 
)

finds matching isoValue for the given value fraction returns the cellStatus

Definition at line 49 of file surfaceIteratorIso.C.

References alpha1, e, f(), UList< T >::first(), forAll, UList< T >::last(), Foam::LUsolve(), M, Foam::mag(), Foam::max(), Foam::min(), Foam::pow(), Foam::pow3(), UList< T >::size(), Foam::sortedOrder(), and Foam::sqr().

Here is the call graph for this function:

◆ subCellCentre()

const point & subCellCentre ( ) const
inline

The centre point of cutted volume.

Definition at line 120 of file surfaceIteratorIso.H.

References cutCellIso::subCellCentre().

Here is the call graph for this function:

◆ subCellVolume()

scalar subCellVolume ( ) const
inline

The volume of cutted volume.

Definition at line 126 of file surfaceIteratorIso.H.

References cutCellIso::subCellVolume().

Here is the call graph for this function:

◆ surfaceCentre()

const point & surfaceCentre ( ) const
inline

The centre of cutting isosurface.

Definition at line 132 of file surfaceIteratorIso.H.

References cutCellIso::faceCentre().

Here is the call graph for this function:

◆ surfaceArea()

const vector & surfaceArea ( ) const
inline

The area vector of cutting isosurface.

Definition at line 138 of file surfaceIteratorIso.H.

References cutCellIso::faceArea().

Here is the call graph for this function:

◆ VolumeOfFluid()

scalar VolumeOfFluid ( ) const
inline

Volume of Fluid for cellI (subCellVolume_/mesh_.V()[cellI])

Definition at line 144 of file surfaceIteratorIso.H.

References cutCellIso::VolumeOfFluid().

Here is the call graph for this function:

◆ cutValue()

scalar cutValue ( ) const
inline

The cutValue.

Definition at line 150 of file surfaceIteratorIso.H.

References cutCellIso::cutValue().

Here is the call graph for this function:

◆ cellStatus()

label cellStatus ( )
inline

The cellStatus.

Definition at line 156 of file surfaceIteratorIso.H.

References cutCellIso::cellStatus().

Here is the call graph for this function:

◆ facePoints()

const DynamicList< point > & facePoints ( )
inline

The points of the cutting isosurface in sorted order.

Definition at line 162 of file surfaceIteratorIso.H.

References cutCellIso::facePoints().

Here is the call graph for this function:

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