View source: R/bandle-function.R
bandle | R Documentation |
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)
)
objectCond1 |
A list of |
objectCond2 |
A list of |
fcol |
The feature meta-data containing marker definitions. Default is
|
hyperLearn |
Algorithm to learn posterior hyperparameters of the Gaussian
processes. Default is |
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 |
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 |
|
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 |
|
PC |
|
pcPrior |
|
nu |
|
propSd |
If MH is used to learn posterior hyperparameters then the proposal standard deviations. A Gaussian random-walk proposal is used. |
numChains |
|
BPPARAM |
BiocParallel parameter. Defaults to machine default backend using bpparam() |
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
bandle
returns an instance of class bandleParams
bandle
returns an instance of class bandleParams
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))
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)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.