--- title: "Introduction to zFPKM Transformation" author: "Ron Ammar" date: "`r Sys.Date()`" output: rmarkdown::html_vignette vignette: > %\VignetteIndexEntry{Introduction to zFPKM Transformation} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} --- ```{r, reproducible_settings, echo=FALSE} rm(list=ls()) set.seed(8675309) options(stringsAsFactors=FALSE) library(printr) knitr::opts_chunk$set(message=FALSE, fig.align="center", fig.width=8, fig.height=6) ``` ## Summary Perform the zFPKM transform on RNA-seq FPKM data. This algorithm is based on the publication by [Hart et al., 2013 (Pubmed ID 24215113)](https://www.ncbi.nlm.nih.gov/pubmed/24215113). The reference recommends using zFPKM > -3 to select expressed genes. Validated with ENCODE open/closed promoter chromatin structure epigenetic data on six of the ENCODE cell lines. It works well for gene level data using FPKM or TPM, but does not appear to calibrate well for transcript level data. ## Identifying *active* genes for subsequent analysis We calculate zFPKM for existing FPKM from [gse94802](https://www.ncbi.nlm.nih.gov/geo/query/acc.cgi?acc=gse94802). ```{r} library(dplyr) library(GEOquery) library(stringr) library(SummarizedExperiment) library(tidyr) getSpecificGEOSupp <- function(url) { temp <- tempfile() download.file(url, temp) out <- read.csv(gzfile(temp), row.names=1, check.names=FALSE) out <- select(out, -MGI_Symbol) return(out) } gse94802_fpkm <- "ftp://ftp.ncbi.nlm.nih.gov/geo/series/GSE94nnn/GSE94802/suppl/GSE94802_Minkina_etal_normalized_FPKM.csv.gz" gse94802_counts <- "ftp://ftp.ncbi.nlm.nih.gov/geo/series/GSE94nnn/GSE94802/suppl/GSE94802_Minkina_etal_raw_counts.csv.gz" if (file.exists("gse94802.rds")) { esetlist <- readRDS("gse94802.rds") } else { esetlist <- getGEO("gse94802") } doe <- pData(esetlist[[1]]) colData <- DataFrame( condition=ifelse(str_detect(doe$title, regex("control", ignore_case=TRUE)), "control", "mutant"), sample_id=str_match(doe$title, "rep\\d_(.+)")[, 2], row.names=str_match(doe$title, "rep\\d_(.+)")[, 2]) se <- SummarizedExperiment(assays=SimpleList(fpkm=getSpecificGEOSupp(gse94802_fpkm), counts=getSpecificGEOSupp(gse94802_counts)), colData=colData) # clear namespace rm(esetlist, gse94802_fpkm, gse94802_counts, doe, colData, getSpecificGEOSupp) ``` We compute zFPKM. ```{r} library(zFPKM) assay(se, "zfpkm") <- zFPKM(se) ``` We can also plot the Guassian fit to the FPKM data for which the z-scores are based. ```{r} zFPKMPlot(se) ``` To determine which genes are *active*, we compute the median expression within each group. ```{r} activeGenes <- assay(se, "zfpkm") %>% mutate(gene=rownames(assay(se, "zfpkm"))) %>% gather(sample_id, zfpkm, -gene) %>% left_join(select(as.data.frame(colData(se)), sample_id, condition), by="sample_id") %>% group_by(gene, condition) %>% summarize(median_zfpkm=median(zfpkm)) %>% ungroup() %>% mutate(active=(median_zfpkm > -3)) %>% filter(active) %>% select(gene) %>% distinct() seActive <- SummarizedExperiment( assays=SimpleList(counts=as.matrix(assay(se, "counts")[activeGenes$gene, ])), colData=colData(se)) ``` In the following DE analysis, we only use genes that were active in either group. ```{r} library(limma) library(edgeR) # Generate normalized log2CPM from counts AFTER we filter for protein-coding # genes that are detectably expressed. dge <- DGEList(counts=assay(seActive, "counts")) dge <- calcNormFactors(dge) design <- model.matrix(~ 0 + condition, data=colData(seActive)) vq <- voomWithQualityWeights(dge, design, plot=TRUE) fit <- lmFit(vq, design) contrastMatrix <- makeContrasts(conditioncontrol - conditionmutant, levels=design) fit <- contrasts.fit(fit, contrastMatrix) fit <- eBayes(fit, robust=TRUE) deGenes <- topTable(fit, number=Inf) ``` ### References Hart T, Komori HK, LaMere S, Podshivalova K, Salomon DR. Finding the active genes in deep RNA-seq gene expression studies. BMC Genomics. 2013 Nov 11;14:778. doi: 10.1186/1471-2164-14-778. ```{r} sessionInfo() ```