Skip to content

Commit

Permalink
Add parallel argument to summary and lrTest.
Browse files Browse the repository at this point in the history
Fixes #169
  • Loading branch information
amcdavid committed Dec 15, 2021
1 parent 0df7fcd commit fb37e77
Show file tree
Hide file tree
Showing 7 changed files with 67 additions and 27 deletions.
6 changes: 3 additions & 3 deletions DESCRIPTION
Original file line number Diff line number Diff line change
@@ -1,8 +1,8 @@
Package: MAST
Type: Package
Title: Model-based Analysis of Single Cell Transcriptomics
Version: 1.21.0
Date: 2021-05-10
Version: 1.21.1
Date: 2021-12-15
Authors@R: c(person("Andrew", "McDavid", email = "[email protected]", role = c("aut", "cre")),
person("Greg", "Finak", email="[email protected]", role='aut'),
person("Masanao", "Yajima", email="[email protected]", role='aut'))
Expand Down Expand Up @@ -88,7 +88,7 @@ Collate:
'thresholdSCRNA.R'
'zeroinf.R'
'zlmHooks.R'
RoxygenNote: 7.1.1
RoxygenNote: 7.1.2
LazyData: true
biocViews: GeneExpression, DifferentialExpression, GeneSetEnrichment,
RNASeq, Transcriptomics, SingleCell
Expand Down
3 changes: 2 additions & 1 deletion R/AllGenerics.R
Original file line number Diff line number Diff line change
Expand Up @@ -58,6 +58,7 @@ setGeneric('se.coef', function(object, ...) standardGeneric('se.coef'))
##' \code{hypothesis} can be one of a \code{character} giving complete factors or terms to be dropped from the model, \code{CoefficientHypothesis} giving names of coefficients to be dropped, \code{Hypothesis} giving contrasts using the symbolically, or a contrast \code{matrix}, with one row for each coefficient in the full model, and one column for each contrast being tested.
##' @param object LMlike or subclass
##' @param hypothesis the hypothesis to be tested. See details.
##' @param ... optional arguments, passed to fitting functions
##' @return array giving test statistics
##' @export
##' @seealso fit
Expand All @@ -68,7 +69,7 @@ setGeneric('se.coef', function(object, ...) standardGeneric('se.coef'))
##' @examples
##' #see ZlmFit-class for examples
##' example('ZlmFit-class')
setGeneric('lrTest', function(object, hypothesis) standardGeneric('lrTest'))
setGeneric('lrTest', function(object, hypothesis, ...) standardGeneric('lrTest'))

##' Run a Wald test
##'
Expand Down
34 changes: 20 additions & 14 deletions R/ZlmFit.R
Original file line number Diff line number Diff line change
Expand Up @@ -25,16 +25,16 @@ collectSummaries <- function(listOfSummaries){

## Refit using new MM and calculate stats
## hString is human-readable hypothesis
.lrtZlmFit <- function(zlmfit, newMM, hString){
.lrtZlmFit <- function(zlmfit, newMM, hString, ...){
o1 <- zlmfit
LMlike <- o1@LMlike
model.matrix(LMlike) <- newMM
message('Refitting on reduced model...')
# may not be necessary, but hopefully will keep backward compatibility with serialized ZlmFit
if('exprs_values' %in% slotNames(zlmfit)) {
o0 <- zlm(sca=o1@sca, LMlike=LMlike, exprs_values = o1@exprs_values)
o0 <- zlm(sca=o1@sca, LMlike=LMlike, exprs_values = o1@exprs_values, ...)
} else{
o0 <- zlm(sca=o1@sca, LMlike=LMlike)
o0 <- zlm(sca=o1@sca, LMlike=LMlike, ...)
}

lambda <- -2*(o0@loglik-o1@loglik)
Expand All @@ -57,8 +57,10 @@ collectSummaries <- function(listOfSummaries){
##' discrete, continuous and hurdle (combined) levels.
##' @param object ZlmFit
##' @param hypothesis See Details
##' @param ... passed to `zlm`
##' @inheritDotParams zlm
##' @return 3D array
setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='character'), function(object, hypothesis){
setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='character'), function(object, hypothesis, ...){
o1 <- object
LMlike <- o1@LMlike
oldMM <- model.matrix(LMlike)
Expand All @@ -73,34 +75,35 @@ setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='character'), funct
newq <- qr.Q(qr(newMM))
diff <- any(abs(newq %*% crossprod(newq, oldMM) - oldMM)>1e-6)
if(!diff) stop('Removing term ', sQuote(hypothesis), " doesn't actually alter the model, maybe due to marginality? Try specifying individual coefficents as a `CoefficientHypothesis`.")
.lrtZlmFit(o1, LMlike@modelMatrix, hypothesis)
.lrtZlmFit(o1, LMlike@modelMatrix, hypothesis, ...)
})

##' @describeIn ZlmFit Returns an array with likelihood-ratio tests on contrasts defined using \code{CoefficientHypothesis()}.
##' @param ... passed to `zlm`
##' @return see "Methods (by generic)"
setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='CoefficientHypothesis'), function(object, hypothesis){
setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='CoefficientHypothesis'), function(object, hypothesis, ...){
h <- generateHypothesis(hypothesis, colnames(object@coefD))
testIdx <- h@index
newMM <- model.matrix(object@LMlike)[,-testIdx, drop=FALSE]
.lrtZlmFit(object, newMM, hypothesis@.Data)
.lrtZlmFit(object, newMM, hypothesis@.Data, ...)
})

##' @describeIn ZlmFit Returns an array with likelihood-ratio tests specified by \code{Hypothesis}, which is a \link{Hypothesis}.
setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='Hypothesis'), function(object, hypothesis){
setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='Hypothesis'), function(object, hypothesis, ...){
## original fit
h <- generateHypothesis(hypothesis, colnames(object@coefD))
## call using coefficient matrix
lrTest(object, h@contrastMatrix)
lrTest(object, h@contrastMatrix, ...)
})

## contrast matrices
##' @describeIn ZlmFit Returns an array with likelihood-ratio tests specified by \code{Hypothesis}, which is a contrast \code{matrix}.
setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='matrix'), function(object, hypothesis){
setMethod('lrTest', signature=c(object='ZlmFit', hypothesis='matrix'), function(object, hypothesis, ...){
## original fit
LMlike <- object@LMlike
MM <- .rotateMM(LMlike, hypothesis)
testIdx <- attr(MM, 'testIdx')
.lrtZlmFit(object, MM[,-testIdx, drop=FALSE], 'Contrast Matrix')
.lrtZlmFit(object, MM[,-testIdx, drop=FALSE], 'Contrast Matrix', ...)
})


Expand Down Expand Up @@ -203,10 +206,11 @@ normalci <- function(center, se, level){
##' @param logFC If TRUE, calculate log-fold changes, or output from a call to \code{getLogFC}.
##' @param doLRT if TRUE, calculate lrTests on each coefficient, or a character vector of such coefficients to consider.
##' @param level what level of confidence coefficient to return. Defaults to 95 percent.
##' @inheritParams zlm
##' @param ... ignored
##' @return \code{data.table}
##' @seealso print.summaryZlmFit
##' @aliases "summary-ZlmFit", "ZlmFit-summary"
##' @aliases "summary-ZlmFit", "ZlmFit-summary", "summary.ZlmFit"
##' @examples
##' data(vbetaFA)
##' z <- zlm(~Stim.Condition, vbetaFA[1:5,])
Expand All @@ -215,8 +219,10 @@ normalci <- function(center, se, level){
##' print(zs)
##' ##Select `datatable` copmonent to get normal print method
##' zs$datatable
##' ## Can use parallel processing for LRT now
##' summary(z, doLRT = TRUE, parallel = TRUE)
##' @export
setMethod('summary', signature=c(object='ZlmFit'), function(object, logFC=TRUE, doLRT=FALSE, level=.95, ...){
setMethod('summary', signature=c(object='ZlmFit'), function(object, logFC=TRUE, doLRT=FALSE, level=.95, parallel = FALSE, ...){
message('Combining coefficients and standard errors')


Expand Down Expand Up @@ -255,7 +261,7 @@ setMethod('summary', signature=c(object='ZlmFit'), function(object, logFC=TRUE,
}
if(!is.logical(doLRT)){
message('Calculating likelihood ratio tests')
llrt <- lapply(doLRT, function(x) lrTest(object, CoefficientHypothesis(x))[,,'Pr(>Chisq)'])
llrt <- lapply(doLRT, function(x) lrTest(object, CoefficientHypothesis(x), parallel = parallel)[,,'Pr(>Chisq)'])
names(llrt) <- doLRT
llrt <- data.table(reshape2::melt(llrt))
setnames(llrt, c('test.type', 'L1', 'value'), c('component', 'contrast', 'Pr(>Chisq)'))
Expand Down
12 changes: 7 additions & 5 deletions man/ZlmFit-class.Rd

Some generated files are not rendered by default. Learn more about how customized files appear on GitHub.

19 changes: 18 additions & 1 deletion man/lrTest-ZlmFit-character-method.Rd

Some generated files are not rendered by default. Learn more about how customized files appear on GitHub.

4 changes: 3 additions & 1 deletion man/lrTest.Rd

Some generated files are not rendered by default. Learn more about how customized files appear on GitHub.

16 changes: 14 additions & 2 deletions man/summary-ZlmFit-method.Rd

Some generated files are not rendered by default. Learn more about how customized files appear on GitHub.

0 comments on commit fb37e77

Please sign in to comment.