isoAdvection Class Reference

Calculates the new VOF (alpha) field after time step dt given the initial VOF field and a velocity field U and face fluxes phi. The fluid transport calculation is based on an idea of using isosurfaces to estimate the internal distribution of fluid in cells and advecting such isosurfaces across the mesh faces with the velocity field interpolated to the isosurfaces. More...

Public Member Functions

 TypeName ("isoAdvection")
 Runtime type information. More...
 
 isoAdvection (volScalarField &alpha1, const surfaceScalarField &phi, const volVectorField &U)
 Constructors. More...
 
virtual ~isoAdvection ()=default
 Destructor. More...
 
void advect ()
 Advect the free surface. Updates alpha field, taking into account. More...
 
void applyBruteForceBounding ()
 Apply the bounding based on user inputs. More...
 
const volScalarFieldalpha () const
 Return alpha field. More...
 
const dictionarydict () const
 Return the controls dictionary. More...
 
void writeSurfaceCells () const
 Return cellSet of surface cells. More...
 
void writeBoundedCells () const
 Return cellSet of bounded cells. More...
 
tmp< surfaceScalarFieldgetRhoPhi (const dimensionedScalar rho1, const dimensionedScalar rho2) const
 Return mass flux. More...
 
scalar advectionTime () const
 
void writeIsoFaces (const DynamicList< List< point >> &isoFacePts) const
 Write isoface points to .obj file. More...
 

Detailed Description

Calculates the new VOF (alpha) field after time step dt given the initial VOF field and a velocity field U and face fluxes phi. The fluid transport calculation is based on an idea of using isosurfaces to estimate the internal distribution of fluid in cells and advecting such isosurfaces across the mesh faces with the velocity field interpolated to the isosurfaces.

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

Original code supplied by Johan Roenby, DHI (2016)

Source files

Definition at line 74 of file isoAdvection.H.

Constructor & Destructor Documentation

◆ isoAdvection()

isoAdvection ( volScalarField alpha1,
const surfaceScalarField phi,
const volVectorField U 
)

Constructors.

Construct given alpha, phi and velocity field. Note: phi should be

divergence free up to a sufficient tolerance

◆ ~isoAdvection()

virtual ~isoAdvection ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "isoAdvection"  )

Runtime type information.

◆ advect()

void advect ( )

Advect the free surface. Updates alpha field, taking into account.

multiple calls within a single time step.

Definition at line 842 of file isoAdvection.C.

References DebugInFunction, Foam::endl(), limitedSurfaceInterpolationScheme< Type >::flux(), Foam::gMax(), Foam::gMin(), Foam::Info, and Foam::fvc::surfaceIntegrate().

Here is the call graph for this function:

◆ applyBruteForceBounding()

void applyBruteForceBounding ( )

Apply the bounding based on user inputs.

Definition at line 888 of file isoAdvection.C.

References Foam::max(), Foam::min(), Foam::neg0(), and Foam::pos0().

Here is the call graph for this function:

◆ alpha()

const volScalarField& alpha ( ) const
inline

Return alpha field.

Definition at line 340 of file isoAdvection.H.

◆ dict()

const dictionary& dict ( ) const
inline

Return the controls dictionary.

Definition at line 346 of file isoAdvection.H.

◆ writeSurfaceCells()

void writeSurfaceCells ( ) const

Return cellSet of surface cells.

Definition at line 917 of file isoAdvection.C.

References HashSet< Key, Hash >::insert(), IOobject::NO_READ, and regIOobject::write().

Here is the call graph for this function:

◆ writeBoundedCells()

void writeBoundedCells ( ) const

Return cellSet of bounded cells.

Definition at line 941 of file isoAdvection.C.

References HashSet< Key, Hash >::insert(), IOobject::NO_READ, and regIOobject::write().

Here is the call graph for this function:

◆ getRhoPhi()

tmp<surfaceScalarField> getRhoPhi ( const dimensionedScalar  rho1,
const dimensionedScalar  rho2 
) const
inline

Return mass flux.

Definition at line 359 of file isoAdvection.H.

References TimeState::deltaT(), rho1, rho2, and fvMesh::time().

Here is the call graph for this function:

◆ advectionTime()

scalar advectionTime ( ) const
inline

Definition at line 374 of file isoAdvection.H.

◆ writeIsoFaces()

void writeIsoFaces ( const DynamicList< List< point >> &  isoFacePts) const

Write isoface points to .obj file.

Definition at line 969 of file isoAdvection.C.

References Foam::endl(), f(), forAll, Pstream::gatherList(), Foam::identity(), Foam::Info, UPstream::master(), Foam::mkDir(), UPstream::myProcNo(), Foam::nl, UPstream::nProcs(), UPstream::parRun(), fileName::path(), and word::printf().

Here is the call graph for this function:

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