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...
 
template<class SpType , class SuType >
void advect (const SpType &Sp, const SuType &Su)
 Advect the free surface. Updates alpha field, taking into account. More...
 
void applyBruteForceBounding ()
 Apply the bounding based on user inputs. More...
 
reconstructionSchemessurf ()
 Return reconstructionSchemes. 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...
 
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
 reference to alphaPhi More...
 
scalar advectionTime () const
 time spend in the advection step More...
 
void writeIsoFaces (const DynamicList< List< point >> &isoFacePts) const
 Write isoface points to .obj file. More...
 
template<class SpType , class SuType >
void boundFlux (const scalarField &alpha1, surfaceScalarField &dVfcorrectionValues, DynamicList< label > &correctedFaces, const SpType &Sp, const SuType &Su)
 

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) Modified Henning Scheufler, DLR

Source files

Definition at line 78 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 366 of file isoAdvectionTemplates.C.

References 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:

◆ applyBruteForceBounding()

void applyBruteForceBounding ( )

Apply the bounding based on user inputs.

Definition at line 519 of file isoAdvection.C.

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

Here is the call graph for this function:

◆ surf()

reconstructionSchemes& surf ( )
inline

Return reconstructionSchemes.

Definition at line 311 of file isoAdvection.H.

◆ alpha()

const volScalarField& alpha ( ) const
inline

Return alpha field.

Definition at line 317 of file isoAdvection.H.

◆ dict()

const dictionary& dict ( ) const
inline

Return the controls dictionary.

Definition at line 323 of file isoAdvection.H.

◆ writeSurfaceCells()

void writeSurfaceCells ( ) const

Return cellSet of surface cells.

Definition at line 548 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 333 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 350 of file isoAdvection.H.

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

Here is the call graph for this function:

◆ alphaPhi()

const surfaceScalarField& alphaPhi ( ) const
inline

reference to alphaPhi

Definition at line 367 of file isoAdvection.H.

◆ advectionTime()

scalar advectionTime ( ) const
inline

time spend in the advection step

Definition at line 373 of file isoAdvection.H.

◆ writeIsoFaces()

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

Write isoface points to .obj file.

Definition at line 572 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:

◆ boundFlux()

void boundFlux ( const scalarField alpha1,
surfaceScalarField dVfcorrectionValues,
DynamicList< label > &  correctedFaces,
const SpType &  Sp,
const SuType &  Su 
)

Definition at line 214 of file isoAdvectionTemplates.C.

References alpha1, DynamicList< T, SizeMin >::append(), DynamicList< T, SizeMin >::clear(), DebugInfo, DebugInFunction, Foam::endl(), forAll, Foam::mag(), Foam::min(), Foam::neg0(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), phi, Foam::pos0(), Foam::sign(), Sp, and Su.

Here is the call graph for this function:

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