Reference
Contents
Index
LevelSetMethods.AbstractMeshLevelSetMethods.AdvectionTermLevelSetMethods.BoundaryConditionLevelSetMethods.CartesianGridLevelSetMethods.CartesianMeshFieldLevelSetMethods.CurvatureTermLevelSetMethods.DirichletBCLevelSetMethods.LevelSetLevelSetMethods.LevelSetEquationLevelSetMethods.LevelSetTermLevelSetMethods.MeshFieldLevelSetMethods.MeshFieldLevelSetMethods.NeumannBCLevelSetMethods.NeumannGradientBCLevelSetMethods.NormalMotionTermLevelSetMethods.PeriodicBCLevelSetMethods.RK2LevelSetMethods.RK3LevelSetMethods.ReinitializationTermLevelSetMethods.TimeIntegratorLevelSetMethods.D2LevelSetMethods.D2⁰LevelSetMethods.D2⁺⁺LevelSetMethods.D2⁻⁻LevelSetMethods.D⁰LevelSetMethods.D⁺LevelSetMethods.D⁻LevelSetMethods.curvatureLevelSetMethods.curvatureLevelSetMethods.export_surface_meshLevelSetMethods.export_surface_meshLevelSetMethods.export_volume_meshLevelSetMethods.export_volume_meshLevelSetMethods.gradientLevelSetMethods.hessianLevelSetMethods.integrate!LevelSetMethods.makie_themeLevelSetMethods.normalLevelSetMethods.normalLevelSetMethods.perimeterLevelSetMethods.set_makie_theme!LevelSetMethods.volume
LevelSetMethods.AbstractMesh — Typeabstract type AbstractMesh{N,T}An abstract mesh structure in dimension N with primite data of type T.
LevelSetMethods.AdvectionTerm — TypeAdvectionTerm(𝐮, scheme = WENO5())Advection term representing 𝐮 ⋅ ∇ϕ. Available schemes are Upwind and WENO5.
LevelSetMethods.BoundaryCondition — Typeabstract type BoundaryConditionTypes used to specify boundary conditions.
LevelSetMethods.CartesianGrid — MethodCartesianGrid(lc, hc, n)Create a uniform cartesian grid with lower corner lc, upper corner hc and and n nodes in each direction.
Examples
using LevelSetMethods
a = (0, 0)
b = (1, 1)
n = (10, 4)
grid = CartesianGrid(a, b, n)
# output
CartesianGrid{2, Int64}([0, 0], [1, 1], (10, 4))LevelSetMethods.CartesianMeshField — Typeconst CartesianMeshField{V,M<:CartesianGrid} = MeshField{V,M}MeshField over a CartesianGrid.
LevelSetMethods.CurvatureTerm — Typestruct CurvatureTerm{V,M} <: LevelSetTermLevel-set curvature term representing bκ|∇ϕ|, where κ = ∇ ⋅ (∇ϕ/|∇ϕ|) is the curvature.
LevelSetMethods.DirichletBC — Typestruct DirichletBC{T} <: BoundaryConditionA Dirichlet boundary condition taking values of f(x) at the boundary.
LevelSetMethods.LevelSet — TypeLevelSetAlias for MeshField with vals as an AbstractArray of Reals.
LevelSetMethods.LevelSetEquation — MethodLevelSetEquation(; terms, levelset, boundary_conditions, t = 0, integrator = RK3())Create a of a level-set equation of the form ϕₜ + sum(terms) = 0, where each t ∈ terms is a LevelSetTerm and levelset is the initial LevelSet.
Calling integrate!(ls, tf) will evolve the level-set equation up to time tf, modifying the current_state(eq) and current_time(eq) of the object eq in the process (and therefore the original levelset).
Boundary conditions can be specified in two ways. If a single BoundaryCondition is provided, it will be applied uniformly to all boundaries of the domain. To apply different boundary conditions to each boundary, pass a tuple of the form (bc_x, bc_y, ...) with as many elements as dimensions in the domain. If bc_x is a BoundaryCondition, it will be applied to both boundaries in the x direction. If bc_x is a tuple of two BoundaryConditions, the first will be applied to the left boundary and the second to the right boundary. The same logic applies to the other dimensions.
The optional parameter t specifies the initial time of the simulation, and integrator is the TimeIntegrator used to evolve the level-set equation.
using LevelSetMethods, StaticArrays
grid = CartesianGrid((-1, -1), (1, 1), (50, 50)) # define the grid
ϕ = LevelSet(x -> x[1]^2 + x[2]^2 - 0.5^2, grid) # initial shape
𝐮 = MeshField(x -> SVector(1, 0), grid) # advection velocity
terms = (AdvectionTerm(𝐮),) # advection and curvature terms
bc = PeriodicBC() # periodic boundary conditions
eq = LevelSetEquation(; terms, levelset = ϕ, bc) # level-set equation
# output
Level-set equation given by
ϕₜ + 𝐮 ⋅ ∇ ϕ = 0
Current time 0.0
LevelSetMethods.LevelSetTerm — Typeabstract type LevelSetTermA typical term in a level-set evolution equation.
LevelSetMethods.MeshField — Typestruct MeshField{V,M,B}A field described by its discrete values on a mesh.
Base.getindex of an MeshField is overloaded to handle indices that lie outside the CartesianIndices of its MeshField by using bcs.
LevelSetMethods.MeshField — MethodMeshField(f::Function, m)Create a MeshField by evaluating a function f on a mesh m.
LevelSetMethods.NeumannBC — Typestruct NeumannBC <: BoundaryConditionHomogenous Neumann boundary condition, i.e. ∂x ϕ = 0.
LevelSetMethods.NeumannGradientBC — Typestruct NeumannGradientBC <: BoundaryConditionHomogenous Neumann gradient boundary condition, i.e. ∂xx ϕ = 0.
LevelSetMethods.NormalMotionTerm — Typestruct NormalMotionTerm{V,M} <: LevelSetTermLevel-set advection term representing v |∇ϕ|. This LevelSetTerm should be used for internally generated velocity fields; for externally generated velocities you may use AdvectionTerm instead.
LevelSetMethods.PeriodicBC — Typestruct PeriodicBC <: BoundaryConditionSingleton type representing periodic boundary conditions.
LevelSetMethods.RK2 — Typestruct RK2Second order total variation dimishing Runge-Kutta scheme, also known as Heun's predictor-corrector method.
LevelSetMethods.RK3 — Typestruct RK3Third order total variation dimishing Runge-Kutta scheme.
LevelSetMethods.ReinitializationTerm — Typestruct ReinitializationTerm <: LevelSetTermLevel-set term representing sign(ϕ) (|∇ϕ| - 1). This LevelSetTerm should be used for reinitializing the level set into a signed distance function: for a sufficiently large number of time steps this term allows one to solve the Eikonal equation |∇ϕ| = 1.
LevelSetMethods.TimeIntegrator — Typeabstract type TimeIntegrator endAbstract type for time integrators. See subtypes(TimeIntegrator) for a list of available time integrators.
LevelSetMethods.D2 — MethodD2(ϕ::CartesianMeshField,I,dims)Finite difference scheme for second order derivative at grid point I along the dimensions dims.
If dims[1] == dims[2], it is more efficient to call D2⁰(ϕ,I,dims[1]).
LevelSetMethods.D2⁰ — MethodD2⁰(ϕ::CartesianMeshField,I,dim)Centered finite difference scheme for second order derivative at grid point I along dimension dim. E.g. if dim=1, this approximates ∂ₓₓ.
LevelSetMethods.D2⁺⁺ — MethodD2⁺⁺(ϕ::CartesianMeshField,I,dim)Upward finite difference scheme for second order derivative at grid point I along dimension dim. E.g. if dim=1, this approximates ∂ₓₓ.
LevelSetMethods.D2⁻⁻ — MethodD2⁻⁺(ϕ::CartesianMeshField,I,dim)Backward finite difference scheme for second order derivative at grid point I along dimension dim. E.g. if dim=1, this approximates ∂ₓₓ.
LevelSetMethods.D⁰ — MethodD⁰(ϕ::CartesianMeshField,I,dim)Centered finite difference scheme for first order derivative at grid point I along dimension dim.
LevelSetMethods.D⁺ — MethodD⁺(ϕ::CartesianMeshField,I,dim)Forward finite difference scheme for first order derivative at grid point I along dimension dim.
LevelSetMethods.D⁻ — MethodD⁻(ϕ::CartesianMeshField,I,dim)Backward finite difference scheme for first order derivative at grid point I along dimension dim.
LevelSetMethods.curvature — Methodcurvature(ϕ::LevelSet, I)Compute the mean curvature of ϕ at I using κ = ∇ ⋅ (∇ϕ / |∇ϕ|). We use the formula κ = (Δϕ |∇ϕ|^2 - ∇ϕ^T Hϕ ∇ϕ) / |∇ϕ|^3 with first order finite differences. https://en.wikipedia.org/wiki/Meancurvature#Implicitformofmean_curvature
LevelSetMethods.curvature — Methodcurvature(ϕ::LevelSet)Compute the mean curvature of ϕ at I using κ = ∇ ⋅ (∇ϕ / |∇ϕ|). See curvature(ϕ::LevelSet, I) for more details.
using LevelSetMethods
N = 50
grid = CartesianGrid((-1, -1), (1, 1), (N, N))
ϕ = LevelSetMethods.star(grid)
using GLMakie
coeff = exp.(-40.0 * values(ϕ) .^ 2)
κ = curvature(ϕ) .* coeff
xs = LevelSetMethods.grid1d(grid, 1)
ys = LevelSetMethods.grid1d(grid, 2)
fig, ax, hm = heatmap(xs, ys, κ)
Colorbar(fig[:, end+1], hm)
contour!(xs, ys, values(ϕ); levels = [0.0])LevelSetMethods.gradient — Methodgradient(ϕ::LevelSet, I)Return the gradient vector ∇ϕ of ϕ at I
LevelSetMethods.hessian — Methodhessian(ϕ::LevelSet, I)Return the Hessian matrix Hϕ of ϕ at I
LevelSetMethods.integrate! — Functionintegrate!(ls::LevelSetEquation,tf,Δt=Inf)Integrate the LevelSetEquation ls up to time tf, mutating the levelset and current_time of the object ls in the process.
An optional parameter Δt can be passed to specify a maximum time-step allowed for the integration. Note that the internal time-steps taken to evolve the level-set up to tf may be smaller than Δt due to stability reasons related to the terms and integrator employed.
LevelSetMethods.makie_theme — Functionmakie_theme()Return a Makie theme for plots of level-set functions.
LevelSetMethods.normal — Methodnormal(ϕ::LevelSet, I)Compute the unit exterior normal vector of ϕ at I using n = ∇ϕ/|∇ϕ|
LevelSetMethods.normal — Methodnormal(ϕ::LevelSet)Compute the unit exterior normal vector of ϕ using n = ∇ϕ/|∇ϕ|
using LevelSetMethods
N = 50
grid = CartesianGrid((-1, -1), (1, 1), (N, N))
ϕ = LevelSetMethods.star(grid)
using GLMakie
n = normal(ϕ)
xs = LevelSetMethods.grid1d(grid, 1)
ys = LevelSetMethods.grid1d(grid, 2)
coeff = exp.(-40.0 * values(ϕ) .^ 2)
us = getindex.(n, 1) .* coeff
vs = getindex.(n, 2) .* coeff
arrows(xs, ys, us, vs; arrowsize = 10 * vec(coeff), lengthscale = 2.0 / (N - 1))
contour!(xs, ys, values(ϕ); levels = [0.0])LevelSetMethods.perimeter — Methodperimeter(ϕ::LevelSet)Compute the perimeter area of the level-set function.
Note: this function does not compute the perimeter on the borders of the domain.
using LevelSetMethods
R = 0.5
S0 = 2π * R
grid = CartesianGrid((-1, -1), (1, 1), (200, 200))
ϕ = LevelSetMethods.circle(grid; center = (0, 0), radius = R)
LevelSetMethods.perimeter(ϕ), S0
# output
(3.1426415491430366, 3.141592653589793)LevelSetMethods.set_makie_theme! — Functionset_makie_theme!()Set the Makie theme to LevelSetMethods.makie_theme().
LevelSetMethods.volume — Methodvolume(ϕ::LevelSet)Compute the volume of the level-set function.
using LevelSetMethods
R = 0.5
V0 = π * R^2
grid = CartesianGrid((-1, -1), (1, 1), (200, 200))
ϕ = LevelSetMethods.circle(grid; center = (0, 0), radius = R)
LevelSetMethods.volume(ϕ), V0
# output
(0.7854362890190668, 0.7853981633974483)LevelSetMethods.export_surface_mesh — Methodexport_surface_mesh(ϕ::LevelSet, output::String;
hgrad = nothing, hmin = nothing, hmax = nothing, hausd = nothing)Compute a mesh of the LevelSet ϕ zero contour using MMGs_O3.
hgrad control the growth ratio between two adjacent edges
hmin and hmax control the edge sizes to be (respectively) greater than the hmin parameter and lower than the hmax one
hausd control the maximal distance between the piecewise linear representation of the boundary and the reconstructed ideal boundary
Only works for 3 dimensional level-set.
LevelSetMethods.export_surface_mesh — Methodexport_surface_mesh(eq::LevelSetEquation, args...; kwargs...)Call export_surface_mesh on current_state(eq).
LevelSetMethods.export_volume_mesh — Methodexport_volume_mesh(ϕ::LevelSet, output::String;
hgrad = nothing, hmin = nothing, hmax = nothing, hausd = nothing)Compute a mesh of the domains associated with LevelSet eq using either MMG2dO3 or MMG3dO3.
hgrad control the growth ratio between two adjacent edges.
hmin and hmax control the edge sizes to be (respectively) greater than the hmin parameter and lower than the hmax one.
hausd control the maximal distance between the piecewise linear representation of the boundary and the reconstructed ideal boundary.
For more information, see the official MMG documentation.
Only works for 2 and 3 dimensional level-set.
LevelSetMethods.export_volume_mesh — Methodexport_volume_mesh(eq::LevelSetEquation, output; kwargs...)Call export_volume_mesh on current_state(eq).