Introduction

Structure

An AgentBasedModel object embeds all the properties of the system to be simulated and maps unique IDs to microbe instances. During the simulation, the model is evolved in discrete time steps, with each microbe's position, velocity and "state" being updated according to specified rules. Standard rules for motion, reorientations and chemotaxis are available by default, but custom behaviors can be implemented via user-defined functions.

The typical workflow to run a simulation in MicrobeAgents.jl goes as follows:

  1. Define the size and properties of the space in which the microbes will move.
  2. Choose an appropriate microbe type to represent the desired behavior, or define a new one.
  3. Initialize an AgentBasedModel object with the desired space, microbe type, integration timestep, and any extra property needed for the simulation.
  4. Populate the ABM with microbe instances.
  5. Choose the observables to collect during production and run the model.

MicrobeAgents.jl re-exports and extends various function from Agents.jl in order to work as a standalone, but it is generally recommended to use it in combination with Agents.jl for extra goodies.

Space

MicrobAgents.jl only supports continuous spaces with dimensions 1, 2 or 3. Spaces can be created with the ContinuousSpace function (reexported from Agents.jl). The extent of the space must be given as a tuple or SVector, and periodicity is set with the periodic kwarg (defaults to true).

# one-dimensional periodic space
extent = (1.0,)
ContinuousSpace(extent)

# two-dimensional non-periodic space
extent = (1.0, 2.0)
ContinuousSpace(extent; periodic=false)

# three-dimensional space periodic only along the x direction
extent = (100.0, 20.0, 20.0)
ContinuousSpace(extent; periodic=(true,false,false))

Microbes

Microbes are represented by subtypes of the AbstractMicrobe type, which is in turn a subtype of AbstractAgent from Agents.jl

MicrobeAgents.AbstractMicrobeType
AbstractMicrobe{D} <: AbstractAgent where {D<:Integer}

All microbe types in MicrobeAgents.jl simulations must be instances of user-defined types that are subtypes of AbstractMicrobe. YourMicrobeType{D} <: AbstractMicrobe{D} The parameter D defines the dimensionality of the space in which the microbe type lives (1, 2 and 3 are supported).

All microbe types must have at least the following fields:

  • id::Int id of the microbe (used internally by Agents.jl)
  • pos::SVectpr{D,Float64} position of the microbe
  • vel::SVector{D,Float64} velocity of the microbe
  • motility::AbstractMotility motile pattern of the microbe
  • turn_rate::Real average reorientation rate of the microbe
  • rotational_diffusivity::Real coefficient of brownian rotational diffusion
  • radius::Real equivalent spherical radius of the microbe
  • state::Real generic variable for a scalar internal state
source

MicrobeAgents provides different AbstractMicrobe subtypes representing different models of bacterial behavior from the scientific literature. The list of implemented models can be obtained with subtypes(AbstractMicrobe).

A basic type, which is typically sufficient for simple motility simulations and does not include chemotaxis, is the Microbe type.

MicrobeAgents.MicrobeType
Microbe{D,N} <: AbstractMicrobe{D,N}

Base microbe type for simple simulations. D is the space dimensionality and N is the number of states of the microbe's motility pattern.

Microbe has the required fields

  • id::Int an identifier used internally
  • pos::SVector{D,Float64} spatial position
  • vel::SVector{D,Float64} unit velocity vector
  • speed::Float64 magnitude of the velocity vector
  • motility::Motility{N} motility pattern

and the default parameters

  • rotational_diffusivity::Float64 = 0.0 coefficient of brownian rotational diffusion
  • radius::Float64 = 0.0 equivalent spherical radius of the microbe
  • state::Float64 = 0.0 generic variable for a scalar internal state
source

Microbe instances should only be created within an AgentBasedModel, the fundamental structure which embeds everything that has to do with the agent-based simulations you want to run. In MicrobeAgents.jl, models are created through the StandardABM function.

Agents.StandardABMType
StandardABM(MicrobeType, space, timestep; kwargs...)

Extension of the Agents.StandardABM method for microbe types. Implementation of AgentBasedModel where agents can be added and removed at any time. If agents removal is not required, it is recommended to use the keyword argument container = Vector for better performance. See Agents.AgentBasedModel for detailed information on the keyword arguments.

Arguments

  • MicrobeType: subtype of AbstractMicrobe{D}, with explicitly specified dimensionality D. A list of available options can be obtained by running subtypes(AbstractMicrobe).
  • space: a ContinuousSpace{D} with the same dimensionality D as MicrobeType which specifies the spatial properties of the simulation domain.
  • timestep: the integration timestep of the simulation.

Keywords

  • properties: additional container of data to specify model-level properties. MicrobeAgents.jl includes a set of default properties (detailed at the end).
  • scheduler = Schedulers.fastest
  • rng = Random.default_rng()
  • warn = true

Default properties

When a model is created, a default set of properties is included in the model (MicrobeAgents.default_ABM_properties):

DEFAULT_ABM_PROPERTIES = Dict(
    :chemoattractant => GenericChemoattractant{D,Float64}()
    :affect! => chemotaxis!
)

By including these default properties, we make sure that all the chemotaxis models will work even without extra user intervention. All these properties can be overwritten by simply passing an equivalent key to the properties dictionary when creating the model.

source

To initialize a model we must provide the microbe type, the simulation space, and the integration timestep of the simulation. All other parameters are optional. To setup a model for Microbes living in a 1-dimensional space we can therefore run

space = ContinuousSpace((100.0,); periodic=false)
dt = 0.1
model = StandardABM(Microbe{1}, space, dt)

Now, with the add_agent! function we will populate the model with microbes of the specified type (Microbe{1}). The only argument we must always specify for add_agent! is the motility of the microbe, via the motility keyword. An overview of the motility interface is given later; For now we will just use a Run-Tumble motility with an average run duration of 1 second and a constant swimming speed of 20 micron / second. The third required argument to the call is the distribution of reorientation angles, but it is irrelevant in 1-dimensional systems so we can just pass an empty array. If unspecified, position, direction and speed of the microbe will be assigned randomly; all the other fields will be assigned default values from the constructor (unless specified). To select a position, it can be passed as the first argument to the add_agent! call, and any other bacterial parameter can be defined via keyword arguments. All of the following are valid calls

motility = RunTumble(1.0, [20.0], [])
# a Microbe with large radius
add_agent!(model; motility, radius=10.0)
# a Microbe with custom position and high coefficient of rotational diffusion
add_agent!((53.2,), model; motility, rotational_diffusivity=0.5)
# a Microbe initialized with velocity to the right
add_agent!(model; motility, vel=(1.0,))

All the other subtypes of AbstractMicrobe work in a similar way, although they will have distinct default values and extra fields. When possible, default values are typically assigned following the original implementation in the literature.

MicrobeAgents.BrownBergType
BrownBerg{D} <: AbstractMicrobe{D}

Model of chemotactic E.coli from 'Brown and Berg (1974) PNAS'

Default parameters:

  • motility = RunTumble(0.67, [30.0], Isotropic(D), 0.1)
  • rotational_diffusivity = 0.035 rad²/s coefficient of brownian rotational diffusion
  • radius = 0.5 μm equivalent spherical radius of the microbe
  • state = 0.0 corresponds to 'weighted dPb/dt' in the paper
  • gain = 660 s
  • receptor_binding_constant = 100 μM
  • memory = 1 s
source
MicrobeAgents.BrumleyType
Brumley{D} <: AbstractMicrobe{D}

Model of chemotactic bacterium from 'Brumley et al. (2019) PNAS'. The model is optimized for simulation of marine bacteria and accounts for the presence of (gaussian) sensing noise in the chemotactic pathway.

Default parameters:

  • motility = RunReverseFlick(0.45, [46.5], 0.45, [46.5])
  • state = 0.0 → 'S'
  • rotational_diffusivity = 0.035 rad²/s
  • memory = 1.3 s → 'τₘ'
  • gain_receptor = 50.0 μM⁻¹ → 'κ'
  • gain = 50.0 → 'Γ'
  • chemotactic_precision = 6.0 → 'Π'
  • radius = 0.5 μm → 'a'
source
MicrobeAgents.CelaniType
Celani{D} <: AbstractMicrobe{D}

Model of chemotactic bacterium using the response kernel from 'Celani and Vergassola (2010) PNAS', extracted from experiments on E. coli.

Sensing noise (not present in the original model) is customarily introduced through the molecular counting noise formula by Berg and Purcell, and can be tuned through a chemotactic_precision factor inspired by 'Brumley et al. (2019) PNAS' (defaults to 0, i.e. no noise).

Default parameters:

  • `motility = RunTumble(0.67, [30.0], Isotropic(D), 0.1)
  • state = 0
  • rotational_diffusivity = 0.26 rad²/s
  • gain = 50.0
  • memory = 1 s
  • radius = 0.5 μm
source
MicrobeAgents.XieType
Xie{D} <: AbstractMicrobe{D}

Model of chemotactic bacterium adapted from 'Xie et al. (2019) Biophys J'. The model is developed based on experimental measurements of the chemotactic response function in the marine bacterium V. alginolyticus. The peculiarity of the model is the presence of distinct parameters for the forward and backward swimming states.

Sensing noise (not present in the original model) is customarily introduced through the molecular counting noise formula by Berg and Purcell, and can be tuned through a chemotactic_precision factor inspired by 'Brumley et al. (2019) PNAS' (defaults to 0, i.e. no noise).

Default parameters:

  • motility
  • turn_rate_forward = 2.3 Hz
  • turn_rate_backward = 1.9 Hz
  • state = 0.0 s
  • state_m = 0.0 s
  • state_z = 0.0 s
  • rotational_diffusivity = 0.26 rad²/s
  • adaptation_time_m = 1.29 s
  • adaptation_time_z = 0.28 s
  • gain_forward = 2.7 1/s
  • gain_backward = 1.6 1/s
  • binding_affinity = 0.39 μM
  • chemotactic_precision = 0.0
  • radius = 0.5 μm
source

More about models

MicrobeAgents.jl exploits the AgentBasedModel interface from Agents.jl. While the standard Agents.jl syntax will always work, it is typically more convenient to use the method extensions provided by MicrobeAgents.jl, which also includes some default parameters required by the simulations. If the simulation requires removal/addition of microbes, it is recommended to call StandardABM with the container=Dict keyword argument, otherwise MicrobeAgents.jl defaults to container=Vector which provides better performance.

In addition to the microbe instances, the model should also wrap all the other information required to perform the simulation.

MicrobeAgents.jl defines default timestepping functions which are used to evolve the microbes and the model, and are accessible through the microbe_step! and model_step! keywords in StandardABM. By default, the microbe_step! function performs, in order:

  • update microbe position according to current velocity
  • randomize the microbe orientation through rotational diffusion (if present)
  • update internal state of the microbe (e.g. chemotaxis or other user-defined behavior)
  • perform reorientation events following Poissonian statistics

The model_step! function instead defaults to a dummy function which does nothing. Any custom behavior can be implemented by simply modifying these two functions.

Any type of external parameter that should be used during the simulation should be passed to StandardABM through the properties dictionary.

Running a model

After the model has been created and populated with the desired number of microbes, we are ready to run the simulation. We just need to specify how many steps we want to simulate and what data to collect during the run:

nsteps = 100
adf, mdf = run!(model, nsteps; adata=[position])

run! will return two dataframes, one for the agent-level data (adf) and one for the model-level data (mdf, which in this case will be empty). This way, we have produced our first random walk. Since adf.position is a vector of static vectors, we first have to unpack the x and y values and then we are ready to plot our trajectory.

using Plots
x = first.(adf.position)
y = last.(adf.position)
plot(x, y)

Motility patterns

In MicrobeAgents.jl, motility patterns are represented through the Motility{N} type. A Motility is composed of N instances of MotileState and of a set of transition probabilities between these N states. A MotileState contains information about the speed distribution, turn angle distributions, and average lifetime of a particular state of motion.

There are two kinds of MotileStates, RunState and TurnState. A RunState is used to represent states associated with translational motion and no angular reorientation (e.g. Run). A TurnState is instead used to represent states where no translational motion happens, but where reorientations may occur (e.g. Tumble, Reverse, Flick, Stop).

Instances of MotileStates can be arbitrarily combined into a Motility.

MicrobeAgents.MotilityType
Motility(motile_states::NTuple{N,MotileState}, rates...)

Construct a Motility from a sequence of MotileStates and their transition rates.

motile_states must be a tuple of MotileState instances. rates must be specified in the form (i => j, p) where i is the index of the starting state (as ordered in the motile_states tuple), j is the index of the state towards which the transition occurs, and p is the probability that, after the duration of i is over, the transition occurs towards the state j. It is not necessary to indicate impossible transitions; any pair i => j which is not explicitly specified is assumed to have a transition probability p=0.

For example, if we want to define a 3-state motility, composed by a slow but long run, an instantaneous tumble with some angle_pdf, and a fast but short run, where the two runs can only transition towards the tumble, but the tumble can transition with probability 30% towards the slow run and 70% towards the fast run, we will call:

Motility((Run(10.0, 2.0), Tumble(0.0, angle_pdf), Run(60.0, 0.5)),
    (1 => 2, 1.0),
    (2 => 1, 0.3), (2 => 3, 0.7),
    (3 => 2, 1.0)
)

Some default constructors for common motility patterns are provided: RunTumble, RunReverse, RunReverseFlick, RunStop

source

The most common motility patterns are, however, pre-implemented.

MicrobeAgents.RunReverseFunction
RunReverse(
    run_speed_forward, run_duration_forward,
    run_speed_backward, run_duration_backward;
    angle=(π,), reverse_duration=0.0,
source
MicrobeAgents.RunReverseFlickFunction
RunReverseFlick(
    run_speed_forward, run_duration_forward,
    run_speed_backward, run_duration_backward;
    reverse_angle=(π,) flick_angle=(-π/2, +π/2),
    reverse_duration=0.0, flick_duration=0.0,
)
source