bandle 1.2.2
Bayesian ANalysis of Differential Localisation Experiments (BANDLE) is an integrative semi-supervised functional mixture model, developed by (Crook et al. 2021), to obtain the probability of a protein being differentially localised between two conditions.
In this vignette we walk users through how to install and use the R (R Development Core Team 2011)
Bioconductor (Gentleman et al. 2004) bandle
package
by simulating a well-defined differential localisation experiment from spatial
proteomics data from the pRolocdata
package (Gatto et al. 2014).
The BANDLE method uses posterior Bayesian computations performed using Markov-chain Monte-Carlo (MCMC) and thus uncertainty estimates are available (Gilks, Richardson, and Spiegelhalter 1995). Throughout this vignette we use the term differentially localised to pertain to proteins which are assigned to different sub-cellular localisations between two conditions. One of the main outputs of BANDLE is the probability that a protein is differentially localised between two conditions.
The package can be installed with the BiocManager
package:
if (!requireNamespace("BiocManager", quietly=TRUE))
install.packages("BiocManager")
BiocManager::install("bandle")
and then loaded,
library("bandle")
For visualisation we also load the packages,
library("pheatmap")
library("viridis")
In this vignette and Crook et al. (2021), the main data source that we use to study
differential protein sub-cellular localisation are data from high-throughput
mass spectrometry-based experiments. The data from these types of experiments
traditionally yield a matrix of measurements wherein we have, for example, PSMs,
peptides or proteins along the rows, and samples/channels/fractions along the
columns. The bandle
package uses the MSnSet
class as implemented in the
Bioconductor MSnbase package and thus requires users to import
and store their data as a MSnSet
instance. For more details on how to create a
MSnSet
see the relevant vignettes in pRoloc. The
pRolocdata experiment data package is a good starting place to
look for test data. This data package contains tens of quantitative proteomics
experiments, stored as MSnSet
s.
To get started with the basics of using bandle
we begin by generating a simple
example dataset which simulates a differential localisation experiment (please
see the second vignette in this package for a full real-life biological use
case). In this example data, the key elements are replicates, and a perturbation
of interest. There is code within the bandle package to simulate
an example experiment.
In the code chunk below we begin by loading the pRolocdata
package to obtain a spatial proteomics dataset. This will be the basis of our
simulation which will use boostrapping to generate new datasets. The dataset we
have chosen to load is a dataset from 2009 (tan2009r1
). This is data from a
early LOPIT experiment performed on Drosophila embryos (Tan et al. 2009). The aim of
this experiment was to apply LOPIT to an organism with heterogeneous cell types.
This experiment used four isotopes across four distinct fractions and thus
yielded four measurements (features) per protein profile. We visualise the
data by using principal components analysis.
library("pRolocdata")
data("tan2009r1")
## Let's set the stock colours of the classes to plot to be transparent
setStockcol(NULL)
setStockcol(paste0(getStockcol(), "90"))
## Plot the data using plot2D from pRoloc
plot2D(tan2009r1,
main = "An example spatial proteomics datasets",
grid = FALSE)
addLegend(tan2009r1, where = "topleft", cex = 0.7, ncol = 2)
The following code chuck simulates a differential localisation experiment. It
will generate numRep/2
of each a control and treatment condition. We will also
simulate relocalisations for numDyn
proteins.
set.seed(1)
tansim <- sim_dynamic(object = tan2009r1,
numRep = 6L,
numDyn = 100L)
## [1] "markers"
The list of the 6 simulated experiments are found in tansim$lopitrep
. Each one
is an MSnSet
instance (the standard data container for proteomics experimental
data). The first 3 are the simulated control experiments (see
tansim$lopitrep[1:3]
), and the following 3 in the list are the treatment
condition simulated experiments (see tansim$lopitrep[4:6]
). We can plot them
using the plot2D
function from pRoloc
.
plot_title <- c(paste0("Replicate ", seq(3), " condition", " A"),
paste0("Replicate ", seq(3), " condition", " B"))
par(mfrow = c(2, 3))
out <- lapply(seq(tansim$lopitrep), function(z)
plot2D(tansim$lopitrep[[z]], grid = FALSE, main = plot_title[z]))
For understanding, exploring and visualizing individual spatial proteomics
experiments, see the vignettes in pRoloc
and MSnbase
packages.
tansim$lopitrep[[1]]
## MSnSet (storageMode: lockedEnvironment)
## assayData: 888 features, 4 samples
## element names: exprs
## protocolData: none
## phenoData
## sampleNames: X114 X115 X116 X117
## varLabels: Fractions
## varMetadata: labelDescription
## featureData
## featureNames: P20353 P53501 ... P07909 (888 total)
## fvarLabels: FBgn Protein.ID ... knn.scores (18 total)
## fvarMetadata: labelDescription
## experimentData: use 'experimentData(object)'
## pubMedIds: 19317464
## Annotation:
## - - - Processing information - - -
## Added markers from 'mrk' marker vector. Thu Jul 16 22:53:44 2015
## Performed knn prediction (k=10) Thu Jan 19 15:55:26 2023
## MSnbase version: 1.17.12
bandle
analysisThe main function of the package is bandle
, this uses a complex model
to analyse the data. Markov-Chain Monte-Carlo (MCMC) is used to sample the
posterior distribution of parameters and latent variables. From which statistics
of interest can be computed. Here we only run a few iterations for brevity but
typically one needs to run thousands of iterations to ensure convergence, as
well as multiple parallel chains.
First, we need to fit non-parametric regression functions to the markers
profiles, upon which we place our analysis. This uses Gaussian processes. The
fitGPmaternPC
function can be used and fits some default penalised complexity
priors (see ?fitGP
), which works well. However, these can be altered, which is
demonstrated in the next code chunk
par(mfrow = c(4, 3))
gpParams <- lapply(tansim$lopitrep, function(x)
fitGPmaternPC(x, hyppar = matrix(c(10, 60, 250), nrow = 1)))
We apply the fitGPmaternPC
function to each datasets by calling lapply
over
the tansim$lopitrep
list of datasets. The output of fitGPmaternPC
returns a
list of posterior predictive means and standard deviations. As well as MAP
hyperparamters for the GP.
Note here we the use the default hyppar = matrix(c(10, 60, 250), nrow = 1)
as
a starting point for fitting the GPs to the marker profile distributions. In the
Crook et al. (2021) manuscript we found that these values worked well for smaller spatial
proteomics datasets. This was visually assessed by passing these values and
visualising the GP fit using the plotGPmatern
function.
The plotGPmatern
function can be used to plot the profiles for each
class in each replicate condition with the posterior predictive distributions
overlayed with the markers protein profiles.
For example, to plot the predictive distributions of the first dataset,
par(mfrow = c(4, 3))
plotGPmatern(tansim$lopitrep[[1]], params = gpParams[[1]])
The prior needs to form a K*3
matrix. K
corresponds to the number of
subcellular classes in the data, and 3 columns for (1) the prior, (2)
length-scale amplitude and (3) standard deviation parameters (see hyppar
in
?fitGPmaternPC
). Increasing these values, increases the shrinkage. For more
details see the manuscript by Crook et al. (2021). We strongly recommend users start with
the recommended hyppar
parameters and change and assess them as necessary for
their dataset by visually evaluating the fit of the GPs using the plotGPmatern
function.
K <- length(getMarkerClasses(tansim$lopitrep[[1]], fcol = "markers"))
pc_prior <- matrix(NA, ncol = 3, K)
pc_prior[seq.int(1:K), ] <- matrix(rep(c(10, 60, 250),
each = K), ncol = 3)
Now we have generated these complexity priors we can pass them as an
argument to the fitGPmaternPC
function. For example,
gpParams <- lapply(tansim$lopitrep,
function(x) fitGPmaternPC(x, hyppar = pc_prior))
By looking at the plot of posterior predictives using the gpParams
we can see
the GP fit looks sensible.
The next step is to set up the matrix Dirichlet prior on the mixing weights. These
weights are defined across datasets so these are slightly different to mixture
weights in usual mixture models. The \((i,j)^{th}\) entry is the prior probability
that a protein localises to organelle \(i\) in the control and \(j\) in the treatment.
This mean that off-diagonal terms have a different interpretation to diagonal terms.
Since we expect re-localisation to be rare, off-diagonal terms should be small.
The following functions help set up the priors and how to interpret them. The
parameter q
allow us to check the prior probability that more than q
differential localisations are expected.
set.seed(1)
dirPrior = diag(rep(1, K)) + matrix(0.001, nrow = K, ncol = K)
predDirPrior <- prior_pred_dir(object = tansim$lopitrep[[1]],
dirPrior = dirPrior,
q = 15)
The mean number of re-localisations is small:
predDirPrior$meannotAlloc
## [1] 0.3457542
The prior probability that more than q
differential localisations are
expected is small
predDirPrior$tailnotAlloc
## [1] 4e-04
The full prior predictive can be visualised as histogram. The prior probability that proteins are allocated to different components between datasets concentrates around 0.
hist(predDirPrior$priornotAlloc, col = getStockcol()[1])
For most use-cases we indeed expect the number of differential
localisations to be small. However, there may be specific cases where one may
expect the either a smaller or larger number of differential localisations.
Users could try testing different values for the dirPrior
for example,
replacing 0.001 with 0.0005 or smaller, for larger datasets to bring the number
of expected re-localisations inline with the biological expectation, and
visa-versa when we expect the number of proteins to have changed to be higher.
bandle
functionWe are now ready to run the main bandle
function. Remember to carefully
select the datasets and replicates that define the control and treatment.
As a reminder, in this introductory vignette we have used a small dataset
and generated theoretical triplicates of each theoretical condition. Please see
the second vignette in this package for a more detailed workflow and real
biological use-case. In the below code chunk we run bandle
for only 50
iterations for the convenience of building the vignette, but typically we’d
recommend you run the number of iterations (numIter
) in the 1000s.
Remember: the first 3 datasets are the first 3 elements of tansim
and the
final 3 elements are the “treatment” triplicate datasets.
control <- tansim$lopitrep[1:3]
treatment <- tansim$lopitrep[4:6]
bandleres <- bandle(objectCond1 = control,
objectCond2 = treatment,
numIter = 20, # usually 10,000
burnin = 5L, # usually 5,000
thin = 1L, # usually 20
gpParams = gpParams,
pcPrior = pc_prior,
numChains = 1, # usually >=4
dirPrior = dirPrior)
The bandle
function generates an object of class bandleParams
. The show
method indicates the number of parallel chains that were run, this should
typically be greater than 4 (here we use 1 just as a demo). Normally, we
should also assess convergence but this omitted for the moment so that we
can move forward with the analysis. Please see the end of the vignette for
convergence plots.
bandleres
## Object of class "bandleParams"
## Method: bandle
## Number of chains: 1
bandle
outputBefore we can begin to extract protein allocation information and a list of
proteins which are differentially localised between conditions, we first need to
populate the bandleres
object by calling the bandleProcess
function.
bandleres
objectCurrently, the summary slots of the bandleres
object are empty. The
summaries
function accesses them.
summaries(bandleres)
## [[1]]
## An object of class "bandleSummary"
## Slot "posteriorEstimates":
## <S4 Type Object>
## attr(,"elementType")
## [1] "ANY"
## attr(,"elementMetadata")
## `\001NULL\001`
## attr(,"metadata")
## list()
##
## Slot "diagnostics":
## <0 x 0 matrix>
##
## Slot "bandle.joint":
## <0 x 0 matrix>
##
##
## [[2]]
## An object of class "bandleSummary"
## Slot "posteriorEstimates":
## <S4 Type Object>
## attr(,"elementType")
## [1] "ANY"
## attr(,"elementMetadata")
## `\001NULL\001`
## attr(,"metadata")
## list()
##
## Slot "diagnostics":
## <0 x 0 matrix>
##
## Slot "bandle.joint":
## <0 x 0 matrix>
These can be populated as follows
bandleres <- bandleProcess(bandleres)
These slots have now been populated
summary(summaries(bandleres))
## Length Class Mode
## [1,] 1 bandleSummary S4
## [2,] 1 bandleSummary S4
bandle
resultsWe can save the results by calling summaries
. We see that it is
of length 2. 1 for control and 1 for treatment.
res <- summaries(bandleres)
length(res)
## [1] 2
There are a number of slots,
str(res[[1]])
## Formal class 'bandleSummary' [package "bandle"] with 3 slots
## ..@ posteriorEstimates:Formal class 'DFrame' [package "S4Vectors"] with 6 slots
## .. .. ..@ rownames : chr [1:677] "P20353" "P53501" "Q7KU78" "P04412" ...
## .. .. ..@ nrows : int 677
## .. .. ..@ elementType : chr "ANY"
## .. .. ..@ elementMetadata: NULL
## .. .. ..@ metadata : list()
## .. .. ..@ listData :List of 7
## .. .. .. ..$ bandle.allocation : Named chr [1:677] "PM" "Golgi" "Ribosome 40S" "PM" ...
## .. .. .. .. ..- attr(*, "names")= chr [1:677] "P20353" "P53501" "Q7KU78" "P04412" ...
## .. .. .. ..$ bandle.probability : Named num [1:677] 1 0.991 1 0.82 0.921 ...
## .. .. .. .. ..- attr(*, "names")= chr [1:677] "P20353" "P53501" "Q7KU78" "P04412" ...
## .. .. .. ..$ bandle.outlier : num [1:677] 0 0 0 0 0 0 0 0 0 0 ...
## .. .. .. ..$ bandle.probability.lowerquantile: Named num [1:677] 1 0.9812 1 0.0657 0.3553 ...
## .. .. .. .. ..- attr(*, "names")= chr [1:677] "P20353" "P53501" "Q7KU78" "P04412" ...
## .. .. .. ..$ bandle.probability.upperquantile: Named num [1:677] 1 0.998 1 0.98 0.996 ...
## .. .. .. .. ..- attr(*, "names")= chr [1:677] "P20353" "P53501" "Q7KU78" "P04412" ...
## .. .. .. ..$ bandle.mean.shannon : Named num [1:677] 0.00 0.00 8.78e-08 0.00 0.00 ...
## .. .. .. .. ..- attr(*, "names")= chr [1:677] "P20353" "P53501" "Q7KU78" "P04412" ...
## .. .. .. ..$ bandle.differential.localisation: Named num [1:677] 0 0 0 0.4 0.2 0 0 0 0 0 ...
## .. .. .. .. ..- attr(*, "names")= chr [1:677] "P20353" "P53501" "Q7KU78" "P04412" ...
## ..@ diagnostics : logi [1, 1] NA
## ..@ bandle.joint : num [1:677, 1:11] 5.31e-148 3.55e-138 1.98e-35 1.14e-123 2.97e-100 ...
## .. ..- attr(*, "dimnames")=List of 2
## .. .. ..$ : chr [1:677] "P20353" "P53501" "Q7KU78" "P04412" ...
## .. .. ..$ : chr [1:11] "Cytoskeleton" "ER" "Golgi" "Lysosome" ...
The main one of interest is the posteriorEstimates
slot,
posteriorEstimates(res[[1]])
## DataFrame with 677 rows and 7 columns
## bandle.allocation bandle.probability bandle.outlier
## <character> <numeric> <numeric>
## P20353 PM 0.999997 0
## P53501 Golgi 0.990784 0
## Q7KU78 Ribosome 40S 0.999995 0
## P04412 PM 0.820140 0
## Q7KJ73 Ribosome 60S 0.921478 0
## ... ... ... ...
## Q95TL8 ER 1.000000 0
## P25007 Proteasome 0.999994 0
## P41374 Cytoskeleton 0.999990 0
## Q8SZM1 Peroxisome 0.756867 0
## P07909 Nucleus 0.999998 0
## bandle.probability.lowerquantile bandle.probability.upperquantile
## <numeric> <numeric>
## P20353 0.9999951 0.999999
## P53501 0.9812102 0.997560
## Q7KU78 0.9999621 1.000000
## P04412 0.0657153 0.980140
## Q7KJ73 0.3553239 0.996129
## ... ... ...
## Q95TL8 1.000000 1.000000
## P25007 0.999980 0.999999
## P41374 0.999950 0.999999
## Q8SZM1 0.301359 0.979029
## P07909 0.999995 1.000000
## bandle.mean.shannon bandle.differential.localisation
## <numeric> <numeric>
## P20353 0.00000e+00 0.0
## P53501 0.00000e+00 0.0
## Q7KU78 8.78429e-08 0.0
## P04412 0.00000e+00 0.4
## Q7KJ73 0.00000e+00 0.2
## ... ... ...
## Q95TL8 0.0000e+00 0.0
## P25007 4.0882e-06 0.0
## P41374 0.0000e+00 0.0
## Q8SZM1 0.0000e+00 0.6
## P07909 0.0000e+00 0.0
This output object is a data.frame
containing the protein allocations and
associated localisation probabilities (including the upper and lower quantiles
of the allocation probability distribution), the mean Shannon entropy and the
bandle.differential.localisation
probability.
We create two new objects pe1
and pe2
in the below code chunk which contain
the output of the posteriorEstimates
slot.
pe1 <- posteriorEstimates(res[[1]])
pe2 <- posteriorEstimates(res[[2]])
One quantity of interest is the protein allocations, which we can plot as a barplot.
alloc1 <- pe1$bandle.allocation
alloc2 <- pe2$bandle.allocation
par(mfrow = c(1, 2), oma = c(6,2,2,2))
barplot(table(alloc1), col = getStockcol()[2],
las = 2, main = "Protein allocation: control")
barplot(table(alloc2), col = getStockcol()[2],
las = 2, main = "Protein allocation: treatment")
The barplot tells us for this example that bandle
has allocated the majority
of unlabelled proteins to the ER, followed by the Golgi (irrespective of the
posterior probability).
The associated posterior estimates are located in the bandle.probability
column.
pe_alloc1 <- pe1$bandle.probability
pe_alloc2 <- pe1$bandle.probability
The full allocation probabilities are stored in the tagm.joint
slot. These can
be visualised in a heatmap
bjoint_control <- bandleJoint(summaries(bandleres)[[1]])
pheatmap(bjoint_control, cluster_cols = FALSE, color = viridis(n = 25))
bjoint_treatment <- bandleJoint(summaries(bandleres)[[2]])
pheatmap(bjoint_treatment, cluster_cols = FALSE, color = viridis(n = 25))
We can append the results to our original MSnSet
datasets using the
bandlePredict
function.
xx <- bandlePredict(control,
treatment,
params = bandleres,
fcol = "markers")
res_control <- xx[[1]]
res_treatment <- xx[[2]]
The output is a list
of MSnSets
. In this example,
we have 3 for the control and 3 for the treatment.
length(res_control)
## [1] 3
length(res_treatment)
## [1] 3
The results are appended to the first MSnSet
feature data slot
for each condition.
fvarLabels(res_control[[1]])
## [1] "FBgn" "Protein.ID"
## [3] "Flybase.Symbol" "AccessionNo"
## [5] "EntryName" "AccessionNoAll"
## [7] "EntryNameAll" "No.peptide.IDs"
## [9] "Mascot.score" "No.peptide.quantified"
## [11] "PLSDA" "pd.2013"
## [13] "pd.markers" "markers.orig"
## [15] "markers" "markers.tl"
## [17] "knn" "knn.scores"
## [19] "bandle.allocation" "bandle.probability"
## [21] "bandle.probability.lowerquantile" "bandle.probability.upperquantile"
## [23] "bandle.mean.shannon" "bandle.differential.localisation"
## [25] "bandle.outlier" "bandle.joint"
To access them use the fData
function
fData(res_control[[1]])$bandle.probability
fData(res_control[[1]])$bandle.allocation
It is common practice in supervised machine learning to set a specific threshold
on which to define new assignments/allocations, below which classifications are
left unassigned/unknown. Indeed, we do not expect the whole subcellular
diversity to be represented by the 11 niches defined here, we expect there to be
many more, many of which will be multiply localised within the cell. It is
important to allow for the possibility of proteins to reside in multiple
locations (this information is available in the bandle.joint
slot - see above
for more details on multiple location).
As we are using a Bayesian model the outputs of the classifier are probabilities. This not only allows us to look at the distribution of probabilities over all subcellular classes but also allows us to extract a probability threshold on which we can define new assignments.
The subcellular allocations are located in the bandle.allocation
column of the
fData
slot and the posteriors are located in the bandle.probability
slot. We
can use the getPredictions
function from the pRoloc
package to return a set
of predicted localisations according to if they meet a probability threshold.
For example, in the below code chunk we set a 1% FDR for assigning proteins a subcellular nice, below which we leave them unlabelled.
res_control[[1]] <- getPredictions(res_control[[1]],
fcol = "bandle.allocation",
scol = "bandle.probability",
mcol = "markers",
t = .99)
## ans
## Cytoskeleton ER Golgi Lysosome Nucleus
## 15 220 135 9 39
## PM Peroxisome Proteasome Ribosome 40S Ribosome 60S
## 133 7 42 73 46
## mitochondrion unknown
## 81 88
res_treatment[[1]] <- getPredictions(res_treatment[[1]],
fcol = "bandle.allocation",
scol = "bandle.probability",
mcol = "markers",
t = .99)
## ans
## Cytoskeleton ER Golgi Lysosome Nucleus
## 19 191 125 10 41
## PM Peroxisome Proteasome Ribosome 40S Ribosome 60S
## 138 14 46 79 46
## mitochondrion unknown
## 84 95
We may also wish to take into account the probability of the protein being an
outlier and thus use the results in the bandle.outlier
column of the feature
data. We could calculate the product of the posterior and the outlier (as they
are both probabilities) to obtain a localisation score that takes into account
the outlier model. More details on this are found in the second vignette of this
package.
As previously mentioned the term “differentially localised” is used to pertain to proteins which are assigned to different sub-cellular localisations between two conditions. For the majority of users this is the main output they are keen to extract using the BANDLE method.
Following on from the above example, after extracting posterior estimates for
each condition using the summaries
function we can also access the
differential localisation probability as it is stored in the
bandle.differential.localisation
column of the data.frames
of pe1
and
pe2
, in the above sections.
The differential localisation probability tells us which proteins are most likely to differentially localise. We can for example, examine how many proteins get a differential probability greater than 0.99 to look for the most confident differentially localised candidates.
diffloc_probs <- pe1$bandle.differential.localisation
head(diffloc_probs, 50)
## P20353 P53501 Q7KU78 P04412 Q7KJ73 Q9VM65 Q9VCK0 B7Z0W3
## 0.0000000 0.0000000 0.0000000 0.4000000 0.2000000 0.0000000 0.0000000 0.0000000
## Q9V415 Q00174 Q9V769 Q27593 Q9V780 P19109 Q95NR4 O77047
## 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000
## Q8SXD0 Q9VBV5 Q9VBU5 Q8IA62 Q9Y105 Q9W2M4 P26308 Q9VP77
## 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 1.0000000 1.0000000
## Q5U0Z2 Q960E2 Q9VLJ6 Q9VJ39 Q9VTX8 Q9VTZ5 B7Z0C1 Q9VRJ4
## 0.0000000 0.0000000 1.0000000 1.0000000 0.0000000 0.6000000 0.0000000 0.0000000
## M9PCB7 P46150 A1ZBH5 Q9W3M7 A8DZ29 Q9VN14 Q9VZL3 M9PC99
## 0.0000000 0.0000000 0.0000000 0.0000000 1.0000000 0.1333333 0.0000000 0.4666667
## Q86BP3 Q9W3N9 Q7K5M6 P16620 P48375 Q9VMB9 Q9VI55 Q9VU58
## 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.2000000
## Q9V498 Q9V496
## 0.0000000 1.0000000
length(which(diffloc_probs[order(diffloc_probs, decreasing = TRUE)] > 0.99))
## [1] 88
We find there are 88 proteins above this threshold.
This can also be seen on a rank plot
plot(diffloc_probs[order(diffloc_probs, decreasing = TRUE)],
col = getStockcol()[3], pch = 19, ylab = "Probability",
xlab = "Rank", main = "Differential localisation rank plot")
In-line with our expectations, the rank plot indicates that most proteins are not differentially localised.
bootstrapdiffLocprob
functionWe can examine the top n
proteins (here we use an example of top = 100
) and
produce bootstrap estimates of the uncertainty (note here the uncertainty is
likely to be underestimated as we did not produce many MCMC samples). These can
be visualised as ranked boxplots.
set.seed(1)
boot_t <- bootstrapdiffLocprob(params = bandleres, top = 100,
Bootsample = 5000, decreasing = TRUE)
boxplot(t(boot_t), col = getStockcol()[5],
las = 2, ylab = "Probability", ylim = c(0, 1),
main = "Differential localisation \nprobability plot (top 100 proteins)")
binomDiffLoc
functionInstead of applying the bootstrapdiffLocprob
we could use the binomDiffLoc
function to obtain credible intervals from the binomial distribution.
bin_t <- binomialDiffLocProb(params = bandleres, top = 100,
nsample = 5000, decreasing = TRUE)
boxplot(t(bin_t), col = getStockcol()[5],
las = 2, ylab = "Probability", ylim = c(0, 1),
main = "Differential localisation \nprobability plot (top 100 proteins)")
There are many ways we could obtain probability estimates from either of the above methods. We could, for example, take the mean of each protein estimate, or compute the cumulative error (there is not really a false discovery rate in Bayesian statistics) or we could threshold on the interval to reduce the number of differential localisations if you feel the model has been overconfident.
# more robust estimate of probabilities
dprobs <- rowMeans(bin_t)
# compute cumulative error, there is not really a false discovery rate in
# Bayesian statistics but you can look at the cumulative error rate
ce <- cumsum(1 - dprobs)
# you could threshold on the interval and this will reduce the number of
# differential localisations
qt <- apply(bin_t, 1, function(x) quantile(x, .025))
Instead of estimating the false discovery rate we can estimate the expected false discovery rate from the posterior probabilities at a particular threshold. This mean that for fixed threshold, we compute the expected proportion of false discoveries. Here is an example below. We can see that setting a probability threshold of 0.9 leads to an expected false discovery rate of less than \(0.5\%\)
EFDR(diffloc_probs, threshold = 0.90)
## [1] 0.002197802
We can visualise the changes in localisation between conditions on an alluvial
plot using the plotTranslocations
function
plotTranslocations(bandleres)
Or alternatively, on a chord (circos) diagram
plotTranslocations(bandleres, type = "chord")
These visualisations are showing the change in class label between the two
conditions (as assigned by bandle
i.e. the result stored in
bandle.allocation
). The results are taken directly from bandleres
and thus
no thresholding on the class label and posterior to allow for proteins to be
left “unassigned” or unknown, is conducted. Furthermore, there is no
thresholding on the bandle.differential.localisation
probability.
It would be better to re-plot these figures with some thresholds on the above
quantities, to get a better representation of what is moving. The easiest way to
do this is to pass the MSnSets
output after performing bandlePredict
and
getPredictions
.
For example, first let’s identify which proteins get a high differential localisation probability,
## identify which proteins get a high differential localisation probability
ind <- which(fData(res_control[[1]])$bandle.differential.localisation > 0.99)
## create two new MSnSets with only these proteins
res_dl_control <- res_control[[1]][ind, ]
res_dl_treatment <- res_treatment[[1]][ind, ]
Now we can plot only these 88
proteins which are deemed to move/differentially localise. We also specify where the prediction results are
located e.g. fcol = "bandle.allocation.pred"
.
## specify colour palette
mycols <- c(getStockcol()[seq(getMarkerClasses(res_control[[1]]))], "grey")
names(mycols) <- c(getMarkerClasses(res_control[[1]]), "unknown")
## Create a list of the datasets for plotTranslocations
res <- list(res_dl_control, res_dl_treatment)
plotTranslocations(res, fcol = "bandle.allocation.pred", col = mycols)
## 88 features in common
## ------------------------------------------------
## If length(fcol) == 1 it is assumed that the
## same fcol is to be used for both datasets
## setting fcol = c(bandle.allocation.pred,bandle.allocation.pred)
## ----------------------------------------------
We can also use the function plotTable
to display a summary table of the
number of proteins that have changed in location between conditions.
(tbl <- plotTable(res, fcol = "bandle.allocation.pred"))
## Condition1 Condition2 value
## 1 ER Golgi 1
## 2 ER Lysosome 2
## 3 ER Nucleus 1
## 4 ER PM 3
## 5 ER Peroxisome 7
## 6 ER Proteasome 4
## 7 ER Ribosome 40S 4
## 8 ER Ribosome 60S 1
## 10 ER unknown 5
## 11 ER Cytoskeleton 2
## 12 Golgi ER 1
## 14 Golgi Nucleus 3
## 19 Golgi Ribosome 60S 1
## 20 Golgi mitochondrion 2
## 21 Golgi unknown 4
## 31 Lysosome mitochondrion 1
## 37 Nucleus PM 1
## 38 Nucleus Peroxisome 1
## 46 PM Golgi 1
## 50 PM Proteasome 2
## 51 PM Ribosome 40S 4
## 55 PM Cytoskeleton 1
## 57 Peroxisome Golgi 1
## 67 Proteasome ER 1
## 71 Proteasome PM 1
## 73 Proteasome Ribosome 40S 2
## 74 Proteasome Ribosome 60S 1
## 75 Proteasome mitochondrion 2
## 76 Proteasome unknown 2
## 79 Ribosome 40S Golgi 1
## 81 Ribosome 40S Nucleus 1
## 82 Ribosome 40S PM 1
## 86 Ribosome 40S mitochondrion 1
## 87 Ribosome 40S unknown 2
## 95 Ribosome 60S Proteasome 1
## 96 Ribosome 60S Ribosome 40S 2
## 106 mitochondrion Proteasome 1
## 109 mitochondrion unknown 1
## 110 mitochondrion Cytoskeleton 1
## 111 unknown ER 1
## 112 unknown Golgi 3
## 115 unknown PM 1
## 117 unknown Proteasome 4
## 118 unknown Ribosome 40S 1
## 120 unknown mitochondrion 1
## 121 unknown Cytoskeleton 1
In this section, we demonstrate how to visually assess convergence of the MCMC algorithm. In the chunk below, we use 4 chains so that we can assess the convergence of the method.
control <- tansim$lopitrep[1:3]
treatment <- tansim$lopitrep[4:6]
bandleres <- bandle(objectCond1 = control,
objectCond2 = treatment,
numIter = 50, # usually 10,000
burnin = 10L, # usually 5,000
thin = 1L, # usually 20
gpParams = gpParams,
pcPrior = pc_prior,
numChains = 4,
dirPrior = dirPrior)
We then use the plotConvergence
function which plots the ranks of the total
number of outliers within in each of the chains. If convergence has been reached
these plots should be uniform. As we can see these plots are not uniform and
are skewed towards the extremes, so they have not converged. Clearly one
of the chains has higher values on average than the other chains producing
the skew. Running the algorithm for more iterations (typically 10,000)
should produce convergence. This is not done here for brevity of the vignette.
par(mfrow = c(2, 2))
out <- plotConvergence(bandleres)
bandle
parametersThe bandle
function has a significant number of parameters to allow flexible
and bespoke analysis. Here, we describe these parameters in more detail to
allow user to make decisions on the level of flexibility they wish to exploit.
objectCond1
. This is a list of MSnSets
containing the first condition.
objectCond2
. This is a list of MSnSets
containing the second condition.
fcol
indicates the feature column in the MSnSets
that indicated the
markers. Proteins that are not markers should be labels unknown
. The default
is markers
.
hyperLearn
is the algorithm used to learn the hyperparameters of the
Gaussian processes. For speed the default is an optimization algorithm called
“LBFGS”, however is users want to perform uncertainty quantification on these
parameters we can use Markov-chain Monte Carlo (MCMC) methods. This is implemented
using the Metropolis-Hastings algorithm. Though this latter methodology provides
more information, it is much more costly. The analysis is expected to take
several days rather than hours.
numIter
is the number of MCMC iterations for the algorithm. We typically
suggest around 10,000 iterations is plenty for convergence. Though some cases
may take longer. If resources are constrained, we suggest 4,000 iterations
as acceptable. A minimum number of iterations is around 1,000 though at this
level we expect the posterior estimates to suffer considerably. If possible
more parallel chains should be run in this case by changing numChains
to,
say, 9. The more chains and iterations the more computationally expensive
the algorithm. The time taken for the algorithm scales roughly linearly
in the number of iterations
burnin
is the number of samples that should be discarded from the
beginning of the chain due to the bias induced by the starting point of the
algorithm. We suggest sensible burnin
values to be roughly \(10-50\%\) of the
number of iterations
thin
reduces auto-correlation in the MCMC samples. The default is \(5\),
which means every 5th sample is taken. If memory requirements are an issue,
we suggest to increase the thinning amount. Though above \(20\), you will see
a decrease in performance.
u
and v
represent the prior hyperparameters of the proportion of outliers.
This is modelled using a Beta(u,v)
with u = 2
and v = 10
a default. This
suggest that roughly \(\frac{u}{u = V} = 16%\) of proteins are believed to be
outliers and that it is quite unlikely that more than \(50%\) of proteins
are outliers. Users can examine the quantiles of the Beta(u,v)
distribution
if they wish to place a more bespoke prior. For example, increasing u
will increase the number of a prior believed outliers.
lambda
is a ridge parameter used for numerical stability and is set to
\(0.01\). If you experience the algorithm fails due to numerical issue then you
can set this value larger. If you require values above \(1\) it is likely that
there are other issues with the analysis. We suggest checking the method
is appropriate for your problem and opening issue detailing the problems.
gpParams
results from fitting Gaussian proccess (Gaussian random fields).
We refer the users to those functions. The default is NULL
which will fit
GPs internally but we recommend setting these outside the bandle function
because it leads to more stable results.
hyperIter
if the hyperparameters of the GP are learnt using MH algorithm
then this is the frequency at which these are updated relative to the bandle
algorithm. By default this is unused, but if hyperLearn
is set to MH
then
this proceed at every 20 iterations.
hyperMean
is the mean of the log normal prior used on the hyperparameters.
Though by default this is not used unless PC
is set to false
hyperSd
is the standard deviation of the log normal prior used on the
hyperparameters. The default is c(1,1,1)
for the 3 hyperparameters, increasing
these values increases the uncertainty in the prior values of the hyperparameters.
seed
is the random-number seed.
pg
indicates whether or not to use the Polya-Gamma (PG) prior. The default
is false and a Dirichlet prior is used instead. If set to true the pg
is used.
In which case a default PG prior is used. This prior attempts to match the
default Dirichlet prior that is used when PG prior is set to false. The PG
prior is more computationally expensive but can provide prior information on
correlations
pgPrior
is by default NULL. We suggest using the pg_prior
function
to help set this parameter and the documentation therein. This function
uses an empirical approach to compute a sensible default.
tau
is a parameter used by the Polya-Gamma prior and we refer to
BANDLE manuscript for details. By default it is only used if pg
prior is true,
when the default becomes 0.2
. At this value the pg
prior is similar to the
Dirichlet prior but with information on correlations.
dirPrior
is the Dirichlet matrix prior on the correlations. This should
be provided as a K by K matrix, where K is the number of subcellular niches.
The diagonal component should represent the prior belief that organelles do
not re-localise (same compartment), where as the off-diagonal terms represent
the prior terms of re-localisation. The prior_pred_dir
can be used to provide
a prior predictive check based on the provided prior. It is recommended
that the off-diagonal terms are at least two orders of magnitude smaller than the
diagonal terms. An example is given in the vignette.
maternCov
is this true the covariance function is the matern covariance,
otherwise a Gaussian covariance is used.
PC
indicates whether a penalised complexity (PC) is used. The default
is true and otherwise log normal priors are used.
pcPrior
is a numeric of length 3 indicating the parameters of the PC prior.
The prior is placed on the parameters of length-scale, amplitude, and variance
in that order. The default values are \(0.5,3,100\), and increasing the value
increases the shrinkage towards straight-lines with zero variance.
nu
which defaults to 2 is the smoothness of the matern covariance. By
increasing nu
you encourage smoother solutions. nu
should be an integer,
though for values of nu
above 3, we have observed numerical instability.
propSd
is the standard deviation of the random-walk update used in the MH
algorithm. We do not recommend changing this unless you are familiar with
Bayesian analysis. The default is c(0.3,0.1,0.05)
for the 3 hyperparameters.
Changing these will alter the efficiency of the underlying samplers.
numChains
is the number of parrallel chains and defaults to 4. We recommend
using as much processing resources as you have and frequently have used 9 in
practise.
BPPARAM
is the BiocParallel back-end which defaults to
BiocParallel::bpparam()
. We refer you to the BiocParallel
package for details on setting this dependent on your computing system.
All software and respective versions used to produce this document are listed below.
sessionInfo()
## R version 4.2.2 (2022-10-31)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] pRolocdata_1.36.0 viridis_0.6.2 viridisLite_0.4.1
## [4] pheatmap_1.0.12 bandle_1.2.2 pRoloc_1.38.2
## [7] BiocParallel_1.32.5 MLInterfaces_1.78.0 cluster_2.1.4
## [10] annotate_1.76.0 XML_3.99-0.13 AnnotationDbi_1.60.0
## [13] IRanges_2.32.0 MSnbase_2.24.2 ProtGenerics_1.30.0
## [16] mzR_2.32.0 Rcpp_1.0.9 Biobase_2.58.0
## [19] S4Vectors_0.36.1 BiocGenerics_0.44.0 BiocStyle_2.26.0
##
## loaded via a namespace (and not attached):
## [1] utf8_1.2.2 tidyselect_1.2.0
## [3] RSQLite_2.2.20 htmlwidgets_1.6.1
## [5] grid_4.2.2 lpSolve_5.6.17
## [7] pROC_1.18.0 munsell_0.5.0
## [9] codetools_0.2-18 preprocessCore_1.60.2
## [11] future_1.30.0 withr_2.5.0
## [13] colorspace_2.0-3 filelock_1.0.2
## [15] highr_0.10 knitr_1.41
## [17] ggalluvial_0.12.3 mzID_1.36.0
## [19] listenv_0.9.0 MatrixGenerics_1.10.0
## [21] GenomeInfoDbData_1.2.9 farver_2.1.1
## [23] bit64_4.0.5 coda_0.19-4
## [25] parallelly_1.34.0 vctrs_0.5.1
## [27] generics_0.1.3 ipred_0.9-13
## [29] xfun_0.36 timechange_0.2.0
## [31] BiocFileCache_2.6.0 randomForest_4.7-1.1
## [33] R6_2.5.1 doParallel_1.0.17
## [35] GenomeInfoDb_1.34.7 clue_0.3-63
## [37] MsCoreUtils_1.10.0 bitops_1.0-7
## [39] cachem_1.0.6 DelayedArray_0.24.0
## [41] assertthat_0.2.1 scales_1.2.1
## [43] nnet_7.3-18 gtable_0.3.1
## [45] affy_1.76.0 globals_0.16.2
## [47] timeDate_4022.108 rlang_1.0.6
## [49] GlobalOptions_0.1.2 splines_4.2.2
## [51] lazyeval_0.2.2 ModelMetrics_1.2.2.2
## [53] impute_1.72.3 hexbin_1.28.2
## [55] BiocManager_1.30.19 yaml_2.3.6
## [57] reshape2_1.4.4 caret_6.0-93
## [59] tools_4.2.2 lava_1.7.1
## [61] bookdown_0.32 ggplot2_3.4.0
## [63] affyio_1.68.0 ellipsis_0.3.2
## [65] jquerylib_0.1.4 RColorBrewer_1.1-3
## [67] proxy_0.4-27 plyr_1.8.8
## [69] progress_1.2.2 zlibbioc_1.44.0
## [71] purrr_1.0.1 RCurl_1.98-1.9
## [73] prettyunits_1.1.1 rpart_4.1.19
## [75] sampling_2.9 SummarizedExperiment_1.28.0
## [77] LaplacesDemon_16.1.6 ggrepel_0.9.2
## [79] magrittr_2.0.3 magick_2.7.3
## [81] data.table_1.14.6 circlize_0.4.15
## [83] pcaMethods_1.90.0 mvtnorm_1.1-3
## [85] matrixStats_0.63.0 hms_1.1.2
## [87] evaluate_0.20 xtable_1.8-4
## [89] mclust_6.0.0 shape_1.4.6
## [91] gridExtra_2.3 compiler_4.2.2
## [93] biomaRt_2.54.0 tibble_3.1.8
## [95] ncdf4_1.21 crayon_1.5.2
## [97] htmltools_0.5.4 segmented_1.6-2
## [99] tidyr_1.2.1 lubridate_1.9.0
## [101] DBI_1.1.3 dbplyr_2.3.0
## [103] MASS_7.3-58.1 rappdirs_0.3.3
## [105] Matrix_1.5-3 cli_3.6.0
## [107] vsn_3.66.0 gdata_2.18.0.1
## [109] parallel_4.2.2 gower_1.0.1
## [111] GenomicRanges_1.50.2 pkgconfig_2.0.3
## [113] plotly_4.10.1 recipes_1.0.4
## [115] MALDIquant_1.22 xml2_1.3.3
## [117] foreach_1.5.2 bslib_0.4.2
## [119] hardhat_1.2.0 XVector_0.38.0
## [121] prodlim_2019.11.13 stringr_1.5.0
## [123] digest_0.6.31 Biostrings_2.66.0
## [125] rmarkdown_2.19 dendextend_1.16.0
## [127] curl_5.0.0 kernlab_0.9-31
## [129] gtools_3.9.4 lifecycle_1.0.3
## [131] nlme_3.1-161 jsonlite_1.8.4
## [133] limma_3.54.0 fansi_1.0.3
## [135] pillar_1.8.1 lattice_0.20-45
## [137] lbfgs_1.2.1.2 KEGGREST_1.38.0
## [139] fastmap_1.1.0 httr_1.4.4
## [141] survival_3.5-0 glue_1.6.2
## [143] FNN_1.1.3.1 png_0.1-8
## [145] iterators_1.0.14 bit_4.0.5
## [147] class_7.3-20 stringi_1.7.12
## [149] sass_0.4.4 mixtools_2.0.0
## [151] blob_1.2.3 memoise_2.0.1
## [153] dplyr_1.0.10 e1071_1.7-12
## [155] future.apply_1.10.0
Crook, Oliver M, Colin TR Davies, Laurent Gatto, Paul DW Kirk, and Kathryn S Lilley. 2021. “Inferring Differential Subcellular Localisation in Comparative Spatial Proteomics Using Bandle.” bioRxiv. doi: https://doi.org/10.1101/2021.01.04.425239.
Gatto, Laurent, Lisa M. Breckels, Samuel Wieczorek, Thomas Burger, and Kathryn S. Lilley. 2014. “Mass-Spectrometry Based Spatial Proteomics Data Analysis Using pRoloc and pRolocdata.” Bioinformatics.
Gentleman, Robert C., Vincent J. Carey, Douglas M. Bates, Ben Bolstad, Marcel Dettling, Sandrine Dudoit, Byron Ellis, et al. 2004. “Bioconductor: Open Software Development for Computational Biology and Bioinformatics.” Genome Biol 5 (10): –80. https://doi.org/10.1186/gb-2004-5-10-r80.
Gilks, Walter R, Sylvia Richardson, and David Spiegelhalter. 1995. Markov Chain Monte Carlo in Practice. CRC press.
R Development Core Team. 2011. R: A Language and Environment for Statistical Computing. Vienna, Austria: R Foundation for Statistical Computing. http://www.R-project.org/.
Tan, Denise JL, Heidi Dvinge, Andrew Christoforou, Paul Bertone, Alfonso Martinez Arias, and Kathryn S Lilley. 2009. “Mapping Organelle Proteins and Protein Complexes in Drosophila Melanogaster.” Journal of Proteome Research 8 (6): 2667–78.