API reference

Main.LagrangianVoronoi.PolygonMultiType
PolygonMulti(; x::RealVector, kwargs...)

Predefined Voronoi Polygon for Navier-Stokes equations with implicit viscosity, optimized for multiphase problems with high density ratios, including surface tension. Mandatory keyword variable x is the position of the generating seed.

source
Main.LagrangianVoronoi.PolygonNSType
PolygonNS(; x::RealVector, kwargs...)

Predefined Voronoi Polygon for Navier-Stokes equations. Mandatory keyword variable x is the position of the generating seed.

source
Main.LagrangianVoronoi.PolygonNSFType
PolygonNSF(; x::RealVector, kwargs...)

Predefined Voronoi Polygon for Navier-Stokes-Fourier equations. Mandatory keyword variable x is the position of the generating seed.

source
Main.LagrangianVoronoi.bdary_friction!Method
bdary_friction!(grid::VoronoiGrid, vDirichlet::Function, dt::Float64)

Applies a viscous drag at boundaries. This is used for Dirichlet condition for velocity when the the velocity is tangent of the surface. It cannot be used for inflows or outflows. The function vDirichlet should accept a single argument x::RealVector and return a RealVector corresponding to the velocity of the boundary.

source
Main.LagrangianVoronoi.find_D!Method
find_D!(grid::VoronoiGrid)

Compute the velocity deformation tensor $D = \frac{1}{2}(\nabla v + \nabla v^T)$ and assign it to every polygon.

source
Main.LagrangianVoronoi.viscous_step!Method
viscous_step!(grid::VoronoiGrid, dt::Float64; artificial_viscosity::Bool = true)

Applies one forward viscous step of size dt to all Voronoi polygons. It assumes that the velocity deformation tensor is already computed using the find_D! function. It is also assumed that every polygon p has its dynamic viscosity p.mu assigned by the initial condition or otherwise. A keyword parameter controls whether Stone-Norman viscosity tensor is used to damp oscillations near shocks (yes by default).

source
Main.LagrangianVoronoi.FastVectorType
FastVector{T} <: AbstractVector{T}

A vector with preallocated memory than can only increase in size. Compared to normal vector, this improves the performance of mesh generation by about 50%. If this is already implemented somewhere, please let me know.

source
Main.LagrangianVoronoi.fourier_step!Method
fourier_step!(grid::VoronoiGrid, dt::Float64)

Update the energy of Voronoi polygon by Fourier heat conduction. This assumes that every polygon p has its value of heat conductivity p.k assigned by the initial condition or otherwise.

source
Main.LagrangianVoronoi.heat_from_bdary!Method
heat_from_bdary!(grid::VoronoiGrid, dt::Float64, T_bc::Function)

Update the internal energy of a Voronoi cell by considering heat flux from boundary. Function T_bc specifies the temperature at that boundary and should return NaN to indicate adiabatic walls. This currently works only for ideal gas.

source
Main.LagrangianVoronoi.EdgeType
Edge(v1::RealVector, v2::RealVector; label::Int = 0)

Edge defined by two endpoints and a label which indicates the neighbor's index. Zero or negative indices are used for domain boundaries. Edges are stack-allocated and immutable.

source
Main.LagrangianVoronoi.RealMatrixType
RealMatrix(x11::Float64, x21::Float64, x12::Float64, x22::Float64)

Static Float64 matrix with 2x2 elements. Be warned that Julia has column-major ordering of elements!

source
Main.LagrangianVoronoi.cross2Method
function cross2(a::RealVector, b::RealVector)::Float64

A cross product in 2d. Returns a scalar equal to the z component of the corresponding 3d cross product

source
Main.LagrangianVoronoi.vertsMethod
verts(r::Rectangle)::NTuple{4, RealVector}

Return the four vertices of a Rectangle in counter-clockwise order, starting with the bottomleft corner.

source
Main.LagrangianVoronoi.export_gridMethod
export_grid(grid::VoronoiGrid, filename::String, vars::Symbol...)

Export grid to a vtk file. Append datasets specified by vars. Currently only numbers and vectors can be exported.

source
Main.LagrangianVoronoi.boundariesMethod
boundaries(p::VoronoiPolygon)

Create an iterator through all edges of 'p' that lie on the boundary of the domain. This can be used to implement boundary forces.

source
Main.LagrangianVoronoi.neighborsMethod
neighbors(p::VoronoiPolygon, grid::VoronoiGrid)

Create an iterator through all Voronoi polygons that neighbor the polygon p. Use it in a for loop to iterate through all triplets (q,e,y) where

  • q = the neighboring polygon
  • e = the edge connecting p and q
  • y = the position of q as a neighbor of p (equivalent to q.x for non-periodic grids)
source
Main.LagrangianVoronoi.move!Method
move!(grid::VoronoiGrid, dt::Float64)

Update the positions of all polygons, moving each by dt*p.v. The function ensures that points will never escape the computational domain (this would lead to undefined behavior). If the grid is peridoic, points will be wrapped around automatically. The grid is remeshed after this update.

source
Main.LagrangianVoronoi.integralMethod
integral(p::VoronoiPolygon, val::Float64, taylor::SVector)::Float64

Compute an integral of a polynomial over polygon p defined by value val at p.x and the Taylor expansion taylor.

source
Main.LagrangianVoronoi.movinglsMethod
movingls(::Type{T}, grid::VoronoiGrid, p::VoronoiPolygon, fun; h,  kernel)

Finds the Taylor expansion of a given function using moving least squares. The first argument T should be one of following:

  • LinearExpansion
  • QuadraticExpansion
  • CubicExpansion

Polygon p is where the Taylor expansion of a function fun is computed. A keyword argument h is the moving radius. This needs to be sufficiently big, otherwise the Taylor expension may be undefined. Kernel is the weighting function used for defining the (weighted) least squared problem.

source
Main.LagrangianVoronoi.point_valueMethod
point_value(grid::VoronoiGrid, x::RealVector, fun::Function)

Obtain the interpolation of a function fun defined on polygons at an arbitrary point x. This is a slow code and should never be used in perfomance-critical areas.

source
Main.LagrangianVoronoi.poly_evalMethod
poly_eval(val::Float64, taylor::T, dx::RealVector)::Float64 where {T <: SVector}

Use this function to interpolate a value from the known Taylor expansion.

source
Main.LagrangianVoronoi.VoronoiPolygonType
VoronoiPolygon

Abstract supertype for Voronoi Polygons. We call them Voronoi Polygons and not Voronoi cells to avoid confusion with cells of a cell list. No functions need to be defined for subtypes of Voronoi Polygon but we require that it has these fields

  • x::RealVector
  • edges::FastVector{Edge}

where x is the generating seed and edges is a list of all edges in no particular order. The edges can be initialized as edges = emptypolygon() and are generated by the remesh! function. The subtypes may contain additional fields depending on the physics involved. See src/celldefs.jl for examples.

source
Main.LagrangianVoronoi.lr_ratioMethod
lr_ratio(dx::RealVector, e::Edge)::Float64

The ratio between the length of a vector dx and length of the edge e. This may sound like a dubious function but is used so very often in this method. It is slightly faster than norm(dx)/len(e) because it computes only one square root.

source
Main.LagrangianVoronoi.tri_areaMethod
tri_area(a::RealVector, b::RealVector, c::RealVector)::Float64

Return the unsigned area of a triangle given by three points. The points can be in arbitrary order.

source
Main.LagrangianVoronoi.populate_circ!Method
populate_circ!(grid::VoronoiGrid{T}; charfun, center, ic!)

Populate computational domain with polygons arranged in concentric circles. Keyword parameters:

  • charfun: the characteristic function; only those areas where charfun(x) = true are populated
  • center: the center of each circle
  • ic!: the initial condition; ic!(p) is called on every polygon after the mesh is generated
source
Main.LagrangianVoronoi.populate_hex!Method
populate_hex!(grid::VoronoiGrid{T}; charfun, ic!)

Populate computational domain with polygons arranged in hexagonal grid. Use this initializing method when you are not sure. Keyword parameters:

  • charfun: the characteristic function; only those areas where charfun(x) = true are populated
  • ic!: the initial condition; ic!(p) is called on every polygon after the mesh is generated
source
Main.LagrangianVoronoi.populate_lloyd!Method
populate_lloyd!(grid::VoronoiGrid{T}; charfun, niterations, ic!)

Populate computational domain with polygons arranged in concentric circles. Keyword parameters:

  • charfun: the characteristic function; only those areas where charfun(x) = true are populated
  • center: the center of each circle
  • ic!: the initial condition; ic!(p) is called on every polygon after the mesh is generated
source
Main.LagrangianVoronoi.populate_rand!Method
populate_rand!(grid::VoronoiGrid{T}; charfun, ic!)

Populate computational domain with polygons arranged randomly. This initializing method is not recommended because the mesh will have very low quality. Keyword parameters:

  • charfun: the characteristic function; only those areas where charfun(x) = true are populated
  • ic!: the initial condition; ic!(p) is called on every polygon after the mesh is generated
source
Main.LagrangianVoronoi.populate_rect!Method
populate_rect!(grid::VoronoiGrid{T}; charfun, ic!)

Populate computational domain with Cartesian grid. Keyword parameters:

  • charfun: the characteristic function; only those areas where charfun(x) = true are populated
  • ic!: the initial condition; ic!(p) is called on every polygon after the mesh is generated
source
Main.LagrangianVoronoi.populate_vogel!Method
populate_vogel!(grid::VoronoiGrid{T}; charfun, center, ic!)

Populate computational domain with polygons arranged on Vogel spiral. Keyword parameters:

  • charfun: the characteristic function; only those areas where charfun(x) = true are populated
  • center: the center of the spiral
  • ic!: the initial condition; ic!(p) is called on every polygon after the mesh is generated
source
Main.LagrangianVoronoi.find_pressure!Function
find_pressure!(solver::PressureSolver, dt::Float64, niter::Int64)

Find the estimated value of pressure field at the next time step. Integer niter is the number fixed point iteratrions.

source
Main.LagrangianVoronoi.ideal_eos!Function
ideal_eos!(grid::VoronoiGrid, gamma = 1.4; Pmin)

Compute pressure and sound speed from internal energy and density using ideal gas equation of state. Number gamma is adiabatic index and Pmin can specify least possible value pressure. In the semi-implicit scheme, this value of pressure is used as an initial condition for an implicit solver.

source
Main.LagrangianVoronoi.stiffened_eos!Function
stiffened_eos!(grid::VoronoiGrid, gamma = 1.4; P0)

Compute pressure and sound speed from internal energy and density using ideal gas equation of state. Number gamma is adiabatic index and P0 is the stiffness constant. Stiffened equation of state fares better for flows with very low Mach number.

source
Main.LagrangianVoronoi.MultiphaseSolverType
MultiphaseSolver(grid::VoronoiGrid, quality_threshold::Float64 = 0.25)

Construct a solver for the multiphase projection step. A keyword argument quality_threshold determines when dv should not be projected because the quality of the cell is too poor. In these problematic cases, stability of the simulation takes priority over physical accuracy.

source
Main.LagrangianVoronoi.find_dv!Function
find_dv!(grid::VoronoiGrid, dt::Float64, alpha::Float64 = 1.0)

Determine dv, the repair velocity. This is the first step of mesh relaxation procedure and is needed to maintain sufficiently high mesh quality. The repair velocity is chosen to be proportional to Frobenius norm of velocity deformation tensor D and inversely proportional to the mesh quality squared. There is also multiplicative dimensionless parameter alpha with default value alpha = 1. Higher alpha causes the Voronoi cells to be rounder (= better) but also less genuinly Lagrangian.

source
Main.LagrangianVoronoi.multiphase_projection!Method
multiphase_projection!(solver::MultiphaseSolver)

Solve the linear system to find an orthogonal projection of dv to a constraint space which gurantees conservation of area for every fluid phase. You only need this for multiphase flows. You can also run multiphase problems without it but expect some artificial density ridges between phases of different densities.

source
Main.LagrangianVoronoi.relaxation_step!Method
relaxation_step!(grid::VoronoiGrid, dt::Float64; rusanov::Bool = true)

Update the variables by taking into account the repair velocity dv. The repair velocity should be known before the relalaxation_step! is called. See find_dv!. The mass, momentum and energy of each cell is updated by solving one step of an advection problem. A keyword boolean argument rusanov controls the use of a Rusanov approximate Riemann solver. The Rusanov flux is recommended to prevent the generation of new extrema (like negative density) and decrease of entropy. This is a second step of the mesh relaxation procedure (or maybe third, if your simulation involves a multiphase projector).

source
Main.LagrangianVoronoi.SimulationWorkspaceType
SimulationWorkspace

Abstract type for containing the computational grid, solvers and global variables. This methods should be defined for an instance sim of SimulationWorkspace:

  • step!(sim, time::Float64)
  • postproc!(sim, time::Float64) [optional]
source
Main.LagrangianVoronoi.run!Method
run!(sim::SimulationWorkspace, dt::Float64, t_end::Float64, step!::Function; kwargs...)

A utility function for time-marching and exporting the results. Specify the time step dt, final time t_end and a function $step!(sim::SimulationWorkspace, time::Float64).$ Keyword arguments:

  • nframes::Bool (how many times should I export the data?)
  • path::String (where to export the data?)
  • save_points::Bool (should I export the point data?)
  • save_grid:::Bool (should I export the grid data?)
  • save_csv::Bool (should I export global variables as a csv data?)
  • vtp_vars (which local variables should I export?)
  • csv_vars (which global variables should I export?)
  • postproc!::Function (If provided, postproc!(sim, dt) is called each time just before the data is saved.

Use this for some (expensive) post proccessing that affects the output files, not the simulation itself.)

Note this assumes constant time step. See examples/sedov.jl for an example with adaptive time step.

source
Main.LagrangianVoronoi.VoronoiGridType
VoronoiGrid{T}(boundary_rect::Rectangle, dr::Float64, kwargs...)

This struct contains all information about geometry, the Voronoi mesh and a cell list. Type variable T specifies the Voronoi Polygon, boundary_rect defines the computational domain and dr is the default resolution (a particle will typically occupy an area of size dr^2).

Keyword arguments:

  • xperiodic::Bool: Is our domain periodic in the horizontal direction?
  • yperiodic::Bool: Is our domain periodic in the vertical direction?
  • h the size of cells in the cell_list
  • r_max the maximum possible distance between neighboring cells
source
Main.LagrangianVoronoi.get_arrowMethod
get_arrow(x::RealVector, y::RealVector, grid::VoronoiGrid)::RealVector

Return a vector from y to x. Equal to x - y on non-periodic domains. For periodic rectangle, the shortest possible vector is returned.

source
Main.LagrangianVoronoi.remesh!Method
remesh!(grid::VoronoiGrid)

Generate all Voronoi cells by performing all Voronoi cuts. This is called automatically when the VoronoiGrid is populated and each the Voronoi generators move in space.

source