Populations

A population is a collection of

  • either a lattice<:RealLattice if the model is spatial (see lattices), or a dummy placeholder NoLattice if it is not.
  • a phylogenetic tree, which is a directed graph where the root(s) represent the wild type(s).
  • metadata about the population and its genetics. See Metadata.
  • the time the state has been evolved for by invoking dynamics on it.

They are of type Population.

It is normally not required to manipulate metadata, phylogeny or lattice directly. Convenient methods to manipulate state are provided.

julia>     state, _ = uniform(CubicLattice, 32; g=0)(CubicLattice{Int64, Array{Int64, 3}}
0	genotypes
0	population, nothing)
julia> # fails because the new genotype is yet unknown state[16,16,16] = 1ERROR: ArgumentError: Genotype 1 is not know. Use push!(::Population, 1) first.
julia> # make it known first # last paremeter is the parental genotype # because there is none yet, set it to `nothing` add_genotype!(state, 1, nothing)1
julia> state[16,16,16] = 1;
julia> state.meta1-element MetaData{Int64}: (genotype = 1, npop = 1, fitness = 1.0, snps = nothing, age = (0, 0.0))
julia> state[15:17,15:17,15:17]3×3×3 Array{Int64, 3}: [:, :, 1] = 0 0 0 0 0 0 0 0 0 [:, :, 2] = 0 0 0 0 1 0 0 0 0 [:, :, 3] = 0 0 0 0 0 0 0 0 0
julia> push!(state, 2)2
julia> state[1:size(state,1), 1, 1] .= 2;
julia> state.meta2-element MetaData{Int64}: (genotype = 1, npop = 1, fitness = 1.0, snps = nothing, age = (0, 0.0)) (genotype = 2, npop = 33, fitness = 1.0, snps = nothing, age = (0, 0.0))

Convenience constructors

A handful of methods to generate common geometries are provided. They return a tuple (population,auxilliary). The latter holds information that might be nice to have. It is often nothing, but for example spherer returns a vector with all indices contained in the ball.

GrowthDynamics.Populations.uniformFunction
uniform(T, L; g=0, a=1.0)

A population on a lattice of type T with linear extension L, filled with genotype g.

Return (population, nothing).

Example

julia> uniform(HexagonalLattice, 128; g=1)
(HexagonalLattice{Int64, Matrix{Int64}}
1	genotypes
16641	population, nothing)
source
GrowthDynamics.Populations.spherefFunction
spheref(T, L::Int; r = 0, g1=0, g2=1) where LT<:Lattices.RealLattice

Fill lattice of type T (e.g CubicLattice) with genotype g1 and put a (L2-)ball with genotype g2 that occupies approx. a fraction f of the lattice at the center.

Return (population, idx_ball) with idx_ball being the indices of sites within the ball.

source
GrowthDynamics.Populations.sphererFunction
spherer(T, L::Int; r = 0, g1=0, g2=1) where LT<:Lattices.RealLattice

Fill lattice of type T (e.g CubicLattice) with genotype g1 and put a (L2-)ball of approx. radius r with genotype g2 at the center.

Return (population, idx_ball) with idx_ball being the indices of sites within the ball.

source
GrowthDynamics.Populations.half_spaceFunction
half_space(T, L; g1=1, g2=2)

A population of genotype g1 on a lattice of type T. Fill the last dimension with g2 up to fraction f

Return (population, fill_to) with fill_to defined as state[:,..., 1:fill_to] == g2 and state[:,..., fill_to+1:end] == g1.

Example

julia> using GrowthDynamics.Populations

julia> state = half_space(CubicLattice, 32, f=1/4, g1=1, g2=2)[1]
CubicLattice{Int64, Array{Int64, 3}}
2	genotypes
35937	population

julia> state.meta
2-element MetaData{Int64}:
 (genotype = 1, npop = 27225, fitness = 1.0, snps = nothing, age = (0, 0.0))
 (genotype = 2, npop = 8712, fitness = 1.0, snps = nothing, age = (0, 0.0))
source
GrowthDynamics.Populations.sphere_with_single_mutant_on_outer_shellFunction
sphere_with_single_mutant_on_outer_shell(T, L::Int; r, s=1.0)

A lattice of type T (e.g CubicLattice) with an (L2-)ball of genotype 1 and radius r at the center. One random site on the outermost shell is populated with genotype 2 that has fitness s.

Return (population, (idx_shell, idx_mutant)) with idx_shell being the indices of the outer shell, and idx_mutant the index of the mutant genotype.

source
GrowthDynamics.Populations.sphere_with_diverse_outer_shellFunction
sphere_with_diverse_outer_shell(T, L; r)

A lattice of type T (e.g CubicLattice) with an (L2-)ball of genotype 1 and radius r at the center. Each site on the outermost shell is populated consecutively with a different genotype, starting from 2.

Return (population, idx_shell) with idx_shell being the indices of the outer ("diverse") shell.

source

An initial population without spatial structure can be constructed with

API

Manipulating genotypes

GrowthDynamics.Populations.PopulationType
Population{G, T}

Represents a population on a lattice of type T with genotypes of data type G.

Fields

  • lattice<:Lattices.AbstractLattice
  • phylogeny: directed graph recording the ancestry of genotypes
  • meta::MetaData: metadata such as fitnesses, mutations, etc. See MetaData.
  • t::Int: age in timesteps
  • treal::Float64: age in "real" time.
  • observables: dictionary to store observables in
source
GrowthDynamics.Populations.add_genotype!Function
add_genotype!(S::Population, G, parent)

Add genotype G to the population and connect it to parent in the phylogeny. G is either a genotype or a full MetaDatum. parent defaults to the first genotype, i.e. the root of the phylogenetic tree. If parent=nothing, the genotype will not be connected to the tree.

See also: MetaDatum, remove_genotype!

source
GrowthDynamics.Populations.remove_genotype!Function
remove_genotype!(S::Population, g; bridge=true)

Remove genotype from the population. Discards it from meta data, prunes it from the phylogeny, and sets all corresponding sites of the lattice to zero.

If bridge=true (default), the gap in the phylogeny is bridged with new edges.

Throw an exception if the requested genotype does not exist.

Return true if successful, else false.

source

Methods to add mutations

GrowthDynamics.Populations.add_snps!Function
add_snps!(state::Population, g, μ)

Randomize and add mutations to genotype g, or replace them. μ is the genome wide mutation rate (Poisson) or mutation count (fixed).

Note: allow_multiple is much faster if the number of mutations already present is large.

Keyword arguments

  • L=10^9: length of the genome
  • allow_multiple=false: allow for a site to mutate more than once.
  • kind=:poisson: either :poisson or :fixed
  • replace=false: replace existing SNPs.
source
add_snps!(state::Population, g, v::Vector{Int})

Add mutations v to genotype g. No checks for duplications are performed.

Return a vector of all mutations.

source
Missing docstring.

Missing docstring for annotate_snps!. Check Documenter's build log for details.

Missing docstring.

Missing docstring for annotate_lineage!. Check Documenter's build log for details.