References

Cells and lineages

SomaticPhylogenies.CellType
Cell

Represents the basic element of a CellLineage.

Fields

  • uid::T : Unique identifier of Cell.
  • generation::T : Position of Cell within a CellLineage.
  • t::F : Time point of the Cells birth.
  • mutations::M : Number of Cell's private mutations.

Base methods

  • == : compares uid

See also: CellLineage, MutationMode, MutationRate

source

Mutations

SomaticPhylogenies.TimeCounterType
TimeCounter <: MutationMode

Describes that the number of newly arising mutations in a cell increases linearly with the cell's age. Together with the already present mutations, these new mutations are passed on to both daughter cells.

See also: Cell, MutationRate

source
SomaticPhylogenies.MutGenProcessType
MutGenProcess

Abstract type that specifies whether or not stochastic fluctuations in the number of newly arising mutations in a cell should be taken into account.

It is either DeterministicProcess or PoissonProcess.

See MutationRate.

source
SomaticPhylogenies.MutationRateType
MutationRate{MM<:MutationMode, MGP<:MutGenProcess}

Determines mutation accumulation along a CellLineage.

The two parameters MutationMode and MutGenProcess specify how the single field μ of MutationRate is interpreted.

  • ReplicationCounter : μ has unit mutations per cell division and counts for example the (average) number replication errors. Upon each cell division, each daughter cell inherits all variants from its mother plus a number of novel, private mutations. For DeterministicProcess this number is given by μ, for PoissonProcess μ serves as expectation value for a poisson distrution.
  • TimeCounter : μ has unit mutations per cell and time and counts for example the (average) number of DNA repair errors in between cell divisions. Since these mutations occur before the mother cell divides, both daughter cells inherit the same mutations, which comprise the mutations the mother cell inherited from its mother plus a number of novel mutations, which depend linearly on the lifespan of the mother cell, Δt, where μ serves as proportionality constant (mutation clock). For DeterministicProcess this number is μ Δt, while for PoissonProcess μ Δt serves as expectation value for a poisson distrution.

Constructor

MutationRate(
    [MM::Type{<:MutationMode} = ReplicationCounter,]
    μ::Real;
    poisson::Bool = false,
)

Construct MutationRate with value μ. The default MutationMode is ReplicationCounter and PoissonProcess may be specified by setting poisson = true.

See also: MutationMode, MutGenProcess, CellLineage

source

Phylogeny

SomaticPhylogenies.PhylogenyType
Phylogeny

Phylogeny of a cell population at time t. Each member of the population is represented by a CellLineage holding ancestral information in terms of Cells. A Cell defines a uid and contains information about the time of its birth, the number of its private mutations, and its position within a CellLineage. Population members descending from the same cell all have that Cell at the same position in their CellLineage.

The Phylogeny-generating dynamics are simulated with AbstractProcesses .

Fields

  • lineages::Vector{CellLineage}
  • t::Real : The time point of observation.
  • μ::MutationRate
  • uid_counter::Int : The global uid counter, which ensures that new mutations in a Cell are uniquely identified (infinte-site hyothesis).

See also: Cell, CellLineage, MutationRate, AbstractProcess


Constructors

Phylogeny(μ::MutationRate[, N0::Integer, t0])

Create an initial Phylogeny with N0 (default 1) somatically unmutated founder cells with initial time t0 (default 0.0).

Phylogeny(μ::MutationRate, N0::Vector{<:Real}[, t0])

In order to equip the founder cells with somatic mutations, a vector N0 may be passed where length(N0) is the number of founder cells and the elements are the number of founder mutations.

See also: MutationRate

source
SomaticPhylogenies.subpopulation_sizesMethod
subpopulation_sizes(phylo::Phylogeny) -> Dict{Cell, Int}

Return the number of cells that derive from distinct founder cells.

subpopulation_sizes(generation::Integer, phylo::Phylogeny) -> Dict{Cell, Int}

Instead of founder cells, a generation may be provided whose Cells define the sub-populations.

subpopulation_sizes(t::AbstractFloat, phylo::Phylogeny) -> Dict{Cell, Int}

Instead of generation, a time point t may be provided where Cells that existed at t define the sub-populations.

source
SomaticPhylogenies.num_mutationsMethod
num_mutations(::CellLineage) -> Real

Return the number of mutations acquired along CellLineage.

num_mutations(::Phylogeny) -> Vector{<:Real}

Return the number of mutations acquired along each CellLineage of Phylogeny.

source
SomaticPhylogenies.cell_vectorMethod
cell_vector(::Phylogeny) -> Vector{Cell}

Return all Cells that ever occurred in Phylogeny.

cell_vector(generation::Integer, ::Phylogeny) -> Vector{Cell}

Return all Cells of generation or greater that ever occurred in Phylogeny.

source
SomaticPhylogenies.cell_countMethod
cell_count(::Phylogeny [, n_min]) -> Int

Return the total number of distinct Cells that ever occurred in Phylogeny.

If n_min is specified, only Cells that occur in at least n_min CellLineages are considered.

source
SomaticPhylogenies.subsamplingFunction
subsampling(::Phylogeny, n::Integer) -> Phylogeny

Return sub-Phylogeny comprising n randomly selected lineages (without replacement).

source
SomaticPhylogenies.MRCAFunction
MRCA(::Phylogeny) -> CellLineage

Return the CellLineage of the Most Recent Common Ancestor of Phylogeny.

If no MRCA exists, an empty CellLineage is returned.

source
SomaticPhylogenies.leading_subclonesFunction
leading_subclones(::Phylogeny, depth::Integer)

Determine the leading subclones of Phylogeny up to depth.

Returns

clone_structure::Vector{Vector{Pair{<:Cell, NamedTuple{(:clonesize, :Δnmut), Tuple{Int64, Float64}}}}} : Vector of length depth where each element represents the subclones of a level, where level = 1:depth. Each level contains 2^(level - 1) clones, which are represented by a Vector of that length. Each clone is represented by a Pair, where first is the last Cell along the sublineage defining the clone, and second is a NamedTuple containing the clonesize and the number of mutations Δnmut accumulated along the sublineage.

The relationship between clones across levels is as follows: The first two clones on level L descent from the first clone on level L-1, the third and fourth clones on level L descent from the second clone on level L-1, and so on. Furthermore, it is possible that a clone on some level L represents a single Cell (clonesize = 1), which is therefore still alive at the time of observation. In order to keep above relationship logic valid,

source

Processes

SomaticPhylogenies.AbstractProcessType
AbstractProcess

Abstract supertype of all population-dynamics processes.

There are three commonalties for the concrete types

  1. A functor is implemented
  2. The functor accepts a Phylogeny instance as input
  3. The functor returns a Phylogeny instance
source
SomaticPhylogenies.BirthDeathProcessType
BirthDeathProcess{F,I} <: AbstractProcess

BirthDeathProcess with division rate λ and loss rate δ.

BirthDeathProcess has three termination conditions: duration Δt, minimal population size N_min, and maximal population size N_max.


Constructors

BirthDeathProcess(
    λ::Real,
    δ::Real;
    Δt::Real = Inf,
    N_min::Real = 0,
    N_max::Real = Inf)
BirthDeathProcess(
    λ::Real = 1.0;
    Δt::Real = Inf,
    N_min::Real = 0,
    N_max::Real = Inf)

Construct BirthDeathProcess. If δ is not provided, it is set to 0.0 (pure birth process).


Functors

(::BirthDeathProcess)(μ::MutationRate [, N0, t0::Real]) -> Phylogeny
(::BirthDeathProcess)(phylo::Phylogeny) -> Phylogeny

Simulate BirthDeathProcess and return Phylogeny.

source
SomaticPhylogenies.MoranProcessType
MoranProcess{F,I} <: AbstractProcess

MoranProcess with turnover rate λ and population size N.

MoranProcess terminates after duration Δt.


Constructors

MoranProcess(λ::Real, N::Integer, Δt::Real)
MoranProcess(N::Integer, Δt::Real)

Construct MoranProcess. If λ is not provided, it is set to 1.0.


Functors

(::MoranProcess)(μ::MutationRate [, t0::Real]) -> Phylogeny
(::MoranProcess)(μ::MutationRate, N0::Vector[, t0::Real]) -> Phylogeny
(::MoranProcess)(phylo::Phylogeny) -> Phylogeny

Simulate MoranProcess and return Phylogeny.

source
SomaticPhylogenies.CompositeProcessType
CompositeProcess <: AbstractProcess

Fields

  • subprocesses::Vector{AbstractProcess} : The subprocesses.

Constructor

CompositeProcess(subprocesses...)

Functor

(::CompositeProcess)(phylo::Phylogeny) -> Phylogeny

Simulate the processes of CompositeProcess successively and return Phylogeny.

source

Whole-genome sequencing

SomaticPhylogenies.BinomialSamplingType
BinomialSampling <: SamplingMode

The simulated read counts for variants are assumed to be binomially distributed while the coverage for each genomic site is the same (= sequencing depth).

source
SomaticPhylogenies.PoissonBinomialSamplingType
PoissonBinomialSampling <: SamplingMode

The simulated read counts for variants are assumed to be binomially distributed and the coverage of a genomic site is assumed to be poisson-distributed (sequencing depth = mean(poisson)).

source
SomaticPhylogenies.BetaSamplingType
BetaSampling <: SamplingMode

Instead of simulating read counts and coverages, the variant allele frequencies are directly determined with a Beta distribution.

  • As the number of variants at true allele frequency x becomes larger, while each individual variant follows BinomialSampling, the measured variant allele frequencies converge against a Beta distribution characterized by α = x * (seq_depth - 1), β = (1 - x) * (seq_depth - 1).
source
SomaticPhylogenies.WholeGenomeSequencingType
WholeGenomeSequencing{SM<:SamplingMode}

Object to simulate whole-genome sequencing.

Fields

  • depth::Int : The (average) coverage of genomic sites.
  • vaf_edges::Vector{Float64} : The bin edges for variant allele frequencies
  • read_min::Int : Only variants with number of reads ≥ read_min are considered in sequencing output.

Constructors

WholeGenomeSequencing(
    ::Type{<:SamplingMode},
    [vaf_edges::AbstractVector{<:AbstractFloat},]
    depth::Integer;
    read_min::Integer=1,
    read_step::Integer=1)

Create WholeGenomeSequencing for SamplingMode and (average) coverage depth. If vaf_edges is not provided, it will be specified according to unique(push!(collect(read_min/depth : read_step/depth : 1), 1)).

SamplingMode is either BetaSampling, BinomialSampling, or PoissonBinomialSampling.

Optional keyword arguments

  • read_min::Integer=1 : The minimal read count to use in sequencing output.
  • read_step::Integer=1 : The number of reads between two adjacent vaf_edges (only relevant if vaf_edges is not supplied).

See also: SamplingMode

Functor

(wgs::WholeGenomeSequencing{<:SamplingMode})(::Phylogeny) -> SequencingResult
(wgs::WholeGenomeSequencing{<:SamplingMode})(::Vector{<:Real}) -> SequencingResult

Simulate WholeGenomeSequencing and return SequencingResult. The input is either a Phylogeny or an already computed site frequency spectrum.

See also: SamplingMode, site_frequency_spectrum

Returns

A SequencingResult comprising two fields:

  • wgs::WholeGenomeSequencing : The WholeGenomeSequencing object used for simulations.
  • vaf_histogram::Vector{<:Real} : The simulated distribution of variant allele frequencies across wgs.vaf_edges.

See also: SequencingResult

source
SomaticPhylogenies.SequencingResultType
SequencingResult

Result of whole-genome sequencing.

Fields

  • wgs::WholeGenomeSequencing : The WholeGenomeSequencing object used for simulations.
  • vaf_histogram::Vector{<:Real} : The simulated number of variant allele frequencies (vafs) across wgs.vaf_edges. Specifically, the k-th value counts vaf_edges[k] ≤ vaf < vaf_edges[k+1], where 1 ≤ k < length(vaf_edges); the number of variants with vaf = 1 is vaf_histogram[end].

See also: WholeGenomeSequencing

source
SomaticPhylogenies.vaf_minMethod
vaf_min(::WholeGenomeSequencing) -> Float64
vaf_min(::SequencingResult) -> Float64

Return the lowest value of vaf_edges.

source
SomaticPhylogenies.vaf_stepMethod
vaf_step(::WholeGenomeSequencing) -> Float64
vaf_step(::SequencingResult) -> Float64

Return the width of the first bin of vaf_edges.

source
SomaticPhylogenies.vafogramFunction
vafogram(::SequencingResult)
vafogram!(::SequencingResult)

Plot the number of variants greater-or-equal than vaf versus inverted vaf.

source