Contents

1 Version Info

## See system.file("LICENSE", package="MotifDb") for use restrictions.

R version: R version 4.4.1 (2024-06-14)
Bioconductor version: 3.20
Package version: 1.30.0

2 Introduction

Eukaryotic gene regulation can be very complex. Transcription factor binding to promoter DNA sequences is a stochastic process, and imperfect matches can be sufficient for binding. Chromatin remodeling, methylation, histone modification, chromosome interaction, distal enhancers, and the cooperative binding of transcription co-factors all play an important role. We avoid most of this complexity in this demonstration workflow in order to examine transcription factor binding sites in a small set of seven broadly co-expressed Saccharomyces cerevisiae genes of related function. These genes exhibit highly correlated mRNA expression across 200 experimental conditions, and are annotated to Nitrogen Catabolite Repression (NCR), the means by which yeast cells switch between using rich and poor nitrogen sources.

We will see, however, that even this small collection of co-regulated genes of similar function exhibits considerable regulatory complexity, with (among other things) activators and repressors competing to bind to the same DNA promoter sequence. Our case study sheds some light on this complexity, and demonstrates how several new Bioconductor packages and methods allow us to

[ Back to top ]

3 Installation and Use

To install the necessary packages and all of their dependencies, evaluate the commands

## try http:// if https:// URLs are not supported
library(BiocManager)
BiocManager::install(c("MotifDb",  "GenomicFeatures", 
           "TxDb.Scerevisiae.UCSC.sacCer3.sgdGene",
           "org.Sc.sgd.db", "BSgenome.Scerevisiae.UCSC.sacCer3",
           "motifStack", "seqLogo"))

Package installation is required only once per R installation. When working with an organism other than S.cerevisiae, substitute the three species-specific packages as needed.

To use these packages in an R session, evaluate these commands:

library(MotifDb)
library(S4Vectors)
library(seqLogo)
library(motifStack)
library(Biostrings)
library(GenomicFeatures)
library(org.Sc.sgd.db)
library(BSgenome.Scerevisiae.UCSC.sacCer3)
library(TxDb.Scerevisiae.UCSC.sacCer3.sgdGene)

These instructions are required once in each R session.

[ Back to top ]

4 Biological Background

The x-y plot below displays expression levels of seven genes across 200 conditions, from a compendium of yeast expression data which accompanies Allocco et al, 2004, “Quantifying the relationship between co-expression, co-regulation and gene function”:

compendium.png


Allocco et al establish that

In S. cerevisiae, two genes have a 50% chance of having a common transcription factor binder if the correlation between their expression profiles is equal to 0.84.

These seven highly-correlated (> 0.85) NCR genes form a connected subnetwork within the complete co-expresson network derived from the compendium data (work not shown). Network edges indicate correlated expression of the two connected genes across all 200 conditions. The edges are colored as a function of that correlation: red for perfect correlation, white indicating correlation of 0.85, and intermediate colors for intermediate values. DAL80 is rendered as an octagon to indicate its special status as a transcription factor. We presume, following Allocco, that such correlation among genes, including one transcription factor, is a plausible place to look for shared transcription factor binding sites.


dal80-subnet.png


Some insight into the co-regulation of these seven genes is obtained from Georis et al, 2009, “The Yeast GATA Factor Gat1 Occupies a Central Position in Nitrogen Catabolite Repression-Sensitive Gene Activation”:

Saccharomyces cerevisiae cells are able to adapt their metabolism according to the quality of the nitrogen sources available in the environment. Nitrogen catabolite repression (NCR) restrains the yeast’s capacity to use poor nitrogen sources when rich ones are available. NCR-sensitive expression is modulated by the synchronized action of four DNA-binding GATA factors. Although the first identified GATA factor, Gln3, was considered the major activator of NCR-sensitive gene expression, our work positions Gat1 as a key factor for the integrated control of NCR in yeast for the following reasons: (i) Gat1 appeared to be the limiting factor for NCR gene expression, (ii) GAT1 expression was regulated by the four GATA factors in response to nitrogen availability, (iii) the two negative GATA factors Dal80 and Gzf3 interfered with Gat1 binding to DNA, and (iv) Gln3 binding to some NCR promoters required Gat1. Our study also provides mechanistic insights into the mode of action of the two negative GATA factors. Gzf3 interfered with Gat1 by nuclear sequestration and by competition at its own promoter. Dal80-dependent repression of NCR-sensitive gene expression occurred at three possible levels: Dal80 represses GAT1 expression, it competes with Gat1 for binding, and it directly represses NCR gene transcription. (emphasis added)

Thus DAL80 is but one of four interacting transcription factors which all bind the GATA motif. We will see below that DAL80 lacks the GATA sequence in its own promoter, but that the motif is well-represented in the promoters of the other six.

In order to demonstrate Bioconductor capabilities for finding binding sites for known transcription factors via sequence matching, we will use the shared DNA-binding GATA sequence as retrieved from one of those factors from MotifDb, DAL80.

[ Back to top ]

6 Minimal Example

Only eight lines of code (excluding library statements) are required to find two matches to the JASPAR DAL80 motif in the promoter of DAL1.

library(MotifDb)
library(seqLogo)
library(motifStack)
library(Biostrings)
library(GenomicFeatures)
library(org.Sc.sgd.db)
library(BSgenome.Scerevisiae.UCSC.sacCer3)
library(TxDb.Scerevisiae.UCSC.sacCer3.sgdGene)

query(MotifDb, "DAL80")   
## MotifDb object of length 6
## | Created from downloaded public sources, last update: 2022-Mar-04
## | 6 position frequency matrices from 6 sources:
## |        JASPAR_2014:    1
## |        JASPAR_CORE:    1
## |             ScerTF:    1
## |         jaspar2016:    1
## |         jaspar2018:    1
## |         jaspar2022:    1
## | 1 organism/s
## |        Scerevisiae:    6
## Scerevisiae-ScerTF-DAL80-harbison 
## Scerevisiae-JASPAR_CORE-DAL80-MA0289.1 
## Scerevisiae-JASPAR_2014-DAL80-MA0289.1 
## Scerevisiae-jaspar2016-DAL80-MA0289.1 
## Scerevisiae-jaspar2018-DAL80-MA0289.1 
## Scerevisiae-jaspar2022-DAL80-MA0289.1
pfm.dal80.jaspar <- query(MotifDb,"DAL80")[[1]]
seqLogo(pfm.dal80.jaspar)

dal1 <- "YIR027C"
chromosomal.loc <- 
  transcriptsBy(TxDb.Scerevisiae.UCSC.sacCer3.sgdGene, by="gene") [dal1]
promoter.dal1 <- 
  getPromoterSeq(chromosomal.loc, Scerevisiae, upstream=1000, downstream=0)
pcm.dal80.jaspar <- round(100 * pfm.dal80.jaspar)
matchPWM(pcm.dal80.jaspar, unlist(promoter.dal1)[[1]], "90%")
## Views on a 1000-letter DNAString subject
## subject: TTGAGGAGTTGTCCACATACACATTAGTGTTGAT...GCAAAAAAAAAGTGAAATACTGCGAAGAACAAAG
## views:
##       start end width
##   [1]   621 625     5 [GATAA]
##   [2]   638 642     5 [GATAA]

[ Back to top ]

7 Sample Workflow: an extended example

We begin by visualizing DAL80’s TF binding motif using either of two Bioconductor packages: seqLogo, and motifStack. First, query MotifDb for the PFM (position frequency matrix):

query(MotifDb,"DAL80")
## MotifDb object of length 6
## | Created from downloaded public sources, last update: 2022-Mar-04
## | 6 position frequency matrices from 6 sources:
## |        JASPAR_2014:    1
## |        JASPAR_CORE:    1
## |             ScerTF:    1
## |         jaspar2016:    1
## |         jaspar2018:    1
## |         jaspar2022:    1
## | 1 organism/s
## |        Scerevisiae:    6
## Scerevisiae-ScerTF-DAL80-harbison 
## Scerevisiae-JASPAR_CORE-DAL80-MA0289.1 
## Scerevisiae-JASPAR_2014-DAL80-MA0289.1 
## Scerevisiae-jaspar2016-DAL80-MA0289.1 
## Scerevisiae-jaspar2018-DAL80-MA0289.1 
## Scerevisiae-jaspar2022-DAL80-MA0289.1

There are two motifs. How do they compare? The seqlogo package has been the standard tool for viewing sequence logos, but can only portray one logo at a time.

dal80.jaspar <- query(MotifDb,"DAL80")[[1]]
dal80.scertf <-query(MotifDb,"DAL80")[[2]]
seqLogo(dal80.jaspar)

seqLogo(dal80.scertf)

With a little preparation, the new (October 2012) package motifStack can plot both motifs together. First, create instances of the pfm class:

pfm.dal80.jaspar <- new("pfm", mat=query(MotifDb, "dal80")[[1]], 
                        name="DAL80-JASPAR")
pfm.dal80.scertf <- new("pfm", mat=query(MotifDb, "dal80")[[2]], 
                        name="DAL80-ScerTF")
plotMotifLogoStack(DNAmotifAlignment(c(pfm.dal80.scertf, pfm.dal80.jaspar)))
## Loading required namespace: Cairo
## Warning in checkValidSVG(doc, warn = warn): This picture may not have been
## generated by Cairo graphics; errors may result
## Warning in checkValidSVG(doc, warn = warn): This picture may not have been
## generated by Cairo graphics; errors may result
## Warning in checkValidSVG(doc, warn = warn): This picture may not have been
## generated by Cairo graphics; errors may result
## Warning in checkValidSVG(doc, warn = warn): This picture may not have been
## generated by Cairo graphics; errors may result

Of these two, the JASPAR motif has more detail, but the ScerTF motif is more recently published. ScerTF has a reputation for careful yeast-specific curation. We will use the ScerTF version.

Georis et al mention that DAL80 “competes with Gat1 for binding” – suggesting that they would have highly similar motifs. MotifDb has 3 entries for GAT1:

query(MotifDb, "gat1")
## MotifDb object of length 7
## | Created from downloaded public sources, last update: 2022-Mar-04
## | 7 position frequency matrices from 7 sources:
## |        JASPAR_2014:    1
## |        JASPAR_CORE:    1
## |             ScerTF:    1
## |           UniPROBE:    1
## |         jaspar2016:    1
## |         jaspar2018:    1
## |         jaspar2022:    1
## | 1 organism/s
## |        Scerevisiae:    7
## Scerevisiae-ScerTF-GAT1-zhu 
## Scerevisiae-JASPAR_CORE-GAT1-MA0300.1 
## Scerevisiae-JASPAR_2014-GAT1-MA0300.1 
## Scerevisiae-jaspar2016-GAT1-MA0300.1 
## Scerevisiae-jaspar2018-GAT1-MA0300.1 
## Scerevisiae-jaspar2022-GAT1-MA0300.1 
## Scerevisiae-UniPROBE-Gat1.UP00287

Plot the three together:

pfm.gat1.jaspar = new("pfm", mat=query(MotifDb, "gat1")[[1]], 
                       name="GAT1-JASPAR")
pfm.gat1.scertf = new("pfm", mat=query(MotifDb, "gat1")[[2]], 
                       name="GAT1-ScerTF")
pfm.gat1.uniprobe = new("pfm", mat=query(MotifDb, "gat1")[[3]], 
                       name="GAT1-UniPROBE")
plotMotifLogoStack(c(pfm.gat1.uniprobe, pfm.gat1.scertf, pfm.gat1.jaspar))
## Warning in checkValidSVG(doc, warn = warn): This picture may not have been
## generated by Cairo graphics; errors may result
## Warning in checkValidSVG(doc, warn = warn): This picture may not have been
## generated by Cairo graphics; errors may result
## Warning in checkValidSVG(doc, warn = warn): This picture may not have been
## generated by Cairo graphics; errors may result
## Warning in checkValidSVG(doc, warn = warn): This picture may not have been
## generated by Cairo graphics; errors may result

The GAT1-JASPAR motif is very similar to DAL80’s GATAA motif, and thus consistent with the Georis claim that GAT1 and DAL80 compete for binding. The GAT1-ScerTF and GAT1-UniPROBE motifs are similar, but differ in length. They are reverse complements of the canonical GATAA motif.

To match motifs in a promoter, these steps are required:

The three search motifs, one for DAL80, and two for GAT1, must be transformed before then can be matched to DNA sequence. MotifDb returns a position frequency matrix (a PFM) with all columns summing to 1.0, but the Biostrings matchPWM method expects a position count matrix (a PCM) with integer values. Transform the frequency matrix into a count matrix using the somewhat arbitrary but functionally reliable scaling factor of 100:

pfm.dal80.scertf <- query(MotifDb, "dal80")[[2]]
pcm.dal80.scertf <- round(100 * pfm.dal80.scertf)

pfm.gat1.jaspar <- query(MotifDb, "gat1")[[1]]
pcm.gat1.jaspar <- round(100 * pfm.gat1.jaspar)

pfm.gat1.scertf <- query(MotifDb, "gat1")[[2]]
pcm.gat1.scertf <- round(100 * pfm.gat1.scertf)

Create a list of the seven genes from the DAL80 co-regulated subnetwork (displayed above). Lookup their systematic names, which will be needed immediately below.

genes <- c("DAL1", "DAL2", "DAL4", "DAL5", "DAL7", "DAL80", "GAP1")
orfs <- as.character(mget(genes, org.Sc.sgdCOMMON2ORF))

Obtain the coordinates of the transcripts for the orfs.

    grl <- transcriptsBy(TxDb.Scerevisiae.UCSC.sacCer3.sgdGene, by="gene") [orfs]

These coordinates, returned in a GRangesList object, specify the start location (chromosome and base pair) of every known transcript for each gene. With this information, GenomicFeatures::getPromoterSeq calculates and returns the DNA sequence of the promoter:

promoter.seqs <- getPromoterSeq(grl, Scerevisiae, upstream=1000,
                                downstream=0)

We will next search for a match of the motif to the first of these sequences, the promoter for DAL1. Note that here, and below, we use a 90% “min.score” when we call matchPWM. This high minimum match score seems reasonable given the relative absence of variability in DAL80’s PFM:

pfm.dal80.scertf
##            1 2 3 4 5          6    7
## A 0.10891089 0 1 0 1 0.90909091 0.03
## C 0.66336634 0 0 0 0 0.01010101 0.19
## G 0.05940594 1 0 0 0 0.01010101 0.75
## T 0.16831683 0 0 1 0 0.07070707 0.03

The GATAA pattern is a very strong signal in this motif.

Note that some restructuring is needed for us to use the results of getPromoterSeqs as an argument to matchPWM. We call the getPromoterSeq method with a GRangesList, which contains a unique set of genomic ranges, representing transcript coordinates, for each of several genes. The corresponding result is a DNAStringSetList in which there is one DNAStringSet (essentially a list of DNAStrings) for each gene in the input list. Both of these variables are therefore lists of lists, in which the outer list is named with genes, and the inner lists are per-transcript coordinates or DNA strings.

Since we need DNA strings without that overarching by-gene-name structure, we call unlist to strip off that structure, leaving us with the desired DNAStringSet:

print (class(promoter.seqs))
## [1] "DNAStringSetList"
## attr(,"package")
## [1] "Biostrings"
promoter.seqs <- unlist(promoter.seqs)
print (class(promoter.seqs))
## [1] "DNAStringSet"
## attr(,"package")
## [1] "Biostrings"
matchPWM(pcm.dal80.scertf, promoter.seqs[[1]], "90%")
## Views on a 1000-letter DNAString subject
## subject: TTGAGGAGTTGTCCACATACACATTAGTGTTGAT...GCAAAAAAAAAGTGAAATACTGCGAAGAACAAAG
## views:
##       start end width
##   [1]   620 626     7 [TGATAAG]
##   [2]   637 643     7 [CGATAAG]

The close proximity of these two GATAA hits suggests that dimeric DAL80, or some other GATAA-binding dimer, may bind DAL1.

All of the matches in the promoters of all seven genes to one binding motif may be found at once with this command:

pwm.hits <- sapply(promoter.seqs, 
                      function(pseq) 
                         matchPWM(pcm.dal80.scertf, pseq, min.score="90%"))

And we can summarize the motif hits for each of the three motifs (dal80.scertf, gat1.jaspar, gat1.scertf) by creating a data.frame of counts, by gene and motif. First, determine the hits:

dal80.scertf.hits <- sapply(promoter.seqs, function(pseq) 
                            matchPWM(pcm.dal80.scertf, pseq, min.score="90%"))
gat1.scertf.hits  <- sapply(promoter.seqs, function(pseq) 
                            matchPWM(pcm.gat1.scertf, pseq, min.score="90%"))
gat1.jaspar.hits  <- sapply(promoter.seqs, function(pseq) 
                            matchPWM(pcm.gat1.jaspar, pseq, min.score="90%"))

Now count their lengths:

dal80.scertf <- sapply(dal80.scertf.hits, length)
gat1.jaspar  <- sapply(gat1.jaspar.hits,  length)
gat1.scertf  <- sapply(gat1.scertf.hits,  length)

Finally, create a data.frame from this information:

tbl.gata     <- data.frame(gene=genes, dal80.scertf, gat1.jaspar, gat1.scertf)

The simple dal80.scertf 5-base motif has the most hits. The more complex 8-base gat1.jaspar mtoif has fewer hits: perhaps it is over-specified. The ‘other’(non-GATAA) motif of GAT1 obtained from ScerTF has fewer matches in the promoters of these genes than do the GATA motifs. The non-GATAA motif hits may in fact, be not much different from chance, as could be revealed by sampling the distribution of motif hits in the promoters of randomly selected genes. Such analyses will be left as an exercise for the reader.

[ Back to top ]

8 Biological Summary

This dataset and our exploration has revealed a number of GATAA binding sites within these tighly co-regulated NCR genes, but leaves unanswered questions, some of which are:

One speculative explanation for the counter-intuitive DAL80 expression is “nuclear sequestration”, a mechanism by which a gene is expressed but the mRNA is held in reserve for later use. See Lavut A, Raveh D 2012.

That GAT1 has multiple binding motifs (we show two, SGD four is yet another indication of the incompletely understood complexity of gene regulation.

The four GATAA-binding regulators, two positive and two negative, and their many downstream targets, some of whose binding sequences we have studied here, can thus be seen to be parts of complex regulatory circuits whose full elucidation has not yet been worked out. Judicious integration of many other kinds of data, careful laboratory work, and the right computational tools, will eventually clarify them.

[ Back to top ]

9 Exploring Package Content

The packages used here have extensive help pages, and include vignettes highlighting common use cases. The help pages and vignettes are available from within an R session. After loading a package, type, for instance:

help(package="MotifDb")
?query

Though it is quite simple, with only a few methods, it will be worthwhile understand the MotifDb package in detail. To access the vignette:

browseVignettes(package="MotifDb")

Finally, you can open a web page containing comprehensive help resources for all installed packages:

help.start()

[ Back to top ]

sessionInfo()
## R version 4.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0
## 
## 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       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] grid      stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] generegulation_1.30.0                      
##  [2] seqLogo_1.72.0                             
##  [3] org.Sc.sgd.db_3.20.0                       
##  [4] motifStack_1.50.0                          
##  [5] TxDb.Scerevisiae.UCSC.sacCer3.sgdGene_3.2.2
##  [6] MotifDb_1.48.0                             
##  [7] GenomicFeatures_1.58.0                     
##  [8] AnnotationDbi_1.68.0                       
##  [9] Biobase_2.66.0                             
## [10] BSgenome.Scerevisiae.UCSC.sacCer3_1.4.0    
## [11] BSgenome_1.74.0                            
## [12] rtracklayer_1.66.0                         
## [13] BiocIO_1.16.0                              
## [14] Biostrings_2.74.0                          
## [15] XVector_0.46.0                             
## [16] GenomicRanges_1.58.0                       
## [17] GenomeInfoDb_1.42.0                        
## [18] IRanges_2.40.0                             
## [19] S4Vectors_0.44.0                           
## [20] BiocGenerics_0.52.0                        
## [21] BiocStyle_2.34.0                           
## 
## loaded via a namespace (and not attached):
##   [1] jsonlite_1.8.9              magrittr_2.0.3             
##   [3] magick_2.8.5                rmarkdown_2.28             
##   [5] zlibbioc_1.52.0             vctrs_0.6.5                
##   [7] memoise_2.0.1               Cairo_1.6-2                
##   [9] Rsamtools_2.22.0            RCurl_1.98-1.16            
##  [11] base64enc_0.1-3             tinytex_0.53               
##  [13] htmltools_0.5.8.1           S4Arrays_1.6.0             
##  [15] curl_5.2.3                  CNEr_1.42.0                
##  [17] SparseArray_1.6.0           sass_0.4.9                 
##  [19] pracma_2.4.4                bslib_0.8.0                
##  [21] htmlwidgets_1.6.4           plyr_1.8.9                 
##  [23] cachem_1.1.0                GenomicAlignments_1.42.0   
##  [25] lifecycle_1.0.4             pkgconfig_2.0.3            
##  [27] Matrix_1.7-1                R6_2.5.1                   
##  [29] fastmap_1.2.0               GenomeInfoDbData_1.2.13    
##  [31] MatrixGenerics_1.18.0       digest_0.6.37              
##  [33] colorspace_2.1-1            TFMPvalue_0.0.9            
##  [35] RSQLite_2.3.7               fansi_1.0.6                
##  [37] httr_1.4.7                  abind_1.4-8                
##  [39] compiler_4.4.1              bit64_4.5.2                
##  [41] BiocParallel_1.40.0         DBI_1.2.3                  
##  [43] highr_0.11                  R.utils_2.12.3             
##  [45] MASS_7.3-61                 poweRlaw_0.80.0            
##  [47] DelayedArray_0.32.0         rjson_0.2.23               
##  [49] gtools_3.9.5                caTools_1.18.3             
##  [51] tools_4.4.1                 splitstackshape_1.4.8      
##  [53] R.oo_1.26.0                 glue_1.8.0                 
##  [55] restfulr_0.0.15             reshape2_1.4.4             
##  [57] ade4_1.7-22                 TFBSTools_1.44.0           
##  [59] generics_0.1.3              gtable_0.3.6               
##  [61] tzdb_0.4.0                  R.methodsS3_1.8.2          
##  [63] data.table_1.16.2           hms_1.1.3                  
##  [65] utf8_1.2.4                  pillar_1.9.0               
##  [67] stringr_1.5.1               dplyr_1.1.4                
##  [69] lattice_0.22-6              bit_4.5.0                  
##  [71] annotate_1.84.0             tidyselect_1.2.1           
##  [73] DirichletMultinomial_1.48.0 GO.db_3.20.0               
##  [75] knitr_1.48                  grImport2_0.3-3            
##  [77] bookdown_0.41               SummarizedExperiment_1.36.0
##  [79] xfun_0.49                   matrixStats_1.4.1          
##  [81] stringi_1.8.4               UCSC.utils_1.2.0           
##  [83] yaml_2.3.10                 evaluate_1.0.1             
##  [85] codetools_0.2-20            tibble_3.2.1               
##  [87] BiocManager_1.30.25         cli_3.6.3                  
##  [89] xtable_1.8-4                munsell_0.5.1              
##  [91] jquerylib_0.1.4             Rcpp_1.0.13                
##  [93] png_0.1-8                   XML_3.99-0.17              
##  [95] parallel_4.4.1              ggplot2_3.5.1              
##  [97] readr_2.1.5                 blob_1.2.4                 
##  [99] jpeg_0.1-10                 bitops_1.0-9               
## [101] pwalign_1.2.0               scales_1.3.0               
## [103] crayon_1.5.3                rlang_1.1.4                
## [105] KEGGREST_1.46.0

[ Back to top ]