API

Agent

CellBasedModels.ABMType
mutable struct ABM

Basic structure which contains the user defined parmeters of the model, the user rules of the agents, both in high level definition and the functions already compiled.

Elements

FieldDescription
dims::IntDimensions of the model.
parameters::OrderedDict{Symbol,UserParameter}Dictionary of parameters of the model and its properties.
declaredUpdates::Dict{Symbol,Expr}Dictionary of updating rules and their user defined content (high level code).
declaredUpdatesCode::Dict{Symbol,Expr}Dictionary of updating rules after wrapping into code (low level code).
declaredUpdatesFunction::Dict{Symbol,Function}Dictionary of updting rules and the compiled functions (compiled code).
agentDEProblemODEProblem or SDEProblem object of Agent
agentAlgAlgorithm for the ODEProblem or SDEProblem of Agent
agentSolveArgsParameters for the ODEProblem or SDEProblem of Agent
modelDEProblemODEProblem or SDEProblem object of Model
modelAlgAlgorithm for the ODEProblem or SDEProblem of Model
modelSolveArgsParameters for the ODEProblem or SDEProblem of Model
mediumDEProblemODEProblem or SDEProblem object of Medium
mediumAlgAlgorithm for the ODEProblem or SDEProblem of Medium
mediumSolveArgsParameters for the ODEProblem or SDEProblem of Medium
neighborsAlgorithm to compute neighbors
platformPlatform in which to run the model
removalOfAgents_::BoolStores the information to check wether agents are removed in the code. Auxiliar parameter for generating the code.

Constructors

function ABM()

Generates an empty instance of ABM to be filled.

function ABM(
    dims;

    agent=OrderedDict{Symbol,DataType}(),
    agentRule::Expr=quote end,
    agentODE::Expr=quote end,
    agentSDE::Expr=quote end,

    model=OrderedDict{Symbol,DataType}(),
    modelRule::Expr=quote end,
    modelODE::Expr=quote end,
    modelSDE::Expr=quote end,

    medium=OrderedDict{Symbol,DataType}(),
    mediumRule::Expr=quote end,
    mediumODE::Expr=quote end,
    mediumSDE::Expr=quote end,

    baseModelInit::Vector{ABM}=ABM[],
    baseModelEnd::Vector{ABM}=ABM[],

    agentAlg::Union{CustomIntegrator,DEAlgorithm} = CBMIntegrators.Euler(),
    agentSolveArgs::Dict{Symbol,Any} = Dict{Symbol,Any}(),

    modelAlg::Union{CustomIntegrator,DEAlgorithm} = CBMIntegrators.Euler(),
    modelSolveArgs::Dict{Symbol,Any} = Dict{Symbol,Any}(),

    mediumAlg::Union{CustomIntegrator,DEAlgorithm} = DifferentialEquations.AutoTsit5(DifferentialEquations.Rosenbrock23()),
    mediumSolveArgs::Dict{Symbol,Any} = Dict{Symbol,Any}(),

    neighborsAlg::Neighbors = CBMNeighbors.Full(),       
    platform::Platform = CPU(),     
)

Generates an agent based model with defined parameters and rules.

ArgumentDescription
ArgsdimsDimensions of the system.
KwArgsagent=OrderedDict{Symbol,DataType}()Agent parameters
agentRule::Expr=quote endAgent rules
agentODE::Expr=quote endAgent Ordinary Differential Equations definition
agentSDE::Expr=quote endAgent Stochastic Differential Equations term definition
model=OrderedDict{Symbol,DataType}()Model parameters
modelRule::Expr=quote endModel rules
modelODE::Expr=quote endModel Ordinary Differential Equations definition
modelSDE::Expr=quote endModel Ordinary Differential Equations definition
medium=OrderedDict{Symbol,DataType}()Medium parameters
mediumRule::Expr=quote endMedium rules
mediumODE::Expr=quote endMedium Ordinary Differential Equations definition
mediumSDE::Expr=quote endMedium Ordinary Differential Equations definition
baseModelInit::Vector{ABM}=ABM[]ABM model whose rules will act before this ABM rules
baseModelEnd::Vector{ABM}=ABM[]ABM model whose rules will act after this ABM rules

For a more extense explanation of how to define rules and parameters, read Usage in the documentation.

source

Macros

These are macros that can be used in the code.

AgentRule

CellBasedModels.@loopOverNeighborsMacro
macro loopOverNeighbors(code)
macro loopOverNeighbors(it1, code)

Macro that creates the loop function to go over all neighbors of the agent.

It can be declared as

@loopOverNeighbors for iterator in ___
    #code
end

or

@loopOverNeighbors iterator begin
    #code
end

for some iterator symbol. It can only be used in agent rules or DEs.

source

Medium

CellBasedModels.@∂Macro
macro ∂(coord,code)

Discretizes the code term of a drift process. e.g. Agent in 2D

∂(1,m) → (m[i1-1,i2]-m[i1-1,i2])/(2*dx)

the coordinate must be 1, 2 or 3, specifing the axis of differentiation.

source
CellBasedModels.@∂2Macro
macro ∂2(coord,code)

Discretizes the code term of a drift process. e.g. Agent in 2D

∂2(1,m) → (m[i1-1,i2]-2*m[i1,i2]+m[i1-1,i2])/dx^2

the coordinate must be 1, 2 or 3, specifing the axis of differentiation.

source
CellBasedModels.@mediumBorderMacro
mediumBorder(coord, border)

Macro that returns true if mesh position is in the lower (border=-1) or upper (border=1) border of the axis coordinate 1, 2 or 3.

source

Community

CellBasedModels.CommunityType

Basic structure keeping the parameters of all the agents in the current simulation of a model.

Parameters to store information essential for Community simulation

SymbolDescription
abmABM model of the community
uuidUnique identifier of the community
loadedIf loaded into the platform CPU or GPU
pastTimes::Array{Community}Store times when called to saveRAM!
parameters::OrderedDict{Symbol,AbstractArray}Dictionary with the User Defined Parameters
agentDEProblemODEProblem or SDEProblem object of Agent
modelDEProblemODEProblem or SDEProblem object of Model
mediumDEProblemODEProblem or SDEProblem object of Medium

Parameters seen in the kernels that may de used directly by the user

SymbolDescription
tTime of the community
dtStepping time of the community
NNumber of agents
NMediumSize of medium mesh
simBoxSize of simulation box
dxSize of the axis 1 in mesh medium
dySize of the axis 2 in mesh medium
dzSize of the axis 3 in mesh medium

Parameters seen in the kernels that are for internal use

SymbolDescription
idIdentification of the agent
nMax_Maximum number of agents when loading to platform
idMax_Maximum id in the Community at all times
NAdd_Number of added agents in the present step
NRemove_Number of removed agents in the present step
NSurvive_Number of agents survived in this step
flagSurvive_0 is agent survived this step, 1 if dead
holeFromRemoveAt_Holes left from agents that are dead
repositionAgentInPos_List the positions of the agents that have to be reubicated to empty spaces
flagRecomputeNeighbors_1 is neighbors have to be recomputed

Constructors

function Community()

Creates fully empty community. Auxiliar method for the following method of declaration for the users.

function Community(
    abm::ABM; 

    dt::Union{Nothing,AbstractFloat} = nothing,
    t::AbstractFloat = 0.,
    N::Int = 1,
    id::AbstractArray{Int} = 1:N,
    NMedium::Union{Nothing,Vector{<:Int}} = nothing,
    simBox::Union{Nothing,Matrix{<:Number}} = nothing,

    args...
)

Function to construct a Community that accepts to provide integration algorithms, neighboring algorithms, the computing platform and setting parameters of the community.

For a more specific indication of the usage see the UserGuide.

source

Platform loading

CellBasedModels.loadToPlatform!Function
function loadToPlatform!(com::Community;preallocateAgents::Int=0)

Function that converts the Community data into the appropiate format to be executed in the corresponding platform. It locks the possibility of accessing and manipulating the data by indexing and propety. If preallocateAgents is provided, it allocates that additional number of agents to the community. Preallocating is necesssary as if more agents will be added during the evolution of the model (Going over the number of preallocated agents will run into an error.).

source
CellBasedModels.bringFromPlatform!Function
function bringFromPlatform!(com::Community)

Return the Community object from the platform where it is being evolved. It locks the possibility of accessing and manipulating the data by indexing and propety.

source

Evolution functions

CellBasedModels.evolve!Function
function evolve!(community;
    steps,saveEach=1,
    saveToFile=false,fileName=nothing,overwrite=false,
    saveCurrentState=false,
    preallocateAgents=0,
    progressMessage=(com)->nothing)

Performs step number of steps on the community, saving each saveEach number of steps the community instance using the saving function provided in saveFunction (See IO). If saveCurrentState is true, the present instance is saved.

preallocateAgents is an integer to sent to the loadToPlatform! function that allocates empty space for agents if the model has to grow. The maximum number of agents in the final simulation has to be specified in here.

progressMessage is a function that is executed every time we save an step and can be used to print some information of progress.

source
CellBasedModels.step!Function
function step!(community)

Executes all the possible step functions and updates the parameters a single time.

source
CellBasedModels.agentStepDE!Function
function agentStepDE!(community)

Function that computes a integration step of the community a time step dt using the defined Integrator defined in Agent.

source
CellBasedModels.modelStepDE!Function
function modelStepDE!(community)

Function that computes a integration step of the community a time step dt using the defined Integrator defined in Agent.

source
CellBasedModels.mediumStepDE!Function
function mediumStepDE!(community)

Function that computes a integration step of the community a time step dt using the defined Integrator defined in Agent.

source
CellBasedModels.update!Function
function update!(community)

Function that updates all the parameters after making a step in the community.

source

IO

CellBasedModels.saveRAM!Function
function saveRAM!(community::Community)

Function that stores the present configuration of community in the field Community.pastTimes.

source
CellBasedModels.saveJLD2Function
saveJLD2(file::String, community::Community; overwrite=false)

Save in file the current instance of the community. If some other community was being saved here, an error will raise. If overwrite=true is specified, it will remove the previous community.

source

Metrics

CellBasedModels.CBMMetrics.@euclideanMacro
macro euclidean(it2)
macro euclidean(it1,it2)

Macro that given the iterator symbos it1 and it2, give the corresponding euclidean distanve in the correct dimentions. If it1 is not provided it asumes the default iteration index of agents (i1_).

source
CellBasedModels.CBMMetrics.@manhattanMacro
macro manhattan(it2)
macro manhattan(it1,it2)

Macro that given the iterator symbos it1 and it2, give the corresponding manhattan distanve in the correct dimentions. If it1 is not provided it asumes the default iteration index of agents (i1_).

source
CellBasedModels.CBMMetrics.cellInMeshFunction
function cellInMesh(edge,x,xMin,xMax,nX)

Give the integer position in a regular discrete mesh with poits separamtions of edge, given a position in x. The simulation domain being (xMin, xMax) and maximum number of mesh points nX.

e.g Grid with 6 points at [0,.1,.2,.3,.4,.5]

>>> cellInMesh(.1,.29,0.,.5,6)
4

which if the closest point in the mesh.

source
CellBasedModels.CBMMetrics.intersection2linesFunction
function intersection2lines(x1,y1,theta1,x2,y2,theta2,inf_eff=100000)

Finds the point of intersection of two lines. You have to provide a point in space and and angle for eachline: (x1,y1,theta1) and (x2,y2,theta2).

If the lines are parallel, it returns a point effectively in the infinite. The effective distance is described by inf_eff`.

Returns the point of intersection.

source
CellBasedModels.CBMMetrics.point2lineFunction
function point2line(x1,y1,x2,y2,theta2)

Given a point (x1,x2), finds the closest point projected over a line described by a point in the line and the angle: (x2,y2,theta2).

Returns the coordinates of the closest point over the line.

source
CellBasedModels.CBMMetrics.pointInsideRodFunction
function pointInsideRod(x1,y1,l1,theta1,pxAux,pyAux,separation)

Given a line segment described by it central point (x1,y1), its angle in the plate theta1 and its length l1; and given a point over the save line (pxAux,pyAux), returns the point if inside the segment or the closes extreme of the segment. If provided a separation (0,1), it moves the point that separation.

Returns the coordinates of the closest point over the segment.

source
CellBasedModels.CBMMetrics.rodIntersectionFunction
function rodIntersection(x1,y1,l1,theta1,x2,y2,l2,theta2;separation=0.99)

Given a two line segment described by it central point (x1,y2), its angle in the plate theta1 and its length l1; finds the closest spheres of both segments.

Returns the coordinates of the closest spheres (x1Aux,y1Aux), (x2Aux,y2Aux).

source

Random

Random number generators from distributions that are compatible both in CPU and GPU.

Platform

CellBasedModels.GPUType

" mutable struct GPU <: Platform

CPU platform structure.

Threads and blocks for each rule method

  • agentThreads
  • agentBlocks
  • modelThreads
  • modelBlocks
  • mediumThreads
  • mediumBlocks
source

Integrators

Integrators can be used the ones from DifferentialEquations for ODE or SDE problems or the custom made solvers provided in ths package.

Neighbor algorithms

CellBasedModels.CBMNeighbors.VerletTimeType
mutable struct VerletTime <: Verlet <: Neighbors

Method that computes VerletList neighbors and updates it at fixed times.

Constructor

function VerletTime(;skin,dtNeighborRecompute,nMaxNeighbors)
ParameterDescription
KwArgsskin is the maximum distance around center to check neighbors.
dtNeighborRecompute is time step at which recompute Verlet Lists.
nMaxNeighbors is the maximum number of neighbors that a cell may have.
source
CellBasedModels.CBMNeighbors.VerletDisplacementType
mutable struct VerletTime <: Verlet <: Neighbors

Method that computes VerletList neighbors and updates it whenever an agent moves too far from the initial position

Constructor

function VerletDisplacement(;skin,nMaxNeighbors)
ParameterDescription
KwArgsskin is the maximum distance around center to check neighbors.
nMaxNeighbors is the maximum number of neighbors that a cell may have.
source
CellBasedModels.CBMNeighbors.CellLinkedType
mutable struct CellLinked <: Neighbors

Method that computes Cell Linked neighbors and updates it whenever an agent moves too far from the initial position

Constructor

function CellLinked(;cellEdge)
ParameterDescription
KwArgscellEdge is the grid size to use to check for neighbors around neighbor cells.
source
CellBasedModels.CBMNeighbors.CLVDType
mutable struct CLVD <: Verlet <: Neighbors

Method that computes Cell Linked and Verlet Displacement neighbors algorithms together to compute neighbors and only recompute when they left very far from center.

Constructor

function CLVD(;skin,nMaxNeighbors,cellEdge)
ParameterDescription
KwArgsskin is the maximum distance around center to check neighbors.
nMaxNeighbors is the maximum number of neighbors that a cell may have.
cellEdge is the grid size to use to check for neighbors around neighbor cells.
source

Fitting

CellBasedModels.CBMFitting.gridSearchFunction
function gridSearch(evalFunction::Function,  
    searchList::Dict{Symbol,<:Vector{<:Number}}; 
    returnAll::Bool=false,
    saveFileName::Union{Nothing,String} = nothing,
    args::Vector{<:Any} = Any[],
    verbose=false)

Function that evaluates a grid of parameter configurations for a model.

ParameterDescription
ArgsevalFunction:: Function** : Function that takes a DataFrame with parameters, generates the simulations and returns a score of the fit.
searchList::Dict{Symbol,<:Tuple{<:Number,<:Number}}}** : Dictionary of parameters and the ranges of exloration the parameters (e.g. :x => (0,1)).
KwArgsreturnAll::Bool = false** : If return the hole list of parameters explored or the just the most fit.
saveFileName::Union{Nothing,String} = nothing** : If given a string, it saves the parameters explored in a file with the corresponding name.
args::Vector{<:Any} = Any[]** : Additional arguments to give to evalFunction.
source
CellBasedModels.CBMFitting.swarmAlgorithmFunction
function swarmAlgorithm(evalFunction::Function,  
    searchList::Dict{Symbol,<:Union{<:Tuple{<:Number,<:Number},Vector{<:Number}}};
    population::Int=100,
    weightInertia::Number = .1,
    weightGlobalBest::Number = .1, 
    weightPopulationBest::Number = .1, 
    stopMaxGenerations::Int = 10,
    initialisation::Union{Nothing,DataFrame} = nothing,
    returnAll::Bool = false,
    saveFileName::Union{Nothing,String} = nothing,
    args::Vector{<:Any} = Any[],
    verbose=false)

Optimization of the parameter space of a model that uses the Swarm Algorithm.

ParameterDescription
ArgsevalFunction:: FunctionFunction that takes a DataFrame with parameters, generates the simulations and returns a score of the fit.
searchList::Dict{Symbol,<:Tuple{<:Number,<:Number}}}Dictionary of parameters and the ranges of exloration the parameters (e.g. :x => (0,1)).
KwArgspopulation::Int = 100Size of the colony used at each generation for the optimization.
weightInertia::Number = .1Hyperparameter of the colony weighting the current velocity.
weightGlobalBest::Number = .1Hyperparameter of the colony weighting the global best solution.
weightPopulationBest::Number = .1Hyperparameter of the colony weighting the local best solution.
stopMaxGenerations::Int = 10How many generations do before stopping the algorithm.
initialisation::Union{Nothing,DataFrame} = nothingDataFrame defining the initial parameters of the population. If nothing, they are set randomly.
returnAll::Bool = falseIf return the hole list of parameters explored or the just the most fit.
saveFileName::Union{Nothing,String} = nothingIf given a string, it saves the parameters explored in a file with the corresponding name.
args::Vector{<:Any} = Any[]Additional arguments to give to evalFunction.
verbose=falseIf true, show progress bar during optimization.
source
CellBasedModels.CBMFitting.beeColonyAlgorithmFunction
function beeColonyAlgorithm(
    evalFunction::Function, 
    searchList::Dict{Symbol,<:Tuple{<:Number,<:Number}}; 
    population::Int=100,
    limitCycles::Int = 10,
    stopMaxGenerations::Int = 100,
    returnAll::Bool = false,
    initialisation::Union{Nothing,DataFrame} = nothing,
    saveFileName::Union{Nothing,String} = nothing,
    args::Vector{<:Any} = Any[],
    verbose=false)

Optimization of the parameter space of a model that uses the Bee Colony Algorithm.

ParameterDescription
ArgsevalFunction:: FunctionFunction that takes a DataFrame with parameters, generates the simulations and returns a score of the fit.
searchList::Dict{Symbol,<:Tuple{<:Number,<:Number}}}Dictionary of parameters and the ranges of exloration the parameters (e.g. :x => (0,1)).
KwArgspopulation::Int=100Size of the colony used at each generation for the optimization.
limitCycles::Int = 10Hyperparameter of the algorithm that says how many generations without update are waited until jump to other position.
stopMaxGenerations::Int = 100How many generations do before stopping the algorithm.
returnAll::Bool = falseIf return the hole list of parameters explored or the just the most fit.
initialisation::Union{Nothing,DataFrame} = nothingDataFrame defining the initial parameters of the population. If nothing, they are set randomly.
saveFileName::Union{Nothing,String} = nothingIf given a string, it saves the parameters explored in a file with the corresponding name.
args::Vector{<:Any} = Any[]Additional arguments to give to evalFunction.
verbose=falseIf true, show progress bar during optimization.
source
CellBasedModels.CBMFitting.geneticAlgorithmFunction
function geneticAlgorithm(
    evalFunction::Function,  
    searchList::Dict{Symbol,<:Union{<:Tuple{<:Number,<:Number},Vector{<:Number}}};
    population::Int=100,
    parentSelectionAlg::String = "weighted", #weighted or random
    parentSelectionP::Number = .1, 
    mutationRate::Number = .1, 
    stopMaxGenerations::Int = 10,
    initialisation::Union{Nothing,DataFrame} = nothing,
    returnAll::Bool = false,
    saveFileName::Union{Nothing,String} = nothing,
    args::Vector{<:Any} = Any[],
    verbose=false)

Optimization of the parameter space of a model that uses the Particle Swarm Algorithm.

ParameterDescription
ArgsevalFunction:: FunctionFunction that takes a DataFrame with parameters, generates the simulations and returns a score of the fit.
searchList::Dict{Symbol,<:Tuple{<:Number,<:Number}}}Dictionary of parameters and the ranges of exloration the parameters (e.g. :x => (0,1)).
KwArgspopulation::Int=100Size of the colony used at each generation for the optimization.
parentSelectionAlg::String = "weighted"Weigthing method of the population ot chose descendants.
parentSelectionP::Number = .1Hyperparameter of the algorithm indicating the proportion of parameters exchanged between parents.
mutationRate::Number = .1Hyperparameter of the algorithm indicating the probability of resampling the parameter with a uniform.
stopMaxGenerations::Int = 10How many generations do before stopping the algorithm.
initialisation::Union{Nothing,DataFrame} = nothingDataFrame defining the initial parameters of the population. If nothing, they are set randomly.
returnAll::Bool = falseIf return the hole list of parameters explored or the just the most fit.
saveFileName::Union{Nothing,String} = nothingIf given a string, it saves the parameters explored in a file with the corresponding name.
args::Vector{<:Any} = Any[]Additional arguments to give to evalFunction.
verbose=falseIf true, show progress bar during optimization.
source

Plotting

CellBasedModels.CBMPlots.plotRods2D!Function
function plotRods2D!(ax, x, y, d, l, angle; kargs...)

Plot rod shape like cells using Makie functions.

ParameterDescription
ArgsaxAxis where to plot the rods
xCoordinates of rods in x
yCoordinates of rods in y
dDiameter of the rod
lLength of the rods
angleAngle of the rods of rods in the XY plane
KwArgsAll arguments that want to be passed to meshscatter! function as color etc...
function plotRods2D!(ax, x, y, xs1, ys1, xs2, ys2, markerSphere, markerCylinder, angle; kargs...)

Plot rod shape like cells using Makie functions. This method is useful to make videos of cells.

ParameterDescription
ArgsaxAxis where to plot the rods
xCoordinates of center of rod in x
yCoordinates of center of rod in y
xs1Coordinates of rod extreme in x
ys1Coordinates of rod extreme in y
xs2Coordinates of other rod extreme in x
ys2Coordinates of other rod extreme in y
markerSpherePoint3f0 describing the radius of the sphere in (x,y,z)
markerCylinderPoint3f0 describing the cylinder sizes in (l,rx,ry)
KwArgsAll arguments that want to be passed to meshscatter! function as color etc...
source