Integrate scATAC-seq and scRNA-seq data to reconstruct regulatory networks
Compare to reconstruct regulatory network only using scRNA-seq data, integrated analysis of scATAC-seq with scRNA-seq overall substantially improve the reconstruction of gene regulatory networks, and have higher precision of identifying known transcription factors. Here, we provide an example of integrating scRNA-seq data and scATAC-seq data to reconstruct regulatory networks for heart regeneration. IReNA analyzes scRNA-seq to infer potential regulatory relationships, and then uses footprints with high footprint occupancy score to refine potential regulatory relationships. Finally performing network analysis through hypergeometric test.
To run this pipeline, Computer or server of linux system and the following software are needed: samtools and fimo.
1. Test data download
Test data used in the following steps can be downloaded from google drive. The raw data of scATAC-seq used to run scATAC-seq preprocessing pipline are available from https://www.ncbi.nlm.nih.gov/geo/query/acc.cgi?acc=GSE142365
2. Workflow
3. Inputs of this pipeline
(1) bam, peak-to-gene links, and footprint files of scATAC-seq data
bam, Peak-to-gene links, and footprint files of scATAC-seq data can be separately generated by the step0, step 2, and step 3 in scATAC-seq analysis pipline.
(2) Seurat object with Pseudotime in the metadata
This seurat object can be generated in scRNA-seq preprocessing pipeline
Hint: seurat object used here should only contain interest genes(e.g. differentially expressed genes and expressed transcription factors)
If bulk RNA-seq data are used to identify basic regulatory relationships, just load raw counts of bulk RNA-seq data, then use the same codes as the scRNA-seq data for further analysis. We suggest to load the expression profile of only differentially expressed genes, otherwise much more time will be used to calculate the correlation of each gene pair. Deseq2 and edgeR can be used to identify differentially expressed genes in bulk RNA-seq data.
(3) Reference genome of the species
Reference genome should be the same as that used for mapping ATAC-seq data, which can be downloaded from UCSC database.
(4) Motif information of transcription factors
IReNA contains DNA motif datasets for four species (Homo sapiens, Mus musculus, Zebrafish and Chicken) derived from TRANSFAC version 201803. Following codes are used to get the motif dataset from TRANSFAC or user-defined motif dataset which should have the same format as these from TRANSFAC database.
library(IReNA)
###call Mus musculus motif database
<- Tranfac201803_Mm_MotifTFsF
motif1 ###call Homo sapiens motif database
<- Tranfac201803_Hs_MotifTFsF
motif1 ###call Zebrafish motif database
<- Tranfac201803_Zf_MotifTFsF
motif1 ###call Chicken motif database
<- Tranfac201803_Ch_MotifTFsF motif1
(5) Motif position weight matirx (PWM)
PWM of all species can be download from https://github.com/jiang-junyao/MEMEmotif
4. Tutorial
IReNA contains four main parts to reconstruct regulatory networks:
Part 1: Analyze scRNA-seq or bulk RNA-seq data to get basic regulatory relationships
Part 2: Analyze ATAC-seq data to refine regulatory relationships (with ATAC-seq data)
Part 3: Regulatory network analysis and visualization
Part 1: Analyze scRNA-seq or bulk RNA-seq data to get basic regulatory relationships
IReNA supports three formats of the inputs:
Input data of this part: Seurat object with Pseudotime in the metadata, which can be generated in scRNA-seq preprocessing pipeline
The test data for this part (‘seurat_with_time.rds’) can be downloaded from https://github.com/jiang-junyao/IReNA-test-data/blob/master/seurat_with_time.rds. The test data are from the study of human spermatogonial stem cell development which is available in https://www.cell.com/cell-stem-cell/fulltext/S1934-5909(17)30370-3.
First, IReNA divide cells into 60 bins across pseudotime (20 bins for each branch). The bin is removed if all genes in this bin have no expression. The gene is filtered if absolute fold change < 0.01 (set by the parameter ‘FC’). Then, genes will be clustered through K-means algorithm (the number of clusters ‘K’ is set by the parameter ‘K1’).
If bulk RNA-seq data are used to identify regulatory relationships, load your expression matrix as expression_profile that generated by get_SmoothByBin_PseudotimeExp(), then run the following codes. I suggest to input expression profile that only contains differentially expressed genes, or you will a huge of time to calculate correlation of each gene pair.
###Load the test data
seurat_with_time <- readRDS('CM_seurat_with_time.rds')
###Get expression profiles ordered by pseudotime
expression_profile <- get_SmoothByBin_PseudotimeExp(seurat_with_time, Bin = 60, method = 'State')
###Filter noise and logFC in expression profile
expression_profile_filter <- filter_expression_profile(expression_profile, FC=0.01)
###K-means clustering
clustering <- clustering_Kmeans(expression_profile_filter, K1=4)
1:5,1:5] clustering[
## KmeansGroup FoldChangeQ95 SmExp1 SmExp2 SmExp3
## Tfcp2 1 1.42074006 0.86973961 1.04984927 1.00683762
## Cux1 1 -0.43935731 1.16796885 0.64085837 0.93599581
## Arnt 1 -0.07020227 0.64529043 1.55948037 0.13001053
## Ppara 1 0.80449683 -0.22638135 0.64400251 1.09123493
## Ddit3 1 1.13501050 0.05743463 -0.04918435 -0.04926762
Visualize the clustering of gene expression profiles through the heatmap
plot_kmeans_pheatmap(clustering,ModuleColor1 = c('#67C7C1','#5BA6DA','#FFBF0F','#C067A9'))
Because of the characteristics of Kmeans algorithm, different results will be obtained each time clustering.
Then, Add Ensmble ID of the genes in the first column.
###Add Ensembl ID as the first column of clustering results
<- add_ENSID(clustering, Spec1='Mm')
Kmeans_clustering_ENS 1:5,1:5] Kmeans_clustering_ENS[
## Symbol KmeansGroup FoldChangeQ95 SmExp1 SmExp2
## ENSMUSG00000009733 Tfcp2 1 1.42074006 0.86973961 1.04984927
## ENSMUSG00000029705 Cux1 1 -0.43935731 1.16796885 0.64085837
## ENSMUSG00000015522 Arnt 1 -0.07020227 0.64529043 1.55948037
## ENSMUSG00000022383 Ppara 1 0.80449683 -0.22638135 0.64400251
## ENSMUSG00000025408 Ddit3 1 1.13501050 0.05743463 -0.04918435
Then IReNA infers the potential regulation for each gene in the expression profile through GENIE3, one of the top performing method for network inference as listed in Beeline benchmark. Futhermore, based on smoothed expression profile, regulation type of each gene pair is defined by Person’s correlation, correlation > 0 reprsents activation, correlation < 0 reprsents repression.
<- GENIE3(as.matrix(seurat_with_time@assays$RNA@data),nCores = 50)
weightMat <- getLinkList(weightMat)
weightMat <- weightMat[weightMat[,3]>0.0004,]
regulation ### add regulation type for each gene pair
<- add_regulation_type(Kmeans_clustering_ENS,regulation)
regulatory_relationships ### check whether source genes are transcription factors
<- c()
motifTF for (i in 1:nrow(motif1)) {
<- strsplit(motif1[i,5],';')[[1]]
TF <- c(motifTF,TF)
motifTF
}<- regulatory_relationships[regulatory_relationships[,1] %in% motifTF,] regulatory_relationships
However, GENIE3 costs a huge of running time and memory compare to Person’s correlation. So, we also provide functions to infer potential regulation through Person’s correlation.
<- get_cor(Kmeans_clustering_ENS, motif = motif1, correlation_filter = 0.6, start_column = 4) regulatory_relationships
Part 2: Analyze scATAC-seq data to refine regulatory relationships (with scATAC-seq data)
For scATAC-seq data is available, IReNA integrate cell-type-specific footprints and to refine regulatory relationships calculates the footprints of transcription factors and footprint occupancy score (FOS) to refine regulatory relationships. First of all, footprints whose distance is less than 4 are merged and then the sequence of each footprint is obtained from the reference genome through the function ‘get_merged_fasta’. The reference genome should be in fasta/fa format which can be downloaded from UCSC or other genome database.
The footprint file and peak-to-gene links file used below is available in https://github.com/jiang-junyao/IReNA-test-data/blob/master/scATAC-seq.zip
The human genome used below can be downloaded from https://hgdownload.soe.ucsc.edu/goldenPath/hg38/bigZips/hg38.fa.gz
###merge footprints whose distance is less than 4
<- read.table('footprints.bed',sep = '\t')
filtered_footprints <- 'Genome/hg38.fa'
fastadir <- get_merged_fasta(filtered_footprints,fastadir)
merged_fasta write.table(merged_fasta,'merged_footprints.fasta',row.names=F,quote=F,col.names=F)
After obtaining the motif sequences, use fimo software to identify binding motifs in the footprints. Because fimo software only supports linux environment, we generate a shell script to run Fimo software.
First, identify differentially expressed genes related motifs through the function ‘motif_select’, which will reduce the running time of the subsequent analysis process.
Then, you should set the following five parameters:
fimodir: path of Fimo software. If Fimo has been set to the global environment variable, just set ‘fimodir <- fimo’.
outputdir1: output path of the shell scripts.
outputdir: output path of Fimo result.
motifdir: path of the motif file, which can be downloaded from https://github.com/jiang-junyao/MEMEmotif or TRANSFAC database.
sequencedir: path of the sequence which is generated by the function ‘get_merged_fasta’.
Please note that, at the end of ‘outputdir’ and ‘sequencedir’ the symbol ’/’should be contained.
### Identify differentially expressed genes related motifs
<- motifs_select(Tranfac201803_Hs_MotifTFsF, rownames(Kmeans_clustering_ENS)) ###Kmeans_clustering_ENS was obtained in part1
motif1 ### run find_motifs()
<- 'fimo'
fimodir <- '/public/home/user/fimo/output/'
outputdir <- '/public/home/user/fimo/output/'
outputdir <- '/public/home/user/fimo/Mememotif/'
motifdir <- '/public/home/user/fimo/merged_footprints.fasta'
sequencedir find_motifs(motif1,step=20,fimodir, outputdir1, outputdir, motifdir, sequencedir)
### run fimo_all script in shell
<- paste0('sh ',outputdir1,'Fimo_All.sh')
shell_code system(shell_code,wait=TRUE)
Then, we combine these Fimo consequence in ‘outputdir’. Notably outputdir folder should only contain fimo result files. Next, Generate TFs-to-genes regulation relationships through integrating peak-to-gene links, combined motifs and seurat object that only contains expressed TFs and DEGs.
<- combine_footprints(outputdir)
combined_motifs <- readRDS('CM_p2g.rds')
p2g <- generate_scATAC_Candid(p2g,combined_motifs,rownames(refined_seurat))
Candid ### output filtered footprints
write.table(Candid[[1]],'filtered_footprints.bed', quote = F, row.names = F, col.names = F, sep = '\t')
Next, because of the size of original BAM file is too large, so we need to use samtools to extract footprints realated regions in BAM to reduce running time of function which analyze bam files in IReNA. (If you use our test data from google drive, just skip this step)
shell_code1 <- 'samtools view -hb -L filtered_footprint.bed possorted_bam.bam > CM_filter.bam'
system(shell_code1,,wait=TRUE)
Optional Locate the center on the Tn5 cut site thorugh R package ATACseqQC. ATACseqQC shifts reads on the forward strand by +4 bp and reads on the reverse strand by -5 bp from the bam-files.
library(ATACseqQC)
library(Rsamtools)
<- 'CM_filter.bam'
bamfilepath1 indexBam(bamfilepath1)
<- readBamFile(bamfilepath1, tag=tags, asMates=TRUE, bigFile=TRUE)
gal <- shiftGAlignmentsList(gal, 'CM_filter_shift.bam') gal1
Furthermore, we count the cuts of each position in footprints by wig_track(), and use these cuts to calculate the FOS of footprints to identify enriched TFs which determine the regulatory relationship. regulatory_relationships used here was calculated in part1. Parameter FOS_threshold used here is the threshold to filter low quality footprints, you can increase it to reduce the number of exported footprints. Hint workers parameter of ‘cal_footprint_cuts’ function represent the number of threads.
### calculate cuts of each each position in footprints
<- 'CM_filter.bam'
bamfilepath1 ### set parameter 'workers' to make this function run in parallel
<- cal_footprint_cuts(bamfilepath = bamfilepath1,bedfile = list2[[1]],workers = 40)
cuts1 ### If you have more than one sample, just put all of them into cut_list
<- list(cuts1)
cut_list ### get related genes of footprints with high FOS
<- Footprints_FOS(cut_list,Candid[[2]], FOS_threshold = 0.5)
scATAC_FOS05 ### Use information of footprints with high FOS to refine regulatory relationships
<- filter_ATAC(scATAC_FOS05,regulatory_relationships) filtered_regulatory
Part 3: Regulatory network analysis and visualization
After we get ‘filtered_regulatory_relationships’ and ‘Kmeans_clustering_ENS’, we can reconstruct regulatory network. Run network_analysis() to get regulatory, this function will generate a list which contain the following 9 elements:
(1)Cor_TFs.txt: list of expressed TFs in the gene networks.
(2)Cor_EnTFs.txt: list of TFs which significantly regulate gene modules (or enriched TFs).
(3)FOSF_RegMTF_Cor_EnTFs.txt: regulatory pairs in which the source gene is enriched TF.
(4)FOSF_RegMTF_Cor_EnTFs.txt: regulatory pairs in which both source gene and target gene are enriched TFs.
(5)FOSF_RegMTF_Cor_EnTFs.txt: regulatory pairs only including regulations within each module but not those between modules, in this step
(6)TF_list: enriched TFs which significantly regulate gene modules
(7)TF_module_regulation: details of enriched TFs which significantly regulate gene modules
(8)TF_network: regulatory network for enriched transcription factors of each module
(9)intramodular_network: intramodular regulatory network
Here, we use refined regulatory relationships from part2 to reconstruct regulatory networks
<- network_analysis(filtered_regulatory_relationships,Kmeans_clustering_ENS,TFFDR1 = 10,TFFDR2 = 10,ModuleFDR = 0.1) TFs_list
We can also make enrichment analysis for differentially expressed genes in each module. Before you run this function, you need to download the org.db of your species through BiocManager.
### Download Mus musculus org.db
#BiocManger::install('org.Mm.eg.db')
library(org.Mm.eg.db)
### Enrichment analysis (KEGG)
<- enrich_module(Kmeans_clustering_ENS, org.Mm.eg.db, enrich.db = 'KEGG',organism = 'mmu',fun_num = 10)
enrichment_KEGG #enrichment_GO <- enrich_module(Kmeans_cluster_ENS, org.Hs.eg.db, 'GO')
Moreover, you can do GO analysis
library(org.Hs.eg.db)
### Enrichment analysis (GO)
enrichment_GO <- enrich_module(clustering_ENS,org.db = org.Mm.eg.db,fun_num = 10,organism = 'mmu',enrich.db = 'GO')
You can visualize regulatory networks for enriched transcription factors of each module through plot_network() function by setting type parameter as ‘TF’. This plot shows regulatory relationships between transcription factors in different modules that significantly regulating other modules. The size of each vertex determine the significance of this transcription factor. Yellow edges are positive regulation, grey edges are negative regulation.
plot_tf_network(TFs_list)
You can visualize intramodular network with enriched function through plot_intramodular_network() function. Before run this function, you can select enriched functions of each module that you want to present in the plot. If you input all enriched functions, this function will automatically select the function with the highest -log10(qvalue) in each module to present in the plot. What’s more transcription factor with the most edge numbers in each module will be presented in the plot too.
### select functions that you want to present in the figure
<- enrichment_KEGG[11,]
enrichment_KEGG_select ### plotting
plot_intramodular_network(TFs_list,enrichment_KEGG_select,layout = 'circle')
It is strongly recommended to use Cytoscape to display the regulatory networks. We provide a function that can provide different Cytoscape styles. You need to install and open Cytoscape before running the function.
###optional: display the network in cytoscape, open cytoscape before running this function
initiate_cy(TFs_list, layout1='degree-circle', type='TF')
initiate_cy(TFs_list, layout1='grid', type='module')
These are the picture we processed through Cytoscape, which can show the regulatory relationship of modularized transcription factors.
Use Cytoscape to visualize regulatory network for enriched transcription factors of each module
Use Cytoscape to visualize intramodular network