Using the delta method, estimate the log-fold change from a state given by a vector contrast0 and the state(s) given by contrast1.
logFC(zlmfit, contrast0, contrast1) getLogFC(zlmfit, contrast0, contrast1)
Hypothesis
. If missing, then the '(Intercept)' is used as baseline.Hypothesis
. If missing, then all non-(Intercept) coefficients are compared.list of matrices `logFC` and `varLogFC`, giving the log-fold-changes for each contrast (columns) and genes (rows) and the estimated sampling variance thereof
The log-fold change is defined as follows. For each gene, let \(u(x)\) be the expected value of the continuous component, given a covariate x and the estimated coefficients coefC
, ie, \(u(x)=\) crossprod(x, coefC)
.
Likewise, Let \(v(x)=\) 1/(1+exp(-crossprod(coefD, x)))
be the expected value of the discrete component.
The log fold change from contrast0 to contrast1 is defined as
$$u(contrast1)v(contrast1)-u(contrast0)v(contrast0).$$
Note that for this to be a log-fold change, then the regression for u must have been fit on the log scale. This is returned in the matrix logFC
.
An approximation of the variance of logFC
(applying the delta method to formula defined above) is provided in varLogFC
.
getLogFC
: Return results as a perhaps friendlier data.table
1. When method='bayesglm'
(the default), it's no longer necessarily true that the log fold change from condition A to B will be the inverse of the log fold change from B to A if the models are fit separately.
This is due to the shrinkage in bayesglm
.
2. The log fold change can be small, but the Hurdle p-value small and significant when the sign of the discrete and continuous model components are discordant so that the marginal log fold change cancels out. The large sample sizes present in many single cell experiments also means that there is substantial power to detect even small changes.
#> #>##log-fold changes in terms of intercept (which is Stim(SEB) and CD154+VbetaResponsive) lfcStim <- logFC(zz) ##If we want to compare against unstim, we can try the following coefnames <- colnames(coef(zz, 'D')) contrast0 <- setNames(rep(0, length(coefnames)), coefnames) contrast0[c('(Intercept)', 'Stim.ConditionUnstim')] <- 1 contrast1 <- diag(length(coefnames)) rownames(contrast1)<-colnames(contrast1)<-coefnames contrast1['(Intercept)',]<-1 lfcUnstim <- logFC(zz, contrast0, contrast1) ##log-fold change with itself is 0 stopifnot(all(lfcUnstim$logFC[,2]==0)) ##inverse of log-fold change with Stim as reference stopifnot(all(lfcStim$logFC[,1]==(-lfcUnstim$logFC[,1]))) ##As a data.table: getLogFC(zz)#> primerid contrast logFC varLogFC #> 1: B3GAT1 Stim.ConditionUnstim 0.02121690 0.05699408 #> 2: B3GAT1 PopulationCD154-VbetaUnresponsive 0.09823225 0.11759250 #> 3: B3GAT1 PopulationCD154+VbetaResponsive 0.97402709 0.26295565 #> 4: B3GAT1 PopulationCD154+VbetaUnresponsive 0.50171226 0.33509425 #> 5: B3GAT1 PopulationVbetaResponsive 0.07638735 0.12164060 #> 6: B3GAT1 PopulationVbetaUnresponsive NaN NaN #> 7: BAX Stim.ConditionUnstim -0.42956548 1.35480633 #> 8: BAX PopulationCD154-VbetaUnresponsive -1.32492734 2.58024639 #> 9: BAX PopulationCD154+VbetaResponsive 0.62690532 2.65894498 #> 10: BAX PopulationCD154+VbetaUnresponsive -0.66503358 3.49255200 #> 11: BAX PopulationVbetaResponsive 1.70757502 3.10585739 #> 12: BAX PopulationVbetaUnresponsive 1.95779243 3.15029439 #> 13: BCL2 Stim.ConditionUnstim -0.45315003 0.10064027 #> 14: BCL2 PopulationCD154-VbetaUnresponsive 1.47536864 0.56464034 #> 15: BCL2 PopulationCD154+VbetaResponsive 1.70251880 0.57023184 #> 16: BCL2 PopulationCD154+VbetaUnresponsive 2.22719080 1.15781350 #> 17: BCL2 PopulationVbetaResponsive 3.13232752 1.14348444 #> 18: BCL2 PopulationVbetaUnresponsive NaN NaN #> 19: CCL2 Stim.ConditionUnstim 0.80636176 1.56877955 #> 20: CCL2 PopulationCD154-VbetaUnresponsive 0.42039146 0.53040577 #> z #> 1: 0.08887244 #> 2: 0.28646019 #> 3: 1.89945839 #> 4: 0.86670485 #> 5: 0.21901918 #> 6: NaN #> 7: -0.36905472 #> 8: -0.82482435 #> 9: 0.38445656 #> 10: -0.35585423 #> 11: 0.96892249 #> 12: 1.10303953 #> 13: -1.42842066 #> 14: 1.96342626 #> 15: 2.25458278 #> 16: 2.06984640 #> 17: 2.92922030 #> 18: NaN #> 19: 0.64379702 #> 20: 0.57723115 #> [ reached getOption("max.print") -- omitted 11 rows ]