TCGAanalyze
: Analyze data from TCGA.You can easily analyze data using following functions:
TCGAanalyze_Preprocessing
: Preprocessing of Gene Expression data (IlluminaHiSeq_RNASeqV2)You can easily search TCGA samples, download and prepare a matrix of gene expression.
# You can define a list of samples to query and download providing relative TCGA barcodes.
listSamples <- c("TCGA-E9-A1NG-11A-52R-A14M-07","TCGA-BH-A1FC-11A-32R-A13Q-07",
"TCGA-A7-A13G-11A-51R-A13Q-07","TCGA-BH-A0DK-11A-13R-A089-07",
"TCGA-E9-A1RH-11A-34R-A169-07","TCGA-BH-A0AU-01A-11R-A12P-07",
"TCGA-C8-A1HJ-01A-11R-A13Q-07","TCGA-A7-A13D-01A-13R-A12P-07",
"TCGA-A2-A0CV-01A-31R-A115-07","TCGA-AQ-A0Y5-01A-11R-A14M-07")
# Query platform Illumina HiSeq with a list of barcode
query <- GDCquery(project = "TCGA-BRCA",
data.category = "Gene expression",
data.type = "Gene expression quantification",
experimental.strategy = "RNA-Seq",
platform = "Illumina HiSeq",
file.type = "results",
barcode = listSamples,
legacy = TRUE)
# Download a list of barcodes with platform IlluminaHiSeq_RNASeqV2
GDCdownload(query)
# Prepare expression matrix with geneID in the rows and samples (barcode) in the columns
# rsem.genes.results as values
BRCARnaseqSE <- GDCprepare(query)
BRCAMatrix <- assay(BRCARnaseqSE,"raw_count") # or BRCAMatrix <- assay(BRCARnaseqSE,"raw_count")
# For gene expression if you need to see a boxplot correlation and AAIC plot to define outliers you can run
BRCARnaseq_CorOutliers <- TCGAanalyze_Preprocessing(BRCARnaseqSE)
The result is shown below:
TCGA-A2-A0CV-01A-31R-A115-07 | TCGA-A7-A13G-11A-51R-A13Q-07 | |
---|---|---|
C20orf46|55321 | 45 | 43 |
DEGS1|8560 | 16697 | 4104 |
INTS3|65123 | 7885 | 6939 |
KLHL28|54813 | 885 | 1131 |
LOC399815|399815 | 7 | 4 |
SLC27A1|376497 | 3094 | 3658 |
ACRBP|84519 | 205 | 47 |
LOC286135|286135 | 0 | 0 |
SUV420H1|51111 | 3832 | 3062 |
FABP6|2172 | 1 | 0 |
The result from TCGAanalyze_Preprocessing is shown below:
TCGAanalyze_DEA
&
TCGAanalyze_LevelTab
: Differential expression analysis (DEA)Perform DEA (Differential expression analysis) to identify differentially expressed genes (DEGs) using the TCGAanalyze_DEA
function.
TCGAanalyze_DEA
performs DEA using following functions from R :
This function receives as arguments:
Next, we filter the output of dataDEGs by abs(LogFC) >=1, and uses the TCGAanalyze_LevelTab
function to create a table with DEGs (differentially expressed genes), log Fold Change (FC), false discovery rate (FDR), the gene expression level for samples in Cond1type, and Cond2type, and Delta value (the difference of gene expression between the two conditions multiplied logFC).
# Downstream analysis using gene expression data
# TCGA samples from IlluminaHiSeq_RNASeqV2 with type rsem.genes.results
# save(dataBRCA, geneInfo , file = "dataGeneExpression.rda")
library(TCGAbiolinks)
# normalization of genes
dataNorm <- TCGAanalyze_Normalization(tabDF = dataBRCA, geneInfo = geneInfo)
# quantile filter of genes
dataFilt <- TCGAanalyze_Filtering(tabDF = dataNorm,
method = "quantile",
qnt.cut = 0.25)
# selection of normal samples "NT"
samplesNT <- TCGAquery_SampleTypes(barcode = colnames(dataFilt),
typesample = c("NT"))
# selection of tumor samples "TP"
samplesTP <- TCGAquery_SampleTypes(barcode = colnames(dataFilt),
typesample = c("TP"))
# Diff.expr.analysis (DEA)
dataDEGs <- TCGAanalyze_DEA(mat1 = dataFilt[,samplesNT],
mat2 = dataFilt[,samplesTP],
Cond1type = "Normal",
Cond2type = "Tumor",
fdr.cut = 0.01 ,
logFC.cut = 1,
method = "glmLRT")
# DEGs table with expression values in normal and tumor samples
dataDEGsFiltLevel <- TCGAanalyze_LevelTab(dataDEGs,"Tumor","Normal",
dataFilt[,samplesTP],dataFilt[,samplesNT])
The result is shown below:
mRNA | logFC | FDR | Tumor | Normal | Delta |
---|---|---|---|---|---|
FN1 | 2.88 | 1.296151e-19 | 347787.48 | 41234.12 | 1001017.3 |
COL1A1 | 1.77 | 1.680844e-08 | 358010.32 | 89293.72 | 633086.3 |
C4orf7 | 5.20 | 2.826474e-50 | 87821.36 | 2132.76 | 456425.4 |
COL1A2 | 1.40 | 9.480478e-06 | 273385.44 | 91241.32 | 383242.9 |
GAPDH | 1.32 | 3.290678e-05 | 179057.44 | 63663.00 | 236255.5 |
CLEC3A | 6.79 | 7.971002e-74 | 27257.16 | 259.60 | 185158.6 |
IGFBP5 | 1.24 | 1.060717e-04 | 128186.88 | 53323.12 | 158674.6 |
CPB1 | 4.27 | 3.044021e-37 | 37001.76 | 2637.72 | 157968.8 |
CARTPT | 6.72 | 1.023371e-72 | 21700.96 | 215.16 | 145872.8 |
DCD | 7.26 | 1.047988e-80 | 19941.20 | 84.80 | 144806.3 |
Other examples are in the next sections.
CancerProject <- "TCGA-BRCA"
DataDirectory <- paste0("../GDC/",gsub("-","_",CancerProject))
FileNameData <- paste0(DataDirectory, "_","HTSeq_Counts",".rda")
query <- GDCquery(project = CancerProject,
data.category = "Transcriptome Profiling",
data.type = "Gene Expression Quantification",
workflow.type = "HTSeq - Counts")
samplesDown <- getResults(query,cols=c("cases"))
dataSmTP <- TCGAquery_SampleTypes(barcode = samplesDown,
typesample = "TP")
dataSmNT <- TCGAquery_SampleTypes(barcode = samplesDown,
typesample = "NT")
dataSmTP_short <- dataSmTP[1:10]
dataSmNT_short <- dataSmNT[1:10]
queryDown <- GDCquery(project = CancerProject,
data.category = "Transcriptome Profiling",
data.type = "Gene Expression Quantification",
workflow.type = "HTSeq - Counts",
barcode = c(dataSmTP_short, dataSmNT_short))
GDCdownload(query = queryDown,
directory = DataDirectory)
dataPrep <- GDCprepare(query = queryDown,
save = TRUE,
directory = DataDirectory,
save.filename = FileNameData)
dataPrep <- TCGAanalyze_Preprocessing(object = dataPrep,
cor.cut = 0.6,
datatype = "HTSeq - Counts")
dataNorm <- TCGAanalyze_Normalization(tabDF = dataPrep,
geneInfo = geneInfoHT,
method = "gcContent")
boxplot(dataPrep, outline = FALSE)
boxplot(dataNorm, outline = FALSE)
dataFilt <- TCGAanalyze_Filtering(tabDF = dataNorm,
method = "quantile",
qnt.cut = 0.25)
dataDEGs <- TCGAanalyze_DEA(mat1 = dataFilt[,dataSmTP_short],
mat2 = dataFilt[,dataSmNT_short],
Cond1type = "Normal",
Cond2type = "Tumor",
fdr.cut = 0.01 ,
logFC.cut = 1,
method = "glmLRT")
require(TCGAbiolinks)
CancerProject <- "TCGA-BRCA"
DataDirectory <- paste0("../GDC/",gsub("-","_",CancerProject))
FileNameData <- paste0(DataDirectory, "_","miRNA_gene_quantification",".rda")
query.miR <- GDCquery(project = CancerProject,
data.category = "Gene expression",
data.type = "miRNA gene quantification",
file.type = "hg19.mirna",
legacy = TRUE)
samplesDown.miR <- getResults(query.miR,cols=c("cases"))
dataSmTP.miR <- TCGAquery_SampleTypes(barcode = samplesDown.miR,
typesample = "TP")
dataSmNT.miR <- TCGAquery_SampleTypes(barcode = samplesDown.miR,
typesample = "NT")
dataSmTP_short.miR <- dataSmTP.miR[1:10]
dataSmNT_short.miR <- dataSmNT.miR[1:10]
queryDown.miR <- GDCquery(project = CancerProject,
data.category = "Gene expression",
data.type = "miRNA gene quantification",
file.type = "hg19.mirna",
legacy = TRUE,
barcode = c(dataSmTP_short.miR, dataSmNT_short.miR))
GDCdownload(query = queryDown.miR,
directory = DataDirectory)
dataAssy.miR <- GDCprepare(query = queryDown.miR,
save = TRUE,
save.filename = FileNameData,
summarizedExperiment = TRUE,
directory =DataDirectory )
rownames(dataAssy.miR) <- dataAssy.miR$miRNA_ID
# using read_count's data
read_countData <- colnames(dataAssy.miR)[grep("count", colnames(dataAssy.miR))]
dataAssy.miR <- dataAssy.miR[,read_countData]
colnames(dataAssy.miR) <- gsub("read_count_","", colnames(dataAssy.miR))
dataFilt <- TCGAanalyze_Filtering(tabDF = dataAssy.miR,
method = "quantile",
qnt.cut = 0.25)
dataDEGs <- TCGAanalyze_DEA(mat1 = dataFilt[,dataSmNT_short.miR],
mat2 = dataFilt[,dataSmTP_short.miR],
Cond1type = "Normal",
Cond2type = "Tumor",
fdr.cut = 0.01 ,
logFC.cut = 1,
method = "glmLRT")
TCGAanalyze_EAcomplete & TCGAvisualize_EAbarplot
: Enrichment AnalysisResearchers, in order to better understand the underlying biological processes, often want to retrieve a functional profile of a set of genes that might have an important role. This can be done by performing an enrichment analysis.
We will perform an enrichment analysis on gene sets using the TCGAanalyze_EAcomplete
function. Given a set of genes that are up-regulated under certain conditions, an enrichment analysis will identify classes of genes or proteins that are over-represented using annotations for that gene set.
To view the results you can use the TCGAvisualize_EAbarplot
function as shown below.
library(TCGAbiolinks)
# Enrichment Analysis EA
# Gene Ontology (GO) and Pathway enrichment by DEGs list
Genelist <- rownames(dataDEGsFiltLevel)
system.time(ansEA <- TCGAanalyze_EAcomplete(TFname="DEA genes Normal Vs Tumor",Genelist))
# Enrichment Analysis EA (TCGAVisualize)
# Gene Ontology (GO) and Pathway enrichment barPlot
TCGAvisualize_EAbarplot(tf = rownames(ansEA$ResBP),
GOBPTab = ansEA$ResBP,
GOCCTab = ansEA$ResCC,
GOMFTab = ansEA$ResMF,
PathTab = ansEA$ResPat,
nRGTab = Genelist,
nBar = 10)
The result is shown below:
The figure shows canonical pathways significantly overrepresented (enriched) by the DEGs (differentially expressed genes). The most statistically significant canonical pathways identified in DEGs list are listed according to their p value corrected FDR (-Log) (colored bars) and the ratio of list genes found in each pathway over the total number of genes in that pathway (Ratio, red line).
TCGAanalyze_survival
: Survival AnalysisWhen analyzing survival, different problems come up than the ones discussed so far. One question is how do we deal with subjects dropping out of a study. For example, assuming that we test a new cancer drug. While some subjects die, others may believe that the new drug is not effective, and decide to drop out of the study before the study is finished. A similar problem would be faced when we investigate how long a machine lasts before it breaks down.
Using the clinical data, it is possible to create a survival plot with the function TCGAanalyze_survival
as follows:
clin.gbm <- GDCquery_clinic("TCGA-GBM", "clinical")
TCGAanalyze_survival(clin.gbm,
"gender",
main = "TCGA Set\n GBM",height = 10, width=10)
The arguments of TCGAanalyze_survival
are:
The result is shown below:
TCGAanalyze_SurvivalKM
: Correlating gene expression and Survival Analysislibrary(TCGAbiolinks)
# Survival Analysis SA
clinical_patient_Cancer <- GDCquery_clinic("TCGA-BRCA","clinical")
dataBRCAcomplete <- log2(BRCA_rnaseqv2)
tokenStop<- 1
tabSurvKMcomplete <- NULL
for( i in 1: round(nrow(dataBRCAcomplete)/100)){
message( paste( i, "of ", round(nrow(dataBRCAcomplete)/100)))
tokenStart <- tokenStop
tokenStop <-100*i
tabSurvKM<-TCGAanalyze_SurvivalKM(clinical_patient_Cancer,
dataBRCAcomplete,
Genelist = rownames(dataBRCAcomplete)[tokenStart:tokenStop],
Survresult = F,
ThreshTop=0.67,
ThreshDown=0.33)
tabSurvKMcomplete <- rbind(tabSurvKMcomplete,tabSurvKM)
}
tabSurvKMcomplete <- tabSurvKMcomplete[tabSurvKMcomplete$pvalue < 0.01,]
tabSurvKMcomplete <- tabSurvKMcomplete[order(tabSurvKMcomplete$pvalue, decreasing=F),]
tabSurvKMcompleteDEGs <- tabSurvKMcomplete[
rownames(tabSurvKMcomplete) %in% dataDEGsFiltLevel$mRNA,
]
The result is shown below:
pvalue | Cancer Deaths | Cancer Deaths with Top | Cancer Deaths with Down | |
---|---|---|---|---|
DCTPP1 | 6.204170e-08 | 66 | 46 | 20 |
APOO | 9.390193e-06 | 65 | 49 | 16 |
LOC387646 | 1.039097e-05 | 69 | 48 | 21 |
PGK1 | 1.198577e-05 | 71 | 49 | 22 |
CCNE2 | 2.100348e-05 | 65 | 48 | 17 |
Mean Tumor Top | Mean Tumor Down | Mean Normal | |
---|---|---|---|
DCTPP1 | 13.31 | 12.01 | 11.74 |
APOO | 11.40 | 10.17 | 10.01 |
LOC387646 | 7.92 | 4.64 | 5.90 |
PGK1 | 15.66 | 14.18 | 14.28 |
CCNE2 | 11.07 | 8.23 | 7.03 |
TCGAanalyze_DMR
: Differentially methylated regions AnalysisWe will search for differentially methylated CpG sites using the TCGAanalyze_DMR
function. In order to find these regions we use the beta-values (methylation values ranging from 0.0 to 1.0) to compare two groups.
First, it calculates the difference between the mean DNA methylation of each group for each probe.
Second, it test for differential expression between two groups using the wilcoxon test adjusting by the Benjamini-Hochberg method. The default arguments was set to require a minimum absolute beta-values difference of 0.2 and an adjusted p-value of < 0.01.
After these tests, we save a volcano plot (x-axis:diff mean methylation, y-axis: statistical significance) that will help the user identify the differentially methylated CpG sites, then the results are saved in a csv file (DMR_results.groupCol.group1.group2.csv) and finally the object is returned with the calculus in the rowRanges.
The arguments of TCGAanalyze_DMR are:
Argument | Description |
---|---|
data | SummarizedExperiment obtained from the TCGAPrepare |
groupCol | Columns with the groups inside the SummarizedExperiment object. (This will be obtained by the function colData(data)) |
group1 | In case our object has more than 2 groups, you should set the name of the group |
group2 | In case our object has more than 2 groups, you should set the name of the group |
calculate.pvalues.probes | In order to get the probes faster the user can select to calculate the pvalues only for the probes with a difference in DNA methylation. The default is to calculate to all probes. Possible values: “all”, “differential”. Default “all” |
plot.filename | Filename. Default: volcano.pdf, volcano.svg, volcano.png. If set to FALSE, there will be no plot. |
ylab | y axis text |
xlab | x axis text |
title | main title. If not specified it will be “Volcano plot (group1 vs group2) |
legend | Legend title |
color | vector of colors to be used in graph |
label | vector of labels to be used in the figure. Example: c(“Not Significant”,“Hypermethylated in group1”, “Hypomethylated in group1”)) |
xlim | x limits to cut image |
ylim | y limits to cut image |
p.cut | p values threshold. Default: 0.01 |
probe.names | is probe.names |
diffmean.cut | diffmean threshold. Default: 0.2 |
paired | Wilcoxon paired parameter. Default: FALSE |
adj.method | Adjusted method for the p-value calculation |
overwrite | Overwrite the pvalues and diffmean values if already in the object for both groups? Default: FALSE |
cores | Number of cores to be used in the non-parametric test Default = groupCol.group1.group2.rda |
save | Save object with results? Default: TRUE |
data <- TCGAanalyze_DMR(data, groupCol = "methylation_subtype",
group1 = "CIMP.H",
group2="CIMP.L",
p.cut = 10^-5,
diffmean.cut = 0.25,
legend = "State",
plot.filename = "coad_CIMPHvsCIMPL_metvolcano.png")
The output will be a plot such as the figure below. The green dots are the probes that are hypomethylated in group 2 compared to group 1, while the red dots are the hypermethylated probes in group 2 compared to group 1
Also, the TCGAanalyze_DMR
function will save the plot as pdf and return the same SummarizedExperiment that was given as input with the values of p-value, p-value adjusted, diffmean and the group it belongs in the graph (non significant, hypomethylated, hypermethylated) in the rowRanges. The collumns will be (where group1 and group2 are the names of the groups):
This values can be view/acessed using the rowRanges
acessesor (rowRanges(data)
).
Observation: Calling the same function again, with the same arguments will only plot the results, as it was already calculated. If you want to have them recalculated, please set overwrite
to TRUE
or remove the calculated collumns.
TCGAvisualize
: Visualize results from analysis functions with TCGA’s data.You can easily visualize results from some following functions:
TCGAvisualize_Heatmap
: Create heatmaps with cluster barsIn order to have a better view of clusters, we normally use heatmaps. TCGAvisualize_Heatmap
will plot a heatmap and add to each sample bars representing different features. This function is a wrapper to the package ComplexHeatmap package,
The arguments of this function are:
For more information please take a look on case study #2.
The result is shown below:
TCGAvisualize_Volcano
: Create volcano plotCreates a volcano plot for DNA methylation or expression
The arguments of this function are:
Argument | Description |
---|---|
x | x-axis data |
y | y-axis data |
filename | Filename. Default: volcano.pdf, volcano.svg, volcano.png |
ylab | y axis text |
xlab | x axis text |
title | main title. If not specified it will be “Volcano plot (group1 vs group2) |
legend | Legend title |
label | vector of labels to be used in the figure. Example: c(“Not Significant”,“Hypermethylated in group1”, “Hypomethylated in group1”))#’ |
xlim | x limits to cut image |
ylim | y limits to cut image |
color | vector of colors to be used in graph |
names | Names to be ploted if significant. Should be the same size of x and y |
names.fill | Names should be filled in a color box? Default: TRUE |
show.names | What names will be showd? Possibilities: “both”, “significant”, “highlighted” |
x.cut | x-axis threshold. Default: 0.0 If you give only one number (e.g. 0.2) the cut-offs will be -0.2 and 0.2. Or you can give diffenrent cutt-ofs as a vector (e.g. c(-0.3,0.4)) |
y.cut | p-values threshold. |
height | Figure height |
width | Figure width |
highlight | List of genes/probes to be highlighted. It should be in the names argument. |
highlight.color | Color of the points highlighted |
names.size | Size of the names text |
dpi | Figure dpi |
For more information please take a look on case study #3.
TCGAvisualize_PCA
: Principal Component Analysis plot for differentially expressed genesIn order to better understand our genes, we can perform a PCA to reduce the number of dimensions of our gene set. The function TCGAvisualize_PCA
will plot the PCA for different groups.
The arguments of this function are:
# normalization of genes
dataNorm <- TCGAbiolinks::TCGAanalyze_Normalization(dataBRCA, geneInfo)
# quantile filter of genes
dataFilt <- TCGAanalyze_Filtering(tabDF = dataNorm,
method = "quantile",
qnt.cut = 0.25)
# selection of normal samples "NT"
group1 <- TCGAquery_SampleTypes(colnames(dataFilt), typesample = c("NT"))
# selection of normal samples "TP"
group2 <- TCGAquery_SampleTypes(colnames(dataFilt), typesample = c("TP"))
# Principal Component Analysis plot for ntop selected DEGs
pca <- TCGAvisualize_PCA(dataFilt,dataDEGsFiltLevel, ntopgenes = 200, group1, group2)
The result is shown below:
TCGAvisualize_meanMethylation
: Mean DNA Methylation AnalysisUsing the data and calculating the mean DNA methylation per group, it is possible to create a mean DNA methylation boxplot with the function TCGAvisualize_meanMethylation
as follows:
query <- GDCquery(project = "TCGA-GBM",
data.category = "DNA methylation",
platform = "Illumina Human Methylation 27",
legacy = TRUE,
barcode = c("TCGA-02-0058-01A-01D-0186-05", "TCGA-12-1597-01B-01D-0915-05",
"TCGA-12-0829-01A-01D-0392-05", "TCGA-06-0155-01B-01D-0521-05",
"TCGA-02-0099-01A-01D-0199-05", "TCGA-19-4068-01A-01D-1228-05",
"TCGA-19-1788-01A-01D-0595-05", "TCGA-16-0848-01A-01D-0392-05"))
GDCdownload(query, method = "api")
data <- GDCprepare(query)
# "shortLetterCode" is a column in the SummarizedExperiment::colData(data) matrix
TCGAvisualize_meanMethylation(data, groupCol = "shortLetterCode",filename = NULL)
## groups Mean Median Max Min
## 1 NT 0.5059958 0.5104387 0.5229913 0.4814882
## 2 TP 0.5103744 0.5052156 0.5305809 0.4893334
## 3 TR 0.5014347 0.5018155 0.5267410 0.4689107
## NT TP TR
## NT NA 0.6108825 0.6511016
## TP 0.6108825 NA 0.4023347
## TR 0.6511016 0.4023347 NA
The arguments of TCGAvisualize_meanMethylation
are:
Argument | Description |
---|---|
data | SummarizedExperiment object obtained from GDCPrepare |
groupCol | Columns in colData(data) that defines the groups. If no columns defined a columns called “Patients” will be used |
subgroupCol | Columns in colData(data) that defines the subgroups. |
shapes | Shape vector of the subgroups. It must have the size of the levels of the subgroups. Example: shapes = c(21,23) if for two levels |
print.pvalue | Print p-value for two groups |
plot.jitter | Plot jitter? Default TRUE |
jitter.size | Plot jitter size? Default 3 |
filename | The name of the pdf that will be saved |
ylab | y axis text in the plot |
xlab | x axis text in the plot |
title | main title in the plot |
labels | Labels of the groups |
group.legend | Name of the group legend. DEFAULT: groupCol |
subgroup.legend | Name of the subgroup legend. DEFAULT: subgroupCol |
color | vector of colors to be used in graph |
y.limits | Change lower/upper y-axis limit |
sort | Sort boxplot by mean or median. Possible values: mean.asc, mean.desc, median.asc, meadian.desc |
order | Order of the boxplots |
legend.position | Legend position (“top”, “right”,“left”,“bottom”) |
legend.title.position | Legend title position (“top”, “right”,“left”,“bottom”) |
legend.ncols | Number of columns of the legend |
add.axis.x.text | Add text to x-axis? Default: FALSE |
width | Plot width default:10 |
height | Plot height default:10 |
dpi | Pdf dpi default:600 |
axis.text.x.angle | Angle of text in the x axis |
Other examples using the parameters:
# setting y limits: lower 0, upper 1
TCGAvisualize_meanMethylation(data,groupCol = "shortLetterCode",
filename = NULL, y.limits = c(0,1))
## groups Mean Median Max Min
## 1 NT 0.5059958 0.5104387 0.5229913 0.4814882
## 2 TP 0.5103744 0.5052156 0.5305809 0.4893334
## 3 TR 0.5014347 0.5018155 0.5267410 0.4689107
## NT TP TR
## NT NA 0.6108825 0.6511016
## TP 0.6108825 NA 0.4023347
## TR 0.6511016 0.4023347 NA
# setting y limits: lower 0
TCGAvisualize_meanMethylation(data,groupCol = "shortLetterCode",
filename = NULL, y.limits = 0)
## groups Mean Median Max Min
## 1 NT 0.5059958 0.5104387 0.5229913 0.4814882
## 2 TP 0.5103744 0.5052156 0.5305809 0.4893334
## 3 TR 0.5014347 0.5018155 0.5267410 0.4689107
## NT TP TR
## NT NA 0.6108825 0.6511016
## TP 0.6108825 NA 0.4023347
## TR 0.6511016 0.4023347 NA
# Changing shapes of jitters to show subgroups
TCGAvisualize_meanMethylation(data,
groupCol = "subtype_Pan.Glioma.DNA.Methylation.Cluster",
subgroupCol ="vital_status", filename = NULL)
## groups Mean Median Max Min
## 1 group1 0.5082970 0.5141680 0.5267410 0.4814882
## 2 group2 0.5007477 0.5052156 0.5279137 0.4689107
## 3 group3 0.5087602 0.5104387 0.5305809 0.4893334
## group1 group2 group3
## group1 NA 0.4669833 0.9601968
## group2 0.4669833 NA 0.4107977
## group3 0.9601968 0.4107977 NA
# Sorting bars by descending mean methylation
TCGAvisualize_meanMethylation(data,
groupCol = "subtype_Pan.Glioma.DNA.Methylation.Cluster",
sort="mean.desc",
filename=NULL)
## groups Mean Median Max Min
## 1 group1 0.5082970 0.5141680 0.5267410 0.4814882
## 2 group2 0.5007477 0.5052156 0.5279137 0.4689107
## 3 group3 0.5087602 0.5104387 0.5305809 0.4893334
## group1 group2 group3
## group1 NA 0.4669833 0.9601968
## group2 0.4669833 NA 0.4107977
## group3 0.9601968 0.4107977 NA
# Sorting bars by asc mean methylation
TCGAvisualize_meanMethylation(data,
groupCol = "subtype_Pan.Glioma.DNA.Methylation.Cluster",
sort = "mean.asc",
filename=NULL)
## groups Mean Median Max Min
## 1 group1 0.5082970 0.5141680 0.5267410 0.4814882
## 2 group2 0.5007477 0.5052156 0.5279137 0.4689107
## 3 group3 0.5087602 0.5104387 0.5305809 0.4893334
## group1 group2 group3
## group1 NA 0.4669833 0.9601968
## group2 0.4669833 NA 0.4107977
## group3 0.9601968 0.4107977 NA
## groups Mean Median Max Min
## 1 ALIVE 0.5103744 0.5052156 0.5305809 0.4893334
## 2 DEAD 0.5037152 0.5100908 0.5267410 0.4689107
## ALIVE DEAD
## ALIVE NA 0.4147255
## DEAD 0.4147255 NA
TCGAvisualize_starburst
: Integration of gene expression and DNA methylation dataThe starburst plot is proposed to combine information from two volcano plots, and is applied for a study of DNA methylation and gene expression. It first introduced in 2010 (Noushmehr et al. 2010). In order to reproduce this plot, we will use the TCGAvisualize_starburst
function.
The function creates Starburst plot for comparison of DNA methylation and gene expression. The log10 (FDR-corrected P value) for DNA methylation is plotted in the x axis, and for gene expression in the y axis, for each gene. The black dashed line shows the FDR-adjusted P value of 0.01.
The arguments of this function are: | Argument | Description |
|————–|———————————————————————————————————————————————————————————————————————————————————-| | exp | Object obtained by DEArnaSEQ function | | group1 | The name of the group 1 Obs: Column p.value.adj.group1.group2 should exist | | group2 | The name of the group 2. Obs: Column p.value.adj.group1.group2 should exist | | exp.p.cut | expression p value cut-off | | met.p.cut | methylation p value cut-off | | diffmean.cut | If set, the probes with diffmean higher than methylation cut-off will be highlighted in the plot. And the data frame return will be subseted. | | logFC.cut | If set, the probes with expression fold change higher than methylation cut-off will be highlighted in the plot. And the data frame return will be subseted. | | met.platform | DNA methylation platform (“27K”,“450K” or “EPIC”) | | genome | Genome of reference (“hg38” or “hg19”) used to identify nearest probes TSS | | names | Add the names of the significant genes? Default: FALSE | | names.fill | Names should be filled in a color box? Default: TRUE | | filename | The filename of the file (it can be pdf, svg, png, etc) | | return.plot | If true only plot object will be returned (pdf will not be created) | | ylab | y axis text | | xlab | x axis text | | title | main title | | legend | legend title | | color | vector of colors to be used in graph | | label | vector of labels to be used in graph | | xlim | x limits to cut image | | ylim | y limits to cut image | | height | Figure height | | width | Figure width | | dpi | Figure dpi |
starburst <- TCGAvisualize_starburst(coad.SummarizeExperiment,
different.experssion.analysis.data,
group1 = "CIMP.H",
group2 = "CIMP.L",
met.platform = "450K",
genome = "hg19",
met.p.cut = 10^-5,
exp.p.cut = 10^-5,
names = TRUE)
As a result, the function will a plot the figure below and return a matrix with the Gene_symbol and it status in relation to expression (up regulated/down regulated) and to methylation (Hyper/Hypo methylated). The case study 3, shows the complete pipeline for creating this figure.
## R version 3.4.3 (2017-11-30)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.6-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.6-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] grid parallel stats4 stats graphics grDevices utils
## [8] datasets methods base
##
## other attached packages:
## [1] png_0.1-7 DT_0.4
## [3] dplyr_0.7.4 SummarizedExperiment_1.8.1
## [5] DelayedArray_0.4.1 matrixStats_0.53.0
## [7] Biobase_2.38.0 GenomicRanges_1.30.1
## [9] GenomeInfoDb_1.14.0 IRanges_2.12.0
## [11] S4Vectors_0.16.0 BiocGenerics_0.24.0
## [13] TCGAbiolinks_2.6.12
##
## loaded via a namespace (and not attached):
## [1] backports_1.1.2 circlize_0.4.3
## [3] aroma.light_3.8.0 plyr_1.8.4
## [5] selectr_0.3-1 ConsensusClusterPlus_1.42.0
## [7] lazyeval_0.2.1 splines_3.4.3
## [9] BiocParallel_1.12.0 ggplot2_2.2.1
## [11] sva_3.26.0 digest_0.6.15
## [13] foreach_1.4.4 htmltools_0.3.6
## [15] magrittr_1.5 memoise_1.1.0
## [17] cluster_2.0.6 doParallel_1.0.11
## [19] limma_3.34.7 ComplexHeatmap_1.17.1
## [21] Biostrings_2.46.0 readr_1.1.1
## [23] annotate_1.56.1 R.utils_2.6.0
## [25] prettyunits_1.0.2 colorspace_1.3-2
## [27] blob_1.1.0 rvest_0.3.2
## [29] ggrepel_0.7.0 RCurl_1.95-4.10
## [31] jsonlite_1.5 roxygen2_6.0.1
## [33] genefilter_1.60.0 bindr_0.1
## [35] survival_2.41-3 zoo_1.8-1
## [37] iterators_1.0.9 glue_1.2.0
## [39] survminer_0.4.2 gtable_0.2.0
## [41] zlibbioc_1.24.0 XVector_0.18.0
## [43] GetoptLong_0.1.6 shape_1.4.3
## [45] scales_0.5.0 DESeq_1.30.0
## [47] DBI_0.7 edgeR_3.20.8
## [49] ggthemes_3.4.0 Rcpp_0.12.15
## [51] xtable_1.8-2 progress_1.1.2
## [53] cmprsk_2.2-7 foreign_0.8-69
## [55] bit_1.1-12 matlab_1.0.2
## [57] km.ci_0.5-2 htmlwidgets_1.0
## [59] httr_1.3.1 RColorBrewer_1.1-2
## [61] pkgconfig_2.0.1 XML_3.98-1.9
## [63] R.methodsS3_1.7.1 locfit_1.5-9.1
## [65] labeling_0.3 rlang_0.1.6
## [67] reshape2_1.4.3 AnnotationDbi_1.40.0
## [69] munsell_0.4.3 tools_3.4.3
## [71] downloader_0.4 RSQLite_2.0
## [73] devtools_1.13.4 broom_0.4.3
## [75] evaluate_0.10.1 stringr_1.2.0
## [77] yaml_2.1.16 knitr_1.19
## [79] bit64_0.9-7 survMisc_0.5.4
## [81] purrr_0.2.4 bindrcpp_0.2
## [83] EDASeq_2.12.0 nlme_3.1-131
## [85] R.oo_1.21.0 xml2_1.2.0
## [87] biomaRt_2.34.2 compiler_3.4.3
## [89] curl_3.1 testthat_2.0.0
## [91] tibble_1.4.2 geneplotter_1.56.0
## [93] stringi_1.1.6 highr_0.6
## [95] GenomicFeatures_1.30.3 lattice_0.20-35
## [97] Matrix_1.2-12 commonmark_1.4
## [99] psych_1.7.8 KMsurv_0.1-5
## [101] pillar_1.1.0 GlobalOptions_0.0.12
## [103] data.table_1.10.4-3 bitops_1.0-6
## [105] rtracklayer_1.38.3 R6_2.2.2
## [107] latticeExtra_0.6-28 hwriter_1.3.2
## [109] RMySQL_0.10.13 ShortRead_1.36.0
## [111] gridExtra_2.3 codetools_0.2-15
## [113] assertthat_0.2.0 rprojroot_1.3-2
## [115] rjson_0.2.15 withr_2.1.1
## [117] GenomicAlignments_1.14.1 Rsamtools_1.30.0
## [119] mnormt_1.5-5 GenomeInfoDbData_1.0.0
## [121] mgcv_1.8-23 hms_0.4.1
## [123] tidyr_0.8.0 rmarkdown_1.8
## [125] ggpubr_0.1.6
Noushmehr, Houtan, Daniel J Weisenberger, Kristin Diefes, Heidi S Phillips, Kanan Pujara, Benjamin P Berman, Fei Pan, et al. 2010. “Identification of a Cpg Island Methylator Phenotype That Defines a Distinct Subgroup of Glioma.” Cancer Cell 17 (5). Elsevier:510–22.