This vignette provides the necessary instructions for performing the Partial Correlation coefficient with Information Theory (PCIT) (Reverter and Chan 2008) and Regulatory Impact Factors (RIF) (Reverter et al. 2010) algorithm.
The PCIT algorithm identifies meaningful correlations to define edges in a weighted network and can be applied to any correlation-based network including but not limited to gene co-expression networks, while the RIF algorithm identify critical transcript factors (TF) from gene expression data.
These two algorithms when combined provide a very relevant layer of information for gene expression studies (Microarray, RNA-seq and single-cell RNA-seq data).
A gene expression data from microarray, RNA-seq or single-cell RNA-seq spanning two biological conditions of interest (e.g. normal/tumor, healthy/disease, malignant/nonmalignant) is subjected to standard normalization techniques and significance analysis to identify the target genes whose expression is differentially expressed (DE) between the two conditions. Then, the regulators (e.g. Transcript Factors genes) are identified in the data. The TF genes can be obtained from the literature (Wang and Nishida 2015)(Vaquerizas et al. 2009). Next, the co-expression correlation between each TF and the DE genes is computed for each of the two conditions. This allows for the computation of the differential wiring (DW) from the difference in co-expression correlation existing between a TF and a DE genes in the two conditions. As a result, RIF analysis assigns an extreme score to those TF that are consistently most differentially co-expressed with the highly abundant and highly DE genes (case of RIF1 score), and to those TF with the most altered ability to act as predictors of the abundance of DE genes (case of RIF2 score). A given TF may not show a change in expression profile between the two conditions to score highly by RIF as long as it shows a big change in co-expression with the DE genes. To this particular, the profile of the TF gene (triangle, solid line) is identical in both conditions (slightly downwards). Instead, the DE gene (circle, dashed line) is clearly over-expressed in condition B. Importantly, the expression of the TF and the DE gene shows a strong positive correlation in condition A, and a strong negative correlation in condition B.
The proposed PCIT algorithm contains two distinct steps as follows:
For every trio of genes in x, y and z, the three first-order partial correlation coefficients are computed by:
\[r_{xy.z} = \frac{r_{xy} - r_{xz} r_{yz}}{\sqrt{(1-r^{2}_{xz})(1-r^{2}_{yz})}}\], and similarly for \(r_{xz.y}\) and \(r_{yz.x}\).
The partial correlation coefficient between x and y given z (here denoted by \(r_{xy.z}\)) indicates the strength of the linear relationship between x and y that is independent of (uncorrelated with) z. Calculating the ordinary (or unconditional or zero-order) correlation coefficient and comparing it with the partial correlation, we might see that the association between the two variables has been sharply reduced after eliminating the effect of the third variable.
We invoke the Data Processing Inequality (DPI) theorem of Information Theory which states that ‘no clever manipulation of the data can improve the inference that can be made from the data’ (Cover and Thomas 2012). For every trio of genes, and in order to obtain the tolerance level (\(\varepsilon\)) to be used as the local threshold for capturing significant associations, the average ratio of partial to direct correlation is computed as follows:
\[\varepsilon = (\frac{r_{xy.z}}{r_{xy}} + \frac{r_{xz.y}}{r_{xz}} + \frac{r_{yz.x}}{r_{yz}})\] In the context of our network reconstruction, a connection between genes x and y is discarded if:
\[|r_{xy}| \le |\varepsilon r_{xz}| and |r_{xy}| \le |\varepsilon r_{yz}|\] Otherwise, the association is defined as significant, and a connection between the pair of genes is established in the reconstruction of the GCN. To ascertain the significance of the association between genes x and y, the above mentioned Steps 1 and 2 are repeated for each of the remaining n−2 genes (denoted here by z).
To install, just type:
for Linux users is necessary to install libcurl4-openssl-dev, libxml2-dev and libssl-dev dependencies.
There are many ways to perform the analysis. The following sections will be splited by steps, and finishing with the complete analysis with visualization. We will use the airway (Himes et al. 2014) dataset in the following sections. This dataset provides a RNA-seq count data from four human ASM cell lines that were treated with dexamenthasone - a potent synthetic glucocorticoid. Briefly, this dataset has 4 samples untreated and other 4 samples with the treatment.
The first option is to perform the PCIT analysis. The output will be a list with 3 elements. The first one contains a dataframe with the pairwise correlation between genes (corr1) and the significant pairwise correlation (corr2 \(\neq\) 0). The second element of the list stores the adjacency matrix with all correlation. And the last element contains the adjacency matrix with only the significant values:
# Loading packages
library(CeTF)
library(airway)
library(kableExtra)
library(knitr)
# Loading airway data
data("airway")
# Creating a variable with annotation data
anno <- as.data.frame(colData(airway))
anno <- anno[order(anno$dex, decreasing = TRUE), ]
anno <- data.frame(cond = anno$dex, 
                   row.names = rownames(anno))
# Creating a variable with count data
counts <- assay(airway)
# Sorting count data samples by conditions (untrt and trt)
counts <- counts[, rownames(anno)]
colnames(counts) <- paste0(colnames(counts), c(rep("_untrt", 4), rep("_trt", 4)))
# Differential Expression analysis to use only informative genes
DEGenes <- expDiff(exp = counts,
                   anno = anno,
                   conditions = c('untrt', 'trt'),
                   lfc = 4.5,
                   padj = 0.05, 
                   diffMethod = "Reverter")
# Selecting only DE genes from counts data
counts <- counts[rownames(DEGenes$DE_unique), ]
# Converting count data to TPM
tpm <- apply(counts, 2, function(x) {
            (1e+06 * x)/sum(x)
        })
# Count normalization
PCIT_input <- normExp(tpm)
# PCIT input for untrt
PCIT_input_untrt <- PCIT_input[,grep("_untrt", colnames(PCIT_input))]
# PCIT input for trt
PCIT_input_trt <- PCIT_input[,grep("_trt", colnames(PCIT_input))]
# Performing PCIT analysis for untrt condition
PCIT_out_untrt <- PCIT(PCIT_input_untrt, tolType = "mean")
# Performing PCIT analysis for trt condition
PCIT_out_trt <- PCIT(PCIT_input_trt, tolType = "mean")
# Printing first 10 rows for untrt condition
kable(PCIT_out_untrt$tab[1:10, ]) %>% 
  kable_styling(bootstrap_options = "striped", full_width = FALSE)| gene1 | gene2 | corr1 | corr2 | 
|---|---|---|---|
| ENSG00000011465 | ENSG00000026025 | 0.01160 | 0.00000 | 
| ENSG00000011465 | ENSG00000035403 | 0.61286 | 0.00000 | 
| ENSG00000011465 | ENSG00000049323 | -0.71685 | 0.00000 | 
| ENSG00000011465 | ENSG00000067225 | -0.99935 | -0.99935 | 
| ENSG00000011465 | ENSG00000071127 | -0.56817 | 0.00000 | 
| ENSG00000011465 | ENSG00000071242 | -0.89063 | 0.00000 | 
| ENSG00000011465 | ENSG00000075624 | 0.09820 | 0.00000 | 
| ENSG00000011465 | ENSG00000077942 | -0.80938 | 0.00000 | 
| ENSG00000011465 | ENSG00000080824 | -0.72639 | 0.00000 | 
| ENSG00000011465 | ENSG00000087086 | -0.91357 | 0.00000 | 
# Printing first 10 rows for trt condition
kable(PCIT_out_trt$tab[1:10, ]) %>% 
  kable_styling(bootstrap_options = "striped", full_width = FALSE)| gene1 | gene2 | corr1 | corr2 | 
|---|---|---|---|
| ENSG00000011465 | ENSG00000026025 | 0.63800 | 0.00000 | 
| ENSG00000011465 | ENSG00000035403 | 0.96857 | 0.96857 | 
| ENSG00000011465 | ENSG00000049323 | 0.14794 | 0.00000 | 
| ENSG00000011465 | ENSG00000067225 | -0.99687 | -0.99687 | 
| ENSG00000011465 | ENSG00000071127 | -0.92734 | -0.92734 | 
| ENSG00000011465 | ENSG00000071242 | -0.75914 | 0.00000 | 
| ENSG00000011465 | ENSG00000075624 | -0.70489 | 0.00000 | 
| ENSG00000011465 | ENSG00000077942 | -0.29712 | 0.00000 | 
| ENSG00000011465 | ENSG00000080824 | -0.93148 | 0.00000 | 
| ENSG00000011465 | ENSG00000087086 | -0.94620 | 0.00000 | 
After performing the PCIT analysis, it is possible to verify the histogram distribution of the clustering coefficient of the adjacency matrix with the significant values:
It’s possible to generate the density plot with the significance values of correlation. We’ll use the raw adjacency matrix and the adjacency matrix with significant values. It is necessary to define a cutoff of the correlation module (values between -1 and 1) that will be considered as significant:
# Example for trt condition
densityPlot(mat1 = PCIT_out_trt$adj_raw, 
           mat2 = PCIT_out_trt$adj_sig, 
           threshold = 0.5)To perform the RIF analysis we will need the count data, an annotation table and a list with the Transcript Factors of specific organism (Homo sapiens in this case) and follow the following steps in order to get the output (dataframe with the average expression, RIF1 and RIF2 metrics for each TF):
# Loading packages
library(CeTF)
library(airway)
library(kableExtra)
library(knitr)
# Loading airway data
data("airway")
# Creating a variable with annotation data
anno <- as.data.frame(colData(airway))
anno <- anno[order(anno$dex, decreasing = TRUE), ]
anno <- data.frame(cond = anno$dex, 
                   row.names = rownames(anno))
# Creating a variable with count data
counts <- assay(airway)
# Sorting count data samples by conditions (untrt and trt)
counts <- counts[, rownames(anno)]
colnames(counts) <- paste0(colnames(counts), c(rep("_untrt", 4), rep("_trt", 4)))
# Differential Expression analysis to use only informative genes
DEGenes <- expDiff(exp = counts,
                   anno = anno,
                   conditions = c('untrt', 'trt'),
                   lfc = 2,
                   padj = 0.05, 
                   diffMethod = "Reverter")
# Selecting only DE genes from counts data
counts <- counts[rownames(DEGenes$DE_unique), ]
# Converting count data to TPM
tpm <- apply(counts, 2, function(x) {
            (1e+06 * x)/sum(x)
        })
# Count normalization
Clean_Dat <- normExp(tpm)
# Loading the Transcript Factors (TFs) character
data("TFs")
# Verifying which TFs are in the subsetted normalized data
TFs <- rownames(Clean_Dat)[rownames(Clean_Dat) %in% TFs]
# Selecting the Target genes
Target <- setdiff(rownames(Clean_Dat), TFs)
# Ordering rows of normalized count data
RIF_input <- Clean_Dat[c(Target, TFs), ]
# Performing RIF analysis
RIF_out <- RIF(input = RIF_input,
               nta = length(Target),
               ntf = length(TFs),
               nSamples1 = 4,
               nSamples2 = 4)
# Printing first 10 rows
kable(RIF_out[1:10, ]) %>% 
  kable_styling(bootstrap_options = "striped", full_width = FALSE)| TF | avgexpr | RIF1 | RIF2 | 
|---|---|---|---|
| ENSG00000008441 | 10.380590 | -0.1386859 | -0.6770894 | 
| ENSG00000102804 | 9.264843 | 0.0147001 | 0.1789910 | 
| ENSG00000103888 | 13.591598 | 0.0985383 | -0.3361795 | 
| ENSG00000103994 | 10.385093 | -1.2748364 | 0.8113721 | 
| ENSG00000112837 | 8.143177 | -0.2246323 | -1.7903370 | 
| ENSG00000116016 | 11.376028 | -1.2315044 | 0.0216295 | 
| ENSG00000116132 | 9.388008 | 1.4038890 | 0.8823288 | 
| ENSG00000118689 | 9.189257 | 0.9971123 | 0.3981473 | 
| ENSG00000123562 | 9.207009 | 0.3381318 | 0.4609342 | 
| ENSG00000157514 | 7.215664 | -1.5628941 | -0.3290858 | 
Finally, it is possible to run the entire analysis all at once. The output will be a CeTF object with all results generated between the different steps. To access the CeTF object is recommended to use the accessors from CeTF class:
# Loading packages
library(airway)
library(CeTF)
# Loading airway data
data("airway")
# Creating a variable with annotation data
anno <- as.data.frame(colData(airway))
# Creating a variable with count data
counts <- assay(airway)
# Sorting count data samples by conditions (untrt and trt)
counts <- counts[, order(anno$dex, decreasing = TRUE)]
# Loading the Transcript Factors (TFs) character
data("TFs")
# Performing the complete analysis
out <- runAnalysis(mat = counts, 
                   conditions=c("untrt", "trt"),
                   lfc = 5,
                   padj = 0.05,
                   TFs = TFs,
                   nSamples1 = 4,
                   nSamples2= 4,
                   tolType = "mean",
                   diffMethod = "Reverter", 
                   data.type = "counts")Based on CeTF class object resulted from runAnalysis function is possible to get a plot for Differentially Expressed (DE) genes that shows the relationship between log(baseMean) and Difference of Expression or log2FoldChange, enabling the visualization of the distribution of DE genes and TF in both conditions. These two first plot provides a initial graphical visualization of results:
# Using the runAnalysis output (CeTF class object)
SmearPlot(object = out, 
          diffMethod = 'Reverter', 
          lfc = 1.5, 
          conditions = c('untrt', 'trt'), 
          type = "DE")Then is possible to generate the same previous plot but now for a single TF. So, if you have a specific TF that you want to visualize, this is the recommended plot:
# Using the runAnalysis output (CeTF class object)
SmearPlot(object = out,
          diffMethod = 'Reverter',
          lfc = 1.5,
          conditions = c('untrt', 'trt'),
          TF = 'ENSG00000185917',
          label = TRUE, 
          type = "TF")The next output is a network for both conditions resulted from :
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.Then we will associate the genes in both conditions with Gene Ontology (GO) and other databases. getGroupGO function will be performed where a set of genes are counted within all possible ontologies, without statistical power. This function is capable of perform groupGO for any organism as long as it’s annotation package exists (i.e. org.Hs.eg.db, org.Bt.eg.db, org.Rn.eg.db, etc). In this example we will perform only for first condition:
# Loading Homo sapiens annotation package
library(org.Hs.eg.db)
# Accessing the network for condition 1
genes <- unique(c(as.character(NetworkData(out, "network1")[, "gene1"]), 
                  as.character(NetworkData(out, "network1")[, "gene2"])))
# Performing getGroupGO analysis
cond1 <- getGroupGO(genes = genes, 
                    ont = "BP", 
                    keyType = "ENSEMBL", 
                    annoPkg = org.Hs.eg.db, 
                    level = 3)Alternatively, it is possible to perform the enrichment analysis with a statistical power. In this case, there are many databases options to perform this analysis (i.e. GO, KEGG, REACTOME, etc). The getEnrich function will perform the enrichment analysis and returns which pathways are enriched in a given set of genes. This analysis requires a reference list of genes with all genes that will be considered to enrich. In this package there is a function, refGenes that has these references genes for ENSEMBL and SYMBOL nomenclatures for 5 organisms: Human (Homo sapiens), Mouse (Mus musculus), Zebrafish (Danio rerio), Cow (Bos taurus) and Rat (Rattus norvegicus). If the user wants to use a different reference set, simply inputs a character vector with the genes.
# Accessing the network for condition 1
genes <- unique(c(as.character(NetworkData(out, "network1")[, "gene1"]), 
                  as.character(NetworkData(out, "network1")[, "gene2"])))
# Performing getEnrich analysis
enrich <- getEnrich(organism='hsapiens', database='geneontology_Biological_Process', 
                   genes=genes, GeneType='ensembl_gene_id', 
                   refGene=refGenes$Homo_sapiens$ENSEMBL, 
                   fdrMethod = "BH", fdrThr = 0.05, minNum = 5, maxNum = 500)In the next steps we will use the getGroupGO results, but it is worth mentioning that you can use the results of the getEnrich function. To run the following steps with the getEnrich results is recommended to change the lfc parameter in aboce runAnalysis function by 3. That way there will be more pathways to enrich. Remembering that this will require more processing time.
After the association of GOs and genes, we are able to plot the network of previously results. It’s possible to choose which variables (“pathways”, “TFs” or “genes”) we want to group. All these variables are inputted by the user, i.e. the user can choose the pathways, TFs or genes returned from analysis from the runAnalysis function and CeTF class object, or can also choose by specific pathways, TFs or genes of interest. The table used to plot the networks are stored in a list resulted from netGOTFPlot function. To get access is needed to use the object$tab. The list will be splitted by groupBy argument choice:
library(dplyr)
# Subsetting only the first 12 Ontologies with more counts
t1 <- head(cond1$results, 12)
# Subsetting the network for the conditions to make available only the 12 nodes subsetted
t2 <- subset(cond1$netGO, cond1$netGO$gene1 %in% as.character(t1[, "ID"]))
# generating the GO plot grouping by pathways
pt <- netGOTFPlot(netCond = NetworkData(out, "network1") %>% 
                    mutate(across(everything(), as.character)),
                  resultsGO = t1,
                  netGO = t2,
                  anno = NetworkData(out, "annotation"),
                  groupBy = 'pathways',
                  type = 'GO')
pt$plot##               gene1    pathway           gene2 class
## 1   ENSG00000004799 GO:0006807 ENSG00000114861  gene
## 2   ENSG00000004799 GO:0006807 ENSG00000119138  gene
## 3   ENSG00000004799 GO:0006807 ENSG00000120093  gene
## 4   ENSG00000004799 GO:0006807 ENSG00000184271  gene
## 49  ENSG00000013293 GO:0019222 ENSG00000114861  gene
## 108 ENSG00000030419 GO:0019222 ENSG00000105989    TF# generating the GO plot grouping by TFs
TFs <- NetworkData(out, "keytfs")$TF[1:4]
pt <- netGOTFPlot(netCond = NetworkData(out, "network1"),
                  resultsGO = t1,
                  netGO = t2,
                  anno = NetworkData(out, "annotation"),
                  groupBy = 'TFs',
                  TFs = TFs, 
                  type = 'GO')
pt$plot## NULL# generating the GO plot grouping by specifics genes
genes <- c('ENSG00000011465', 'ENSG00000026025', 'ENSG00000075624', 'ENSG00000077942')
pt <- netGOTFPlot(netCond = NetworkData(out, "network1"),
                  resultsGO = t1,
                  netGO = t2,
                  anno = NetworkData(out, "annotation"),
                  groupBy = 'genes',
                  genes = genes,
                  type = 'GO')
pt$plot##              gene1           genes gene2 class
## 1  ENSG00000004799 ENSG00000011465  <NA>  gene
## 5  ENSG00000013293 ENSG00000011465  <NA>  gene
## 9  ENSG00000030419 ENSG00000011465  <NA>    TF
## 13 ENSG00000053254 ENSG00000011465  <NA>    TF
## 17 ENSG00000064961 ENSG00000011465  <NA>    TF
## 21 ENSG00000067082 ENSG00000011465  <NA>    TFFinally, we will merge the results from the network of condition 1 (Gene-Gene and TF-Gene interaction) with groupGO results (GO-Gene or GO-TF interaction). The final results has the objective to integrate all these data in order to results in a complete view of the data:
pt <- netGOTFPlot(netCond = NetworkData(out, "network1") %>% 
                   mutate(across(everything(), as.character)),
                  netGO = t2,
                  keyTFs = NetworkData(out, "keytfs"), 
                  type = 'Integrated')
pt$plot##              gene1           gene2
## 22 ENSG00000004799 ENSG00000113761
## 24 ENSG00000004799 ENSG00000114861
## 26 ENSG00000004799 ENSG00000119138
## 27 ENSG00000004799 ENSG00000120093
## 47 ENSG00000004799 ENSG00000143127
## 77 ENSG00000004799 ENSG00000184271Is possible to use the accessors from CeTF class object to access the outputs generated by the runAnalysis function. These outputs can be used as input in Cytoscape (Shannon et al. 2003). The accessors are:
All accessors can be further explored by looking in more detail at the documentation.
This package has some additional features to plot the results. Let’s use the gene set from airway data previously used. The features are:
Network propagation is an important and widely used algorithm in systems biology, with applications in protein function prediction, disease gene prioritization, and patient stratification. Propagation provides a robust estimate of network distance between sets of nodes. Network propagation uses the network of interactions to find new genes that are most relevant to a set of well understood genes. To perform this analysis is necessary to install the latest Cytoscape software version and know the path that will be installed the software. After running the diffusion analysis is possible to perform the enrichment for the subnetwork resulted. Finally, ih this vignette we’ll not perform the diffusion analysis because this requires Cytoscape to be installed.
This plot makes it possible to visualize the targets of specific TFs or genes. The main idea is to identify which chromosomes are linked to the targets of a given gene or TF to infer whether there are cis (same chromosome) or trans (different chromosomes) links between them. The black links are between different chromosomes while the red links are between the same chromosome.
# Loading the CeTF object class demo
data("CeTFdemo")
CircosTargets(object = CeTFdemo, 
              file = "/path/to/gtf/file/specie.gtf", 
              nomenclature = "ENSEMBL", 
              selection = "ENSG00000185591", 
              cond = "condition1")To visualize the relationship between RIF results (RIF1 and RIF2) is possible to generate the following plot:
Then to get the relationship between RIF1/RIF2 and differential expressed genes:
We provide some option to plot the results from getEnrich function. First of all we need to perform the enrichment analysis for a set of genes:
# Selecting the genes related to the network for condition 1
genes <- unique(c(as.character(NetworkData(out, "network1")[, "gene1"]), 
                  as.character(NetworkData(out, "network1")[, "gene2"])))
# Performing enrichment analysis
cond1 <- getEnrich(organism='hsapiens', database='geneontology_Biological_Process', 
                   genes=genes, GeneType='ensembl_gene_id',
                   refGene=refGenes$Homo_sapiens$ENSEMBL,
                   fdrMethod = "BH", fdrThr = 0.05, minNum = 5, maxNum = 500)## Loading the functional categories...
## Loading the ID list...
## Loading the reference list...
## Performing the enrichment analysis...After performing the enrichment analysis we have 4 options for viewing this data:
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.13-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.13-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              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] parallel  stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] dplyr_1.0.7                 org.Hs.eg.db_3.13.0        
##  [3] AnnotationDbi_1.54.1        knitr_1.34                 
##  [5] kableExtra_1.3.4            airway_1.12.0              
##  [7] SummarizedExperiment_1.22.0 Biobase_2.52.0             
##  [9] GenomicRanges_1.44.0        GenomeInfoDb_1.28.4        
## [11] IRanges_2.26.0              S4Vectors_0.30.0           
## [13] BiocGenerics_0.38.0         MatrixGenerics_1.4.3       
## [15] matrixStats_0.60.1          CeTF_1.4.5                 
## 
## loaded via a namespace (and not attached):
##   [1] pbdZMQ_0.3-5                     GGally_2.1.2                    
##   [3] R.methodsS3_1.8.1                clinfun_1.0.15                  
##   [5] coda_0.19-4                      tidyr_1.1.3                     
##   [7] ggplot2_3.3.5                    bit64_4.0.5                     
##   [9] DelayedArray_0.18.0              R.utils_2.10.1                  
##  [11] data.table_1.14.0                KEGGREST_1.32.0                 
##  [13] RCurl_1.98-1.4                   doParallel_1.0.16               
##  [15] generics_0.1.0                   cowplot_1.1.1                   
##  [17] RSQLite_2.2.8                    shadowtext_0.0.8                
##  [19] bit_4.0.4                        tzdb_0.1.2                      
##  [21] enrichplot_1.12.2                base64url_1.4                   
##  [23] webshot_0.5.2                    xml2_1.3.2                      
##  [25] assertthat_0.2.1                 viridis_0.6.1                   
##  [27] xfun_0.25                        hms_1.1.0                       
##  [29] jquerylib_0.1.4                  evaluate_0.14                   
##  [31] fansi_0.5.0                      readxl_1.3.1                    
##  [33] igraph_1.2.6                     DBI_1.1.1                       
##  [35] geneplotter_1.70.0               reshape_0.8.8                   
##  [37] GenomicTools.fileHandler_0.1.5.9 apcluster_1.4.8                 
##  [39] purrr_0.3.4                      ellipsis_0.3.2                  
##  [41] ggpubr_0.4.0                     backports_1.2.1                 
##  [43] signal_0.7-7                     annotate_1.70.0                 
##  [45] vctrs_0.3.8                      Cairo_1.5-12.2                  
##  [47] abind_1.4-5                      withr_2.4.2                     
##  [49] cachem_1.0.6                     ggforce_0.3.3                   
##  [51] vroom_1.5.4                      sna_2.6                         
##  [53] treeio_1.16.2                    svglite_2.0.0                   
##  [55] cluster_2.1.2                    DOSE_3.18.2                     
##  [57] ape_5.5                          IRdisplay_1.0                   
##  [59] lazyeval_0.2.2                   crayon_1.4.1                    
##  [61] uchardet_1.1.0                   genefilter_1.74.0               
##  [63] pkgconfig_2.0.3                  labeling_0.4.2                  
##  [65] tweenr_1.0.2                     nlme_3.1-153                    
##  [67] rlang_0.4.11                     RJSONIO_1.3-1.5                 
##  [69] lifecycle_1.0.0                  downloader_0.4                  
##  [71] cellranger_1.1.0                 polyclip_1.10-0                 
##  [73] graph_1.70.0                     rngtools_1.5                    
##  [75] Matrix_1.3-4                     aplot_0.1.0                     
##  [77] IRkernel_1.2                     carData_3.0-4                   
##  [79] base64enc_0.1-3                  whisker_0.4                     
##  [81] GlobalOptions_0.1.2              WebGestaltR_0.4.4               
##  [83] png_0.1-7                        viridisLite_0.4.0               
##  [85] rjson_0.2.20                     bitops_1.0-7                    
##  [87] R.oo_1.24.0                      ggnetwork_0.5.10                
##  [89] dplR_1.7.2                       Biostrings_2.60.2               
##  [91] blob_1.2.2                       doRNG_1.8.2                     
##  [93] shape_1.4.6                      stringr_1.4.0                   
##  [95] qvalue_2.24.0                    readr_2.0.1                     
##  [97] rstatix_0.7.0                    gridGraphics_0.5-1              
##  [99] gMWT_1.1.1                       ggsignif_0.6.3                  
## [101] scales_1.1.1                     memoise_2.0.0                   
## [103] magrittr_2.0.1                   plyr_1.8.6                      
## [105] zlibbioc_1.38.0                  compiler_4.1.1                  
## [107] scatterpie_0.1.7                 RColorBrewer_1.1-2              
## [109] clue_0.3-59                      DESeq2_1.32.0                   
## [111] XVector_0.32.0                   patchwork_1.1.1                 
## [113] MASS_7.3-54                      tidyselect_1.1.1                
## [115] stringi_1.7.4                    forcats_0.5.1                   
## [117] highr_0.9                        yaml_2.2.1                      
## [119] GOSemSim_2.18.1                  locfit_1.5-9.4                  
## [121] ggrepel_0.9.1                    grid_4.1.1                      
## [123] sass_0.4.0                       fastmatch_1.1-3                 
## [125] tools_4.1.1                      rio_0.5.27                      
## [127] circlize_0.4.13                  rstudioapi_0.13                 
## [129] uuid_0.1-4                       foreach_1.5.1                   
## [131] foreign_0.8-81                   gridExtra_2.3                   
## [133] farver_2.1.0                     ggraph_2.0.5                    
## [135] digest_0.6.27                    snpStats_1.42.0                 
## [137] Rcpp_1.0.7                       car_3.0-11                      
## [139] broom_0.7.9                      httr_1.4.2                      
## [141] ComplexHeatmap_2.8.0             colorspace_2.0-2                
## [143] rvest_1.0.1                      XML_3.99-0.7                    
## [145] splines_4.1.1                    yulab.utils_0.0.2               
## [147] tidytree_0.3.5                   graphlayouts_0.7.1              
## [149] ggplotify_0.1.0                  systemfonts_1.0.2               
## [151] xtable_1.8-4                     jsonlite_1.7.2                  
## [153] ggtree_3.0.4                     tidygraph_1.2.0                 
## [155] ggfun_0.0.3                      R6_2.5.1                        
## [157] pillar_1.6.2                     htmltools_0.5.2                 
## [159] glue_1.4.2                       fastmap_1.1.0                   
## [161] clusterProfiler_4.0.5            BiocParallel_1.26.2             
## [163] codetools_0.2-18                 fgsea_1.18.0                    
## [165] mvtnorm_1.1-2                    utf8_1.2.2                      
## [167] lattice_0.20-44                  bslib_0.3.0                     
## [169] tibble_3.1.4                     network_1.17.1                  
## [171] curl_4.3.2                       magick_2.7.3                    
## [173] zip_2.2.0                        GO.db_3.13.0                    
## [175] openxlsx_4.2.4                   survival_3.2-13                 
## [177] rmarkdown_2.10                   statnet.common_4.5.0            
## [179] GenomicTools_0.2.9.7             repr_1.1.3                      
## [181] munsell_0.5.0                    DO.db_2.9                       
## [183] GetoptLong_1.0.5                 GenomeInfoDbData_1.2.6          
## [185] iterators_1.0.13                 RCy3_2.12.4                     
## [187] haven_2.4.3                      reshape2_1.4.4                  
## [189] gtable_0.3.0Cover, Thomas M, and Joy A Thomas. 2012. Elements of Information Theory. John Wiley & Sons.
Himes, Blanca E, Xiaofeng Jiang, Peter Wagner, Ruoxi Hu, Qiyu Wang, Barbara Klanderman, Reid M Whitaker, et al. 2014. “RNA-Seq Transcriptome Profiling Identifies Crispld2 as a Glucocorticoid Responsive Gene That Modulates Cytokine Function in Airway Smooth Muscle Cells.” PloS One 9 (6): e99625.
Reverter, Antonio, and Eva KF Chan. 2008. “Combining Partial Correlation and an Information Theory Approach to the Reversed Engineering of Gene Co-Expression Networks.” Bioinformatics 24 (21): 2491–7.
Reverter, Antonio, Nicholas J Hudson, Shivashankar H Nagaraj, Miguel Pérez-Enciso, and Brian P Dalrymple. 2010. “Regulatory Impact Factors: Unraveling the Transcriptional Regulation of Complex Traits from Expression Data.” Bioinformatics 26 (7): 896–904.
Shannon, Paul, Andrew Markiel, Owen Ozier, Nitin S Baliga, Jonathan T Wang, Daniel Ramage, Nada Amin, Benno Schwikowski, and Trey Ideker. 2003. “Cytoscape: A Software Environment for Integrated Models of Biomolecular Interaction Networks.” Genome Research 13 (11): 2498–2504.
Vaquerizas, Juan M, Sarah K Kummerfeld, Sarah A Teichmann, and Nicholas M Luscombe. 2009. “A Census of Human Transcription Factors: Function, Expression and Evolution.” Nature Reviews Genetics 10 (4): 252.
Wang, Kai, and Hiroki Nishida. 2015. “REGULATOR: A Database of Metazoan Transcription Factors and Maternal Factors for Developmental Studies.” BMC Bioinformatics 16 (1): 114.