Title: | Inference of Transmission Tree from a Dated Phylogeny |
---|---|
Description: | Inference of transmission tree from a dated phylogeny. Includes methods to simulate and analyse outbreaks. The methodology is described in Didelot et al. (2014) <doi:10.1093/molbev/msu121>, Didelot et al. (2017) <doi:10.1093/molbev/msw275>. |
Authors: | Xavier Didelot [aut, cre] , Yuanwei Xu [ctb] |
Maintainer: | Xavier Didelot <[email protected]> |
License: | GPL (>=2) |
Version: | 1.4.10 |
Built: | 2024-11-09 04:39:38 UTC |
Source: | https://github.com/xavierdidelot/TransPhylo |
Inference of transmission tree from a dated phylogeny. Includes methods to simulate and analyse outbreaks.
Xavier Didelot [email protected]
Didelot et al. (2014) <doi:10.1093/molbev/msu121> Didelot et al. (2017) <doi:10.1093/molbev/msw275>.
https://github.com/xavierdidelot/TransPhylo
Convert to coda mcmc format
as.mcmc.resTransPhylo(x, burnin = 0.5)
as.mcmc.resTransPhylo(x, burnin = 0.5)
x |
Output from inferTTree |
burnin |
Proportion of the MCMC output to be discarded as burnin |
mcmc object from coda package
Build a matrix indicating for each pairs of individuals how many intermediates there are in the transmission chain
computeMatTDist(record, burnin = 0.5)
computeMatTDist(record, burnin = 0.5)
record |
Output from inferTTree function |
burnin |
Proportion of the MCMC output to be discarded as burnin |
Matrix of intermediates in transmission chains between pairs of hosts
Build a matrix of probability of who infected whom from a MCMC output
computeMatWIW(record, burnin = 0.5)
computeMatWIW(record, burnin = 0.5)
record |
Output from inferTTree function |
burnin |
Proportion of the MCMC output to be discarded as burnin |
Matrix of probability of who infected whom
Build a consensus transmission tree from a MCMC output
consTTree(record, burnin = 0.5, minimum = 0.2, debug = F)
consTTree(record, burnin = 0.5, minimum = 0.2, debug = F)
record |
Output from inferTTree function |
burnin |
Proportion of the MCMC output to be discarded as burnin |
minimum |
Minimum probability for inclusion of a partition in the consensus |
debug |
Used for debugging |
The consensus transmission tree
Convert to coda mcmc format
convertToCoda(record, burnin = 0.5)
convertToCoda(record, burnin = 0.5)
record |
Output from inferTTree function |
burnin |
Proportion of the MCMC output to be discarded as burnin |
Object of class mcmc from coda package
Return the date of last sample from a ttree or ctree or ptree
dateLastSample(x)
dateLastSample(x)
x |
A transmission tree or colored tree or phylogenetic tree |
date of the last sample
Return the combined tree corresponding to a given iteration of the TransPhylo results
extractCTree(res, iteration)
extractCTree(res, iteration)
res |
Output from inferTTree command |
iteration |
Number of the iteration to be extracted |
The colored tree at the specified iteeatino
Extracts phylogenetic tree from a combined phylogenetic/transmission tree
extractPTree(ctree)
extractPTree(ctree)
ctree |
Combined tree |
phylogenetic tree
extractPTree(simulateOutbreak())
extractPTree(simulateOutbreak())
Extracts transmission tree from a combined phylogenetic/transmission tree
extractTTree(ctree)
extractTTree(ctree)
ctree |
Combined tree |
transmission tree
extractTTree(simulateOutbreak())
extractTTree(simulateOutbreak())
Extract and return realised generation time distribution
getGenerationTimeDist( record, burnin = 0.5, maxi = 2, numBins = 20, show.plot = F )
getGenerationTimeDist( record, burnin = 0.5, maxi = 2, numBins = 20, show.plot = F )
record |
MCMC output produced by inferTTree |
burnin |
Proportion of the MCMC output to be discarded as burnin |
maxi |
Maximum generation time to consider |
numBins |
Number of time bins to compute and display distribution |
show.plot |
Show a barplot of the distribution |
Vector of times between becoming infected and infecting others (generation times) in the posterior
Returns and/or plot numbers of sampled and unsampled cases over time
getIncidentCases( record, burnin = 0.5, numBins = 10, dateT = NA, show.plot = FALSE )
getIncidentCases( record, burnin = 0.5, numBins = 10, dateT = NA, show.plot = FALSE )
record |
Output from inferTTree function |
burnin |
Proportion of the MCMC output to be discarded as burnin |
numBins |
Number of time bins to compute and display incident cases |
dateT |
Date when process stops (this can be Inf for fully resolved outbreaks) |
show.plot |
Show a plot of incident cases over time with stacked bars |
List with four entries. Time is a vector of the time points. allCases is the average number of cases at each time in the posterior. sampledCases: average number of sampled cases. unsampCases: average number of unsampled cases.
Extract and return distribution of infection time of given sampled case(s)
getInfectionTimeDist(record, burnin = 0.5, k, numBins = 10, show.plot = F)
getInfectionTimeDist(record, burnin = 0.5, k, numBins = 10, show.plot = F)
record |
MCMC output produced by inferTTree |
burnin |
Proportion of the MCMC output to be discarded as burnin |
k |
Case(s) whose posterior infection times are to be extracted. Either a string matching one of the case names in the data, or a vector of such strings |
numBins |
Number of bins to use for plot |
show.plot |
Show a barplot of the distribution |
Posterior infection times for the case(s) in k. If length(k)==1 then a vector is returned, otherwise a matrix
Extract and return offspring distribution of given sampled case(s)
getOffspringDist(record, burnin = 0.5, k, show.plot = F)
getOffspringDist(record, burnin = 0.5, k, show.plot = F)
record |
MCMC output produced by inferTTree |
burnin |
Proportion of the MCMC output to be discarded as burnin |
k |
Case(s) whose offspring distribution are to be extracted. Either a string matching one of the case names in the data, or a vector of such strings |
show.plot |
Show a barplot of the distribution |
Posterior offspring distribution for the case(s) in k. If length(k)==1 then a vector is returned, otherwise a matrix
Extract and return realised sampling time distribution
getSamplingTimeDist( record, burnin = 0.5, maxi = 2, numBins = 20, show.plot = F )
getSamplingTimeDist( record, burnin = 0.5, maxi = 2, numBins = 20, show.plot = F )
record |
MCMC output produced by inferTTree |
burnin |
Proportion of the MCMC output to be discarded as burnin |
maxi |
Maximum generation time to consider |
numBins |
Number of time bins to compute and display distribution |
show.plot |
Show a barplot of the distribution |
Vector of times between becoming infected and becoming sampled in the posterior
Infer transmission tree given a phylogenetic tree
inferTTree( ptree, w.shape = 2, w.scale = 1, ws.shape = NA, ws.scale = NA, w.mean = NA, w.std = NA, ws.mean = NA, ws.std = NA, mcmcIterations = 1000, thinning = 1, startNeg = 100/365, startOff.r = 1, startOff.p = 0.5, startPi = 0.5, updateNeg = TRUE, updateOff.r = TRUE, updateOff.p = FALSE, updatePi = TRUE, qNeg = NA, qOff.r = NA, qOff.p = NA, qPi = NA, startCTree = NA, updateTTree = TRUE, optiStart = 2, dateT = Inf, delta_t = NA, verbose = F )
inferTTree( ptree, w.shape = 2, w.scale = 1, ws.shape = NA, ws.scale = NA, w.mean = NA, w.std = NA, ws.mean = NA, ws.std = NA, mcmcIterations = 1000, thinning = 1, startNeg = 100/365, startOff.r = 1, startOff.p = 0.5, startPi = 0.5, updateNeg = TRUE, updateOff.r = TRUE, updateOff.p = FALSE, updatePi = TRUE, qNeg = NA, qOff.r = NA, qOff.p = NA, qPi = NA, startCTree = NA, updateTTree = TRUE, optiStart = 2, dateT = Inf, delta_t = NA, verbose = F )
ptree |
Phylogenetic tree |
w.shape |
Shape parameter of the Gamma distribution representing the generation time |
w.scale |
Scale parameter of the Gamma distribution representing the generation time |
ws.shape |
Shape parameter of the Gamma distribution representing the sampling time |
ws.scale |
Scale parameter of the Gamma distribution representing the sampling time |
w.mean |
Mean of the Gamma distribution representing the generation time |
w.std |
Std of the Gamma distribution representing the generation time |
ws.mean |
Mean of the Gamma distribution representing the sampling time |
ws.std |
Std of the Gamma distribution representing the sampling time |
mcmcIterations |
Number of MCMC iterations to run the algorithm for |
thinning |
MCMC thinning interval between two sampled iterations |
startNeg |
Starting value of within-host coalescent parameter Ne*g |
startOff.r |
Starting value of parameter off.r |
startOff.p |
Starting value of parameter off.p |
startPi |
Starting value of sampling proportion pi |
updateNeg |
Whether of not to update the parameter Ne*g |
updateOff.r |
Whether or not to update the parameter off.r |
updateOff.p |
Whether or not to update the parameter off.p |
updatePi |
Whether or not to update the parameter pi |
qNeg |
Proposal kernel range for parameter Ne*g |
qOff.r |
Proposal kernel range for parameter Off.r |
qOff.p |
Proposal kernel range for parameter Off.p |
qPi |
Proposal kernel range for parameter pi |
startCTree |
Optional combined tree to start from |
updateTTree |
Whether or not to update the transmission tree |
optiStart |
Type of optimisation to apply to MCMC start point (0=none, 1=slow, 2=fast) |
dateT |
Date when process stops (this can be Inf for fully simulated outbreaks) |
delta_t |
Grid precision (smaller is better but slower) |
verbose |
Whether or not to use verbose mode (default is false) |
posterior sample set of transmission trees
inferTTree(ptreeFromPhylo(ape::rtree(5),2020),mcmcIterations=100)
inferTTree(ptreeFromPhylo(ape::rtree(5),2020),mcmcIterations=100)
Create a transmission tree compatible with the provided phylogenetic tree
makeCTreeFromPTree( ptree, off.r = NA, off.p = NA, neg = NA, pi = NA, w.shape = NA, w.scale = NA, ws.shape = NA, ws.scale = NA, T = NA, optiStart = 0 )
makeCTreeFromPTree( ptree, off.r = NA, off.p = NA, neg = NA, pi = NA, w.shape = NA, w.scale = NA, ws.shape = NA, ws.scale = NA, T = NA, optiStart = 0 )
ptree |
Phylogenetic tree |
off.r |
First parameter of the negative binomial distribution for offspring number |
off.p |
Second parameter of the negative binomial distribution for offspring number |
neg |
the within-host effective population size (Ne) timesgeneration duration (g) |
pi |
probability of sampling an infected individual |
w.shape |
Shape parameter of the Gamma probability density function representing the generation time |
w.scale |
Scale parameter of the Gamma probability density function representing the generation time |
ws.shape |
Shape parameter of the Gamma probability density function representing the sampling time |
ws.scale |
Scale parameter of the Gamma probability density function representing the sampling time |
T |
Date when process stops (this can be Inf for fully simulated outbreaks) |
optiStart |
Method used to optimised colored tree (0=none, 1=slow, 2=fast) |
A minimal non-zero probability phylogenetic+transmission tree, or an optimised version if parameters are provided
Simulate a transmission tree
makeTTree( off.r, off.p, pi, w.shape, w.scale, ws.shape = w.shape, ws.scale = w.scale, maxTime = Inf, nSampled = NA )
makeTTree( off.r, off.p, pi, w.shape, w.scale, ws.shape = w.shape, ws.scale = w.scale, maxTime = Inf, nSampled = NA )
off.r |
First parameter of the negative binomial distribution for offspring number |
off.p |
Second parameter of the negative binomial distribution for offspring number |
pi |
probability of sampling an infected individual |
w.shape |
Shape parameter of the Gamma probability density function representing the generation time |
w.scale |
Scale parameter of the Gamma probability density function representing the generation time |
ws.shape |
Shape parameter of the Gamma probability density function representing the sampling time |
ws.scale |
Scale parameter of the Gamma probability density function representing the sampling time |
maxTime |
Duration of simulation (can be Inf) |
nSampled |
Number of sampled individuals (can be NA for any) |
A N*3 matrix in the following format with one row per infected host, first column is time of infection, second column is time of sampling, third column is infector
Return the medoid from a MCMC output
medTTree(record, burnin = 0.5)
medTTree(record, burnin = 0.5)
record |
Output from inferTTree function |
burnin |
Proportion of the MCMC output to be discarded as burnin |
The index of the medoid
Converts a phylogenetic tree into an ape phylo object
phyloFromPTree(ptree)
phyloFromPTree(ptree)
ptree |
phylogenetic tree |
phylo object
phyloFromPTree(extractPTree(simulateOutbreak()))
phyloFromPTree(extractPTree(simulateOutbreak()))
Plotting for ctree
## S3 method for class 'ctree' plot(x, ...)
## S3 method for class 'ctree' plot(x, ...)
x |
Object of class ctree, ie a colored phylogenetic tree |
... |
Additional parameters are passed on |
Plot of ctree
plot(simulateOutbreak())
plot(simulateOutbreak())
Plotting for ptree
## S3 method for class 'ptree' plot(x, ...)
## S3 method for class 'ptree' plot(x, ...)
x |
Object of class ptree, ie a phylogenetic tree |
... |
Additional parameters are passed on to ape::plot.phylo |
Plot of ptree
plot(ptreeFromPhylo(ape::rtree(5),2020))
plot(ptreeFromPhylo(ape::rtree(5),2020))
Plotting for resTransPhylo
## S3 method for class 'resTransPhylo' plot(x, ...)
## S3 method for class 'resTransPhylo' plot(x, ...)
x |
Output from inferTTree |
... |
Additional parameters are passed on |
Plot of TransPhylo results
Plotting for ttree
## S3 method for class 'ttree' plot(x, type = "summarised", w.shape = NA, w.scale = NA, ...)
## S3 method for class 'ttree' plot(x, type = "summarised", w.shape = NA, w.scale = NA, ...)
x |
Object of class ttree, ie a transmission tree |
type |
Type of plot to display, can be 'detailed' or 'summarised' (default) |
w.shape |
Shape parameter of the generation time, needed for detailed plot only |
w.scale |
Scale parameter of the generation time, needed for detailed plot only |
... |
Additional parameters are passed on |
Plot of ttree
plot(extractTTree(simulateOutbreak()))
plot(extractTTree(simulateOutbreak()))
Plot both phylogenetic and transmission trees using colors on the phylogeny
plotCTree( tree, showLabels = TRUE, showStars = TRUE, cols = NA, maxTime = NA, cex = 1 )
plotCTree( tree, showLabels = TRUE, showStars = TRUE, cols = NA, maxTime = NA, cex = 1 )
tree |
Combined phylogenetic/transmission tree |
showLabels |
Whether or not to show the labels |
showStars |
Whether or not to show stars representing transmission events |
cols |
Colors to use for hosts |
maxTime |
Maximum time to show on the x axis |
cex |
Expansion factor |
Returns invisibly the first parameter
plotCTree(simulateOutbreak())
plotCTree(simulateOutbreak())
Plot MCMC traces
plotTraces(record, burnin = 0, extend = F)
plotTraces(record, burnin = 0, extend = F)
record |
Output from inferTTree function |
burnin |
Proportion of the MCMC output to be discarded as burnin |
extend |
Whether to also show traces of off.r and off.p |
Returns invisibly the first parameter
Plot a transmission tree in a detailed format
plotTTree(ttree, w.shape, w.scale, showLabels = TRUE, maxTime = NA, cex = 1)
plotTTree(ttree, w.shape, w.scale, showLabels = TRUE, maxTime = NA, cex = 1)
ttree |
Transmission tree |
w.shape |
Shape parameter of the Gamma probability density function representing the generation time |
w.scale |
Scale parameter of the Gamma probability density function representing the generation time |
showLabels |
Whether or not to show the labels |
maxTime |
Maximum value of time to show on x axis |
cex |
Expansion factor |
Returns invisibly the first parameter
plotTTree(extractTTree(simulateOutbreak()),2,1)
plotTTree(extractTTree(simulateOutbreak()),2,1)
Plot a transmission tree in an economic format
plotTTree2( ttree, showLabels = TRUE, showMissingLinks = 0, maxTime = NA, cex = 1 )
plotTTree2( ttree, showLabels = TRUE, showMissingLinks = 0, maxTime = NA, cex = 1 )
ttree |
Transmission tree |
showLabels |
Boolean for whether or not to show the labels |
showMissingLinks |
Option for how to show missing links: (0) as dots, (1) as several gray levels, (2) as a single gray level |
maxTime |
Maximum value of time to show on x axis |
cex |
Expansion factor |
Returns invisibly the first parameter
plotTTree2(extractTTree(simulateOutbreak()))
plotTTree2(extractTTree(simulateOutbreak()))
Print function for ctree objects
## S3 method for class 'ctree' print(x, ...)
## S3 method for class 'ctree' print(x, ...)
x |
Object of class ctree, ie a colored phylogenetic tree |
... |
Additional parameters are passed on |
Print out details of the ctree
print(simulateOutbreak())
print(simulateOutbreak())
Print function for ptree objects
## S3 method for class 'ptree' print(x, ...)
## S3 method for class 'ptree' print(x, ...)
x |
Object of class ptree, ie a phylogenetic tree |
... |
Additional parameters are passed on |
Print out details of the ptree
print(extractPTree(simulateOutbreak()))
print(extractPTree(simulateOutbreak()))
Print function for resTransPhylo objects
## S3 method for class 'resTransPhylo' print(x, ...)
## S3 method for class 'resTransPhylo' print(x, ...)
x |
output from inferTTree |
... |
Additional parameters are passed on |
Print out details of TransPhylo results
Print function for ttree objects
## S3 method for class 'ttree' print(x, ...)
## S3 method for class 'ttree' print(x, ...)
x |
Object of class ttree, ie a transmission tree |
... |
Additional parameters are passed on |
Print out details of the ttree
print(extractTTree(simulateOutbreak()))
print(extractTTree(simulateOutbreak()))
Calculate the probability of a phylogenetic tree given a transmission tree
probPTreeGivenTTree(ctree, neg, w = integer(0))
probPTreeGivenTTree(ctree, neg, w = integer(0))
ctree |
Combined phylogenetic/transmission tree |
neg |
Within-host coalescent rate |
w |
Vector of hosts for which to calculate the probability, or nothing for all |
Probability of phylogeny given transmission tree
Calculate the probability of a phylogenetic tree given a transmission tree
probPTreeGivenTTreeR(ctree, neg, w = NULL)
probPTreeGivenTTreeR(ctree, neg, w = NULL)
ctree |
Combined phylogenetic/transmission tree |
neg |
Within-host coalescent rate |
w |
Vector of hosts for which to calculate the probability, or NULL for all |
Probability of phylogeny given transmission tree
Calculates the log-probability of a transmission tree
probTTree( ttree, rOff, pOff, pi, shGen, scGen, shSam, scSam, dateT, delta_t = 0.01 )
probTTree( ttree, rOff, pOff, pi, shGen, scGen, shSam, scSam, dateT, delta_t = 0.01 )
ttree |
Transmission tree |
rOff |
First parameter of the negative binomial distribution for offspring number |
pOff |
Second parameter of the negative binomial distribution for offspring number |
pi |
probability of sampling an infected individual |
shGen |
Shape parameter of the Gamma probability density function representing the generation time |
scGen |
Scale parameter of the Gamma probability density function representing the generation time |
shSam |
Shape parameter of the Gamma probability density function representing the sampling time |
scSam |
Scale parameter of the Gamma probability density function representing the sampling time |
dateT |
Date when process stops (this can be Inf for fully simulated outbreaks) |
delta_t |
Grid precision |
Probability of the transmission tree
Calculates the log-probability of a transmission tree
probTTreeR( ttree, off.r, off.p, pi, w.shape, w.scale, ws.shape, ws.scale, dateT )
probTTreeR( ttree, off.r, off.p, pi, w.shape, w.scale, ws.shape, ws.scale, dateT )
ttree |
Transmission tree |
off.r |
First parameter of the negative binomial distribution for offspring number |
off.p |
Second parameter of the negative binomial distribution for offspring number |
pi |
probability of sampling an infected individual |
w.shape |
Shape parameter of the Gamma probability density function representing the generation time |
w.scale |
Scale parameter of the Gamma probability density function representing the generation time |
ws.shape |
Shape parameter of the Gamma probability density function representing the sampling time |
ws.scale |
Scale parameter of the Gamma probability density function representing the sampling time |
dateT |
Date when process stops (this can be Inf for fully simulated outbreaks) |
Probability of the transmission tree
Converts an ape phylo object into a phylogenetic tree
ptreeFromPhylo(tr, dateLastSample)
ptreeFromPhylo(tr, dateLastSample)
tr |
phylo object |
dateLastSample |
date of the last sample |
phylogenetic tree
ptreeFromPhylo(ape::rtree(5),2020)
ptreeFromPhylo(ape::rtree(5),2020)
Select the most representative transmission tree from a MCMC output
selectTTree(record, burnin = 0.5)
selectTTree(record, burnin = 0.5)
record |
Output from inferTTree function |
burnin |
Proportion of the MCMC output to be discarded as burnin |
The index of the selected transmission tree
Simulate an outbreak
simulateOutbreak( off.r = 1, off.p = 0.5, neg = 0.25, nSampled = NA, pi = 0.5, w.shape = 2, w.scale = 1, ws.shape = NA, ws.scale = NA, w.mean = NA, w.std = NA, ws.mean = NA, ws.std = NA, dateStartOutbreak = 2000, dateT = Inf )
simulateOutbreak( off.r = 1, off.p = 0.5, neg = 0.25, nSampled = NA, pi = 0.5, w.shape = 2, w.scale = 1, ws.shape = NA, ws.scale = NA, w.mean = NA, w.std = NA, ws.mean = NA, ws.std = NA, dateStartOutbreak = 2000, dateT = Inf )
off.r |
First parameter of the negative binomial distribution for offspring number |
off.p |
Second parameter of the negative binomial distribution for offspring number |
neg |
the within-host effective population size (Ne) timesgeneration duration (g) |
nSampled |
number of sampled infected individuals, or NA for any |
pi |
probability of sampling an infected individual |
w.shape |
Shape parameter of the Gamma probability density function representing the generation time |
w.scale |
Scale parameter of the Gamma probability density function representing the generation time |
ws.shape |
Shape parameter of the Gamma probability density function representing the sampling time |
ws.scale |
Scale parameter of the Gamma probability density function representing the sampling time |
w.mean |
Mean of the Gamma distribution representing the generation time |
w.std |
Std of the Gamma distribution representing the generation time |
ws.mean |
Mean of the Gamma distribution representing the sampling time |
ws.std |
Std of the Gamma distribution representing the sampling time |
dateStartOutbreak |
Date when index case becomes infected |
dateT |
Date when process stops (this can be Inf for fully simulated outbreaks) |
Combined phylogenetic and transmission tree
simulateOutbreak() simulateOutbreak(off.r=2,dateStartOutbreak=2010,dateT=2015)
simulateOutbreak() simulateOutbreak(off.r=2,dateStartOutbreak=2010,dateT=2015)
Summary function for resTransPhylo objects
## S3 method for class 'resTransPhylo' summary(object, ...)
## S3 method for class 'resTransPhylo' summary(object, ...)
object |
output from inferTTree |
... |
Passed on to print.phylo |
Print out details of TransPhylo results