Loading...
Searching...
No Matches
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...

#include <isoAdvection.H>

Public Member Functions

 TypeName ("isoAdvection")
 Runtime type information.
 isoAdvection (volScalarField &alpha1, const surfaceScalarField &phi, const volVectorField &U)
 Constructors.
virtual ~isoAdvection ()=default
 Destructor.
template<class SpType, class SuType>
void advect (const SpType &Sp, const SuType &Su)
 Advect the free surface. Updates alpha field, taking into account.
reconstructionSchemessurf () noexcept
 Return reconstructionSchemes.
const volScalarFieldalpha () const noexcept
 Return alpha field.
const dictionarydict () const noexcept
 Return the controls dictionary.
void writeSurfaceCells () const
 Return cellSet of surface cells.
tmp< surfaceScalarFieldgetRhoPhi (const dimensionedScalar &rho1, const dimensionedScalar &rho2) const
 Return mass flux.
tmp< surfaceScalarFieldgetRhoPhi (const volScalarField &rho1, const volScalarField &rho2)
 Return mass flux.
const surfaceScalarFieldalphaPhi () const noexcept
 reference to alphaPhi
scalar advectionTime () const noexcept
 time spend in the advection step
void writeIsoFaces (const UList< List< point > > &isoFacePts) const
 Write isoface points to .obj file.
template<class SpType, class SuType>
void boundFlux (const bitSet &nextToInterface, surfaceScalarField &dVfcorrectionValues, DynamicList< label > &correctedFaces, const SpType &Sp, const SuType &Su)

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 87 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

Definition at line 48 of file isoAdvection.C.

References alpha1, boundary, cutFace::debug, Foam::dimTime, Foam::dimVol, e, Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::gMinMax(), limits, mesh, Foam::name(), Foam::New(), IOobjectOption::NO_WRITE, UPstream::parRun(), phi, porosityProperties(), primitiveFieldRef(), IOobjectOption::READ_IF_PRESENT, timeName, U, and Foam::Zero.

Here is the call graph for this function:

◆ ~isoAdvection()

virtual ~isoAdvection ( )
virtualdefault

Destructor.

References Foam::Sp(), and Foam::Su().

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

TypeName ( "isoAdvection" )

Runtime type information.

References alpha1, phi, and U.

◆ advect()

template<class SpType, class SuType>
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 423 of file isoAdvectionTemplates.C.

References addProfilingInFunction, DebugInfo, DebugInFunction, Foam::endl(), limitedSurfaceInterpolationScheme< Type >::flux(), Foam::gMinMax(), Foam::Info, limits, Foam::Sp(), Foam::Su(), Foam::fvc::surfaceIntegrate(), and writeSurfaceCells().

Here is the call graph for this function:

◆ surf()

reconstructionSchemes & surf ( )
inlinenoexcept

Return reconstructionSchemes.

Definition at line 426 of file isoAdvection.H.

References Foam::noexcept.

◆ alpha()

const volScalarField & alpha ( ) const
inlinenoexcept

Return alpha field.

Definition at line 434 of file isoAdvection.H.

References Foam::noexcept.

◆ dict()

const dictionary & dict ( ) const
inlinenoexcept

Return the controls dictionary.

Definition at line 442 of file isoAdvection.H.

References Foam::noexcept.

◆ writeSurfaceCells()

void writeSurfaceCells ( ) const

Return cellSet of surface cells.

Definition at line 648 of file isoAdvection.C.

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

Referenced by advect().

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

◆ getRhoPhi() [1/2]

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

Return mass flux.

Definition at line 455 of file isoAdvection.H.

References GeometricField< scalar, fvsPatchField, surfaceMesh >::New(), IOobjectOption::NO_REGISTER, rho1, and rho2.

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 474 of file isoAdvection.H.

References Foam::fvc::interpolate(), GeometricField< scalar, fvsPatchField, surfaceMesh >::New(), IOobjectOption::NO_REGISTER, rho1, and rho2.

Here is the call graph for this function:

◆ alphaPhi()

const surfaceScalarField & alphaPhi ( ) const
inlinenoexcept

reference to alphaPhi

Definition at line 494 of file isoAdvection.H.

References Foam::noexcept.

◆ advectionTime()

scalar advectionTime ( ) const
inlinenoexcept

time spend in the advection step

Definition at line 502 of file isoAdvection.H.

References Foam::noexcept.

◆ writeIsoFaces()

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

Write isoface points to .obj file.

Definition at line 672 of file isoAdvection.C.

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

Here is the call graph for this function:

◆ boundFlux()

template<class SpType, class SuType>
void boundFlux ( const bitSet & nextToInterface,
surfaceScalarField & dVfcorrectionValues,
DynamicList< label > & correctedFaces,
const SpType & Sp,
const SuType & Su )

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