API
Agent
CellBasedModels.ABM
— Typemutable 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
Field | Description |
---|---|
dims::Int | Dimensions 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). |
agentDEProblem | ODEProblem or SDEProblem object of Agent |
agentAlg | Algorithm for the ODEProblem or SDEProblem of Agent |
agentSolveArgs | Parameters for the ODEProblem or SDEProblem of Agent |
modelDEProblem | ODEProblem or SDEProblem object of Model |
modelAlg | Algorithm for the ODEProblem or SDEProblem of Model |
modelSolveArgs | Parameters for the ODEProblem or SDEProblem of Model |
mediumDEProblem | ODEProblem or SDEProblem object of Medium |
mediumAlg | Algorithm for the ODEProblem or SDEProblem of Medium |
mediumSolveArgs | Parameters for the ODEProblem or SDEProblem of Medium |
neighbors | Algorithm to compute neighbors |
platform | Platform in which to run the model |
removalOfAgents_::Bool | Stores 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.
Argument | Description | |
---|---|---|
Args | dims | Dimensions of the system. |
KwArgs | agent=OrderedDict{Symbol,DataType}() | Agent parameters |
agentRule::Expr=quote end | Agent rules | |
agentODE::Expr=quote end | Agent Ordinary Differential Equations definition | |
agentSDE::Expr=quote end | Agent Stochastic Differential Equations term definition | |
model=OrderedDict{Symbol,DataType}() | Model parameters | |
modelRule::Expr=quote end | Model rules | |
modelODE::Expr=quote end | Model Ordinary Differential Equations definition | |
modelSDE::Expr=quote end | Model Ordinary Differential Equations definition | |
medium=OrderedDict{Symbol,DataType}() | Medium parameters | |
mediumRule::Expr=quote end | Medium rules | |
mediumODE::Expr=quote end | Medium Ordinary Differential Equations definition | |
mediumSDE::Expr=quote end | Medium 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.
Macros
These are macros that can be used in the code.
AgentRule
CellBasedModels.@addAgent
— Macromacro addAgentCode(arguments)
Macro that returns the special code based on the arguments provided to addAgent
.
CellBasedModels.@removeAgent
— Macromacro removeAgent()
Macro that returns the special code based on the arguments provided to removeAgent()
.
CellBasedModels.@loopOverNeighbors
— Macromacro 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.
Medium
CellBasedModels.@∂
— Macromacro ∂(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.
CellBasedModels.@∂2
— Macromacro ∂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.
CellBasedModels.@mediumInside
— Macromacro mediumInside()
Macro that returns true if mesh position is not in the border of the region.
CellBasedModels.@mediumBorder
— MacromediumBorder(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.
Community
CellBasedModels.Community
— TypeBasic structure keeping the parameters of all the agents in the current simulation of a model.
Parameters to store information essential for Community simulation
Symbol | Description |
---|---|
abm | ABM model of the community |
uuid | Unique identifier of the community |
loaded | If 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 |
agentDEProblem | ODEProblem or SDEProblem object of Agent |
modelDEProblem | ODEProblem or SDEProblem object of Model |
mediumDEProblem | ODEProblem or SDEProblem object of Medium |
Parameters seen in the kernels that may de used directly by the user
Symbol | Description |
---|---|
t | Time of the community |
dt | Stepping time of the community |
N | Number of agents |
NMedium | Size of medium mesh |
simBox | Size of simulation box |
dx | Size of the axis 1 in mesh medium |
dy | Size of the axis 2 in mesh medium |
dz | Size of the axis 3 in mesh medium |
Parameters seen in the kernels that are for internal use
Symbol | Description |
---|---|
id | Identification 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.
Platform loading
CellBasedModels.loadToPlatform!
— Functionfunction 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.).
CellBasedModels.bringFromPlatform!
— Functionfunction 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.
Evolution functions
CellBasedModels.evolve!
— Functionfunction 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.
CellBasedModels.step!
— Functionfunction step!(community)
Executes all the possible step functions and updates the parameters a single time.
CellBasedModels.agentStepRule!
— Functionfunction agentStepRule!(community)
Function that computes a local step of the community a time step dt
.
CellBasedModels.agentStepDE!
— Functionfunction agentStepDE!(community)
Function that computes a integration step of the community a time step dt
using the defined Integrator defined in Agent.
CellBasedModels.modelStepRule!
— Functionfunction agentStepRule!(community)
Function that computes a local step of the community a time step dt
.
CellBasedModels.modelStepDE!
— Functionfunction modelStepDE!(community)
Function that computes a integration step of the community a time step dt
using the defined Integrator defined in Agent.
CellBasedModels.mediumStepRule!
— Functionfunction agentStepRule!(community)
Function that computes a local step of the community a time step dt
.
CellBasedModels.mediumStepDE!
— Functionfunction mediumStepDE!(community)
Function that computes a integration step of the community a time step dt
using the defined Integrator defined in Agent.
CellBasedModels.update!
— Functionfunction update!(community)
Function that updates all the parameters after making a step in the community.
IO
CellBasedModels.saveRAM!
— Functionfunction saveRAM!(community::Community)
Function that stores the present configuration of community in the field Community.pastTimes
.
CellBasedModels.saveJLD2
— FunctionsaveJLD2(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.
CellBasedModels.loadJLD2
— Functionfunction loadJLD2(file::String)
Load the Community structure saved in file.
Metrics
CellBasedModels.CBMMetrics.euclidean
— Functioneuclidean(x1,x2)
euclidean(x1,x2,y1,y2)
euclidean(x1,x2,y1,y2,z1,z2)
Euclidean distance metric between two positions.
d = (x₁-x₂)²
CellBasedModels.CBMMetrics.@euclidean
— Macromacro 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_).
CellBasedModels.CBMMetrics.manhattan
— Functionmanhattan(x1,x2)
manhattan(x1,x2,y1,y2)
manhattan(x1,x2,y1,y2,z1,z2)
Manhattan distance metric between two positions.
d = |x₁-x₂|
CellBasedModels.CBMMetrics.@manhattan
— Macromacro 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_).
CellBasedModels.CBMMetrics.cellInMesh
— Functionfunction 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.
CellBasedModels.CBMMetrics.intersection2lines
— Functionfunction 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.
CellBasedModels.CBMMetrics.point2line
— Functionfunction 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.
CellBasedModels.CBMMetrics.pointInsideRod
— Functionfunction 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.
CellBasedModels.CBMMetrics.rodIntersection
— Functionfunction 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).
Random
Random number generators from distributions that are compatible both in CPU and GPU.
CellBasedModels.CBMDistributions.normal
— Functionnormal(μ,σ)
Normal distribution radom generation. μ
mean, σ
std.
CellBasedModels.CBMDistributions.uniform
— Functionuniform(l0,l1)
Uniform distribution radom generation. l0
min, l1
max.
CellBasedModels.CBMDistributions.exponential
— Functionexponential(θ)
Exponential distribution radom generation. θ
mean.
Platform
CellBasedModels.CPU
— Type" mutable struct CPU <: Platform
CPU platform structure.
CellBasedModels.GPU
— Type" mutable struct GPU <: Platform
CPU platform structure.
Threads and blocks for each rule method
- agentThreads
- agentBlocks
- modelThreads
- modelBlocks
- mediumThreads
- mediumBlocks
Integrators
Integrators can be used the ones from DifferentialEquations
for ODE or SDE problems or the custom made solvers provided in ths package.
CellBasedModels.CBMIntegrators.Euler
— Typemutable struct Euler <: CustomIntegrator
Euler integrator for ODE problems.
CellBasedModels.CBMIntegrators.Heun
— Typemutable struct Heun <: CustomIntegrator
Heun integrator for ODE problems.
CellBasedModels.CBMIntegrators.RungeKutta4
— Typemutable struct RungeKutta4 <: CustomIntegrator
RungeKutta4 for ODE integrators
CellBasedModels.CBMIntegrators.EM
— Typemutable struct EM <: CustomIntegrator
Euler-Majurana integrator for SDE poblems.
CellBasedModels.CBMIntegrators.EulerHeun
— Typemutable struct EulerHeun <: CustomIntegrator
Euler-Heun method for SDE integration.
Neighbor algorithms
CellBasedModels.CBMNeighbors.Full
— Typemutable struct Full <: Neighbors
Method that computes all against all neighbors.
CellBasedModels.CBMNeighbors.VerletTime
— Typemutable struct VerletTime <: Verlet <: Neighbors
Method that computes VerletList neighbors and updates it at fixed times.
Constructor
function VerletTime(;skin,dtNeighborRecompute,nMaxNeighbors)
Parameter | Description | |
---|---|---|
KwArgs | skin 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. |
CellBasedModels.CBMNeighbors.VerletDisplacement
— Typemutable 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)
Parameter | Description | |
---|---|---|
KwArgs | skin is the maximum distance around center to check neighbors. | |
nMaxNeighbors is the maximum number of neighbors that a cell may have. |
CellBasedModels.CBMNeighbors.CellLinked
— Typemutable 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)
Parameter | Description | |
---|---|---|
KwArgs | cellEdge is the grid size to use to check for neighbors around neighbor cells. |
CellBasedModels.CBMNeighbors.CLVD
— Typemutable 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)
Parameter | Description | |
---|---|---|
KwArgs | skin 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. |
Fitting
CellBasedModels.CBMFitting.gridSearch
— Functionfunction 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.
Parameter | Description | |
---|---|---|
Args | evalFunction:: 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)). | ||
KwArgs | returnAll::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 . |
CellBasedModels.CBMFitting.swarmAlgorithm
— Functionfunction 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.
Parameter | Description | |
---|---|---|
Args | evalFunction:: 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)). | |
KwArgs | population::Int = 100 | Size of the colony used at each generation for the optimization. |
weightInertia::Number = .1 | Hyperparameter of the colony weighting the current velocity. | |
weightGlobalBest::Number = .1 | Hyperparameter of the colony weighting the global best solution. | |
weightPopulationBest::Number = .1 | Hyperparameter of the colony weighting the local best solution. | |
stopMaxGenerations::Int = 10 | How many generations do before stopping the algorithm. | |
initialisation::Union{Nothing,DataFrame} = nothing | DataFrame defining the initial parameters of the population. If nothing, they are set randomly. | |
returnAll::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 . | |
verbose=false | If true , show progress bar during optimization. |
CellBasedModels.CBMFitting.beeColonyAlgorithm
— Functionfunction 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.
Parameter | Description | |
---|---|---|
Args | evalFunction:: 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)). | |
KwArgs | population::Int=100 | Size of the colony used at each generation for the optimization. |
limitCycles::Int = 10 | Hyperparameter of the algorithm that says how many generations without update are waited until jump to other position. | |
stopMaxGenerations::Int = 100 | How many generations do before stopping the algorithm. | |
returnAll::Bool = false | If return the hole list of parameters explored or the just the most fit. | |
initialisation::Union{Nothing,DataFrame} = nothing | DataFrame defining the initial parameters of the population. If nothing, they are set randomly. | |
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 . | |
verbose=false | If true , show progress bar during optimization. |
CellBasedModels.CBMFitting.geneticAlgorithm
— Functionfunction 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.
Parameter | Description | |
---|---|---|
Args | evalFunction:: 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)). | |
KwArgs | population::Int=100 | Size of the colony used at each generation for the optimization. |
parentSelectionAlg::String = "weighted" | Weigthing method of the population ot chose descendants. | |
parentSelectionP::Number = .1 | Hyperparameter of the algorithm indicating the proportion of parameters exchanged between parents. | |
mutationRate::Number = .1 | Hyperparameter of the algorithm indicating the probability of resampling the parameter with a uniform. | |
stopMaxGenerations::Int = 10 | How many generations do before stopping the algorithm. | |
initialisation::Union{Nothing,DataFrame} = nothing | DataFrame defining the initial parameters of the population. If nothing, they are set randomly. | |
returnAll::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 . | |
verbose=false | If true , show progress bar during optimization. |
Plotting
CellBasedModels.CBMPlots.plotRods2D!
— Functionfunction plotRods2D!(ax, x, y, d, l, angle; kargs...)
Plot rod shape like cells using Makie functions.
Parameter | Description | |
---|---|---|
Args | ax | Axis where to plot the rods |
x | Coordinates of rods in x | |
y | Coordinates of rods in y | |
d | Diameter of the rod | |
l | Length of the rods | |
angle | Angle of the rods of rods in the XY plane | |
KwArgs | All 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.
Parameter | Description | |
---|---|---|
Args | ax | Axis where to plot the rods |
x | Coordinates of center of rod in x | |
y | Coordinates of center of rod in y | |
xs1 | Coordinates of rod extreme in x | |
ys1 | Coordinates of rod extreme in y | |
xs2 | Coordinates of other rod extreme in x | |
ys2 | Coordinates of other rod extreme in y | |
markerSphere | Point3f0 describing the radius of the sphere in (x,y,z) | |
markerCylinder | Point3f0 describing the cylinder sizes in (l,rx,ry) | |
KwArgs | All arguments that want to be passed to meshscatter! function as color etc... |