These function implement the bandle model for dynamic mass spectrometry based spatial proteomics datasets using MCMC for inference

These functions implement the bandle model for dynamic mass spectrometry based spatial proteomics datasets using MCMC for inference, this is an internal sampling function

bandle(
  objectCond1,
  objectCond2,
  fcol = "markers",
  hyperLearn = "LBFGS",
  numIter = 1000,
  burnin = 100L,
  thin = 5L,
  u = 2,
  v = 10,
  lambda = 1,
  gpParams = NULL,
  hyperIter = 20,
  hyperMean = c(0, 0, 0),
  hyperSd = c(1, 1, 1),
  seed = NULL,
  pg = FALSE,
  pgPrior = NULL,
  tau = 0.2,
  dirPrior = NULL,
  maternCov = TRUE,
  PC = TRUE,
  pcPrior = matrix(c(0.5, 3, 100), nrow = 1),
  nu = 2,
  propSd = c(0.3, 0.1, 0.05),
  numChains = 4L,
  BPPARAM = BiocParallel::bpparam()
)

diffLoc(
  objectCond1,
  objectCond2,
  fcol = "markers",
  hyperLearn = "MH",
  numIter = 1000,
  burnin = 100L,
  thin = 5L,
  u = 2,
  v = 10,
  lambda = 1,
  gpParams = NULL,
  hyperIter = 20,
  hyperMean = c(0, 0, 0),
  hyperSd = c(1, 1, 1),
  seed = NULL,
  pg = TRUE,
  pgPrior = NULL,
  tau = 0.2,
  dirPrior = NULL,
  maternCov = TRUE,
  PC = TRUE,
  nu = 2,
  pcPrior = NULL,
  propSd = c(0.3, 0.1, 0.05)
)

Arguments

objectCond1

A list of MSnbase::MSnSets where each is an experimental replicate for the first condition, usually a control

objectCond2

A list of MSnbase::MSnSets where each is an experimental replicate for the second condition, usually a treatment

fcol

The feature meta-data containing marker definitions. Default is markers

hyperLearn

Algorithm to learn posterior hyperparameters of the Gaussian processes. Default is LBFGS and MH for metropolis-hastings is also implemented.

numIter

The number of iterations of the MCMC algorithm. Default is 1000. Though usually much larger numbers are used

burnin

The number of samples to be discarded from the begining of the chain. Default is 100.

thin

The thinning frequency to be applied to the MCMC chain. Default is 5.

u

The prior shape parameter for Beta(u, v). Default is 2

v

The prior shape parameter for Beta(u, v). Default is 10.

lambda

Controls the variance of the outlier component. Default is 1.

gpParams

Object of class gpParams. parameters from prior fitting of GPs to each niche to accelerate inference. Default is NULL.

hyperIter

The frequency of MCMC interation to update the hyper-parameters default is 20

hyperMean

The prior mean of the log normal prior of the GP parameters. Default is 0 for each. Order is length-scale, amplitude and noise variance

hyperSd

The prior standard deviation of the log normal prior fo the GP parameters. Default is 1 for each. Order is length-scale, ampliture and noise variance.

seed

The random number seed.

pg

logical indicating whether to use polya-gamma prior. Default is FALSE.

pgPrior

A matrix generated by pgPrior function. If param pg is TRUE but pgPrior is NULL then a pgPrior is generated on the fly.

tau

The tau parameter for the polya-Gamma prior (is used). Defaults to 0.2

dirPrior

A matrix generated by dirPrior function. Default is NULL and dirPrior is generated on the fly.

maternCov

logical indicated whether to use a matern or gaussian covariance. Default is True and matern covariance is used

PC

logical indicating whether to use a penalised complexity prior. Default is TRUE.

pcPrior

matrix with 3 columns indicating the lambda paramters for the penalised complexity prior. Default is null which internally sets the penalised complexity prior to c(0.5, 3, 100) for each organelle and the order is length-scale, amplitude and variance. See vignette for more details.

nu

integer indicating the smoothness of the matern prior. Default is 2.

propSd

If MH is used to learn posterior hyperparameters then the proposal standard deviations. A Gaussian random-walk proposal is used.

numChains

integer indicating the number of parallel chains to run. Defaults to 4.

BPPARAM

BiocParallel parameter. Defaults to machine default backend using bpparam()

Value

bandle returns an instance of class bandleParams

bandle returns an instance of class bandleParams

Details

The bandle function generate the sample from the posterior distributions (object or class bandleParams) based on an annotated quantitative spatial proteomics datasets (object of class MSnbase::MSnSet). Both are then passed to the bandlePredict function to predict the sub-cellular localisation and compute the differential localisation probability of proteins. See the vignette for examples

The diffloc function generate the sample from the posterior distributions (object or class bandleParam) based on an annotated quantitative spatial proteomics datasets (object of class MSnbase::MSnSet). Both are then passed to the bandlePredict function to predict the sub-cellular localisation and compute the differential localisation probability of proteins. See the vignette for examples

Examples

library(pRolocdata)
data("tan2009r1")
set.seed(1)
tansim <- sim_dynamic(object = tan2009r1, 
                      numRep = 6L,
                      numDyn = 100L)
gpParams <- lapply(tansim$lopitrep, function(x) 
fitGPmaternPC(x, hyppar = matrix(c(0.5, 1, 100), nrow = 1)))
d1 <- tansim$lopitrep
control1 <- d1[1:3]
treatment1 <- d1[4:6]
mcmc1 <- bandle(objectCond1 = control1,
 objectCond2 = treatment1, gpParams = gpParams,
 fcol = "markers", numIter = 5L, burnin = 1L, thin = 2L,
 numChains = 1, BPPARAM = SerialParam(RNGseed = 1))
#> You haven't provided a seed, you may wish to provide a seed
#> 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |==============                                                        |  20%
  |                                                                            
  |============================                                          |  40%
  |                                                                            
  |==========================================                            |  60%

library(pRolocdata)
data("tan2009r1")
set.seed(1)
tansim <- sim_dynamic(object = tan2009r1, 
                    numRep = 6L,
                   numDyn = 100L)
gpParams <- lapply(tansim$lopitrep, 
function(x) fitGPmaternPC(x, hyppar = matrix(c(0.5, 1, 100), nrow = 1)))
d1 <- tansim$lopitrep
control1 <- d1[1:3]
treatment1 <- d1[4:6]
mcmc1 <- diffLoc(objectCond1 = control1, objectCond2 = treatment1, gpParams = gpParams,
                                     fcol = "markers", numIter = 5L, burnin = 1L, thin = 2L)
#> You haven't provided a seed, you may wish to provide a seed
#> 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |==============                                                        |  20%
  |                                                                            
  |============================                                          |  40%
  |                                                                            
  |==========================================                            |  60%