isoAdvection Class Reference

An implementation of the isoAdvector geometric Volume-of-Fluid method advancing the provided volume fraction field (alpha1) in time using the given velocity field, U, and corresponding face fluxes, phi. 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...
 
template<class SpType , class SuType >
void advect (const SpType &Sp, const SuType &Su)
 Advect the free surface. Updates alpha field, taking into account. More...
 
reconstructionSchemessurf () noexcept
 Return reconstructionSchemes. More...
 
const volScalarFieldalpha () const noexcept
 Return alpha field. More...
 
const dictionarydict () const noexcept
 Return the controls dictionary. More...
 
void writeSurfaceCells () const
 Return cellSet of surface cells. More...
 
tmp< surfaceScalarFieldgetRhoPhi (const dimensionedScalar rho1, const dimensionedScalar rho2) const
 Return mass flux. More...
 
tmp< surfaceScalarFieldgetRhoPhi (const volScalarField &rho1, const volScalarField &rho2)
 Return mass flux. More...
 
const surfaceScalarFieldalphaPhi () const noexcept
 reference to alphaPhi More...
 
scalar advectionTime () const noexcept
 time spend in the advection step More...
 
void writeIsoFaces (const DynamicList< List< point >> &isoFacePts) const
 Write isoface points to .obj file. More...
 

Detailed Description

An implementation of the isoAdvector geometric Volume-of-Fluid method advancing the provided volume fraction field (alpha1) in time using the given velocity field, U, and corresponding face fluxes, phi.

References:

Main isoAdvector idea:

    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

Calculation of rhoPhi:

    Roenby, J., Bredmose, H., & Jasak, H. (2019).
    IsoAdvector: Geometric VOF on general meshes.
    OpenFOAMĀ® (pp. 281-296). Springer, Cham.

Extension to porous media flows:

    Missios, K., Jacobsen, N. G., Moeller, K., & Roenby, J.
    Using the isoAdvector Geometric VOF Method for Interfacial Flows 
    Through Porous Media. MARINE 2021.

Original code supplied by Johan Roenby, DHI (2016) Modified Henning Scheufler, DLR

Source files

Definition at line 91 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 ( const SpType &  Sp,
const SuType &  Su 
)

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

multiple calls within a single time step.

Definition at line 392 of file isoAdvectionTemplates.C.

References addProfilingInFunction, DebugInfo, DebugInFunction, Foam::endl(), zeroField::field(), limitedSurfaceInterpolationScheme< Type >::flux(), Foam::gMax(), Foam::gMin(), Foam::Info, Sp, Su, and Foam::fvc::surfaceIntegrate().

Here is the call graph for this function:

◆ surf()

reconstructionSchemes& surf ( )
inlinenoexcept

Return reconstructionSchemes.

Definition at line 334 of file isoAdvection.H.

◆ alpha()

const volScalarField& alpha ( ) const
inlinenoexcept

Return alpha field.

Definition at line 340 of file isoAdvection.H.

◆ dict()

const dictionary& dict ( ) const
inlinenoexcept

Return the controls dictionary.

Definition at line 346 of file isoAdvection.H.

◆ writeSurfaceCells()

void writeSurfaceCells ( ) const

Return cellSet of surface cells.

Definition at line 640 of file isoAdvection.C.

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

Here is the call graph for this function:

◆ getRhoPhi() [1/2]

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

Return mass flux.

Definition at line 356 of file isoAdvection.H.

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

Here is the call graph for this function:

◆ getRhoPhi() [2/2]

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

Return mass flux.

Definition at line 373 of file isoAdvection.H.

References Foam::fvc::interpolate(), rho1, and rho2.

Here is the call graph for this function:

◆ alphaPhi()

const surfaceScalarField& alphaPhi ( ) const
inlinenoexcept

reference to alphaPhi

Definition at line 390 of file isoAdvection.H.

◆ advectionTime()

scalar advectionTime ( ) const
inlinenoexcept

time spend in the advection step

Definition at line 396 of file isoAdvection.H.

◆ writeIsoFaces()

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

Write isoface points to .obj file.

Definition at line 665 of file isoAdvection.C.

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

Here is the call graph for this function:

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