Title: | Identifies differentially expressed genes with respect to other local genes |
---|---|
Description: | The goal of DELocal is to identify DE genes compared to their neighboring genes from the same chromosomal location. It has been shown that genes of related functions are generally very far from each other in the chromosome. DELocal utilzes this information to identify DE genes comparing with their neighbouring genes. |
Authors: | Rishi Das Roy [aut, cre] |
Maintainer: | Rishi Das Roy <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.7.0 |
Built: | 2024-11-29 05:11:04 UTC |
Source: | https://github.com/bioc/DELocal |
Finds differentially expressed genes by comparing neighboring genes
DELocal( pSmrExpt, nearest_neighbours, pDesign, pValue_cut = 0.05, pLogFold_cut = 0 )
DELocal( pSmrExpt, nearest_neighbours, pDesign, pValue_cut = 0.05, pLogFold_cut = 0 )
pSmrExpt |
SummarizedExperiment object |
nearest_neighbours |
How many nearest neighbours within 1 Mb window to evaluate? |
pDesign |
design formula |
pValue_cut |
cut off value for adjusted p-value |
pLogFold_cut |
cut off value for relative log fold change compared to neighbouring genes |
A data.frame with top significant genes with the following columns:
relative.logFC: relative logFC compared to neighbouring genes
P.Value: raw p-value
adj.P.Value: adjusted p-value
B: log-odds that the gene is differentially expressed
count_matrix <- as.matrix(read.table(file = system.file("extdata", "tooth_RNASeq_counts.txt", package = "DELocal"))) colData <- data.frame(condition=gsub("\\..*",x=colnames(count_matrix), replacement = "")) gene_location <- read.table(file = system.file("extdata", "gene_location.txt", package = "DELocal")) smrExpt <- SummarizedExperiment::SummarizedExperiment( assays=list(counts=count_matrix), rowData = gene_location, colData=colData) contrast= c("condition","ME13","ME14") require(dplyr) x_genes <- SummarizedExperiment::rowData(smrExpt) %>% as.data.frame() %>% filter(chromosome_name=="X") %>% rownames() DELocal_result <- DELocal(pSmrExpt = smrExpt[x_genes,], nearest_neighbours = 5, pDesign = ~ condition, pValue_cut = 0.05, pLogFold_cut = 0)
count_matrix <- as.matrix(read.table(file = system.file("extdata", "tooth_RNASeq_counts.txt", package = "DELocal"))) colData <- data.frame(condition=gsub("\\..*",x=colnames(count_matrix), replacement = "")) gene_location <- read.table(file = system.file("extdata", "gene_location.txt", package = "DELocal")) smrExpt <- SummarizedExperiment::SummarizedExperiment( assays=list(counts=count_matrix), rowData = gene_location, colData=colData) contrast= c("condition","ME13","ME14") require(dplyr) x_genes <- SummarizedExperiment::rowData(smrExpt) %>% as.data.frame() %>% filter(chromosome_name=="X") %>% rownames() DELocal_result <- DELocal(pSmrExpt = smrExpt[x_genes,], nearest_neighbours = 5, pDesign = ~ condition, pValue_cut = 0.05, pLogFold_cut = 0)
Returns median expression from different conditions of genes from a neighbourhood of a gene of interest
plotNeighbourhood( pSmrExpt, pNearest_neighbours = 5, pDesign = ~condition, colorFactor = "condition", pGene_id )
plotNeighbourhood( pSmrExpt, pNearest_neighbours = 5, pDesign = ~condition, colorFactor = "condition", pGene_id )
pSmrExpt |
SummarizedExperiment object |
pNearest_neighbours |
How many nearest neighbours within 1 Mb window to plot |
pDesign |
design formula |
colorFactor |
The coloring factor |
pGene_id |
The gene of interest |
a list which contains both the data from the neighbourhood and a ggplot object
count_matrix <- as.matrix(read.table(file = system.file("extdata", "tooth_RNASeq_counts.txt", package = "DELocal"))) colData <- data.frame(condition=gsub("\\..*",x=colnames(count_matrix), replacement = "")) gene_location <- read.table(file = system.file("extdata", "gene_location.txt", package = "DELocal")) smrExpt <- SummarizedExperiment::SummarizedExperiment(assays=list(counts=count_matrix), rowData = gene_location, colData = colData) contrast= c("condition","ME13","ME14") require(dplyr) x_genes <- SummarizedExperiment::rowData(smrExpt) %>% as.data.frame() %>% filter(chromosome_name=="X") %>% rownames() DELocal::plotNeighbourhood(pSmrExpt = smrExpt, pGene_id = "ENSMUSG00000059401")
count_matrix <- as.matrix(read.table(file = system.file("extdata", "tooth_RNASeq_counts.txt", package = "DELocal"))) colData <- data.frame(condition=gsub("\\..*",x=colnames(count_matrix), replacement = "")) gene_location <- read.table(file = system.file("extdata", "gene_location.txt", package = "DELocal")) smrExpt <- SummarizedExperiment::SummarizedExperiment(assays=list(counts=count_matrix), rowData = gene_location, colData = colData) contrast= c("condition","ME13","ME14") require(dplyr) x_genes <- SummarizedExperiment::rowData(smrExpt) %>% as.data.frame() %>% filter(chromosome_name=="X") %>% rownames() DELocal::plotNeighbourhood(pSmrExpt = smrExpt, pGene_id = "ENSMUSG00000059401")