## Related to https://stat.ethz.ch/pipermail/bioc-devel/2020-April/016656.html
knitr::opts_chunk$set(
    collapse = TRUE,
    comment = "#>",
    crop = NULL
)

BugSigDB: a comprehensive database of published microbial signatures

BugSigDB is a manually curated database of microbial signatures from the published literature of differential abundance studies of human and other host microbiomes.

BugSigDB provides:

The bugsigdbr package implements convenient access to BugSigDB from within R/Bioconductor. The goal of the package is to facilitate import of BugSigDB data into R/Bioconductor, provide utilities for extracting microbe signatures, and enable export of the extracted signatures to plain text files in standard file formats such as GMT.

The bugsigdbr package is primarily a data package. For descriptive statistics and comprehensive analysis of BugSigDB contents, please see the BugSigDBStats package and analysis vignette.

We start by loading the package.

library(bugsigdbr)

Obtaining published microbial signatures from BugSigDB

The function importBugSigDB can be used to import the complete collection of curated signatures from BugSigDB. The dataset is downloaded once and subsequently cached. Use cache = FALSE to force a fresh download of BugSigDB and overwrite the local copy in your cache.

bsdb <- importBugSigDB()
dim(bsdb)
colnames(bsdb)

Each row of the resulting data.frame corresponds to a microbe signature from differential abundance analysis, i.e. a set of microbes that has been found with increased or decreased abundance in one sample group when compared to another sample group (eg. in a case-vs.-control setup). The curated signatures are richly annotated with additional metadata columns providing information on study design, antibiotics exclusion criteria, sample size, and experimental and statistical procedures, among others.

Subsetting the full dataset to certain conditions, body sites, or other metadata columns of interest can be done along the usual lines for subsetting data.frames.

For example, the following subset command restricts the dataset to signatures obtained from microbiome studies on obesity, based on fecal samples from participants in the US.

us.obesity.feces <- subset(bsdb,
                           `Location of subjects` == "United States of America" &
                           Condition == "obesity" &
                           `Body site` == "feces")

Extracting microbe signatures

Given the full BugSigDB collection (or a subset of interest), the function getSignatures can be used to obtain the microbes annotated to each signature.

Microbes annotated to a signature are returned following the NCBI Taxonomy nomenclature per default.

sigs <- getSignatures(bsdb)
length(sigs)
sigs[1:3]

It is also possible obtain signatures based on the full taxonomic classification in MetaPhlAn format ...

mp.sigs <- getSignatures(bsdb, tax.id.type = "metaphlan")
mp.sigs[1:3]

... or using the taxonomic name only:

tn.sigs <- getSignatures(bsdb, tax.id.type = "taxname")
tn.sigs[1:3]

As metagenomic profiling with 16S RNA sequencing or whole-metagenome shotgun sequencing is typically conducted on a certain taxonomic level, it is also possible to obtain signatures restricted to eg. the genus level ...

gn.sigs <- getSignatures(bsdb, 
                         tax.id.type = "taxname",
                         tax.level = "genus")
gn.sigs[1:3]

... or the species level:

gn.sigs <- getSignatures(bsdb, 
                         tax.id.type = "taxname",
                         tax.level = "species")
gn.sigs[1:3]

Note that restricting signatures to microbes given at the genus level, will per default exclude microbes given at a more specific taxonomic rank such as species or strain.

For certain applications, it might be desirable to not exclude microbes given at a more specific taxonomic rank, but rather extract the more general tax.level for microbes given at a more specific taxonomic level.

This can be achieved by setting the argument exact.tax.level to FALSE, which will here extract genus level taxon names, for taxa given at the species or strain level.

gn.sigs <- getSignatures(bsdb, 
                         tax.id.type = "taxname",
                         tax.level = "genus",
                         exact.tax.level = FALSE)
gn.sigs[1:3]

Writing microbe signatures to file in GMT format

Once signatures have been extracted using a taxonomic identifier type of choice, the function writeGMT allows to write the signatures to plain text files in GMT format.

writeGMT(sigs, gmt.file = "bugsigdb_signatures.gmt")

This is the standard file format for gene sets used by MSigDB and GeneSigDB and is compatible with most enrichment analysis software.

file.remove("bugsigdb_signatures.gmt")

Displaying BugSigDB signature and taxon pages

Leveraging BugSigDB's semantic MediaWiki web interface, we can also programmatically access annotations for individual microbes and microbe signatures.

The browseSignature function can be used to display BugSigDB signature pages in an interactive session. For programmatic access in a non-interactive setting, the URL of the signature page is returned.

browseSignature(names(sigs)[1])

Analogously, the browseTaxon function displays BugSigDB taxon pages in an interactive session, or the URL of the corresponding taxon page otherwise.

browseTaxon(sigs[[1]][1])

Ontology-based queries for experimental factors and body sites

The Semantic MediaWiki curation interface at bugsigdb.org enforces metadata annotation of signatures to follow established ontologies such as the Experimental Factor Ontology (EFO) for condition, and the Uber-Anatomy Ontology (UBERON) for body site.

The getOntology function can be used to import both ontologies into R. The result is an object of class ontology_index from the ontologyIndex package.

efo <- getOntology("efo")
efo
uberon <- getOntology("uberon")
uberon

As demonstrated above, subsets of BugSigDB signatures can be obtained for signatures associated with certain experimental factors or specific body sites of interest. Higher-level queries can be performed with the subsetByOntology function, which implements subsetting by more general ontology terms. This facilitates grouping of signatures that semantically belong together.

More specifically, subsetting BugSigDB signatures by an EFO term then involves subsetting the Condition column to the term itself and all descendants of that term in the EFO ontology and that are present in the Condition column. Here, we demonstrate the usage by subsetting to signatures associated with cancer.

sdf <- subsetByOntology(bsdb,
                        column = "Condition",
                        term = "cancer",
                        ontology = efo)
dim(sdf)
table(sdf[,"Condition"])

And analogously, subsetting by an UBERON term involves subsetting the Body site column to the term itself and all descendants of that term in the UBERON ontology and that are present in the Body site column. For example, we can use subsetByOntology to subset to signatures for which microbiome samples have been obtained from parts of the digestive system.

sdf <- subsetByOntology(bsdb,
                        column = "Body site",
                        term = "digestive system",
                        ontology = uberon)
dim(sdf)
table(sdf[,"Body site"])

Session info

sessionInfo()


waldronlab/bugsigdbr documentation built on May 2, 2024, 1:21 p.m.