Level Set

The LevelSet module implements a fire spread simulator using the level set method. A signed distance function φ(x, y, t) implicitly tracks the fire front:

Region Condition Meaning
Burned φ < 0 Fire has passed
Front φ = 0 Active fire perimeter
Unburned φ > 0 Fuel ahead of fire

The front evolves according to the Hamilton-Jacobi equation:

\[\frac{\partial \phi}{\partial t} + F |\nabla \phi| = 0\]

where \(F(x, y, t) \geq 0\) is the local spread rate (e.g. from the Rothermel model).

LevelSetGrid

grid = LevelSetGrid(100, 100, dx=30.0)
LevelSetGrid{Float64} 100×100 (t=0.0, burned=0/10000)

Coordinates

xs = xcoords(grid)
println("x: $(first(xs)) to $(last(xs)) m  ($(length(xs)) cells)")

ys = ycoords(grid)
println("y: $(first(ys)) to $(last(ys)) m  ($(length(ys)) cells)")
x: 15.0 to 2985.0 m  (100 cells)
y: 15.0 to 2985.0 m  (100 cells)

Ignition

ignite!(grid, 1500.0, 1500.0, 100.0)
grid
LevelSetGrid{Float64} 100×100 (t=0.0, burned=32/10000)

Queries

burn_area(grid)  # m²
28800.0
sum(burned(grid))  # number of burned cells
32

Advancing the Front

The spread rate field F has the same dimensions as the grid. Each cell specifies the local spread rate in m/min.

F = fill(10.0, size(grid))  # uniform 10 m/min everywhere

for _ in 1:20
    advance!(grid, F, 0.5)
end

grid
LevelSetGrid{Float64} 100×100 (t=10.0, burned=264/10000)

CFL Condition

The time step dt must satisfy the Courant–Friedrichs–Lewy (CFL) condition for numerical stability:

\[\text{dt} < \frac{\text{dx}}{\max(F)}\]

If dt is too large relative to dx and the maximum spread rate, the upwind scheme becomes unstable and φ values will diverge. For example, with dx = 30 m and max(F) = 31 m/min, you need dt < 0.97 min.

The simulate! function in the SpreadModel module does not automatically enforce the CFL condition — you are responsible for choosing a safe dt.

Reinitialization

Over many time steps, φ can drift away from a true signed distance function (gradients become too flat or too steep). Periodic reinitialization restores |∇φ| ≈ 1, which keeps the upwind scheme accurate.

reinitialize!(grid)
grid
LevelSetGrid{Float64} 100×100 (t=10.0, burned=264/10000)

Plotting

Standard Makie Plots

The LevelSetGrid integrates with Makie’s heatmap, contour, contourf, and surface via convert_arguments:

fig = Figure(size=(700, 300))
ax1 = Axis(fig[1, 1], title="heatmap", aspect=DataAspect())
heatmap!(ax1, grid, colormap=:RdYlGn)

ax2 = Axis(fig[1, 2], title="contourf", aspect=DataAspect())
contourf!(ax2, grid, colormap=:RdYlGn)
fig

fireplot

The fireplot recipe combines a heatmap of φ with a contour line at the fire front:

fig = Figure()
ax = Axis(fig[1, 1], title="t = $(grid.t) min", aspect=DataAspect(),
    xlabel="x (m)", ylabel="y (m)")
fireplot!(ax, grid)
fig

Full Example

Simulate fire spread with a spatially varying spread rate and visualize snapshots:

# Fresh grid with ignition
g = LevelSetGrid(200, 200, dx=15.0)
ignite!(g, 1500.0, 1500.0, 50.0)

# Spread rate: faster on the right side
F = [5.0 + 10.0 * (j / 200)
    for i in 1:200, j in 1:200]

fig = Figure(size=(700, 250))
times = [0, 50, 100]
for (col, target_t) in enumerate(times)
    while g.t < target_t
        advance!(g, F, 0.25)
    end
    col > 1 && reinitialize!(g)
    t = round(g.t, digits=1)
    ax = Axis(fig[1, col],
        title="t = $t min",
        aspect=DataAspect())
    fireplot!(ax, g)
    hidedecorations!(ax)
end
fig

References

  • Osher, S. & Sethian, J.A. (1988). Fronts propagating with curvature-dependent speed. J. Computational Physics, 79(1), 12–49.
  • Mallet, V., Keyes, D.E., & Fendell, F.E. (2009). Modeling wildland fire propagation with level set methods. Computers & Mathematics with Applications.