1. Basics of specialized workflows
  2. Basics of agent-based modeling: spatial epidemic dynamics with Agents.jl
  • (Just enough) Julia for scientific informatics, modeling, and reasoning
  • Introduction
  • Basic frameworks and mechanisms
    • Orientation
    • Basics of setting up and running Julia
    • Basics of visualizing mathematical models
    • Basics of working with randomness and probabilities
    • Basics of working with data tables
  • Basics of specialized workflows
    • Basics of paleobiological fossil collection analyses
    • Basics of agent-based modeling: spatial epidemic dynamics with Agents.jl
      • Basics of agent-based modeling: spatial epidemic dynamics with Agents.jl
    • Basics of species distribution modeling
  • Primers
    • Bernoulli trial
    • Pathogen fitness as a function of virulence (Frank, 1996)
    • Virulence-transmission trade-off (Frank, 1996)
    • Julia – Environments – Global vs project
    • Julia: Functions, methods, and signatures
    • Markov property
    • Probabilty distributions–Essential concepts
    • Pseudo-random number generators
    • Pseudo-random number generators: best practices
    • Pseudo-random number generators: continuous values from discrete machines

On this page

  • 1 What is an agent-based model?
    • 1.1 The components of an ABM
    • 1.2 The SIR framework as an ABM
  • 2 Setting up the environment
  • 3 Defining the agent
  • 4 Model parameters and initialization
  • 5 The stepping function
  • 6 Running the model
    • 6.1 Single step
    • 6.2 Running for many steps and collecting data
  • 7 Visualizing the epidemic time series
  • 8 Visualizing the spatial state
    • 8.1 Producing a multi-panel snapshot sequence
  • 9 Parameter sweeps with paramscan
    • 9.1 Visualizing the parameter sweep
  • 10 Extending the model: heterogeneous transmission
  • 11 Exercises
  1. Basics of specialized workflows
  2. Basics of agent-based modeling: spatial epidemic dynamics with Agents.jl

Basics of agent-based modeling: spatial epidemic dynamics with Agents.jl

Author

Jeet Sukumaran

Proficiencies
  • Background concepts:
    • What an agent-based model (ABM) is and how it differs from equation-based models
    • The relationship between individual-level rules and population-level emergent dynamics
    • The SIR compartmental framework and what it assumes
    • Why space matters: local contact structure vs. well-mixed assumptions
    • Stochastic vs. deterministic dynamics in small vs. large populations
    • What a scheduler is and how agent update order affects outcomes
  • Core skills:
    • Defining an agent type with @agent and the fields it carries
    • Creating a GridSpace and an AgentBasedModel (ABM)
    • Writing agent stepping functions that implement individual behavioral rules
    • Initializing a model with a defined starting condition
    • Running a model with step! and run!
    • Collecting agent-level and model-level data during a run using adata and mdata
    • Visualizing time series of population-level dynamics with CairoMakie
    • Visualizing the spatial state of the model with abmplot
    • Sweeping over parameter combinations with paramscan

1 What is an agent-based model?

An agent-based model (ABM) — also called an individual-based model (IBM) — is a computational model in which a population of discrete, autonomous agents each follow a defined set of local rules, and the system-level dynamics emerge from the aggregate of all individual interactions over time.

This is fundamentally different from the equation-based models we encountered earlier in this course. In the virulence-transmission trade-off model of Frank (1996), the dynamics of the entire host–pathogen system were summarized as a small number of state variables (the reproductive number \(R_0\), the exploitation severity \(z\)) governed by deterministic equations. These models describe the average behavior of infinitely large, well-mixed populations, and they do so efficiently: a handful of equations can characterize system behavior across the full parameter space.

But real biological populations are not infinitely large, well-mixed, or average. Individual organisms move through space, encounter neighbors, make decisions, vary in their traits, are born, and die. When we want to ask how spatial structure shapes epidemic spread, how rare mutants invade populations where stochastic loss is a real possibility, or how individual behavioral heterogeneity alters collective outcomes, we need a framework that operates at the level of the individual. ABMs provide that framework.

TipEmergent dynamics

Emergence is the appearance of higher-level patterns that are not explicitly encoded in any individual agent’s rules. A classic example: no individual ant is following a rule that says “build a spiral tunnel structure.” Each ant follows simple local rules about where to deposit or pick up soil, and the global spiral architecture emerges from the aggregate of those interactions. In epidemic ABMs, the characteristic spatial wave of infection spreading from a focal source, or the patchy, fragmented pattern of infection under high host heterogeneity, are emergent: no agent is instructed to produce a wave.

1.1 The components of an ABM

Every ABM consists of three essential ingredients:

  1. Agents: discrete entities with internal state (their attributes or fields) and rules governing their behavior (the stepping function).
  2. Space: the environment the agents occupy. This may be a grid, a continuous plane, a network, or no space at all. Space determines which agents can interact with which.
  3. Scheduler: the rule that determines the order in which agents are updated at each time step. For stochastic models, random activation order matters.

At each time step, the model iterates through all agents (in scheduler-determined order), applies the stepping function to each, and optionally collects data.

1.2 The SIR framework as an ABM

We will build a spatial stochastic SIR model: a classic framework in infectious disease epidemiology that partitions the host population into three compartments:

  • S (Susceptible): not yet infected, capable of contracting the disease.
  • I (Infected): currently infected and capable of transmitting to susceptible neighbors.
  • R (Recovered): no longer infectious; immune to reinfection (within this model).

The deterministic well-mixed SIR model — which you may have encountered as a set of ordinary differential equations — describes the dynamics of these three compartments in an average sense, assuming every individual can potentially contact every other.

The ABM version relaxes the well-mixed assumption: agents occupy positions on a spatial grid, and transmission can only occur between neighboring agents. This single change produces qualitatively different dynamics: infection spreads as a spatial wave, local depletion of susceptibles can extinguish epidemics before they burn through the full population, and heterogeneous spatial distributions of hosts produce ragged, patchy outbreak patterns that no well-mixed equation can reproduce.

NoteConnecting to earlier material

The virulence-transmission fitness model from Frank (1996),

\[ w = z(1 - \alpha z), \]

describes the fitness of a pathogen strategy in a well-mixed population at evolutionary equilibrium. The ABM we build here operates at a shorter time scale — ecological rather than evolutionary — and introduces the spatial dimension that is absent from that analysis. A natural extension, explored in the exercises, is to allow pathogens to vary in their exploitation strategy \(z\), making the ABM the mechanistic substrate from which evolutionary outcomes like those Frank analyzes can emerge.

2 Setting up the environment

using Pkg
Pkg.add(["Agents", "CairoMakie", "DataFrames", "Random", "StatsBase"])
using Agents
using CairoMakie
using DataFrames
using Random
using StatsBase
NoteAgents.jl version

This chapter targets Agents.jl v6.x. The package has undergone significant API changes across major versions. Check your installed version with Pkg.status("Agents"). In v6, agent types are defined with the @agent macro, and the ABM constructor takes the agent type, space, and a properties NamedTuple or Dict.

3 Defining the agent

In Agents.jl, every agent type is defined using the @agent macro, which generates a Julia struct with a mandatory unique integer id and, for grid models, a mandatory pos field storing the agent’s grid coordinates. Additional fields carry the agent’s internal state.

For our SIR model, each agent is a host individual. Its internal state is just its epidemiological compartment:

@agent struct Host(GridAgent{2})
    status::Symbol    # :S, :I, or :R
    days_infected::Int
end

Breaking this down:

  • GridAgent{2} is the base type for agents living on a 2-dimensional grid. It provides the id::Int and pos::Tuple{Int,Int} fields automatically.
  • status::Symbol holds the agent’s compartment: :S, :I, or :R.
  • days_infected::Int tracks how many time steps the agent has been infected, which we use to implement recovery after a fixed infectious period.
TipWhy Symbol for status?

We use a Symbol (:S, :I, :R) rather than, say, an Int (0, 1, 2) for readability and because Julia dispatches on types efficiently. An alternative design would use an @enum, which is more robust for larger state spaces:

@enum EpiStatus S I R
@agent struct Host(GridAgent{2})
    status::EpiStatus
    days_infected::Int
end

Either works; the Symbol approach is more concise for a three-state model.

4 Model parameters and initialization

We collect the model’s parameters into a NamedTuple passed to the ABM constructor as properties. This makes parameters accessible inside stepping functions via the model object, and allows systematic parameter sweeps later.

function sir_step!(agent::Host, model)
    if agent.status == :S
        # Check neighbors for infection
        for neighbor in nearby_agents(agent, model, 1)
            if neighbor.status == :I
                if rand(model.rng) < model.infection_prob
                    agent.status = :I
                    agent.days_infected = 1
                    break   # one successful transmission event per step is sufficient
                end
            end
        end

    elseif agent.status == :I
        agent.days_infected += 1
        if agent.days_infected >= model.recovery_time
            agent.status = :R
        end

    end
    # :R agents require no action
    return nothing
end

function initialize_sir(;
    grid_size        = (50, 50),     # dimensions of the spatial grid
    n_initial_infected = 5,          # number of agents seeded as infected at t = 0
    infection_prob   = 0.06,         # probability of transmission per infected neighbor per step
    recovery_time    = 14,           # time steps until an infected agent recovers
    seed             = 42
)
    rng   = Random.Xoshiro(seed)
    space = GridSpace(grid_size; periodic = false, metric = :chebyshev)

    # Model properties accessible as model.infection_prob etc.
    properties = (
        infection_prob  = infection_prob,
        recovery_time   = recovery_time,
        rng             = rng
    )

    model = StandardABM(Host, space; properties = properties, rng = rng, scheduler = Schedulers.Randomly(), agent_step! = sir_step!)

    # Populate the grid: one agent per cell, all initially susceptible
    for pos in positions(model)
        add_agent!(pos, model, :S, 0)
    end

    # Seed a small number of infected agents at random positions
    initial_infected = sample(rng, collect(positions(model)), n_initial_infected; replace = false)
    for pos in initial_infected
        agent = model[ids_in_position(pos, model)[1]]
        agent.status = :I
        agent.days_infected = 1
    end

    return model
end
initialize_sir (generic function with 1 method)
TipSchedulers.Randomly()

The scheduler determines the order in which agents are stepped within each time step. Randomly shuffles the agent order randomly at every step. This is the standard choice for models where the order of interactions matters. Schedulers.Fastest() activates agents in arbitrary (fast, non-random) order and is useful when agent order is irrelevant, but is generally not appropriate for epidemic models where who moves first affects who infects whom.

Notemetric = :chebyshev

metric = :chebyshev defines neighborhood as the Moore neighborhood: the 8 immediately adjacent cells (including diagonals), as opposed to metric = :euclidean which gives the 4-cell Von Neumann neighborhood. The choice of neighborhood metric affects the effective transmission radius and the shape of the spreading wave.

5 The stepping function

The stepping function defines what each agent does during a single time step. It is called once per agent per step, in scheduler-determined order, and is the heart of the ABM.

Our SIR stepping function implements three rules:

  1. A susceptible agent checks its neighbors; for each infected neighbor, it becomes infected with probability infection_prob.
  2. An infected agent increments its days_infected counter; when this reaches recovery_time, it transitions to recovered.
  3. A recovered agent does nothing (permanent immunity in this model).
function sir_step!(agent::Host, model)
    if agent.status == :S
        # Check neighbors for infection
        for neighbor in nearby_agents(agent, model, 1)
            if neighbor.status == :I
                if rand(model.rng) < model.infection_prob
                    agent.status = :I
                    agent.days_infected = 1
                    break   # one successful transmission event per step is sufficient
                end
            end
        end

    elseif agent.status == :I
        agent.days_infected += 1
        if agent.days_infected >= model.recovery_time
            agent.status = :R
        end

    end
    # :R agents require no action
    return nothing
end
sir_step! (generic function with 1 method)
Tipnearby_agents and neighborhood radius

nearby_agents(agent, model, r) returns an iterator over all agents within Chebyshev distance r of the focal agent’s position. With r = 1 this gives the 8 immediate neighbors. Increasing r to 2 or 3 implements longer-range transmission, mimicking airborne or movement-mediated spread. Notably, increasing r does not change the stepping function logic at all — only the argument to nearby_agents — which illustrates one of the advantages of ABMs: spatial transmission range is a single, clearly interpretable parameter.

6 Running the model

6.1 Single step

model = initialize_sir(seed = 42)

# Advance the model by a single time step, applying sir_step! to every agent
step!(model, 1)
StandardABM with 2500 agents of type Host
 agents container: Dict
 space: GridSpace with size (50, 50), metric=chebyshev, periodic=false
 scheduler: Agents.Schedulers.Randomly
 properties: infection_prob, recovery_time, rng

6.2 Running for many steps and collecting data

run! is the workhorse function for multi-step simulation. It calls step! repeatedly and, at each step, records the quantities specified in adata (agent-level data) and mdata (model-level data).

Rather than recording every agent’s status at every step (which would produce a very large table), we aggregate immediately: we count the number of agents in each compartment.

# Aggregation functions: count agents with each status at each step
n_susceptible(agent)  = agent.status == :S
n_infected(agent)     = agent.status == :I
n_recovered(agent)    = agent.status == :R

# adata specifies which agent attributes to collect and how to aggregate them
adata = [
    (n_susceptible, count),
    (n_infected,    count),
    (n_recovered,   count)
]

# Re-initialize and run for 200 steps
model = initialize_sir(seed = 42)
agent_data, _ = run!(model, 200; adata)

first(agent_data, 8)
8×4 DataFrame
Row time count_n_susceptible count_n_infected count_n_recovered
Int64 Int64 Int64 Int64
1 0 2495 5 0
2 1 2492 8 0
3 2 2489 11 0
4 3 2488 12 0
5 4 2481 19 0
6 5 2476 24 0
7 6 2468 32 0
8 7 2457 43 0

The returned agent_data is a DataFrame with columns time, count_n_susceptible, count_n_infected, and count_n_recovered — one row per time step.

7 Visualizing the epidemic time series

fig = Figure(size = (800, 380))
ax  = Axis(
    fig[1, 1];
    xlabel    = "time step",
    ylabel    = "number of agents",
    title     = "SIR epidemic dynamics on a 50 × 50 spatial grid"
)

n_total = maximum(agent_data.count_n_susceptible .+
                  agent_data.count_n_infected    .+
                  agent_data.count_n_recovered)

lines!(ax, agent_data.time, agent_data.count_n_susceptible;
    color = :steelblue3, linewidth = 2, label = "Susceptible (S)")
lines!(ax, agent_data.time, agent_data.count_n_infected;
    color = :firebrick3, linewidth = 2, label = "Infected (I)")
lines!(ax, agent_data.time, agent_data.count_n_recovered;
    color = :seagreen,   linewidth = 2, label = "Recovered (R)")

axislegend(ax; position = :rc)
fig
NoteComparing to the well-mixed ODE SIR

A well-mixed ODE SIR model — one in which every susceptible is equally likely to contact every infected — typically produces a smooth, symmetric epidemic curve with a single sharp peak. The spatial ABM produces a qualitatively different shape: a slower initial rise (because infection can only spread to immediate neighbors), a broader, flatter peak (because local depletion of susceptibles impedes spread before the global susceptible pool is exhausted), and often a residual pool of susceptibles that escape infection entirely because they happen to be surrounded by recovered agents before the epidemic wave reaches them. These are not artifacts — they are biologically real phenomena.

8 Visualizing the spatial state

Agents.jl provides abmplot for visualizing the spatial state of a model at a single point in time. We color agents by their epidemiological status.

# Re-initialize model, step to a moment near peak infection for an interesting snapshot
model_snap = initialize_sir(seed = 42)
step!(model_snap, 40)   # step to approximately peak infection

# Color map: blue = S, red = I, green = R
status_color(agent) = agent.status == :S ? :steelblue3 :
                      agent.status == :I ? :firebrick3 :
                      :seagreen

fig_snap, ax_snap, abmobs = abmplot(
    model_snap;
    agent_color = status_color,
    agent_size  = 6,
    figure      = (; size = (520, 520)),
    axis        = (; title = "Spatial state at step 40 (near peak infection)")
)

fig_snap
TipSpatial wave structure

In the spatial snapshot, you will typically see a visually coherent wavefront of infected agents (red) separating a core of recovered agents (green, in the interior where the infection originated) from the still-susceptible outer region (blue). This wave structure — absent from any well-mixed model — arises because local transmission creates spatial autocorrelation in infection status: infected agents are surrounded by agents they have already infected or will soon infect. The wave speed is a function of infection_prob and recovery_time, and can be derived analytically for simple spatial SIR models.

8.1 Producing a multi-panel snapshot sequence

snapshots = [5, 20, 40, 80, 120, 200]

fig_seq = Figure(size = (1050, 700))

for (i, t) in enumerate(snapshots)
    row = div(i - 1, 3) + 1
    col = mod(i - 1, 3) + 1

    m = initialize_sir(seed = 42)
    step!(m, t)

    ax = Axis(fig_seq[row, col]; title = "t = $t", aspect = DataAspect())
    hidedecorations!(ax)

    # Collect positions and colors for all agents
    xs = [a.pos[1] for a in allagents(m)]
    ys = [a.pos[2] for a in allagents(m)]
    cs = [status_color(a) for a in allagents(m)]

    scatter!(ax, xs, ys; color = cs, markersize = 4)
end

fig_seq

9 Parameter sweeps with paramscan

One of the central questions in epidemic modeling is how the final size of an epidemic — the total fraction of the population that is ever infected — depends on the transmission probability \(\beta\) (here infection_prob) and the infectious period (here recovery_time).

Agents.jl provides paramscan to automate multi-run simulations across parameter grids.

# Define the parameter grid
parameters = Dict(
    :infection_prob  => [0.02, 0.04, 0.06, 0.10],
    :recovery_time   => [7, 14, 21],
    :seed            => [1, 2, 3]   # three replicate runs per parameter combination
)

scan_data, _ = paramscan(
    parameters, initialize_sir;
    n           = 200,
    adata       = [(n_recovered, count)]
)

# Extract the final epidemic size (recovered count at time 200) per parameter combination
final_sizes = combine(
    groupby(
        filter(:time => t -> t == 200, scan_data),
        [:infection_prob, :recovery_time]
    ),
    :count_n_recovered => mean => :mean_recovered,
    :count_n_recovered => std  => :std_recovered
)

final_sizes
12×4 DataFrame
Row infection_prob recovery_time mean_recovered std_recovered
Float64 Int64 Float64 Float64
1 0.02 7 14.3333 4.93288
2 0.02 14 371.667 138.857
3 0.02 21 1473.0 525.963
4 0.04 7 277.333 209.381
5 0.04 14 2405.0 60.6218
6 0.04 21 2472.0 31.225
7 0.06 7 2205.33 17.6163
8 0.06 14 2492.0 3.60555
9 0.06 21 2500.0 0.0
10 0.1 7 2473.33 0.57735
11 0.1 14 2499.33 1.1547
12 0.1 21 2500.0 0.0

9.1 Visualizing the parameter sweep

fig_sweep = Figure(size = (720, 420))
ax_sweep  = Axis(
    fig_sweep[1, 1];
    xlabel    = "transmission probability (β)",
    ylabel    = "mean final epidemic size (recovered agents)",
    title     = "Final epidemic size across transmission and recovery parameter space"
)

recovery_times = sort(unique(final_sizes.recovery_time))
colors = [:steelblue3, :darkorange2, :seagreen]

for (color, rt) in zip(colors, recovery_times)
    subset = filter(:recovery_time => r -> r == rt, final_sizes)
    sort!(subset, :infection_prob)
    errorbars!(ax_sweep,
        subset.infection_prob, subset.mean_recovered, subset.std_recovered;
        color = (color, 0.5), whiskerwidth = 8)
    lines!(ax_sweep, subset.infection_prob, subset.mean_recovered;
        color, linewidth = 2, label = "recovery time = $rt steps")
    scatter!(ax_sweep, subset.infection_prob, subset.mean_recovered;
        color, markersize = 8)
end

axislegend(ax_sweep; position = :rb)
fig_sweep
NoteThe epidemic threshold

In a well-mixed SIR model, the basic reproduction number is

\[ R_0 = \frac{\beta}{\gamma}, \]

where \(\beta\) is the transmission rate and \(\gamma = 1/\tau\) is the recovery rate (\(\tau\) is the infectious period). An epidemic can only grow when \(R_0 > 1\). In the spatial ABM, the effective threshold differs because local depletion of susceptibles around infected agents reduces the effective transmission rate below what a well-mixed calculation would predict. The epidemic threshold in spatial models is therefore higher than in the equivalent well-mixed ODE — spatial structure is generally epidemic-inhibiting.

10 Extending the model: heterogeneous transmission

A powerful feature of the ABM framework is that individual variation is trivially representable: we simply add fields to the agent struct and let the stepping function use them.

Here we extend the model to allow each host to carry a heritable susceptibility trait \(s_i \in [0, 1]\), which scales its probability of becoming infected on contact with an infected neighbor. This models variation in immune competence, nutritional status, age, or any other factor that modifies individual susceptibility.

@agent struct VariableHost(GridAgent{2})
    status::Symbol
    days_infected::Int
    susceptibility::Float64    # individual-level modifier on infection probability
end

function variable_sir_step!(agent::VariableHost, model)
    if agent.status == :S
        for neighbor in nearby_agents(agent, model, 1)
            if neighbor.status == :I
                # Effective transmission probability is scaled by individual susceptibility
                if rand(model.rng) < model.infection_prob * agent.susceptibility
                    agent.status = :I
                    agent.days_infected = 1
                    break
                end
            end
        end
    elseif agent.status == :I
        agent.days_infected += 1
        if agent.days_infected >= model.recovery_time
            agent.status = :R
        end
    end
    return nothing
end

function initialize_variable_sir(;
    grid_size          = (50, 50),
    n_initial_infected = 5,
    infection_prob     = 0.08,
    recovery_time      = 14,
    susceptibility_mean = 0.8,
    susceptibility_sd  = 0.2,
    seed               = 42
)
    rng   = Random.Xoshiro(seed)
    space = GridSpace(grid_size; periodic = false, metric = :chebyshev)
    properties = (infection_prob = infection_prob, recovery_time = recovery_time, rng = rng)

    model = StandardABM(VariableHost, space; properties = properties, rng = rng, scheduler = Schedulers.Randomly(), agent_step! = variable_sir_step!)

    for pos in positions(model)
        # Draw susceptibility from a truncated normal: clamp to [0, 1]
        s = clamp(susceptibility_mean + susceptibility_sd * randn(rng), 0.0, 1.0)
        add_agent!(pos, model, :S, 0, s)
    end

    initial_infected = sample(rng, collect(positions(model)), n_initial_infected; replace = false)
    for pos in initial_infected
        agent = model[ids_in_position(pos, model)[1]]
        agent.status = :I
        agent.days_infected = 1
    end

    return model
end
initialize_variable_sir (generic function with 1 method)
TipWhy individual variation matters

Models with heterogeneous susceptibility routinely produce lower final epidemic sizes than equivalent homogeneous models with the same mean transmission probability, because the most susceptible individuals tend to be infected (and subsequently removed from the susceptible pool) early, leaving a residual population that is disproportionately resistant. This effect — sometimes called susceptibility herd immunity — is one reason that epidemic forecasts based on homogeneous mean-field models consistently overestimate final epidemic sizes in real outbreaks.

11 Exercises


Exercise 1

Exploring neighborhood radius and transmission range

The nearby_agents(agent, model, r) call in sir_step! uses radius r = 1, giving the 8 immediately adjacent cells.

  1. Modify sir_step! to accept the neighborhood radius as a model property (model.transmission_radius), and update initialize_sir to include this property.
  2. Run the model for radii \(r \in \{1, 2, 3, 5\}\) with the same infection_prob = 0.04 and recovery_time = 14.
  3. For each radius, plot the epidemic time series and record the peak infected count and the final recovered count.
  4. How does increasing transmission radius change the shape of the epidemic curve and the final epidemic size? Explain the mechanism in terms of the spatial depletion argument.

Exercise 2

Demographic stochasticity: small populations and epidemic extinction

In very small populations, stochasticity can cause an epidemic to go extinct before it spreads widely — even when \(R_0 > 1\) under the mean-field calculation.

  1. Re-run the model on grids of size \(10 \times 10\), \(20 \times 20\), and \(50 \times 50\), using infection_prob = 0.08, recovery_time = 14, and n_initial_infected = 1.
  2. For each grid size, run 50 replicate simulations (varying seed from 1 to 50).
  3. For each replicate, record whether the epidemic “takes off” (final recovered count \(> 5\%\) of the population) or goes extinct.
  4. Plot the probability of epidemic establishment as a function of grid size.
  5. How does this relate to the concept of the basic reproduction number \(R_0\) and why does a threshold \(R_0 > 1\) not guarantee epidemic spread in finite populations?

Exercise 3

Spatial heterogeneity: a landscape with barriers

Real landscapes are not uniform grids. Physical barriers — rivers, mountain ranges, roads — interrupt dispersal and contact networks, fragmenting the epidemic wave.

  1. Modify initialize_sir to introduce a barrier: a vertical strip of cells spanning the full height of the grid at the horizontal midpoint, in which no agents are placed.
  2. Seed initial infections on the left half of the grid only.
  3. Run the model and produce a spatial snapshot sequence (as in the multi-panel plot above), observing how the barrier affects wave propagation.
  4. Compute and compare the final epidemic sizes on the left and right halves of the grid separately.
  5. Under what conditions does the epidemic cross the barrier? Modify the barrier width (number of empty columns) to find the approximate minimum barrier width that reliably prevents the epidemic from crossing in at least 80% of 20 replicate runs.

Exercise 4

Connecting to the Frank (1996) virulence model

The simple single-lineage pathogen fitness model (described in Frank (1996)), predicts that, in a well-mixed population, parasites evolve toward an intermediate exploitation severity \(z^* = \frac{1}{2\alpha}\) that maximizes pathogen fitness \(w = z(1 - \alpha z)\).

  1. Extend the VariableHost model so that each infected agent also carries a pathogen trait exploitation::Float64 (representing \(z\)), drawn at initialization from a uniform distribution over \([0.1, 0.9]\).
  2. Modify the stepping function so that infection_prob is replaced by a function of the pathogen’s exploitation severity: use \(\beta(z) = z\) (transmission increases with exploitation) and \(\gamma(z) = \alpha z\) (recovery rate, equivalently virulence, also increases with \(z\); use \(\alpha = 1.5\)), so that the infectious period becomes \(\tau(z) = \lceil 1 / (\alpha z) \rceil\) steps.
  3. When an agent is infected by a neighbor, it inherits the neighbor’s exploitation value (with optional small mutation: add Gaussian noise \(\mathcal{N}(0, 0.05)\), clamped to \([0.01, 1.0]\)).
  4. Run the model for 500 steps and record the mean exploitation value among infected agents over time.
  5. Does the mean exploitation value in the population converge toward the analytical prediction \(z^* = \frac{1}{2\alpha}\)? Plot the trajectory and compare to the predicted optimum.
Back to top

References

Frank, Steven A. 1996. “Models of Parasite Virulence.” The Quarterly Review of Biology 71 (1): 37–78. https://doi.org/10.1086/419267.
Basics of paleobiological fossil collection analyses
Basics of agent-based modeling: spatial epidemic dynamics with Agents.jl
  • © Jeet Sukumaran

Please share or adapt under the Creative Commons Attribution-NonCommercial-ShareAlike 4.0 International License (CC BY-NC-SA 4.0).