Cellular Automata

Tip

fireplot and fireplot! work with both LevelSetGrid and CAGrid. Load a Makie backend (e.g. using CairoMakie) to enable them.

The CellularAutomata module implements a deterministic, grid-based fire spread model using a travel-time approach. Each cell occupies one of four discrete states and transitions are governed by directional spread rates from a RothermelModel.

How It Works

Cell States

Each cell is in exactly one of four states:

State Value Meaning
UNBURNED 0 Fuel available, not yet ignited
BURNING 1 Actively on fire
BURNED 2 Fuel consumed, fire passed
UNBURNABLE 3 Cannot ignite (water, road, fuel break)
(UNBURNED, BURNING, BURNED, UNBURNABLE)
(UNBURNED, BURNING, BURNED, UNBURNABLE)

State Transitions

At each time step, three phases occur:

  1. Spread: For each BURNING cell, compute directional spread rate R(θ) to each UNBURNED neighbor. The travel time is distance / R(θ), and the neighbor’s arrival time is updated.

  2. Ignition: Cells whose minimum arrival time has been reached transition to BURNING.

  3. Burnout: BURNING cells whose elapsed time exceeds the residence_time transition to BURNED.

UNBURNED ──(arrival time reached)──▶ BURNING ──(burnout)──▶ BURNED
                                        ▲
UNBURNABLE: never transitions ──────────╳

Neighborhoods

Two neighborhood types control which cells can spread fire:

Type Neighbors Description
Moore() 8 Cardinal + diagonal (default)
VonNeumann() 4 Cardinal only

CAGrid

The simulation grid stores per-cell states and ignition times:

grid = CAGrid(100, 100, dx=30.0)
CAGrid{Float64} 100×100 (t=0.0, burning=0, burned=0/10000)
grid_vn = CAGrid(100, 100, dx=30.0, neighborhood=VonNeumann())
CAGrid{Float64} 100×100 (t=0.0, burning=0, burned=0/10000, neighborhood=VonNeumann)

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
CAGrid{Float64} 100×100 (t=0.0, burning=32, burned=0/10000)

Queries

burn_area(grid)  # m²
28800.0
count(burning(grid))  # number of actively burning cells
32

Simulation with RothermelModel

The same RothermelModel used with LevelSetGrid drives the CA simulation. The CA uses directional_speed to compute direction-dependent spread rates to each neighbor, supporting both CosineBlending and EllipticalBlending.

M = FuelClasses(d1=0.06, d10=0.07, d100=0.08, herb=0.0, wood=0.0)
model = RothermelModel(SHORT_GRASS, UniformWind(speed=8.0), UniformMoisture(M), FlatTerrain())

g = CAGrid(100, 100, dx=30.0)
ignite!(g, 1500.0, 1500.0, 100.0)
simulate!(g, model, steps=50)
g
CAGrid{Float64} 100×100 (t=24.100944217000432, burning=572, burned=0/10000)

Plotting Cell States

fireplot(g)

Full Example

M = FuelClasses(d1=0.06, d10=0.07, d100=0.08, herb=0.0, wood=0.0)

fig = Figure(size=(900, 280))
steps_list = [0, 30, 80]
for (col, nsteps) in enumerate(steps_list)
    g = CAGrid(100, 100, dx=30.0)
    ignite!(g, 1500.0, 1500.0, 100.0)
    m = RothermelModel(SHORT_GRASS, UniformWind(speed=8.0), UniformMoisture(M), FlatTerrain())
    nsteps > 0 && simulate!(g, m, steps=nsteps)
    ax = Axis(fig[1, col], title="steps = $nsteps (t = $(round(g.t, digits=1)) min)", aspect=DataAspect())
    fireplot!(ax, g)
    hidedecorations!(ax)
end
fig

Animation

using Wildfires.SpreadModels: Trace

M = FuelClasses(d1=0.06, d10=0.07, d100=0.08, herb=0.0, wood=0.0)
g = CAGrid(100, 100, dx=30.0)
ignite!(g, 1500.0, 1500.0, 100.0)
model = RothermelModel(SHORT_GRASS, UniformWind(speed=8.0), UniformMoisture(M), FlatTerrain())

trace = Trace(g, 2)
simulate!(g, model, steps=80, trace=trace)
firegif(joinpath(@__DIR__, "ca_spread.gif"), trace, g)
"/home/runner/work/Wildfires.jl/Wildfires.jl/docs/ca_spread.gif"

Blending Modes

The blending_mode controls fire shape. CosineBlending (default) produces wider fires, while EllipticalBlending gives more realistic elongated shapes:

M = FuelClasses(d1=0.06, d10=0.07, d100=0.08, herb=0.0, wood=0.0)

fig = Figure(size=(700, 300))
for (col, (label, blend)) in enumerate([
    ("Cosine (default)", CosineBlending()),
    ("Elliptical", EllipticalBlending()),
])
    g = CAGrid(100, 100, dx=30.0)
    ignite!(g, 1500.0, 1500.0, 100.0)
    m = RothermelModel(SHORT_GRASS, UniformWind(speed=8.0), UniformMoisture(M), FlatTerrain(), blend)
    simulate!(g, m, steps=80)
    ax = Axis(fig[1, col], title=label, aspect=DataAspect())
    fireplot!(ax, g)
    hidedecorations!(ax)
end
fig

Fuel Breaks

Mark cells as UNBURNABLE to create barriers that block fire spread:

M = FuelClasses(d1=0.06, d10=0.07, d100=0.08, herb=0.0, wood=0.0)
model = RothermelModel(SHORT_GRASS, UniformWind(speed=8.0), UniformMoisture(M), FlatTerrain())

fig = Figure(size=(700, 300))

# Without fuel break
g1 = CAGrid(100, 100, dx=30.0)
ignite!(g1, 1500.0, 1500.0, 60.0)
simulate!(g1, model, steps=80)
ax1 = Axis(fig[1, 1], title="No fuel break", aspect=DataAspect())
fireplot!(ax1, g1)
hidedecorations!(ax1)

# With fuel break (vertical strip)
g2 = CAGrid(100, 100, dx=30.0)
ignite!(g2, 1500.0, 1500.0, 60.0)
xs = collect(xcoords(g2))
for j in eachindex(xs)
    if 2000.0 <= xs[j] <= 2100.0
        for i in 1:size(g2, 1)
            set_unburnable!(g2, xs[j], ycoords(g2)[i], 1.0)
        end
    end
end
simulate!(g2, model, steps=80)
ax2 = Axis(fig[1, 2], title="With fuel break", aspect=DataAspect())
fireplot!(ax2, g2)
hidedecorations!(ax2)
fig

Burnout

BURNING cells transition to BURNED after residence_time minutes and stop spreading fire. Shorter residence times mean faster burnout:

M = FuelClasses(d1=0.06, d10=0.07, d100=0.08, herb=0.0, wood=0.0)

fig = Figure(size=(900, 280))
for (col, rt) in enumerate([5.0, 15.0, Inf])
    g = CAGrid(100, 100, dx=30.0)
    ignite!(g, 1500.0, 1500.0, 60.0)
    m = RothermelModel(SHORT_GRASS, UniformWind(speed=8.0), UniformMoisture(M), FlatTerrain())
    simulate!(g, m, steps=80, residence_time=rt)
    nb = count(burning(g))
    nd = count(burned(g))
    label = rt == Inf ? "∞" : string(rt)
    ax = Axis(fig[1, col],
        title="residence_time = $label\nburning=$nb, burned=$nd",
        aspect=DataAspect())
    fireplot!(ax, g)
    hidedecorations!(ax)
end
fig

Animation

M = FuelClasses(d1=0.06, d10=0.07, d100=0.08, herb=0.0, wood=0.0)
g = CAGrid(100, 100, dx=30.0)
ignite!(g, 1500.0, 1500.0, 60.0)
model = RothermelModel(SHORT_GRASS, UniformWind(speed=8.0), UniformMoisture(M), FlatTerrain())

trace = Trace(g, 2)
simulate!(g, model, steps=80, residence_time=5.0, trace=trace)
firegif(joinpath(@__DIR__, "ca_burnout.gif"), trace, g)
"/home/runner/work/Wildfires.jl/Wildfires.jl/docs/ca_burnout.gif"

References

  • Carrasco, J. et al. (2021). Cell2Fire: A cell-based forest fire growth model. Frontiers in Forests and Global Change, 4, 692706.