Skip to content

Commit

Permalink
quiet user supplied jabber and a fix
Browse files Browse the repository at this point in the history
  • Loading branch information
mikelove committed Jul 5, 2018
1 parent ea9df73 commit 911c45b
Show file tree
Hide file tree
Showing 4 changed files with 22 additions and 7 deletions.
2 changes: 1 addition & 1 deletion DESCRIPTION
Original file line number Diff line number Diff line change
Expand Up @@ -2,7 +2,7 @@ Package: DESeq2
Type: Package
Title: Differential gene expression analysis based on the negative
binomial distribution
Version: 1.21.7
Version: 1.21.8
Author: Michael Love, Simon Anders, Wolfgang Huber
Maintainer: Michael Love <[email protected]>
Description: Estimate variance-mean dependence in count data from
Expand Down
9 changes: 9 additions & 0 deletions NEWS
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,15 @@ CHANGES IN VERSION 1.22.0
o No replicate designs no longer supported (previous
version began deprecation with a warning).

CHANGES IN VERSION 1.21.8
-------------------------

o DESeq() now only says one time 'using supplied model matrix',
previously this was repeated three times from sub-functions.
Sub-functions therefore no longer print this message.
o Fixed bug when lfcShrink run directly after LRT
with supplied model matrices.

CHANGES IN VERSION 1.20.0
--------------------------

Expand Down
11 changes: 5 additions & 6 deletions R/core.R
Original file line number Diff line number Diff line change
Expand Up @@ -311,6 +311,8 @@ DESeq <- function(object, test=c("Wald","LRT"),
}
modelMatrix <- NULL
} else {
# model not as formula, so DESeq() is using supplied model matrix
if (!quiet) message("using supplied model matrix")
if (betaPrior == TRUE) {
stop("betaPrior=TRUE is not supported for user-provided model matrices")
}
Expand Down Expand Up @@ -602,8 +604,6 @@ estimateDispersionsGeneEst <- function(object, minDisp=1e-8, kappa_0=1,

if (is.null(modelMatrix)) {
modelMatrix <- getModelMatrix(object)
} else {
if (!quiet) message("using supplied model matrix")
}
checkFullRank(modelMatrix)
if (nrow(modelMatrix) == ncol(modelMatrix)) {
Expand Down Expand Up @@ -818,8 +818,6 @@ estimateDispersionsMAP <- function(object, outlierSD=2, dispPriorVar,

if (is.null(modelMatrix)) {
modelMatrix <- getModelMatrix(object)
} else {
if (!quiet) message("using supplied model matrix")
}

# fill in the calculated dispersion prior variance
Expand Down Expand Up @@ -1198,7 +1196,6 @@ nbinomWaldTest <- function(object,
if (betaPrior) {
if (missing(betaPriorVar)) stop("user-supplied model matrix with betaPrior=TRUE requires supplying betaPriorVar")
}
if (!quiet) message("using supplied model matrix")
modelAsFormula <- FALSE
attr(object, "modelMatrixType") <- "user-supplied"
renameCols <- FALSE
Expand Down Expand Up @@ -1613,7 +1610,6 @@ nbinomLRT <- function(object, full=design(object), reduced,
data=as.data.frame(colData(object)))
df <- ncol(fullModelMatrix) - ncol(reducedModelMatrix)
} else {
if (!quiet) message("using supplied model matrix")
df <- ncol(full) - ncol(reduced)
}

Expand Down Expand Up @@ -1691,6 +1687,9 @@ nbinomLRT <- function(object, full=design(object), reduced,
# record maximum of Cook's
maxCooks <- recordMaxCooks(design(object), colData(object), dispModelMatrix, cooks, nrow(objectNZ))

# store hat matrix diagonals
assays(object)[["H"]] <- buildMatrixWithNARows(H, mcols(object)$allZero)

# store Cook's distance for each sample
assays(object)[["cooks"]] <- buildMatrixWithNARows(cooks, mcols(object)$allZero)

Expand Down
7 changes: 7 additions & 0 deletions tests/testthat/test_lfcShrink.R
Original file line number Diff line number Diff line change
Expand Up @@ -86,4 +86,11 @@ test_that("LFC shrinkage works", {
res.normal <- lfcShrink(dds=dds, coef=2, res=res, type="normal")
res.ape <- lfcShrink(dds=dds, coef=2, res=res, type="apeglm")

# only running LRT upstream
dds <- makeExampleDESeqDataSet(m=4)
mm <- model.matrix(~condition, colData(dds))
mm0 <- model.matrix(~1, colData(dds))
dds <- DESeq(dds, full=mm, reduced=mm0, test="LRT")
res.normal <- lfcShrink(dds, coef="conditionB", type="normal")

})

0 comments on commit 911c45b

Please sign in to comment.