Title: | Infer Ligand-Receptor Interactions from bulk expression (transcriptomics/proteomics) data, or spatial transcriptomics |
---|---|
Description: | Inference of ligand-receptor (LR) interactions from bulk expression (transcriptomics/proteomics) data, or spatial transcriptomics. BulkSignalR bases its inferences on the LRdb database included in our other package, SingleCellSignalR available from Bioconductor. It relies on a statistical model that is specific to bulk data sets. Different visualization and data summary functions are proposed to help navigating prediction results. |
Authors: | Jacques Colinge [aut] , Jean-Philippe Villemin [cre] |
Maintainer: | Jean-Philippe Villemin <[email protected]> |
License: | CeCILL | file LICENSE |
Version: | 0.99.22 |
Built: | 2025-01-11 03:00:20 UTC |
Source: | https://github.com/bioc/BulkSignalR |
We note discrepancy between format available over internet.
.formatPathwaysFromGmt(file, resourceName = NULL)
.formatPathwaysFromGmt(file, resourceName = NULL)
file |
Path to GMT file |
resourceName |
Two options "GO-BP" or "REACTOME" |
Here we consider a valid gmt file format defined on each lines as follows : First is Pathway name, Then comes the ID, Finally you will find genes symbols according to the pathway defined on the line.
You can find an example here. - For Reactome. (Directly from their website) https://reactome.org/download/current/ReactomePathways.gmt.zip Note that you need to unzip the file to read the content. The code is inspired from read.gmt function from the gsa R package.
Dataframe with pathwayID, geneName and pathwayName
Format dataframe according to json input
.formatPathwaysFromJson(file, resourceName = NULL)
.formatPathwaysFromJson(file, resourceName = NULL)
file |
Path to file. |
resourceName |
Two options "GO-BP" or "REACTOME". |
Dataframe with pathwayID, geneName and pathwayName
Read dataframe from txt file
.formatPathwaysFromTxt(file, resourceName = NULL)
.formatPathwaysFromTxt(file, resourceName = NULL)
file |
Path to a tabular file. |
resourceName |
Two options "GO-BP" "REACTOME". |
Dataframe with pathwayID, geneName and pathwayName
Add a comparison to a BSRDataModelComp object.
## S4 method for signature 'BSRDataModelComp' addClusterComp(obj, cmp, cmp.name)
## S4 method for signature 'BSRDataModelComp' addClusterComp(obj, cmp, cmp.name)
obj |
A BSRDataModelComp object output by
|
cmp |
A BSRClusterComp object to add. |
cmp.name |
The name of the comparison to add. |
Add cmp
to the list of comparisons contained in
obj
.
A BSRDataModelComp object.
# prepare data data(sdc, package = "BulkSignalR") normal <- grep("^N", names(sdc)) bsrdm <- BSRDataModel(sdc[, -normal]) # define the comparison bsrdm.comp <- as(bsrdm, "BSRDataModelComp") colA <- as.integer(1:3) colB <- as.integer(12:15) n <- nrow(ncounts(bsrdm.comp)) stats <- data.frame( pval = runif(n), logFC = rnorm(n, 0, 2), expr = runif(n, 0, 10) ) rownames(stats) <- rownames(ncounts(bsrdm.comp)) bsrcc <- BSRClusterComp(bsrdm.comp, colA, colB, stats) bsrdm.comp <- addClusterComp(bsrdm.comp, bsrcc, "random.example")
# prepare data data(sdc, package = "BulkSignalR") normal <- grep("^N", names(sdc)) bsrdm <- BSRDataModel(sdc[, -normal]) # define the comparison bsrdm.comp <- as(bsrdm, "BSRDataModelComp") colA <- as.integer(1:3) colB <- as.integer(12:15) n <- nrow(ncounts(bsrdm.comp)) stats <- data.frame( pval = runif(n), logFC = rnorm(n, 0, 2), expr = runif(n, 0, 10) ) rownames(stats) <- rownames(ncounts(bsrdm.comp)) bsrcc <- BSRClusterComp(bsrdm.comp, colA, colB, stats) bsrdm.comp <- addClusterComp(bsrdm.comp, bsrcc, "random.example")
Representation of the links between Ligands,Receptors and Pathways.
alluvialPlot(bsrinf, keywords, type = c("L", "R", "pw.id"), qval.thres = 0.01)
alluvialPlot(bsrinf, keywords, type = c("L", "R", "pw.id"), qval.thres = 0.01)
bsrinf |
object bsrinf inference. |
keywords |
vector of pathways. |
type |
filter on Ligand, Receptor or pathway id. |
qval.thres |
threshold over Q-value. |
NULL
This is a convenience function that relies on the ggalluvial
package to propose a simple way
of representing Ligands, Receptors
data(bsrinf, package = "BulkSignalR") alluvialPlot(bsrinf, keywords = c("LAMC1"), type = "L", qval.thres = 0.01)
data(bsrinf, package = "BulkSignalR") alluvialPlot(bsrinf, keywords = c("LAMC1"), type = "L", qval.thres = 0.01)
Dataframe subset describing the spatial spots
data(annotation.spa)
data(annotation.spa)
Dataframe that contains the following columns : barcode_id,sample_id, in_tissue,array_row array_col,ground_truth,reference,cell_count,idSpatial
barcode_id is the id of the spot idSpatial is the spatial id of the spot(array_rowXarray_col) ground_truth is the label (Layer1/2 were only kept)
They are the mandatory informations in order to make plots for the spatial workflow.
http://spatial.libd.org/spatialLIBD/
Generate a data.frame linking interactions to cell types.
assignCellTypesToInteractions( bsrdm, bsrinf, ct.scores, normalize.scores = TRUE, min.weight = 0.1, min.r2 = 0.25, min.r2.after = 0.35, lasso = TRUE, qval.thres = 0.001 )
assignCellTypesToInteractions( bsrdm, bsrinf, ct.scores, normalize.scores = TRUE, min.weight = 0.1, min.r2 = 0.25, min.r2.after = 0.35, lasso = TRUE, qval.thres = 0.001 )
bsrdm |
A BSRDataModel object. |
bsrinf |
A BSRInference object. |
ct.scores |
A matrix of cell type signature scores. |
normalize.scores |
A logical indicating whether scores should be normalized before assigning cell types. |
min.weight |
Minimum weight to keep in the linear model (cell types
with lower weights will be discarded) if |
min.r2 |
Minimum r2 between a candidate cell type and a L-R gene signature score. |
min.r2.after |
Minimum r2 between the proposed linear model and a L-R gene signature score to retain the model. |
lasso |
Logical indicating that the LASSO (or linear regression if
only one cell type satisfies the |
qval.thres |
Maximum Q-value of the L-R pairs to be considered. |
A data.frame containing the cell type assignments for each
L-R interaction. Unique interactions are considered only (thanks to
"reduceToBestPathway"
that is applied
internally). An interaction can be associated with several cell types
or none. In case it is associated with a single cell type, it is labelled
autocrine (indicative only).
Cell type signature scores must be provided. They can be computed
with BulkSignalR utility function "scoreSignatures"
, but
also any other external tool such as CIBERSORT or BisqueRNA. In case
such a tool would score cell types in a nonlinear fashion, we
recommend to transform the score matrix to restore a linear relationship
cell type abundance/score. By default, cell type (and L-R gene
signature) scores are normalized between 0 and 1 to make the weights of each
cel type in the linear models as comparable as possible.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") data(tme.signatures, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores)
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") data(tme.signatures, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores)
A dataset containing rpkm values of brain and liver.
data(bodyMap.mouse)
data(bodyMap.mouse)
A data frame with 24543 rows and 8 variables.
Bin Li & al., Scientific Reports, 2017;
Define the columns of the expression matrix that belong to each cluster, and store the result of the cluster differences statistical analysis obtained by an external tool such as edgeR, DESeq2, etc.
BSRClusterComp(obj, col.clusterA, col.clusterB, differential.stats)
BSRClusterComp(obj, col.clusterA, col.clusterB, differential.stats)
obj |
A BSRDataModelComp object output by
|
col.clusterA |
Cluster A column indices. |
col.clusterB |
Cluster B column indices. |
differential.stats |
A data.frame containing statistics about
the differential
analysis cluster A versus B. |
Create a BSRClusterComp object describing a comparison
of two clusters of columns taken from the expression matrix
in the BSRDataModelComp object obj
. Such a cluster comparison
description is the basis for inferring LRIs from differential
expression P-values instead of correlation analysis.
The rows of differentialStats
must be in the same order
as those of the count
matrix in obj
. Alternatively, differentialStats
rows can be named and a 1-1 correspondence must exist between these names
and those of the count matrix.
A BSRClusterComp object.
# prepare data data(sdc, package = "BulkSignalR") normal <- grep("^N", names(sdc)) bsrdm <- BSRDataModel(sdc[, -normal]) # define the comparison bsrdm.comp <- as(bsrdm, "BSRDataModelComp") colA <- as.integer(1:3) colB <- as.integer(12:15) n <- nrow(ncounts(bsrdm.comp)) stats <- data.frame( pval = runif(n), logFC = rnorm(n, 0, 2), expr = runif(n, 0, 10) ) rownames(stats) <- rownames(ncounts(bsrdm.comp)) bsrcc <- BSRClusterComp(bsrdm.comp, colA, colB, stats)
# prepare data data(sdc, package = "BulkSignalR") normal <- grep("^N", names(sdc)) bsrdm <- BSRDataModel(sdc[, -normal]) # define the comparison bsrdm.comp <- as(bsrdm, "BSRDataModelComp") colA <- as.integer(1:3) colB <- as.integer(12:15) n <- nrow(ncounts(bsrdm.comp)) stats <- data.frame( pval = runif(n), logFC = rnorm(n, 0, 2), expr = runif(n, 0, 10) ) rownames(stats) <- rownames(ncounts(bsrdm.comp)) bsrcc <- BSRClusterComp(bsrdm.comp, colA, colB, stats)
An S4 class to represent the comparison of two clusters of samples to infer LR interactions based on the resulting P-values, log-fold-changes (logFC), and expression values.
col.clusterA
Column indices for the samples in cluster A.
col.clusterB
Column indices for the samples in cluster B.
differential.stats
Comparison statistics A versus B as a data.frame and containing at least two columns named 'pval', 'logFC', and 'expr'.
new("BSRClusterComp")
new("BSRClusterComp")
Take a matrix or data frame containing RNA sequencing, microarray, or expression proteomics data and return a BSRDataModel object ready for subsequent training. Normally, BSRDataModel objects are not instantiated directly, but through this function.
BSRDataModel( counts, normalize = TRUE, symbol.col = NULL, min.count = 10, prop = 0.1, method = c("UQ", "TC"), log.transformed = FALSE, min.LR.found = 80, species = "hsapiens", conversion.dict = NULL, UQ.pc = 0.75, x.col = NULL, y.col = NULL, barcodeID.col = NULL )
BSRDataModel( counts, normalize = TRUE, symbol.col = NULL, min.count = 10, prop = 0.1, method = c("UQ", "TC"), log.transformed = FALSE, min.LR.found = 80, species = "hsapiens", conversion.dict = NULL, UQ.pc = 0.75, x.col = NULL, y.col = NULL, barcodeID.col = NULL )
counts |
A table or matrix of read counts. |
normalize |
A logical indicating whether |
symbol.col |
The index of the column containing the gene symbols in case
those are not the row names of |
min.count |
The minimum read count of a gene to be considered expressed in a sample. |
prop |
The minimum proportion of samples where a gene must be
expressed higher than |
method |
The normalization method ('UQ' for upper quartile or 'TC'
for total count). If |
log.transformed |
A logical indicating whether expression data were already log2-transformed, e.g., some microarray data. |
min.LR.found |
The minimum number of ligands or receptors found in
|
species |
Data were obtained for this organism. |
conversion.dict |
Correspondence table of HUGO gene symbols human/nonhuman. Not used unless the organism is different from human. |
UQ.pc |
Percentile for upper-quartile normalization, number between 0 and 1 (in case the default 0.75 - hence the name - is not appropriate). |
x.col |
In a SpatialExperiment object, the index of the column containing x coordinates in the dafaframe returned by rowData(), usually named array_row |
y.col |
In a SpatialExperiment object, the index of the column containing y coordinates in the dafaframe returned by rowData(), usually named array_col |
barcodeID.col |
In a SpatialExperiment object, the index of the column containing barcodeID in the dafaframe returned by colData(), usually named barcode_id |
The counts
matrix or table should be provided with expression
levels of protein coding genes in each samples (column) and
rownames(counts)
set to HUGO official gene symbols.
For commodity, it is also possible
to provide counts
with the
gene symbols stored in one of its columns. This column must be specified
with symbol.col
. In such a case, BSRDataModel
will extract
this column and use it to set the row names. Because row names must be
unique, BSRDataModel
will eliminate rows with duplicated gene
symbols by keeping the rows with maximum average expression. Gene symbol
duplication may occur in protein coding genes after genome alignment
due to errors in genome feature annotation files (GTF/GFF), where a handful
of deprecated gene annotations might remain, or
some genes are not given their fully specific symbols. If your read count
extraction pipeline does not take care of this phenomenon, the maximum mean
expression selection strategy implemented here should solve this difficulty
for the sake of inferring ligand-receptor interactions.
If normalize
is TRUE
then normalization is performed
according to method
. If those two simple methods are not satisfying,
then it is possible to provide a pre-normalized matrix setting
normalize
to FALSE
. In such a case, the parameter
method
must be used to document the name of the normalization
algorithm used.
In case proteomic or microarray data are provided, min.count
must be
understood as its equivalent with respect to those data.
A BSRModelData object with empty model parameters.
data(sdc, package = "BulkSignalR") idx <- sample(nrow(sdc), 4000) bsrdm <- BSRDataModel(sdc[idx, c("N22","SDC17")], normalize = FALSE,method="UQ")
data(sdc, package = "BulkSignalR") idx <- sample(nrow(sdc), 4000) bsrdm <- BSRDataModel(sdc[idx, c("N22","SDC17")], normalize = FALSE,method="UQ")
An S4 class to represent the expression data used for inferring ligand-receptor interactions.
ncounts
Normalized read count matrix. Row names must be set to HUGO official gene symbols.
log.transformed
Logical indicating whether values in
ncounts
were log2-transformed.
normalization
Name of the normalization method.
param
List containing the statistical model parameters.
initial.organism
Organism for which the data were obtained.
initial.orthologs
List of genes for which human orthologs exist.
new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" )
new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" )
An S4 class to represent the expression data used for inferring ligand-receptor interactions based on sample cluster comparisons.
comp
A named list of BSRClusterComp objects, one per comparison.
mu
A number representing the average value in the normalized and lop1p-transformed gene expression matrix. This value is used to compute the LR-score (cf. SingleCellSignalR paper, Cabello-Aguilar, et al., Nucleic Acids Res, 2020)
new("BSRDataModelComp")
new("BSRDataModelComp")
Output from the 'learnParameters' function to get BulkSignalR statistical model parameters.
data(bsrdm)
data(bsrdm)
An example of an object created by 'BSRDataModel' applied to an sdc subset (Patients N20,N22,SDC17,SDC25) and 10 000 genes sampled (seed set to 123) 'learnParameters' was also called to get statistical model parameters.
See Vignette BulkSignalR-Differential.
data(bsrdm.comp)
data(bsrdm.comp)
An example of an BSR-dataModelComp object
obtained from STexampleData::Visium_humanDLPFC. A single sample (sample 151673) of human brain dorsolateral prefrontal cortex (DLPFC) in the human brain, measured using the 10x Genomics Visium platform. This is a subset of the full dataset published by Maynard and Collado-Torres et al. (2021). The subset is reproduced in the vignette. name.idx <- c("10x32","3x47","4x50", "17x111","5x59","0x20","8x100", "8x108","14x30","11x39")
data(bsrdm.spa)
data(bsrdm.spa)
An example of an object created by 'BSRDataModel' applied to a subset of a spatial dataset. 'learnParameters' was also called to get statistical model parameters.
Output from the 'learnParameters' function to get BulkSignalR statistical model parameters for a subset of a spatial dataset.
http://spatial.libd.org/spatialLIBD/
From the previous object 'bsrdm', you can generate inferences by calling its method 'BSRInference'. The resulting BSR-Inference object is 'bsrinf', It contains all the inferred L-R interactions with their associated pathways and corrected p-values.
data(bsrinf)
data(bsrinf)
An example of an object created by inference function
See Vignette BulkSignalR-Differential.
data(bsrinf.comp)
data(bsrinf.comp)
An example of an BSR-InferenceComp object
see related workflow for non human organism in the vignette
data(bsrinf.mouse)
data(bsrinf.mouse)
An example of an object created by inference function
Output from the 'learnParameters' function to get BulkSignalR statistical model parameters.
data(bsrinf.spa)
data(bsrinf.spa)
From the previous object 'bsrdm.spa', you can generate inferences by calling its method 'BSRInference'. The resulting BSR-Inference object is 'bsrinf.spa', It contains all the inferred L-R interactions with their associated pathways and corrected p-values. 'learnParameters' was also called to get statistical model parameters.
http://spatial.libd.org/spatialLIBD/
Computes putative LR interactions along with their statistical confidence. In this initial inference, all the relevant pathways are reported, see reduction functions to reduce this list.
BSRInference( obj, rank.p = 0.55, min.cor = 0.25, restrict.genes = NULL, reference = c("REACTOME-GOBP", "REACTOME", "GOBP"), max.pw.size = NULL, min.pw.size = NULL, min.positive = NULL, use.full.network = FALSE, restrict.pw = NULL, with.complex = NULL, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
BSRInference( obj, rank.p = 0.55, min.cor = 0.25, restrict.genes = NULL, reference = c("REACTOME-GOBP", "REACTOME", "GOBP"), max.pw.size = NULL, min.pw.size = NULL, min.positive = NULL, use.full.network = FALSE, restrict.pw = NULL, with.complex = NULL, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
obj |
A BSRDataModel output by method. |
rank.p |
A number between 0 and 1 defining the rank of the last considered target genes. |
min.cor |
The minimum Spearman correlation required between the ligand and the receptor. |
restrict.genes |
A list of gene symbols that restricts ligands and receptors. |
reference |
Which pathway reference should be used ("REACTOME" for Reactome, "GOBP" for GO Biological Process, or "REACTOME-GOBP" for both). |
max.pw.size |
Maximum pathway size to consider from the pathway reference. |
min.pw.size |
Minimum pathway size to consider from the pathway reference. |
min.positive |
Minimum number of target genes to be found in a given pathway. |
use.full.network |
A logical to avoid limiting the reference network to the detected genes and use the whole reference network. |
restrict.pw |
A list of pathway IDs to restrict the application of the function. |
with.complex |
A logical indicating whether receptor co-complex members should be included in the target genes. |
fdr.proc |
The procedure for adjusting P-values according to
|
Perform the initial ligand-receptor inference. Initial means that no reduction is applied. All the (ligand, receptor, downstream pathway) triples are reported, i.e., a given LR pair may appear multiple times with different pathways downstream the receptor. Specific reduction functions are available from the package to operate subsequent simplifications based on the BSRInference object created by the initial inference.
Parameters defining minimum/maximum pathway sizes, etc. are set to NULL
by default, meaning that their values will be taken from what was set
during the training of the statistical model with
"learnParameters"
To use different values at the time of inference sounds like a bad idea, although this could be used to explore without retraining the underlying model. Retraining of the model with adjusted parameters is advised following such an exploration.
A BSRInference object with initial inferences set.
data(bsrdm, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") # We use a subset of the reference to speed up # inference in the context of the example. reactSubset <- getResource(resourceName = "Reactome", cache = FALSE) subset <- c("REACTOME_BASIGIN_INTERACTIONS", "REACTOME_SYNDECAN_INTERACTIONS", "REACTOME_ECM_PROTEOGLYCANS", "REACTOME_CELL_JUNCTION_ORGANIZATION") reactSubset <- reactSubset[ reactSubset$`Reactome name` %in% subset,] resetPathways(dataframe = reactSubset, resourceName = "Reactome") bsrinf <- BSRInference(bsrdm, min.cor = 0.2,restrict.genes=immune.signatures$gene, reference="REACTOME")
data(bsrdm, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") # We use a subset of the reference to speed up # inference in the context of the example. reactSubset <- getResource(resourceName = "Reactome", cache = FALSE) subset <- c("REACTOME_BASIGIN_INTERACTIONS", "REACTOME_SYNDECAN_INTERACTIONS", "REACTOME_ECM_PROTEOGLYCANS", "REACTOME_CELL_JUNCTION_ORGANIZATION") reactSubset <- reactSubset[ reactSubset$`Reactome name` %in% subset,] resetPathways(dataframe = reactSubset, resourceName = "Reactome") bsrinf <- BSRInference(bsrdm, min.cor = 0.2,restrict.genes=immune.signatures$gene, reference="REACTOME")
An S4 class to represent inferred ligand-receptor interactions.
This class contains inferred LR interactions along with
their statistical significance. Data representation supports subsequent
reductions to pathways, etc. See reduction functions
"reduceToBestPathway"
,
"reduceToLigand"
,
"reduceToReceptor"
and
"reduceToPathway"
.
LRinter
A data frame describing the (ligand,receptor,pathway) triples with P- and Q-values.
ligands
A list of ligands, one entry per LR interaction.
receptors
A list of receptors, one entry per LR interaction.
tg.genes
A list of target genes, one entry per LR interaction.
tg.corr
A list of target gene correlations to the receptor, one entry per interaction
inf.param
The parameters used for the inference.
new("BSRInference")
new("BSRInference")
This method supports two configurations that we refer to as paracrine and autocrine.
BSRInferenceComp( obj, cmp.name, src.cmp.name = NULL, rank.p = 0.55, max.pval = 0.01, min.logFC = 1, neg.receptors = FALSE, pos.targets = FALSE, neg.targets = FALSE, min.t.logFC = 0.5, restrict.genes = NULL, use.full.network = FALSE, reference = c("REACTOME-GOBP", "REACTOME", "GOBP"), max.pw.size = 400, min.pw.size = 5, min.positive = 2, restrict.pw = NULL, with.complex = TRUE, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
BSRInferenceComp( obj, cmp.name, src.cmp.name = NULL, rank.p = 0.55, max.pval = 0.01, min.logFC = 1, neg.receptors = FALSE, pos.targets = FALSE, neg.targets = FALSE, min.t.logFC = 0.5, restrict.genes = NULL, use.full.network = FALSE, reference = c("REACTOME-GOBP", "REACTOME", "GOBP"), max.pw.size = 400, min.pw.size = 5, min.positive = 2, restrict.pw = NULL, with.complex = TRUE, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
obj |
A BSRDataModelComp object. |
cmp.name |
The name of the cluster comparison that should be used for the inference. Autocrine interactions if only this comparison name is provided, paracrine if a source comparison name is provided as well. |
src.cmp.name |
The name of the source cluster comparison that should be used for paracrine interaction inferences. |
rank.p |
A number between 0 and 1 defining the rank of the last considered target genes. |
max.pval |
The maximum P-value imposed to both the ligand and the receptor. |
min.logFC |
The minimum log2 fold-change allowed for both the receptor and the ligand. |
neg.receptors |
A logical indicating whether receptors are only allowed to be upregulated (FALSE), or up- and downregulated (TRUE). |
pos.targets |
A logical imposing that all the network targets must display positive logFC, i.e. logFC >= min.t.logFC. |
neg.targets |
A logical imposing that all the network targets must display negative logFC, i.e. logFC <= - min.t.logFC. |
min.t.logFC |
The minimum log2 fold-change allowed for targets in case pos.targets or neg.targets are used. |
restrict.genes |
A list of gene symbols that restricts ligands and receptors. |
use.full.network |
A logical to avoid limiting the reference network to the detected genes and use the whole reference network. |
reference |
Which pathway reference should be used ("REACTOME" for Reactome, "GOBP" for GO Biological Process, or "REACTOME-GOBP" for both). |
max.pw.size |
Maximum pathway size to consider from the pathway reference. |
min.pw.size |
Minimum pathway size to consider from the pathway reference. |
min.positive |
Minimum number of target genes to be found in a given pathway. |
restrict.pw |
A list of pathway IDs to restrict the application of the function. |
with.complex |
A logical indicating whether receptor co-complex members should be included in the target genes. |
fdr.proc |
The procedure for adjusting P-values according to
|
In the autocrine case, a single cluster comparison name is provided. In the corresponding cluster comparison, a group of samples A was compared to a group of samples B to determine fold-changes and associated P-values. The inferred ligand-receptor interactions take place in the samples of group A. They are paracrine interactions in the case of single-cell data or they take place in the same tissue represented by cluster A. A typical single-cell example would be a population of macrophages (group A) compared to all the other populations (group B) to represent specific increased or decreased expression in macrophages. The resulting ligand-receptor interactions will be autocrine interactions that are exacerbated (or reduced depending on the chosen parameters) in macrophages.
In the paracrine case, two cluster comparison names must be provided. For instance, a first comparison coul involved macrophages versus all the other cell populations as above. The second comparison could be B-cells against all the other populations. Now, calling BSRInferenceComp() with comparison macrophages vs. the rest and, as source comparison, B-cells vs. the rest, will result in inferring interactions between B-cells (ligands) and macrophages (receptors and downstream pathways). To obtain macrophages to B-cells paracrine interactions, it is necessary to call the method a second time with permuted cluster comparison names. Another example in spatial transcriptomics could be two thin bands at the boundary of two tissue regions, one emitting the ligand and the other one expressing the receptor.
In this initial inference, all the receptor-containing pathways are reported, see reduction functions to reduce this list.
Perform the initial ligand-receptor inference. Initial means that no reduction is applied. All the (ligand, receptor, downstream pathway) triples are reported, i.e., a given LR pair may appear multiple times with different pathways downstream the receptor. Specific reduction functions are available from the package to operate subsequent simplifications based on the BSRInferenceComp object created by this method.
Here, ligand-receptor interactions are inferred based on gene or protein
regulation-associated P-values when comparing two clusters of samples. Since
a BSRDataModelComp object can contain several such comparisons, the name
of the comparison to use must be specified (parameter cmp.name
).
Note that since the introduction of the use.full.network parameter (April 29, 2024), the pathway sizes are always computed before potential intersection with the observed data (use.full.network set to FALSE) for consistency. Accordingly, the minimum and maximum pathway default values have been raised from 5 & 200 to 5 & 400 respectively. By default, use.full.network is set to FALSE.
In addition to statistical significance estimated according to BulkSignalR statistical model, we compute SingleCellSignalR original LR-score, based on L and R cluster average expression. In the paracrine case, L average expression is taken from the source cluster.
A BSRInferenceComp object with initial inferences set.
data(bsrdm.comp, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") # infer ligand-receptor interactions from the comparison bsrinf.comp <- BSRInferenceComp(bsrdm.comp, max.pval = 1, reference="REACTOME", "random.example")
data(bsrdm.comp, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") # infer ligand-receptor interactions from the comparison bsrinf.comp <- BSRInferenceComp(bsrdm.comp, max.pval = 1, reference="REACTOME", "random.example")
An S4 class to represent ligand-receptor interactions inferred from a comparison between two clusters of samples. This class inherits from BSRInference.
This class is contains inferred LR interactions along with
their statistical significance. Data representation supports subsequent
reductions to pathways, etc. See reduction functions
"reduceToBestPathway"
,
"reduceToLigand"
,
"reduceToReceptor"
and
"reduceToPathway"
.
cmp.name
The name of the BSRClusterComp object in a BSRDataModelComp object comp list.
src.cmp.name
The name of an optional BSRClusterComp object in a BSRDataModelComp object comp list in case paracrine inferences were performed.
tg.pval
A list of target gene P-values, one entry per interaction
tg.logFC
A list of target gene logFC, one entry per interaction
tg.expr
A list of target gene expression, one entry per interaction
new("BSRInferenceComp")
new("BSRInferenceComp")
Obtains gene signatures reflecting ligand-receptor as well as
receptor downstream activity to
score ligand-receptor pairs across samples subsequently with
"scoreLRGeneSignatures"
BSRSignature(obj, pval.thres = NULL, qval.thres = NULL, with.pw.id = FALSE)
BSRSignature(obj, pval.thres = NULL, qval.thres = NULL, with.pw.id = FALSE)
obj |
BSRinference object. |
pval.thres |
P-value threshold. |
qval.thres |
Q-value threshold. |
with.pw.id |
A logical indicating whether the ID of a pathway should be concatenated to its name. |
A BSRSignature object containing a gene signature for each triple
ligand-receptor pair. A reduction to the best pathway
for each pair is automatically performed and the gene signature is
comprised of the ligand, the receptor,
and all the target genes with rank equal or superior to pairs$rank
.
data(bsrinf, package = "BulkSignalR") bsrinf.redP <- reduceToPathway(bsrinf) bsrsig.redP <- BSRSignature(bsrinf, qval.thres = 0.001)
data(bsrinf, package = "BulkSignalR") bsrinf.redP <- reduceToPathway(bsrinf) bsrsig.redP <- BSRSignature(bsrinf, qval.thres = 0.001)
S4 class to represent gene signatures of inferred ligand-receptor interactions, including their reduced versions.
ligands
A list of ligands, one entry per LR interaction.
receptors
A list of receptors, one entry per LR interaction.
tg.genes
A list of target genes, one entry per LR interaction.
pathways
An atomic vector of pathway names, one per interaction.
tg.corr
A list of target genes correlation.
new("BSRSignature")
new("BSRSignature")
Obtains gene signatures reflecting ligand-receptor as well as
receptor downstream activity to
score ligand-receptor pairs across samples subsequently with
"scoreLRGeneSignatures"
BSRSignatureComp(obj, pval.thres = NULL, qval.thres = NULL, with.pw.id = FALSE)
BSRSignatureComp(obj, pval.thres = NULL, qval.thres = NULL, with.pw.id = FALSE)
obj |
BSRInferenceComp object. |
pval.thres |
P-value threshold. |
qval.thres |
Q-value threshold. |
with.pw.id |
A logical indicating whether the ID of a pathway should be concatenated to its name. |
A BSRSignatureComp object containing a gene signature for each triple
ligand-receptor pair. A reduction to the best pathway
for each pair is automatically performed and the gene signature is
comprised of the ligand, the receptor,
and all the target genes with rank equal or superior to pairs$rank
.
data(bsrinf.comp, package = "BulkSignalR") bsrinf.redP <- reduceToPathway(bsrinf.comp) bsrsig.redP <- BSRSignatureComp(bsrinf.redP, qval.thres = 0.001)
data(bsrinf.comp, package = "BulkSignalR") bsrinf.redP <- reduceToPathway(bsrinf.comp) bsrsig.redP <- BSRSignatureComp(bsrinf.redP, qval.thres = 0.001)
S4 class to represent gene signatures associated with ligand-receptor interactions that were inferred from the comparison of two clusters of samples. This class inherits from BSRSignature.
cmp.name
The name of the comparison.
tg.pval
A list of target genes P-values.
tg.logFC
A list of target genes logFC.
tg.expr
A list of target genes expression
new("BSRSignatureComp")
new("BSRSignatureComp")
Quick check to observe LR - Pathways association with their respective correlation and Q-values.
bubblePlotPathwaysLR( bsrinf, pathways, qval.thres = 1, filter.L = NULL, filter.R = NULL, color = "#16a647", pointsize = 6 )
bubblePlotPathwaysLR( bsrinf, pathways, qval.thres = 1, filter.L = NULL, filter.R = NULL, color = "#16a647", pointsize = 6 )
bsrinf |
BulkSignalR inference object. |
pathways |
Vector of pathway names to keep. |
qval.thres |
Maximum Q-value. |
filter.L |
Vector of ligands to keep. |
filter.R |
Vector of receptors to keep. |
color |
Main color used for the gradient. |
pointsize |
Global pointsize. |
A bubble plot displayed in the current viewport or in a file in case a filename was provided.
This is a convenience function to propose a simple way of representing LR - Pathways association with their respective correlation and Q-values.
data(bsrinf, package = "BulkSignalR") pathways <- LRinter(bsrinf)[1,c("pw.name")] bubblePlotPathwaysLR(bsrinf, pathways = pathways, qval.thres = 0.1, color = "red", pointsize = 8 )
data(bsrinf, package = "BulkSignalR") pathways <- LRinter(bsrinf)[1,c("pw.name")] bubblePlotPathwaysLR(bsrinf, pathways = pathways, qval.thres = 0.1, color = "red", pointsize = 8 )
Delete the content of cache directory.
cacheClear(dir = c("both", "resources", "database"))
cacheClear(dir = c("both", "resources", "database"))
dir |
Directory to remove. Can be only 'resources' or 'database'. |
Returns 'NULL', invisibly.
cacheClear(dir="database") # need to recreate database in order to run examples well createDatabase(verbose=TRUE)
cacheClear(dir="database") # need to recreate database in order to run examples well createDatabase(verbose=TRUE)
Get cache content informations for specific cache dir.
cacheInfo(dir = c("both", "resources", "database"))
cacheInfo(dir = c("both", "resources", "database"))
dir |
Directory to remove in order to clean the cache. Can be only 'resources', 'database' or 'both'. |
Returns 'NULL', invisibly.
cacheInfo()
cacheInfo()
Check to see if some ressources has has been updated.
cacheVersion(dir = c("both", "resources", "database"))
cacheVersion(dir = c("both", "resources", "database"))
dir |
Directory for which you want to check Version. Can be only 'resources', 'database' or 'both'. |
Returns 'NULL', invisibly.
cacheVersion()
cacheVersion()
Count how many times and with which weights cell types were involved in the (L,R,pathway) triples that targeted genes in a gene set.
cellTypeFrequency(rel, lr, min.n.genes = 1)
cellTypeFrequency(rel, lr, min.n.genes = 1)
rel |
The data.frame output by
|
lr |
The data.frame output by
|
min.n.genes |
Minimum number of genes in the gene set for one (L,R,pathway) triple. |
A list of two slots: t for counting how many times each cell type is involved; s for summing the weights of each involved cell type.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") data(tme.signatures, package = "BulkSignalR") data(p.EMT, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores) # relate to p-EMT (should be done in HNSCC normally, not in SDC) p.EMT <- p.EMT$gene triggers <- relateToGeneSet(bsrinf, p.EMT) # counts cf <- cellTypeFrequency(triggers, lr2ct)
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") data(tme.signatures, package = "BulkSignalR") data(p.EMT, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores) # relate to p-EMT (should be done in HNSCC normally, not in SDC) p.EMT <- p.EMT$gene triggers <- relateToGeneSet(bsrinf, p.EMT) # counts cf <- cellTypeFrequency(triggers, lr2ct)
Generate a igraph object including all the links between cell types.
cellularNetwork(tab)
cellularNetwork(tab)
tab |
The data.frame output by |
A igraph object containing all the links in the cellular network.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data("tme.signatures", package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores) # cellular network g.table <- cellularNetworkTable(lr2ct) gCN <- cellularNetwork(g.table) #plot(gCN, edge.width=5*E(gCN)$score)
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data("tme.signatures", package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores) # cellular network g.table <- cellularNetworkTable(lr2ct) gCN <- cellularNetwork(g.table) #plot(gCN, edge.width=5*E(gCN)$score)
Generate a data.frame including all the links between cell types mediated by L-R interactions with their respective weights.
cellularNetworkTable(lr, autocrine = FALSE)
cellularNetworkTable(lr, autocrine = FALSE)
lr |
The data.frame output by
|
autocrine |
A logical indicating whether autocrine interactions should be included. |
A data.frame containing all the links in the cellular network. A link is created between two cell types as soon as there was a L-R interaction that was associated with both cell types. The link is given a score equal to the geometric mean of each cell type assignment r2.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") data(tme.signatures, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores) # cellular network g.table <- cellularNetworkTable(lr2ct)
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") data(tme.signatures, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores) # cellular network g.table <- cellularNetworkTable(lr2ct)
Chord diagram.
chordDiagramLR( bsrinf, pw.id.filter = NULL, qval.thres = 1, ligand = NULL, receptor = NULL, limit = 20 )
chordDiagramLR( bsrinf, pw.id.filter = NULL, qval.thres = 1, ligand = NULL, receptor = NULL, limit = 20 )
bsrinf |
bsrinf object |
pw.id.filter |
One Pathway ID accepted only to |
qval.thres |
threshold over Q-value. |
ligand |
Ligand of the LR pair that you want to highlight in the chord diagram. |
receptor |
Receptor of the LR pair that you want to highlight in the chord diagram. |
limit |
Number of interactions you can visualize. |
Circos Plot on the screen or a file
data(bsrinf, package = "BulkSignalR") chordDiagramLR(bsrinf, pw.id.filter = "R-HSA-3000178", limit = 20, ligand="ADAM15", receptor="ITGAV" )
data(bsrinf, package = "BulkSignalR") chordDiagramLR(bsrinf, pw.id.filter = "R-HSA-3000178", limit = 20, ligand="ADAM15", receptor="ITGAV" )
Convert BSRDataModel to BSRDataModelComp
from |
BSRDataModel object |
A BSRDataModelComp object
bsrdm <- new("BSRDataModel") bsrdm.comp <- as(bsrdm, "BSRDataModelComp")
bsrdm <- new("BSRDataModel") bsrdm.comp <- as(bsrdm, "BSRDataModelComp")
Cluster A columns accessor
## S4 method for signature 'BSRClusterComp' colClusterA(x)
## S4 method for signature 'BSRClusterComp' colClusterA(x)
x |
object BSRClusterComp |
col.clusterA
bsrcc <- new("BSRClusterComp") colClusterA(bsrcc)
bsrcc <- new("BSRClusterComp") colClusterA(bsrcc)
Cluster B columns accessor
## S4 method for signature 'BSRClusterComp' colClusterB(x)
## S4 method for signature 'BSRClusterComp' colClusterB(x)
x |
object BSRClusterComp |
col.clusterB
bsrcc <- new("BSRClusterComp") colClusterB(bsrcc)
bsrcc <- new("BSRClusterComp") colClusterB(bsrcc)
Comparisons list accessor
## S4 method for signature 'BSRDataModelComp' comparison(x)
## S4 method for signature 'BSRDataModelComp' comparison(x)
x |
object BSRDataModelComp |
comp
bsrdm.comp <- new("BSRDataModelComp") comparison(bsrdm.comp)
bsrdm.comp <- new("BSRDataModelComp") comparison(bsrdm.comp)
Comparison name accessor
Comparison name accessor
## S4 method for signature 'BSRInferenceComp' comparisonName(x) ## S4 method for signature 'BSRSignatureComp' comparisonName(x)
## S4 method for signature 'BSRInferenceComp' comparisonName(x) ## S4 method for signature 'BSRSignatureComp' comparisonName(x)
x |
BSRSignatureComp object |
cmp.name
cmp.name
bsrinf <- new("BSRInferenceComp") comparisonName(bsrinf)
bsrinf <- new("BSRInferenceComp") comparisonName(bsrinf)
By default, BulkSignalR is designed to work with Homo sapiens. In order to work with other organisms, gene names need to be first converted to human by orthology.
convertToHuman(counts, dictionary)
convertToHuman(counts, dictionary)
counts |
A table or matrix of read counts. |
dictionary |
A data frame where the first column is made of gene symbols for the actual organism and row names are the ortholog human gene symbols. |
Return a counts matrix transposed for Human.
data(bodyMap.mouse) idx <- sample(nrow(bodyMap.mouse), 500) bodyMap.mouse <- bodyMap.mouse[idx,] ortholog.dict <- findOrthoGenes( from_organism = "mmusculus", from_values = rownames(bodyMap.mouse) ) matrix.expression.human <- convertToHuman( counts = bodyMap.mouse, dictionary = ortholog.dict )
data(bodyMap.mouse) idx <- sample(nrow(bodyMap.mouse), 500) bodyMap.mouse <- bodyMap.mouse[idx,] ortholog.dict <- findOrthoGenes( from_organism = "mmusculus", from_values = rownames(bodyMap.mouse) ) matrix.expression.human <- convertToHuman( counts = bodyMap.mouse, dictionary = ortholog.dict )
Fetch LR database from remote location.
createDatabase(onRequest = TRUE, verbose = FALSE)
createDatabase(onRequest = TRUE, verbose = FALSE)
onRequest |
logical True if you force download again. This will overwrite pre-existing database. Default is True. |
verbose |
Logical TRUE/FALSE |
Returns 'NULL', invisibly.
print("Function already called elsewhere by cacheClear()") # createDatabase(onRequest = FALSE)
print("Function already called elsewhere by cacheClear()") # createDatabase(onRequest = FALSE)
Create cache for all resources (pathways, or PWC network) downloaded from the web when library is first loaded. This part is handled with BiocFileCache. Otherwise datatabase, is handled by another process not relying on BiocFileCache instance.
createResources(onRequest = TRUE, verbose = FALSE)
createResources(onRequest = TRUE, verbose = FALSE)
onRequest |
logical True if you force download again. This will overwrite pre-existing database. Default is True. |
verbose |
Default is FALSE |
Returns 'NULL', invisibly.
createResources(onRequest=FALSE)
createResources(onRequest=FALSE)
Cluster comparison statistics accessor
## S4 method for signature 'BSRClusterComp' differentialStats(x)
## S4 method for signature 'BSRClusterComp' differentialStats(x)
x |
BSRClusterComp object |
diffferential.stats
bsrcc <- new("BSRClusterComp") differentialStats(bsrcc)
bsrcc <- new("BSRClusterComp") differentialStats(bsrcc)
By default, BulkSignalR is designed to work with Homo sapiens. In order to work with other organisms, gene names need to be first converted to human following an orthology mapping process.
findOrthoGenes( from_organism, from_values, method = c("gprofiler", "homologene", "babelgene") )
findOrthoGenes( from_organism, from_values, method = c("gprofiler", "homologene", "babelgene") )
from_organism |
An organism defined as in Ensembl: drerio, mmusculus, celegans, dmelanogaster, etc. This is the source organism from which you want to convert the gene names to Homo sapiens. |
from_values |
A vector of gene names from the current species studied. |
method |
Ortholog mapping method. |
Return a data frame with 2 columns containing the gene names for the two species. First column is the gene name from the source organism and the second column corresponds to the homologous gene name in Homo sapiens.
data(bodyMap.mouse) idx <- sample(nrow(bodyMap.mouse), 20) bodyMap.mouse <- bodyMap.mouse[idx,] ortholog.dict <- findOrthoGenes( from_organism = "mmusculus", from_values = rownames(bodyMap.mouse) )
data(bodyMap.mouse) idx <- sample(nrow(bodyMap.mouse), 20) bodyMap.mouse <- bodyMap.mouse[idx,] ortholog.dict <- findOrthoGenes( from_organism = "mmusculus", from_values = rownames(bodyMap.mouse) )
Generate a series of individual spatial score plots in a folder. Not limited to BulkSignalR gene signature scores.
generateSpatialPlots( scores, areas, plot.folder, width = 5, height = 3, pointsize = 8, rev.y = TRUE, ref.plot = TRUE, image.raster = NULL, x.col = "array_col", y.col = "array_row", label.col = "label", idSpatial.col = "idSpatial", cut.p = 0.01, low.color = "royalblue3", mid.color = "white", high.color = "orange", title.fs = 12, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 0.5, ref.colors = NULL )
generateSpatialPlots( scores, areas, plot.folder, width = 5, height = 3, pointsize = 8, rev.y = TRUE, ref.plot = TRUE, image.raster = NULL, x.col = "array_col", y.col = "array_row", label.col = "label", idSpatial.col = "idSpatial", cut.p = 0.01, low.color = "royalblue3", mid.color = "white", high.color = "orange", title.fs = 12, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 0.5, ref.colors = NULL )
scores |
A matrix of scores, one L-R interaction per row and
spatial locations in the columns. This matrix is typically obtained
from BulkSignalR functions |
areas |
A data.frame containing at least the x and y
coordinates of the locations as well as the unique IDs of spatial locations.
In case |
plot.folder |
The folder name in which the plot files will be written. |
width |
The width of each individual plot. |
height |
The height of each individual plot. |
pointsize |
PDF font point size. |
rev.y |
A Boolean indicating whether low y coordinates should be at the top of the plot. |
ref.plot |
A Boolean indicating whether a reference map of the tissue with area labels should be plot aside. |
image.raster |
Raster object image to plot raw tissue image as reference. |
x.col |
Column name in |
y.col |
Column name in |
label.col |
Column name in |
idSpatial.col |
Column name in |
cut.p |
Proportion of top and bottom values for thresholding. |
low.color |
Color for low score values. |
mid.color |
Color for score = 0. |
high.color |
Color for high score values. |
title.fs |
Title font size. |
legend.fs |
Legend items font size. |
axis.fs |
Axis ticks font size. |
label.fs |
Legend titles and axis names font size. |
dot.size |
Dot size. |
ref.colors |
A vector of colors to bypass those automatically chosen by ggplot2 for the tissue areas in the reference plot. |
A set of PDF files are created in the provided folder.
Create PDF file and returns 'NULL', invisibly.
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) generateSpatialPlots(scores.red[1:2,], annotation.spa, ".", label.col = "ground_truth")
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) generateSpatialPlots(scores.red[1:2,], annotation.spa, ".", label.col = "ground_truth")
Fetch LR complexes from database and and return a dataframe
getComplexes(idRelease = NULL)
getComplexes(idRelease = NULL)
idRelease |
integer id version Release Default is NULL so last version is selected. |
Returns dataframe Complexex, invisibly.
getComplexes(idRelease=1)
getComplexes(idRelease=1)
Fetch LR interactions from database and and return a dataframe
getInteractions(idRelease = NULL)
getInteractions(idRelease = NULL)
idRelease |
integer id version Release Default is NULL so last version is selected. |
Returns dataframe LR interactions, invisibly.
getInteractions(idRelease=1)
getInteractions(idRelease=1)
Generate a ligand-receptor network from a BSRInference object and add the shortest paths from the receptors to correlated target genes following Reactome and KEGG pathways.
getLRIntracellNetwork( bsrinf, pval.thres = NULL, qval.thres = NULL, min.cor = 0.25, max.pval = NULL, min.logFC = NULL, pos.targets = FALSE, neg.targets = FALSE, restrict.pw = NULL, node.size = 5 )
getLRIntracellNetwork( bsrinf, pval.thres = NULL, qval.thres = NULL, min.cor = 0.25, max.pval = NULL, min.logFC = NULL, pos.targets = FALSE, neg.targets = FALSE, restrict.pw = NULL, node.size = 5 )
bsrinf |
A BSRInference or BSRInference Comp object. |
pval.thres |
P-value LR interaction threshold. |
qval.thres |
Q-value LR interaction threshold. |
min.cor |
Minimum correlation required for the target genes. |
max.pval |
Maximum P-value required for the target genes in case a BSRInferenceComp object is provided. |
min.logFC |
Minimum logFC required for the target genes in case a BSRInferenceComp object is provided. |
pos.targets |
A logical imposing that all the network targets must display positive correlation or logFC in case of a BSRInferenceComp object. |
neg.targets |
A logical imposing that all the network targets must display negative correlation or logFC in case of a BSRInferenceComp object. Correlations must be <= -min.cor or logFC <= - min.logFC with this option activated. |
restrict.pw |
A vector of pathway IDs to which receptor downstream signaling is restricted. |
node.size |
Default node size in the network. |
An igraph
object featuring the ligand-receptor-downstream
signaling network. Default colors and node sizes are assigned,
which can be changed afterwards if necessary.
The target genes to which the min.cor
correlation is imposed are
those listed in tgGenes(bsrinf)
, correlations are in
tgCorr(bsrinf)
.
The construction of shortest paths from the receptors to those selected
targets adds other genes, which were either some targets with too low
correlation or genes along the shortest paths to reach the selected targets.
data(bsrinf, package = "BulkSignalR") bsrinf.redBP <- reduceToBestPathway(bsrinf) pairs <- LRinter(bsrinf.redBP) top <- unique(pairs[pairs$pval < 1e-20, c("pw.id", "pw.name")]) gLRintra.res <- getLRIntracellNetwork(bsrinf.redBP, qval.thres = 0.01, restrict.pw = top[1,]$pw.id ) # write.graph(gLRintra, file="SDC-LR-intracellular-network.reduced.graphml", # format="graphml")
data(bsrinf, package = "BulkSignalR") bsrinf.redBP <- reduceToBestPathway(bsrinf) pairs <- LRinter(bsrinf.redBP) top <- unique(pairs[pairs$pval < 1e-20, c("pw.id", "pw.name")]) gLRintra.res <- getLRIntracellNetwork(bsrinf.redBP, qval.thres = 0.01, restrict.pw = top[1,]$pw.id ) # write.graph(gLRintra, file="SDC-LR-intracellular-network.reduced.graphml", # format="graphml")
Generate a ligand-receptor network from a ligand-receptor table.
getLRNetwork( bsrinf, pval.thres = NULL, qval.thres = NULL, node.size = 5, red.pairs = NULL )
getLRNetwork( bsrinf, pval.thres = NULL, qval.thres = NULL, node.size = 5, red.pairs = NULL )
bsrinf |
A BSRInference object. |
pval.thres |
P-value threshold. |
qval.thres |
Q-value threshold. |
node.size |
Default node size in the network. |
red.pairs |
A data frame with columns L (ligands) and R (receptors) that restrict LR pairs to those listed. |
An igraph
object featuring the ligand-receptor network.
Default colors and node sizes are assigned,
which can be changed afterwards if necessary.
data(bsrinf, package = "BulkSignalR") gLR <- getLRNetwork(bsrinf, qval.thres = 1e-4) # plot(gLR) # write.graph(gLR, file="SDC-LR-network.graphml", format="graphml")
data(bsrinf, package = "BulkSignalR") gLR <- getLRNetwork(bsrinf, qval.thres = 1e-4) # plot(gLR) # write.graph(gLR, file="SDC-LR-network.graphml", format="graphml")
Basic statistics about hit pathways
## S4 method for signature 'BSRInference' getPathwayStats(obj, pval.thres = NULL, qval.thres = NULL)
## S4 method for signature 'BSRInference' getPathwayStats(obj, pval.thres = NULL, qval.thres = NULL)
obj |
BSRinf object. |
pval.thres |
P-value threshold. |
qval.thres |
Q-value threshold. |
A table with the pathways selected after the chosen threshold was
applied to rows in LRinter(obj)
.
Each pathway is reported along with various statistics:
the number of selected receptors in this pathway, the total number of
receptors described this pathway,
the number of selected ligand-receptor pairs hitting this pathway,
and the total number of
ligand-receptor pairs described that could hit this pathway.
Obviously, one could imagine computing enrichment in receptors or ligand-receptor pairs based on such statistics, but the actual meaning of such an analysis would be ambiguous since the pathways were already selected as significantly regulated by the receptor. We thus did not implement this (hypergeometric test) computation.
data(bsrinf, package = "BulkSignalR") pw.stat <- getPathwayStats(bsrinf)
data(bsrinf, package = "BulkSignalR") pw.stat <- getPathwayStats(bsrinf)
Get resources (pathways, or PathwayCommons network from https://www.pathwaycommons.org/) stored in the cache.
getResource(resourceName = NULL, cache = FALSE)
getResource(resourceName = NULL, cache = FALSE)
resourceName |
Ressource name. |
cache |
True/False. Defautlt is False If True, you will use environment variables. |
Returns a dataframe of the requested resource.
reactome <- getResource(resourceName = "Reactome",cache=TRUE)
reactome <- getResource(resourceName = "Reactome",cache=TRUE)
A dataset containing gene signatures for general immune cell populations.
data(immune.signatures)
data(immune.signatures)
A data frame with 1541 rows and 2 variables:
HUGO gene symbol
cell population name
PanglaoDB (Franzén et al., Database, 2019).
Inference parameters accessor
## S4 method for signature 'BSRInference' inferenceParameters(x)
## S4 method for signature 'BSRInference' inferenceParameters(x)
x |
BRSInference object. |
inf.param
bsrinf <- new ("BSRInference") inferenceParameters(bsrinf)
bsrinf <- new ("BSRInference") inferenceParameters(bsrinf)
organism accessor
## S4 method for signature 'BSRDataModel' initialOrganism(x)
## S4 method for signature 'BSRDataModel' initialOrganism(x)
x |
Object BSRDataModel |
initialOrganism
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) initialOrganism(bsrdm)
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) initialOrganism(bsrdm)
Model parameter accessor
## S4 method for signature 'BSRDataModel' initialOrthologs(x)
## S4 method for signature 'BSRDataModel' initialOrthologs(x)
x |
Object BSRDataModel |
initialOrthologs
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) initialOrthologs(bsrdm)
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) initialOrthologs(bsrdm)
Unique entry point for training the parameters behind BulkSignalR statistical models.
## S4 method for signature 'BSRDataModel' learnParameters( obj, plot.folder = NULL, verbose = FALSE, n.rand.LR = 5L, n.rand.RT = 2L, with.complex = TRUE, max.pw.size = 400, min.pw.size = 5, min.positive = 4, quick = FALSE, null.model = c("automatic", "mixedNormal", "normal", "kernelEmpirical", "empirical", "stable"), filename = NULL, min.corr.LR = -1 )
## S4 method for signature 'BSRDataModel' learnParameters( obj, plot.folder = NULL, verbose = FALSE, n.rand.LR = 5L, n.rand.RT = 2L, with.complex = TRUE, max.pw.size = 400, min.pw.size = 5, min.positive = 4, quick = FALSE, null.model = c("automatic", "mixedNormal", "normal", "kernelEmpirical", "empirical", "stable"), filename = NULL, min.corr.LR = -1 )
obj |
A BSRDatamodel without learned paramaters. |
plot.folder |
A folder name for generating control plots. |
verbose |
A logical activating progress messages for the user. |
n.rand.LR |
The number of random expression matrices to use for learning the ligand-receptor correlation distribution. |
n.rand.RT |
The number of random expression matrices to use for learning the receptor-target genes correlation distribution. |
with.complex |
A logical indicating whether receptor co-complex members should be included in the target genes. |
max.pw.size |
Maximum pathway size to consider from the pathway reference. |
min.pw.size |
Minimum pathway size to consider from the pathway reference. |
min.positive |
Minimum number of target genes to be found in a given pathway. |
quick |
A logical indicating whether approximate parameters for the receptor-target correlations should be used. |
null.model |
The null model to use for Spearman correlation null distributions. |
filename |
Name of the output plot. |
min.corr.LR |
The minimum ligand-receptor correlation required. |
Estimates the model parameters that are stored in the
slot param
.
In a reference pathway, i.e., a Reactome pathway or the genes of a
GOBP term, the target genes are the genes coding for proteins forming a
complex with the receptor and the genes in the pathway downstream the
receptor, which are given as regulated by the pathway. If
with.complex
is set to FALSE
, then only the regulated genes
are considered. Participation to a complex, being regulated, and
pathway directed topologies are defined by Reactome and KEGG pathways
as provided by PathwayCommons.
The min.pw.size
, max.pw.size
, and min.positive
parameters should be identical to the
values intended when searching for ligand-receptor pairs with
.getCorrelatedLR
) and .checkReceptorSignaling
)
Although the statistical distributions are rather robust, it is not
advisable to use different parameters that could introduce unanticipated
biases, but for saving compute time and exploring.
The maximum pathway size is used to limit the redundancy inherent to GOBP and Reactome. The minimum pathway size is used to avoid overspecific, noninformative results.
BulkSignalR approach relies on modeling (Spearman) correlations and
different models of null distributions are available for this purpose
(parameter null.model
). By default, the "automatic" option is
selected meaning that censored normal and mixed normal as well as
an empirical model based on Gaussian kernels (R density()
function)
are compared to pick the one closest to the data. Preference is given
to normal and then mixture of normal over the empirical version for
comparable quality of fit. It is also to bypass the automatic selection.
Fitting of an alpha-stable distribution is quite time consuming as the
computation of its PDF is compute-intensive. Finally, in the automaic
selection mode, the choice of the actual model will be done based on
the L-R null assuming a similar shape for the R-T null (with
different parameters though, unless quick
was set to TRUE
).
Note that since the introduction of the use.full.network parameter (April 29, 2024) in the BSRInference method parameters, the pathway sizes are always computed before potential intersection with the observed data (use.full.network set to FALSE) for consistency. Accordingly, the minimum and maximum pathway default values have been raised from 5 & 200 to 5 & 400 respectively. By default, use.full.network is set to TRUE, meaning no intersection and hence larger pathways.
A BSRDataModel object with trained model parameters
data(sdc, package = "BulkSignalR") idx <- sample(nrow(sdc), 4000) bsrdm <- BSRDataModel(sdc[idx, c("N22","SDC17")],min.LR.found = 20) bsrdm <- learnParameters(bsrdm, n.rand.LR = 1L, verbose=FALSE,quick=TRUE)
data(sdc, package = "BulkSignalR") idx <- sample(nrow(sdc), 4000) bsrdm <- BSRDataModel(sdc[idx, c("N22","SDC17")],min.LR.found = 20) bsrdm <- learnParameters(bsrdm, n.rand.LR = 1L, verbose=FALSE,quick=TRUE)
ligands accessor
ligands accessor
## S4 method for signature 'BSRInference' ligands(x) ## S4 method for signature 'BSRSignature' ligands(x)
## S4 method for signature 'BSRInference' ligands(x) ## S4 method for signature 'BSRSignature' ligands(x)
x |
BSRSignature |
ligands
ligands
bsr.sig <- new("BSRSignature") ligands(bsr.sig)
bsr.sig <- new("BSRSignature") ligands(bsr.sig)
log.transformed accessor
## S4 method for signature 'BSRDataModel' logTransformed(x)
## S4 method for signature 'BSRDataModel' logTransformed(x)
x |
Object BRSDataModel |
logTransformed
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) logTransformed(bsrdm)
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) logTransformed(bsrdm)
LRinter accessor
## S4 method for signature 'BSRInference' LRinter(x)
## S4 method for signature 'BSRInference' LRinter(x)
x |
BSRInference object |
LRinter
bsrinf <- new ("BSRInference") LRinter(bsrinf)
bsrinf <- new ("BSRInference") LRinter(bsrinf)
Simplified LRinter accessor with focus on the LR-score
## S4 method for signature 'BSRInferenceComp' LRinterScore(x)
## S4 method for signature 'BSRInferenceComp' LRinterScore(x)
x |
BSRInferenceComp object |
LRinterScore
data(bsrinf.comp, package = "BulkSignalR") LRinterScore(bsrinf.comp)[5,]
data(bsrinf.comp, package = "BulkSignalR") LRinterScore(bsrinf.comp)[5,]
Simplified LRinter accessor reporting the essential columns
Simplified LRinter accessor reporting the essential columns
## S4 method for signature 'BSRInference' LRinterShort(x) ## S4 method for signature 'BSRInferenceComp' LRinterShort(x)
## S4 method for signature 'BSRInference' LRinterShort(x) ## S4 method for signature 'BSRInferenceComp' LRinterShort(x)
x |
BSRInferenceComp object |
LRinterShort
LRinterShort
data(bsrinf.comp, package = "BulkSignalR") LRinterShort(bsrinf.comp)[5,]
data(bsrinf.comp, package = "BulkSignalR") LRinterShort(bsrinf.comp)[5,]
Get maximal ligand expression at nearby locations
maxLigandSpatialCounts( bsrdm, areas, nnn = 4, radius = NULL, x.col = "array_col", y.col = "array_row" )
maxLigandSpatialCounts( bsrdm, areas, nnn = 4, radius = NULL, x.col = "array_col", y.col = "array_row" )
bsrdm |
A BSRDataModel object containing the expression data to smooth. |
areas |
A data.frame containing at least the x and y coordinates of the locations. |
nnn |
Number of nearest-neighbor locations to use for smoothing
each location. In case |
radius |
A maximal distance to include neighbors in the smoothing. |
x.col |
Column name in |
y.col |
Column name in |
Ligand expression data contained in a BSRDataModel object are modified to consider the possibility that the ligand of a L-R interaction might be expressed at nearby locations. This is achieved replacing each ligand expression by its maximum over the central location and its neighbors. Since ligands and receptors are never used as gene targets in computing the receptor downstream signal correlations, this substitution is compatible with our statistical model. Moreover, the reciprocal configuration where the ligand is expressed at the central location and hits a receptors at a neighbor location is covered when the same ligand maximization scheme is applied to the neighbor. L-R localization and gene signature scoring is defined by the location at which the receptor is expressed after applying this function.
Two strategies are available to identify the neighbors. It is possible to
simply set the number of nearest-neighbors (parameter nnn
). An
alternative consists in providing a distance radius (radius
) along
with a a maximum number of nearest-neighbors within the radius
(nnn.radius
). To properly define the radius, the user must know the
location coordinates. The strategy with the radius enables having corner
locations with two neighbors only and border locations with three
neighbors only, whereas to simply set a maximum of four neighbors for
instance would retrieve the four closest neighbors in every case.
A BSRDataModel object containing the maximized ligand expressions.
data(bsrdm.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") max.bsrdm <- maxLigandSpatialCounts(bsrdm.spa, annotation.spa, radius = 1.2, nnn = 4)
data(bsrdm.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") max.bsrdm <- maxLigandSpatialCounts(bsrdm.spa, annotation.spa, radius = 1.2, nnn = 4)
Mu accessor
## S4 method for signature 'BSRDataModelComp' mu(x)
## S4 method for signature 'BSRDataModelComp' mu(x)
x |
object BSRDataModelComp |
mu
bsrdm.comp <- new("BSRDataModelComp") mu(bsrdm.comp)
bsrdm.comp <- new("BSRDataModelComp") mu(bsrdm.comp)
Normalized count matrix accessor
## S4 method for signature 'BSRDataModel' ncounts(x)
## S4 method for signature 'BSRDataModel' ncounts(x)
x |
object BSRDataModel |
ncounts
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) ncounts(bsrdm)
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) ncounts(bsrdm)
Normalization accessor
## S4 method for signature 'BSRDataModel' normalization(x)
## S4 method for signature 'BSRDataModel' normalization(x)
x |
object BSRDatamModel |
normalization
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) normalization(bsrdm)
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) normalization(bsrdm)
Synthetic object used during the call to the function 'resetToInitialOrganism“
data(ortholog.dict)
data(ortholog.dict)
An example of a dataframe created by findOrthoGenes
A dataset containing a partial EMT gene signature.
data(p.EMT)
data(p.EMT)
A data frame with 100 rows and 1 variables:
HUGO gene symbol
Puram, SV & al., Cell, 2017.
Model parameter accessor
## S4 method for signature 'BSRDataModel' parameters(x)
## S4 method for signature 'BSRDataModel' parameters(x)
x |
BSRDataModel oject |
param
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) parameters(bsrdm)
bsrdm <- new("BSRDataModel", ncounts = matrix(1.5, nrow = 2, ncol = 2, dimnames = list(c("A", "B"), c("C", "D")) ), log.transformed = TRUE, normalization = "TC" ) parameters(bsrdm)
pathways accessor
## S4 method for signature 'BSRSignature' pathways(x)
## S4 method for signature 'BSRSignature' pathways(x)
x |
BSRSignature |
pathways
bsr.sig <- new("BSRSignature") pathways(bsr.sig)
bsr.sig <- new("BSRSignature") pathways(bsr.sig)
receptors accessor
receptors accessor
## S4 method for signature 'BSRInference' receptors(x) ## S4 method for signature 'BSRSignature' receptors(x)
## S4 method for signature 'BSRInference' receptors(x) ## S4 method for signature 'BSRSignature' receptors(x)
x |
BSRSignature |
receptors
receptors
bsr.sig <- new("BSRSignature") ligands(bsr.sig)
bsr.sig <- new("BSRSignature") ligands(bsr.sig)
Keep one pathway per ligand-receptor pair
Keep one pathway per ligand-receptor pair
## S4 method for signature 'BSRInference' reduceToBestPathway(obj) ## S4 method for signature 'BSRInferenceComp' reduceToBestPathway(obj)
## S4 method for signature 'BSRInference' reduceToBestPathway(obj) ## S4 method for signature 'BSRInferenceComp' reduceToBestPathway(obj)
obj |
BSRInferenceComp object |
Ligand-receptor pairs are evaluated in relation with pathways that allow checking receptor downstream correlations. It is thus possible that several pathways are reported for a same LR pair.
Ligand-receptor pairs are evaluated in relation with pathways that allow checking receptor downstream correlations. It is thus possible that several pathways are reported for a same LR pair.
A BSRInference object reduced to only report one pathway per ligand-receptor pair. The pathway with the smallest P-value is selected.
A BSRInferenceComp object reduced to only report one pathway per ligand-receptor pair. The pathway with the smallest P-value is selected.
data(bsrinf, package = "BulkSignalR") bsrinf.redBP <- reduceToBestPathway(bsrinf) data(bsrinf.comp, package = "BulkSignalR") reduceToBestPathway(bsrinf.comp)
data(bsrinf, package = "BulkSignalR") bsrinf.redBP <- reduceToBestPathway(bsrinf) data(bsrinf.comp, package = "BulkSignalR") reduceToBestPathway(bsrinf.comp)
Simplifies a ligand-receptor table to focus on the ligands.
Simplifies a ligand-receptor table to focus on the ligands.
## S4 method for signature 'BSRInference' reduceToLigand(obj) ## S4 method for signature 'BSRInferenceComp' reduceToLigand(obj)
## S4 method for signature 'BSRInference' reduceToLigand(obj) ## S4 method for signature 'BSRInferenceComp' reduceToLigand(obj)
obj |
BSRInferenceComp object |
A BSRInference object but reduced to one row per ligand.
All the receptors are combined in a
semi-colon-separated list surrounded by curly brackets in the tabular
slot LRinter
, and in vectors in the ligands
(list) slot.
The reported P-value and target genes are those from the pathway with the smallest P-value.
A BSRInferenceComp object but reduced to one row per ligand.
All the receptors are combined in a
semi-colon-separated list surrounded by curly brackets in the tabular
slot LRinter
, and in vectors in the ligands
(list) slot.
The reported P-value and target genes are those from the pathway with the smallest P-value. The same logic applies to the LR-score, and the receptor expression.
data(bsrinf, package = "BulkSignalR") bsrinf.redL <- reduceToLigand(bsrinf) data(bsrinf.comp, package = "BulkSignalR") bsrinf.redL <- reduceToLigand(bsrinf.comp)
data(bsrinf, package = "BulkSignalR") bsrinf.redL <- reduceToLigand(bsrinf) data(bsrinf.comp, package = "BulkSignalR") bsrinf.redL <- reduceToLigand(bsrinf.comp)
Simplifies a ligand-receptor inference object to focus on the pathways.
Simplifies a ligand-receptor inference object to focus on the pathways.
## S4 method for signature 'BSRInference' reduceToPathway(obj) ## S4 method for signature 'BSRInferenceComp' reduceToPathway(obj)
## S4 method for signature 'BSRInference' reduceToPathway(obj) ## S4 method for signature 'BSRInferenceComp' reduceToPathway(obj)
obj |
BSRInferenceComp object |
A BSRInference object reduced to only report one row per pathway.
The information of which ligand interacted with which receptor is lost as
all the ligands and all the receptors forming pairs related to a certain
pathway are combined.
For a given pathway, the reported P-values and target genes are those of
the best ligand-receptor pair that
was in this pathway.
Receptors and ligands are combined in two semi-colon-separated
lists surrounded by curly brackets in the tabular slot LRinter
,
while the list representation slots (ligands
and
receptors
) are update accordingly.
A BSRInferenceComp object reduced to only report one row per pathway.
The information of which ligand interacted with which receptor is lost as
all the ligands and all the receptors forming pairs related to a certain
pathway are combined.
For a given pathway, the reported P-values and target genes are those of
the best ligand-receptor pair that was in this pathway.
The same logic applies to the LR-score, and the ligand and receptor
expression.
Receptors and ligands are combined in two semi-colon-separated
lists surrounded by curly brackets in the tabular slot LRinter
,
while the list representation slots (ligands
and
receptors
) are update accordingly.
data(bsrinf, package = "BulkSignalR") bsrinf.redP <- reduceToPathway(bsrinf) data(bsrinf.comp, package = "BulkSignalR") bsrinf.redP <- reduceToPathway(bsrinf.comp)
data(bsrinf, package = "BulkSignalR") bsrinf.redP <- reduceToPathway(bsrinf) data(bsrinf.comp, package = "BulkSignalR") bsrinf.redP <- reduceToPathway(bsrinf.comp)
Simplifies a ligand-receptor table to focus on the receptors.
Simplifies a ligand-receptor table to focus on the receptors.
## S4 method for signature 'BSRInference' reduceToReceptor(obj) ## S4 method for signature 'BSRInferenceComp' reduceToReceptor(obj)
## S4 method for signature 'BSRInference' reduceToReceptor(obj) ## S4 method for signature 'BSRInferenceComp' reduceToReceptor(obj)
obj |
BRSInferenceComp object |
BSRInference object reduced to one row per receptor.
All the ligands are combined in a
semi-colon-separated list surrounded by curly brackets in the tabular
slot LRinter
, and in vectors in the ligands
(list) slot.
The reported P-value and target genes are those from the line with the pathway featuring the smallest P-value.
BSRInferenceComp object reduced to one row per receptor.
All the ligands are combined in a
semi-colon-separated list surrounded by curly brackets in the tabular
slot LRinter
, and in vectors in the ligands
(list) slot.
The reported P-value and target genes are those from the line with the pathway featuring the smallest P-value. The same logic applies to the LR-score, and the ligand expression.
data(bsrinf, package = "BulkSignalR") bsrinf.redR <- reduceToReceptor(bsrinf) data(bsrinf.comp, package = "BulkSignalR") # reduction bsrinf.redR <- reduceToReceptor(bsrinf.comp)
data(bsrinf, package = "BulkSignalR") bsrinf.redR <- reduceToReceptor(bsrinf) data(bsrinf.comp, package = "BulkSignalR") # reduction bsrinf.redR <- reduceToReceptor(bsrinf.comp)
Finds ligands related to a gene set by following receptor, and receptor downstream pathway targets.
relateToGeneSet(bsrinf, gs, min.cor = 0.25, qval.thres = 0.001)
relateToGeneSet(bsrinf, gs, min.cor = 0.25, qval.thres = 0.001)
bsrinf |
BSRInference object. |
gs |
The gene set. |
min.cor |
Minimum Spearman correlation between the receptor of a triple (L,R,pw) and a gene of the gene set. |
qval.thres |
Maximum Q-value imposed to the (L,R,pw) triples to be considered. |
A data.frame listing all the (L,R,pathway) triples that lead
to at least one gene in the gene set. The number of genes found
by each triple is indicated in the column n.genes
.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data(p.EMT, package = "BulkSignalR") p.EMT <- p.EMT$gene triggers <- relateToGeneSet(bsrinf, p.EMT)
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data(p.EMT, package = "BulkSignalR") p.EMT <- p.EMT$gene triggers <- relateToGeneSet(bsrinf, p.EMT)
Remove a comparison from a BSRDataModelComp object.
## S4 method for signature 'BSRDataModelComp' removeClusterComp(obj, cmp.name)
## S4 method for signature 'BSRDataModelComp' removeClusterComp(obj, cmp.name)
obj |
A BSRDataModelComp object output by
|
cmp.name |
The name of the comparison to remove. |
Remove the comparison with cmp.name
from the list of
comparisons contained in obj
.
A BSRDataModelComp object.
# prepare data data(sdc, package = "BulkSignalR") normal <- grep("^N", names(sdc)) bsrdm <- BSRDataModel(sdc[, -normal]) # define the comparison bsrdm.comp <- as(bsrdm, "BSRDataModelComp") colA <- as.integer(1:3) colB <- as.integer(12:15) n <- nrow(ncounts(bsrdm.comp)) stats <- data.frame( pval = runif(n), logFC = rnorm(n, 0, 2), expr = runif(n, 0, 10) ) rownames(stats) <- rownames(ncounts(bsrdm.comp)) bsrcc <- BSRClusterComp(bsrdm.comp, colA, colB, stats) bsrdm.comp <- addClusterComp(bsrdm.comp, bsrcc, "random.example") bsrdm.comp <- removeClusterComp(bsrdm.comp, "random.example")
# prepare data data(sdc, package = "BulkSignalR") normal <- grep("^N", names(sdc)) bsrdm <- BSRDataModel(sdc[, -normal]) # define the comparison bsrdm.comp <- as(bsrdm, "BSRDataModelComp") colA <- as.integer(1:3) colB <- as.integer(12:15) n <- nrow(ncounts(bsrdm.comp)) stats <- data.frame( pval = runif(n), logFC = rnorm(n, 0, 2), expr = runif(n, 0, 10) ) rownames(stats) <- rownames(ncounts(bsrdm.comp)) bsrcc <- BSRClusterComp(bsrdm.comp, colA, colB, stats) bsrdm.comp <- addClusterComp(bsrdm.comp, bsrcc, "random.example") bsrdm.comp <- removeClusterComp(bsrdm.comp, "random.example")
A method to re-score an existing BSRInference object (P- and Q-value estimations).
A method to re-score an existing BSRInferenceComp object (P- and Q-value estimations).
## S4 method for signature 'BSRInference' rescoreInference( obj, param, rank.p = 0.55, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") ) ## S4 method for signature 'BSRInferenceComp' rescoreInference( obj, param = NULL, rank.p = 0.55, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
## S4 method for signature 'BSRInference' rescoreInference( obj, param, rank.p = 0.55, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") ) ## S4 method for signature 'BSRInferenceComp' rescoreInference( obj, param = NULL, rank.p = 0.55, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
obj |
BSRInferecenceComp object. |
param |
NULL by default |
rank.p |
A number between 0 and 1 defining the rank of the last considered target genes. |
fdr.proc |
The procedure for adjusting P-values according to
|
A BSRInference object should be created by calling
"BSRInference"
Parameters controlling the estimation
of the statistical significance of the ligand/receptor/pathway triples
(param
) are provided at the time of calling the latter method.
Nonetheless, it might be useful to change the initially-provided parameters, in which case this method should not be called.
A BSRInferenceComp object should be created by calling
"BSRInferenceComp"
A BSRInference object.
A BSRInferenceComp object.
data(bsrinf, package = "BulkSignalR") data(bsrdm, package = "BulkSignalR") bsrinf.new <- rescoreInference(bsrinf, param = parameters(bsrdm)) data(bsrinf.comp, package = "BulkSignalR") bsrinf.less <- rescoreInference(bsrinf.comp, rank.p = 0.75)
data(bsrinf, package = "BulkSignalR") data(bsrdm, package = "BulkSignalR") bsrinf.new <- rescoreInference(bsrinf, param = parameters(bsrdm)) data(bsrinf.comp, package = "BulkSignalR") bsrinf.less <- rescoreInference(bsrinf.comp, rank.p = 0.75)
User can provide a data frame with 2 columns named ligand and receptor. This can be used to extend or replace the existing LRdb.
resetLRdb(db, switch = FALSE)
resetLRdb(db, switch = FALSE)
db |
A dataframe with 2 columns named ligand and receptor. |
switch |
A logical indicating whether LRdb should be extended only (FALSE, default) or completely replaced (TRUE). |
Returns 'NULL', invisibly.
resetLRdb(db = data.frame(ligand = "A2M", receptor = "LRP1"), switch = FALSE)
resetLRdb(db = data.frame(ligand = "A2M", receptor = "LRP1"), switch = FALSE)
Network is a dataframe that gives relation between genes. It's composed of 3 columns annoted as follows :
resetNetwork(network)
resetNetwork(network)
network |
Network dataframe is defined with 3 columns a.gn, b.gn & type. 'a.gn' & 'b.gn' should be gene symbols of gene interactions. 'type' should be set as 'controls-expression-of' when user provide his own file. |
a.gn : Gene Symbol 1 type : controls-expression-of b.gn : Gene Symbol 2
When the user provide his own network 'type' should be set to 'controls-expression-of'.
Returns 'NULL', invisibly.
BulkSignalR_Network <- getResource(resourceName = "Network", cache = FALSE) resetNetwork(BulkSignalR_Network)
BulkSignalR_Network <- getResource(resourceName = "Network", cache = FALSE) resetNetwork(BulkSignalR_Network)
resetPathways
is a function
we provide to user to refresh REACTOME
and GO-BP content included in BulkSignalR.
resetPathways( dataframe = NULL, file = NULL, fileType = c("json", "gmt", "txt"), resourceName = NULL )
resetPathways( dataframe = NULL, file = NULL, fileType = c("json", "gmt", "txt"), resourceName = NULL )
dataframe |
Dataframe formated as
When |
file |
Path to file. |
fileType |
Default is Json. Other options are gmt or txt files. |
resourceName |
Two options "GO-BP" or "Reactome". |
Pathways are defined in Reactome and GoBP databases. Those can be updated using json files from the Human Molecular Signatures Database (MSigDB) at https://www.gsea-msigdb.org/ Gmt file format also can be imported. A dataframe can be used directly also.
Returns 'NULL', invisibly.
reactSubset <- getResource(resourceName = "Reactome", cache = TRUE) subset <- c("REACTOME_BASIGIN_INTERACTIONS", "REACTOME_SYNDECAN_INTERACTIONS", "REACTOME_ECM_PROTEOGLYCANS", "REACTOME_CELL_JUNCTION_ORGANIZATION") reactSubset <- reactSubset[ reactSubset$`Reactome name` %in% subset,] resetPathways(dataframe = reactSubset, resourceName = "Reactome")
reactSubset <- getResource(resourceName = "Reactome", cache = TRUE) subset <- c("REACTOME_BASIGIN_INTERACTIONS", "REACTOME_SYNDECAN_INTERACTIONS", "REACTOME_ECM_PROTEOGLYCANS", "REACTOME_CELL_JUNCTION_ORGANIZATION") reactSubset <- reactSubset[ reactSubset$`Reactome name` %in% subset,] resetPathways(dataframe = reactSubset, resourceName = "Reactome")
Reset gene names to initial organism providen in first instance
## S4 method for signature 'BSRInference' resetToInitialOrganism(obj, conversion.dict)
## S4 method for signature 'BSRInference' resetToInitialOrganism(obj, conversion.dict)
obj |
BSRInference object |
conversion.dict |
A dictionnary |
An BSRInference object updated for gene names. The gene names are replaced by the ones from the organism providen in first instance.
data(bodyMap.mouse, package = "BulkSignalR") data(bsrinf.mouse, package = "BulkSignalR") data(ortholog.dict, package = "BulkSignalR") #idx <- sample(nrow(bodyMap.mouse), 7500) #bodyMap.mouse <- bodyMap.mouse[idx,1:3] #ortholog.dict <- findOrthoGenes( # from_organism = "mmusculus", # from_values = rownames(bodyMap.mouse) #) #matrix.expression.human <- convertToHuman( # counts = bodyMap.mouse, # dictionary = ortholog.dict #) #bsrdm <- BSRDataModel( # counts = matrix.expression.human, # species = "mmusculus", # conversion.dict = ortholog.dict #) #bsrdm <- learnParameters(bsrdm, # quick = TRUE #) #reactSubset <- getResource(resourceName = "Reactome", #cache = TRUE) #subset <- c("REACTOME_BASIGIN_INTERACTIONS", #"REACTOME_SYNDECAN_INTERACTIONS", #"REACTOME_ECM_PROTEOGLYCANS", #"REACTOME_CELL_JUNCTION_ORGANIZATION") #reactSubset <- reactSubset[ #reactSubset$`Reactome name` %in% subset,] #bsrinf.mouse <- BSRInference(bsrdm,reference="REACTOME") bsrinf <- resetToInitialOrganism(bsrinf.mouse, conversion.dict = ortholog.dict)
data(bodyMap.mouse, package = "BulkSignalR") data(bsrinf.mouse, package = "BulkSignalR") data(ortholog.dict, package = "BulkSignalR") #idx <- sample(nrow(bodyMap.mouse), 7500) #bodyMap.mouse <- bodyMap.mouse[idx,1:3] #ortholog.dict <- findOrthoGenes( # from_organism = "mmusculus", # from_values = rownames(bodyMap.mouse) #) #matrix.expression.human <- convertToHuman( # counts = bodyMap.mouse, # dictionary = ortholog.dict #) #bsrdm <- BSRDataModel( # counts = matrix.expression.human, # species = "mmusculus", # conversion.dict = ortholog.dict #) #bsrdm <- learnParameters(bsrdm, # quick = TRUE #) #reactSubset <- getResource(resourceName = "Reactome", #cache = TRUE) #subset <- c("REACTOME_BASIGIN_INTERACTIONS", #"REACTOME_SYNDECAN_INTERACTIONS", #"REACTOME_ECM_PROTEOGLYCANS", #"REACTOME_CELL_JUNCTION_ORGANIZATION") #reactSubset <- reactSubset[ #reactSubset$`Reactome name` %in% subset,] #bsrinf.mouse <- BSRInference(bsrdm,reference="REACTOME") bsrinf <- resetToInitialOrganism(bsrinf.mouse, conversion.dict = ortholog.dict)
Compute ligand-receptor gene signature scores over a BSRDataModel.
Compute ligand-receptor gene signature scores over a BSRDataModelComp specific comparison.
## S4 method for signature 'BSRDataModel' scoreLRGeneSignatures( obj, sig, LR.weight = 0.5, robust = FALSE, name.by.pathway = FALSE, abs.z.score = FALSE, rownames.LRP = FALSE ) ## S4 method for signature 'BSRDataModelComp' scoreLRGeneSignatures( obj, sig, LR.weight = 0.5, robust = FALSE, name.by.pathway = FALSE, abs.z.score = FALSE, rownames.LRP = FALSE )
## S4 method for signature 'BSRDataModel' scoreLRGeneSignatures( obj, sig, LR.weight = 0.5, robust = FALSE, name.by.pathway = FALSE, abs.z.score = FALSE, rownames.LRP = FALSE ) ## S4 method for signature 'BSRDataModelComp' scoreLRGeneSignatures( obj, sig, LR.weight = 0.5, robust = FALSE, name.by.pathway = FALSE, abs.z.score = FALSE, rownames.LRP = FALSE )
obj |
A BSRDataModelComp object. |
sig |
A BSRSignatureComp object. |
LR.weight |
A number between 0 and 1 defining the relative weight of the ligand and the receptor in the signature. |
robust |
A logical indicating that z-scores should be computed with median and MAD instead of mean and standard deviation. |
name.by.pathway |
A logical indicating whether row names of the resulting score matrix should be pathway names. |
abs.z.score |
A logical to use absolute z-scores (useful if the activity of a paythway is reported by a mixture of up- and down-genes whose z-score averages might hide actual activity). |
rownames.LRP |
A logical indicating, in case |
A matrix containing the scores of each ligand-receptor gene signature in each sample.
A matrix containing the scores of each ligand-receptor gene signature in each sample.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") bsrinf.redBP <- reduceToBestPathway(bsrinf) bsrsig.redBP <- BSRSignature(bsrinf.redBP, qval.thres = 0.001) res <-scoreLRGeneSignatures(bsrdm, bsrsig.redBP, name.by.pathway = FALSE ) # prepare data data(bsrdm.comp, package = "BulkSignalR") data(bsrinf.comp, package = "BulkSignalR") # reduction bsrinf.red <- reduceToBestPathway(bsrinf.comp) # signature extraction and scoring bsrsig.red <- BSRSignatureComp(bsrinf.red, qval.thres = 1e-6) scores.red <- scoreLRGeneSignatures(bsrdm.comp, bsrsig.red, name.by.pathway = TRUE, rownames.LRP = TRUE )
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") bsrinf.redBP <- reduceToBestPathway(bsrinf) bsrsig.redBP <- BSRSignature(bsrinf.redBP, qval.thres = 0.001) res <-scoreLRGeneSignatures(bsrdm, bsrsig.redBP, name.by.pathway = FALSE ) # prepare data data(bsrdm.comp, package = "BulkSignalR") data(bsrinf.comp, package = "BulkSignalR") # reduction bsrinf.red <- reduceToBestPathway(bsrinf.comp) # signature extraction and scoring bsrsig.red <- BSRSignatureComp(bsrinf.red, qval.thres = 1e-6) scores.red <- scoreLRGeneSignatures(bsrdm.comp, bsrsig.red, name.by.pathway = TRUE, rownames.LRP = TRUE )
Scores generic gene signatures over the samples of a BSRDataModel object.
scoreSignatures(ds, ref.signatures, robust = FALSE)
scoreSignatures(ds, ref.signatures, robust = FALSE)
ds |
A BSRDataModel object. |
ref.signatures |
Gene signatures. |
robust |
A logical indicating that z-scores should be computed with median and MAD instead of mean and standard deviation. |
This function relies on a simple average of gene z-scores over each signature. It is no replacement for mode advanced methods such as CIBERSORT or BisqueRNA. It is provided for convenience.
A matrix containing the scores of
each gene signature in each sample.
Note that ligand-receptor gene
signature scores should be computed with
"scoreLRGeneSignatures"
instead.
data(sdc, package = "BulkSignalR") data(bsrdm, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") imm.scores <- scoreSignatures(bsrdm, immune.signatures)
data(sdc, package = "BulkSignalR") data(bsrdm, package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") imm.scores <- scoreSignatures(bsrdm, immune.signatures)
A dataset containing the read counts of salivary duct carcinomas and adjacent normal tissues.
data(sdc)
data(sdc)
A data frame with 19764 rows and 26 variables.
https://www.ncbi.nlm.nih.gov/geo/query/acc.cgi?acc=GSE138581
Generate a detailed view related to a chosen interaction made of series of small individual spatial plots: tissue organization (optional), gene signature score, ligand and receptor expression.
separatedLRPlot( v, L, R, ncounts, areas, inter.name = NULL, rev.y = TRUE, ref.plot = TRUE, image.raster = NULL, x.col = "array_col", y.col = "array_row", label.col = "label", idSpatial.col = "idSpatial", cut.p = 0.01, low.color = "royalblue3", mid.color = "white", high.color = "orange", title.fs = 12, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 0.5, legend.dot.factor = 10, ref.colors = NULL )
separatedLRPlot( v, L, R, ncounts, areas, inter.name = NULL, rev.y = TRUE, ref.plot = TRUE, image.raster = NULL, x.col = "array_col", y.col = "array_row", label.col = "label", idSpatial.col = "idSpatial", cut.p = 0.01, low.color = "royalblue3", mid.color = "white", high.color = "orange", title.fs = 12, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 0.5, legend.dot.factor = 10, ref.colors = NULL )
v |
A named vector containing the gene signature scores for the
L-R interaction including the contribution of the pathway, names must be
the IDs of each location. Alternatively, v can be a gene signature score
matrix such as those returned by |
L |
The name of the ligand. |
R |
The name of the receptor. |
ncounts |
The (normalized) expression matrix with column names equal to the IDs of each location. |
areas |
A data.frame containing at leastcluster_columns the x and y
coordinates of the locations as well as the unique IDs of spatial locations.
In case |
inter.name |
Interaction name to display as plot title,
equal to " |
rev.y |
A Boolean indicating whether low y coordinates should be at the top of the plot. |
ref.plot |
A Boolean indicating whether a reference map of the tissue with area labels should be plot aside. |
image.raster |
Raster object image to plot raw tissue image as reference. |
x.col |
Column name in |
y.col |
Column name in |
label.col |
Column name in |
idSpatial.col |
Column name in |
cut.p |
Proportion of top and bottom values for thresholding. |
low.color |
Color for low score values. |
mid.color |
Color for score = 0. |
high.color |
Color for high score values. |
title.fs |
Title font size. |
legend.fs |
Legend items font size. |
axis.fs |
Axis ticks font size. |
label.fs |
Legend titles and axis names font size. |
dot.size |
Dot size. |
legend.dot.factor |
A factor applied to obtain the legend dot size. |
ref.colors |
A vector of colors to bypass those automatically chosen by ggplot2 for the tissue areas in the reference plot. |
A set of spatial plots are generated including an optional reference tissue plot (image or areas represented), the gene signature scores, the ligand expression values, and the receptor expression values.
A set of spatial plots.
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) separatedLRPlot(scores.red, "SLIT2", "GPC1", ncounts(bsrdm.spa), annotation.spa, label.col = "ground_truth")
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) separatedLRPlot(scores.red, "SLIT2", "GPC1", ncounts(bsrdm.spa), annotation.spa, label.col = "ground_truth")
Generate a list of heatmaps for ligand, receptor and target genes for a specific pathway
signatureHeatmaps( pathway, bsrdm, bsrsig, h.width = 6, h.height = 9, fontsize = 6, show_column_names = FALSE )
signatureHeatmaps( pathway, bsrdm, bsrsig, h.width = 6, h.height = 9, fontsize = 6, show_column_names = FALSE )
pathway |
Pathway name |
bsrdm |
BulkSignalR data model object. |
bsrsig |
BulkSignalR signature object. to display on screen. |
h.width |
Heatmap width in cm. |
h.height |
Heatmap height in cm. |
fontsize |
Fontsize. |
show_column_names |
Add column names on heatmap. |
A plot is created.
This is a convenience function to propose a simple way of representing expression of genes involved in a specific pathway.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") if(FALSE){ bsrinf.redP <- reduceToPathway(bsrinf) bsrinf.redPBP <- reduceToBestPathway(bsrinf) bsrsig.redPBP <- BSRSignature(bsrinf, qval.thres = 1) pathway1 <- pathways(bsrsig.redPBP)[1] signatureHeatmaps( pathway = pathway1, bsrdm = bsrdm, bsrsig = bsrsig.redPBP, h.width = 3, h.height = 4, fontsize = 1, show_column_names = TRUE) }
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") if(FALSE){ bsrinf.redP <- reduceToPathway(bsrinf) bsrinf.redPBP <- reduceToBestPathway(bsrinf) bsrsig.redPBP <- BSRSignature(bsrinf, qval.thres = 1) pathway1 <- pathways(bsrsig.redPBP)[1] signatureHeatmaps( pathway = pathway1, bsrdm = bsrdm, bsrsig = bsrsig.redPBP, h.width = 3, h.height = 4, fontsize = 1, show_column_names = TRUE) }
Generate a heatmap representing ligand-receptor gene signature scores.
simpleHeatmap( mat.c, width = 4, height = 3, dend.row = NULL, dend.spl = NULL, cols = NULL, pointsize = 4, bottom.annotation = NULL, n.col.clust = 0, n.row.clust = 0, gap.size = 0.5, cut.p = 0.01, row.names = TRUE, column.names = TRUE, hcl.palette = NULL, reverse = FALSE )
simpleHeatmap( mat.c, width = 4, height = 3, dend.row = NULL, dend.spl = NULL, cols = NULL, pointsize = 4, bottom.annotation = NULL, n.col.clust = 0, n.row.clust = 0, gap.size = 0.5, cut.p = 0.01, row.names = TRUE, column.names = TRUE, hcl.palette = NULL, reverse = FALSE )
mat.c |
A matrix with the signature scores such as output by
|
width |
Heatmap width. |
height |
Heatmap height. |
dend.row |
A precomputed row dendrogram. |
dend.spl |
A precompute sample (column) dendrogram. |
cols |
A vector of colors to use for the heatmap. |
pointsize |
Heatmap fontsize |
bottom.annotation |
|
n.col.clust |
Number of column clusters. |
n.row.clust |
Number of row clusters. |
gap.size |
Gap size between clusters. |
cut.p |
Proportion of top and bottom values for thresholding. |
row.names |
A logical to turn on/off the display of row names. |
column.names |
A logical to turn on/off the display of column (sample) names. |
hcl.palette |
support for HCL colormaps in ComplexHeatmap using color mapping function with circlize::colorRamp2(). palettes are listed in grDevides::hcl.pals(). of row (gene) names. |
reverse |
A logicial to reverse or not colors in hcl.palette. |
A heatmap. Since heatmap plotting tend to be slow on the screen, it is advisable to provide a PDF file name and plot in a file (much faster).
If hcl.palette is set, the colors parameter won't be used.
Extreme values (top and bottom) can be replaced by global quantiles at
cut.p
and 1-cut.p
to avoid color scales shrunk by a few outliers.
This is a convenience function that relies on the ComplexHeatmap
package to propose a simple way
of representing signature scores. If more advanced features are needed or
more graphic parameters
should be controlled, users should implement their own function.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") bsrinf.redBP <- reduceToBestPathway(bsrinf) bsrsig.redBP <- BSRSignature(bsrinf, qval.thres = 0.001 ) scoresLR <- scoreLRGeneSignatures(bsrdm, bsrsig.redBP, name.by.pathway = FALSE ) simpleHeatmap(scoresLR[1:3, ], column.names = TRUE, hcl.palette = "Cividis", width=2, height=1.5)
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") bsrinf.redBP <- reduceToBestPathway(bsrinf) bsrsig.redBP <- BSRSignature(bsrinf, qval.thres = 0.001 ) scoresLR <- scoreLRGeneSignatures(bsrdm, bsrsig.redBP, name.by.pathway = FALSE ) simpleHeatmap(scoresLR[1:3, ], column.names = TRUE, hcl.palette = "Cividis", width=2, height=1.5)
Smooth spatial expression data
smoothSpatialCounts( bsrdm, areas, nnn = 4, radius = NULL, weight.ratio = 0.5, x.col = "array_col", y.col = "array_row" )
smoothSpatialCounts( bsrdm, areas, nnn = 4, radius = NULL, weight.ratio = 0.5, x.col = "array_col", y.col = "array_row" )
bsrdm |
A BSRDataModel object containing the expression data to smooth. |
areas |
A data.frame containing at least the x and y coordinates of the locations. |
nnn |
Number of nearest-neighbor locations to use for smoothing
each location. In case |
radius |
A maximal distance to include neighbors in the smoothing. |
weight.ratio |
The weight given to the central location. |
x.col |
Column name in |
y.col |
Column name in |
The expression data contained in a BSRDataModel object are smoothed using a weighted average of nearby locations.
Two strategies are available to identify the neighbors. It is possible to
simply set the number of nearest-neighbors (parameter nnn
). An
alternative consists in providing a distance radius (radius
) along
with a a maximum number of nearest-neighbors within the radius
(nnn.radius
). To properly define the radius, the user must know the
location coordinates. The strategy with the radius enables having corner
locations with two neighbors only and border locations with three
neighbors only, whereas to simply set a maximum of four neighbors for
instance would retrieve the four closest neighbors in every case.
For each location, its nearest-neighbors are found and a weighted average
computed with weight.ratio
given to the central location itself
and a total weight of 1-weight.ratio
shared within the neighbors
based on the inverse of their distances. In case radius
is set,
some locations may have less than nnn
neighbors (see above).
At such locations, the weight given to
the central location is augmented according to
1-(1-weight.ratio
)*(number of neighbors)/nnn
.
A BSRDataModel object containing the smoothed ncounts.
data(bsrdm.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") sm.bsrdm <- smoothSpatialCounts(bsrdm.spa, annotation.spa, radius = 1.2, nnn = 4)
data(bsrdm.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") sm.bsrdm <- smoothSpatialCounts(bsrdm.spa, annotation.spa, radius = 1.2, nnn = 4)
Source comparison name accessor
## S4 method for signature 'BSRInferenceComp' sourceComparisonName(x)
## S4 method for signature 'BSRInferenceComp' sourceComparisonName(x)
x |
BSRInferenceComp object |
src.comp.name
bsrinf <- new("BSRInferenceComp") sourceComparisonName(bsrinf)
bsrinf <- new("BSRInferenceComp") sourceComparisonName(bsrinf)
Compute the statistical association of L-R interaction score spatial distributions with tissue area labels. Not limited to BulkSignalR gene signature scores.
spatialAssociation( scores, areas, test = c("Kruskal-Wallis", "ANOVA", "Spearman", "r2"), label.col = "label", idSpatial.col = "idSpatial", fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
spatialAssociation( scores, areas, test = c("Kruskal-Wallis", "ANOVA", "Spearman", "r2"), label.col = "label", idSpatial.col = "idSpatial", fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
scores |
A matrix of scores, one L-R interaction per row and
spatial locations in the columns. This matrix is typically obtained
from BulkSignalR functions |
areas |
A data.frame containing at least the x and y coordinates of the locations, the unique IDs of spatial locations, and a label column. |
test |
The chosen statistical test or statistics (see details below). |
label.col |
Column name in |
idSpatial.col |
Column name in |
fdr.proc |
Multiple hypothesis correction procedure, see
|
In case the nonparametric Kruskal-Wallis test is chosen, additional columns are provided testing each label for significantly larger scores (Kruskal-Wallis is global and only says whether one or several labels show a bias). Individual labels are tested with Wilcoxon and two columns are added *per* label, one for the statistics and one for a Bonferroni-corrected P-value over all the labels.
In case an actual statistical test is chosen, a parametric test (ANOVA) and a non-parametric test (Kruskal-Wallis) are available for the global analysis. Individual labels are tested with T-tests or Wilcoxon (Bonferroni-corrected) accordingly.
In case a statistics is preferred, Spearman correlation or explained variance
(r2 or coefficient of determination, through linear models)
are available. They mesure the relationship
between each individual area and scores
. For the explained variance,
a global value (R2) is also computed from a multi-linear model (the same as
what is used for the ANOVA).
A data.frame with the names of the interactions, the value of the chosen statistics, and the corresponding Q-value.
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 #bsrinf.red <- reduceToBestPathway(bsrinf.spa) #s.red <- BSRSignature(bsrinf.red, qval.thres=thres) #scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) # Run in other examples no need to be run again # spatialAssociation(scores.red[c(1:2),], areas = annotation.spa, # label.col = "ground_truth")
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 #bsrinf.red <- reduceToBestPathway(bsrinf.spa) #s.red <- BSRSignature(bsrinf.red, qval.thres=thres) #scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) # Run in other examples no need to be run again # spatialAssociation(scores.red[c(1:2),], areas = annotation.spa, # label.col = "ground_truth")
Plot a heatmap featuring Q-values or values of statistical association between L-R interaction score spatial distributions and tissue area labels.
spatialAssociationPlot( associations, qval.thres = 0.01, absval.thres = 0, colors = NULL )
spatialAssociationPlot( associations, qval.thres = 0.01, absval.thres = 0, colors = NULL )
associations |
A statistical association data.frame generated
by the function |
qval.thres |
The maximum Q-value to consider in the plot (a L-R interaction must associate with one label at least with a Q-value smaller or equal to this threshold). |
absval.thres |
The minimum value to consider in the plot (a L-R interaction must associate with one label at least with an absolute value larger or equal to this threshold). |
colors |
A function returning a color for a given value such as
generated by |
Display a heatmap linking L-R interactions to labels.
ComplexHeatmap::Heatmap object
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) # statistical association with tissue areas based on correlations assoc.bsr.corr <- spatialAssociation(scores.red[c(1:10), ], areas = annotation.spa, label.col = "ground_truth",test = "Spearman") spatialAssociationPlot(assoc.bsr.corr)
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) # statistical association with tissue areas based on correlations assoc.bsr.corr <- spatialAssociation(scores.red[c(1:10), ], areas = annotation.spa, label.col = "ground_truth",test = "Spearman") spatialAssociationPlot(assoc.bsr.corr)
Use PCA or t-SNE to obtain a 2D-projection of a set of spatial scores
or associations.
This plot summarizes the diversity of patterns occuring in a spatial
dataset. Use the function spatialIndexPlot
to create a large
visual index of many spatial distributions.
Not limited to BulkSignalR gene signature scores.
spatialDiversityPlot( scores, associations, proj = c("PCA", "tSNE"), score.based = FALSE, qval.thres = 0.01, val.thres = 0, with.names = FALSE, text.fs = 2.5, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 1, perplexity = 10 )
spatialDiversityPlot( scores, associations, proj = c("PCA", "tSNE"), score.based = FALSE, qval.thres = 0.01, val.thres = 0, with.names = FALSE, text.fs = 2.5, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 1, perplexity = 10 )
scores |
A matrix of scores, one L-R interaction per row and
spatial locations in the columns. This matrix is typically obtained
from BulkSignalR functions |
associations |
A statistical association data.frame generated
by the function |
proj |
Projection method : 'PCA' or 'tSNE' are available arguements. |
score.based |
A logical indicating whether the plot should be based on scores or the associations directly. |
qval.thres |
The maximum Q-value to consider in the plot (a
L-R interaction must associate with one label at least with a Q-value
smaller or equal to this threshold). Relevant for Kruskal-Wallis and
ANOVA tests in |
val.thres |
The minimum value to consider in the plot (a
L-R interaction must associate with one label at least with a value
larger or equal to this threshold). Relevant for Spearman and r2
associations in |
with.names |
A logical indicating whether L-R names should be plotted. |
text.fs |
Point label font size in case |
legend.fs |
Legend items font size. |
axis.fs |
Axis ticks font size. |
label.fs |
Legend titles and axis names font size. |
dot.size |
Dot size. |
perplexity |
Perplexity parameter for t-SNE. |
Display a 2D-projection of the score spatial distributions.
Display a 2D-projection of the score spatial distributions.
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) # statistical association with tissue areas based on correlations # For display purpose, we only use a subset here assoc.bsr.corr <- spatialAssociation(scores.red[c(1:3), ], annotation.spa, label.col = "ground_truth",test = "Spearman") spatialDiversityPlot(scores.red[c(1:3),],assoc.bsr.corr)
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) # statistical association with tissue areas based on correlations # For display purpose, we only use a subset here assoc.bsr.corr <- spatialAssociation(scores.red[c(1:3), ], annotation.spa, label.col = "ground_truth",test = "Spearman") spatialDiversityPlot(scores.red[c(1:3),],assoc.bsr.corr)
Generate an index made of series of small individual spatial score plots in a PDF. Not limited to BulkSignalR gene signature scores.
spatialIndexPlot( scores, areas, out.file, ref.plot = TRUE, image.raster = NULL, x.col = "array_col", y.col = "array_row", label.col = "label", idSpatial.col = "idSpatial", cut.p = 0.01, low.color = "royalblue3", mid.color = "white", high.color = "orange", title.fs = 12, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 0.25, ratio = 1.25, base.v = 2.5, base.h = 3, ref.colors = NULL )
spatialIndexPlot( scores, areas, out.file, ref.plot = TRUE, image.raster = NULL, x.col = "array_col", y.col = "array_row", label.col = "label", idSpatial.col = "idSpatial", cut.p = 0.01, low.color = "royalblue3", mid.color = "white", high.color = "orange", title.fs = 12, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 0.25, ratio = 1.25, base.v = 2.5, base.h = 3, ref.colors = NULL )
scores |
A matrix of scores, one L-R interaction per row and
spatial locations in the columns. This matrix is typically obtained
from BulkSignalR functions |
areas |
A data.frame containing at least the x and y coordinates of the locations, the unique IDs of spatial locations, and a tissue label column. |
out.file |
File name for the output PDF. |
ref.plot |
A Boolean indicating whether a reference map of the tissue with area labels should be plot first. |
image.raster |
Raster object image to plot raw tissue image as reference. |
x.col |
Column name in |
y.col |
Column name in |
label.col |
Column name in |
idSpatial.col |
Column name in |
cut.p |
Proportion of top and bottom values for thresholding. |
low.color |
Color for low score values. |
mid.color |
Color for score = 0. |
high.color |
Color for high score values. |
title.fs |
Title font size. |
legend.fs |
Legend items font size. |
axis.fs |
Axis ticks font size. |
label.fs |
Legend titles and axis names font size. |
dot.size |
Dot size. |
ratio |
the vertical/horizontal ratio. |
base.v |
Height of each plot. |
base.h |
Width of each plot. |
ref.colors |
A vector of colors to bypass those automatically chosen by ggplot2 for the tissue areas in the reference plot. |
A PDF file is created that contains the index.
Create PDF file and returns 'NULL', invisibly.
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) # generate visual index on disk in pdf file spatialIndexPlot(scores.red[1:2,], annotation.spa, label.col = "ground_truth", out.file = "spatialIndexPlot")
data(bsrdm.spa, package = "BulkSignalR") data(bsrinf.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) # generate visual index on disk in pdf file spatialIndexPlot(scores.red[1:2,], annotation.spa, label.col = "ground_truth", out.file = "spatialIndexPlot")
Generate a plot with scores at the spatial coordinates of the corresponding sample locations. Not limited to BulkSignalR gene signature scores.
spatialPlot( v, areas, inter.name, rev.y = TRUE, ref.plot = FALSE, ref.plot.only = FALSE, image.raster = NULL, x.col = "array_col", y.col = "array_row", label.col = "label", idSpatial.col = "idSpatial", cut.p = 0.01, low.color = "royalblue3", mid.color = "white", high.color = "orange", title.fs = 12, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 0.5, legend.dot.factor = 10, ref.colors = NULL )
spatialPlot( v, areas, inter.name, rev.y = TRUE, ref.plot = FALSE, ref.plot.only = FALSE, image.raster = NULL, x.col = "array_col", y.col = "array_row", label.col = "label", idSpatial.col = "idSpatial", cut.p = 0.01, low.color = "royalblue3", mid.color = "white", high.color = "orange", title.fs = 12, legend.fs = 10, axis.fs = 10, label.fs = 12, dot.size = 0.5, legend.dot.factor = 10, ref.colors = NULL )
v |
A named vector containing the scores, names must be the IDs of each location. |
areas |
A data.frame containing at least the x and y
coordinates of the locations as well as the unique IDs of spatial locations.
In case |
inter.name |
Interaction name to display as plot title. |
rev.y |
A Boolean indicating whether low y coordinates should be at the top of the plot. |
ref.plot |
A Boolean indicating whether a reference map of the tissue with area labels should be plot aside. |
ref.plot.only |
A Boolean indicating that only the reference plot should be output. |
image.raster |
Raster object image to plot raw tissue image as reference. |
x.col |
Column name in |
y.col |
Column name in |
label.col |
Column name in |
idSpatial.col |
Column name in |
cut.p |
Proportion of top and bottom values for thresholding. |
low.color |
Color for low score values. |
mid.color |
Color for score = 0. |
high.color |
Color for high score values. |
title.fs |
Title font size. |
legend.fs |
Legend items font size. |
axis.fs |
Axis ticks font size. |
label.fs |
Legend titles and axis names font size. |
dot.size |
Dot size. |
legend.dot.factor |
A factor applied to obtain the legend dot size. |
ref.colors |
A vector of colors to bypass those automatically chosen by ggplot2 for the tissue areas in the reference plot. |
A single (scores) or side-by-side (reference tissue & scores) plot is generated.
A spatial plot
data(bsrinf.spa, package = "BulkSignalR") data(bsrdm.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) inter <- "{SLIT2} / {GPC1}" spatialPlot(scores.red[inter, ], annotation.spa, inter, ref.plot = TRUE, ref.plot.only = FALSE, image.raster = NULL, dot.size = 1, label.col = "ground_truth")
data(bsrinf.spa, package = "BulkSignalR") data(bsrdm.spa, package = "BulkSignalR") data(annotation.spa, package = "BulkSignalR") thres <- 0.01 bsrinf.red <- reduceToBestPathway(bsrinf.spa) s.red <- BSRSignature(bsrinf.red, qval.thres=thres) scores.red <- scoreLRGeneSignatures(bsrdm.spa,s.red) inter <- "{SLIT2} / {GPC1}" spatialPlot(scores.red[inter, ], annotation.spa, inter, ref.plot = TRUE, ref.plot.only = FALSE, image.raster = NULL, dot.size = 1, label.col = "ground_truth")
Generate a igraph object with one link between each cell type.
summarizedCellularNetwork(tab)
summarizedCellularNetwork(tab)
tab |
The data.frame output by |
A igraph object containing a summary cellular network with edge weights proportional to the sum of individual link scores. Edge weight are normalized to a total of one.
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data("tme.signatures", package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores) # cellular network g.table <- cellularNetworkTable(lr2ct[c(1:25),]) gSummary <- summarizedCellularNetwork(g.table) # plot(gSummary, edge.width=1+30*E(gSummary)$score)
data(bsrdm, package = "BulkSignalR") data(bsrinf, package = "BulkSignalR") data("tme.signatures", package = "BulkSignalR") data(immune.signatures, package = "BulkSignalR") immune.signatures <- immune.signatures[immune.signatures$signature %in% c("T cells"), ] signatures <- rbind(immune.signatures, tme.signatures[ tme.signatures$signature %in% c("Fibroblasts"), ]) tme.scores <- scoreSignatures(bsrdm, signatures) # assignment lr2ct <- assignCellTypesToInteractions(bsrdm, bsrinf, tme.scores) # cellular network g.table <- cellularNetworkTable(lr2ct[c(1:25),]) gSummary <- summarizedCellularNetwork(g.table) # plot(gSummary, edge.width=1+30*E(gSummary)$score)
Target gene correlations accessor
Target gene correlations accessor
## S4 method for signature 'BSRInference' tgCorr(x) ## S4 method for signature 'BSRSignature' tgCorr(x)
## S4 method for signature 'BSRInference' tgCorr(x) ## S4 method for signature 'BSRSignature' tgCorr(x)
x |
BSRSignature |
tgCorr
bsr.sig <- new("BSRSignature") tgCorr(bsr.sig)
bsr.sig <- new("BSRSignature") tgCorr(bsr.sig)
Target gene expression accessor
Target gene expression accessor
## S4 method for signature 'BSRInferenceComp' tgExpr(x) ## S4 method for signature 'BSRSignatureComp' tgExpr(x)
## S4 method for signature 'BSRInferenceComp' tgExpr(x) ## S4 method for signature 'BSRSignatureComp' tgExpr(x)
x |
BSRSignatureComp object |
tgExpr
tg.expr
bsrinf <- new("BSRInferenceComp") tgExpr(bsrinf)
bsrinf <- new("BSRInferenceComp") tgExpr(bsrinf)
Target genes accessor
Target genes accessor
## S4 method for signature 'BSRInference' tgGenes(x) ## S4 method for signature 'BSRSignature' tgGenes(x)
## S4 method for signature 'BSRInference' tgGenes(x) ## S4 method for signature 'BSRSignature' tgGenes(x)
x |
BSRSignature |
tgGenes
tgGenes
bsr.sig <- new("BSRSignature") tgGenes(bsr.sig)
bsr.sig <- new("BSRSignature") tgGenes(bsr.sig)
Target gene logFC accessor
Target gene logFC accessor
## S4 method for signature 'BSRInferenceComp' tgLogFC(x) ## S4 method for signature 'BSRSignatureComp' tgLogFC(x)
## S4 method for signature 'BSRInferenceComp' tgLogFC(x) ## S4 method for signature 'BSRSignatureComp' tgLogFC(x)
x |
BSRSignatureComp object |
tgLogFC
tg.logFC
bsrinf <- new("BSRInferenceComp") tgLogFC(bsrinf)
bsrinf <- new("BSRInferenceComp") tgLogFC(bsrinf)
Target gene P-values accessor
Target gene P-values accessor
## S4 method for signature 'BSRInferenceComp' tgPval(x) ## S4 method for signature 'BSRSignatureComp' tgPval(x)
## S4 method for signature 'BSRInferenceComp' tgPval(x) ## S4 method for signature 'BSRSignatureComp' tgPval(x)
x |
BSRSignatureComp object |
tgPval
tg.pval
bsrinf <- new("BSRInferenceComp") tgPval(bsrinf)
bsrinf <- new("BSRInferenceComp") tgPval(bsrinf)
A dataset containing gene signatures for some immune and stromal cell populations that are present in the microenvironment of a tumor.
data(tme.signatures)
data(tme.signatures)
A data frame with 209 rows and 2 variables:
HUGO gene symbol
cell population name
Becht & al., Genome Biol, 2016; Angelova et al., Genome Biol, 2015.
A method to update the data underlying statistical significance estimations prior to rescoring for an existing BSRInferenceComp object (P- and Q-value estimations as well as LR-score).
## S4 method for signature 'BSRInferenceComp' updateInference( obj, bsrcc, ncounts, src.bsrcc = NULL, rank.p = 0.55, max.pval = 0.01, min.logFC = 1, min.LR.score = 0, neg.receptors = FALSE, pos.targets = FALSE, neg.targets = FALSE, min.t.logFC = 0.5, min.positive = 2, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
## S4 method for signature 'BSRInferenceComp' updateInference( obj, bsrcc, ncounts, src.bsrcc = NULL, rank.p = 0.55, max.pval = 0.01, min.logFC = 1, min.LR.score = 0, neg.receptors = FALSE, pos.targets = FALSE, neg.targets = FALSE, min.t.logFC = 0.5, min.positive = 2, fdr.proc = c("BH", "Bonferroni", "Holm", "Hochberg", "SidakSS", "SidakSD", "BY", "ABH", "TSBH") )
obj |
BSRInferenceComp object. |
bsrcc |
BSRClusterComp object relative to target cells. |
ncounts |
Matrix counts normalized. |
src.bsrcc |
BSRClusterComp object relative to source cells. |
rank.p |
A number between 0 and 1 defining the rank of the last considered target genes. |
max.pval |
The maximum P-value imposed to both the ligand and the receptor. |
min.logFC |
The minimum log2 fold-change allowed for both the receptor and the ligand. |
min.LR.score |
The minimum LR-score allowed for the interaction. |
neg.receptors |
A logical indicating whether receptors are only allowed to be upregulated (FALSE), or up- and downregulated (TRUE). |
pos.targets |
A logical imposing that all the network targets must display positive logFC, i.e. logFC >= min.t.logFC. |
neg.targets |
A logical imposing that all the network targets must display negative logFC, i.e. logFC <= - min.t.logFC. |
min.t.logFC |
The minimum log2 fold-change allowed for targets in case pos.targets or neg.targets are used. |
min.positive |
Minimum number of target genes to be found in a given pathway. |
fdr.proc |
The procedure for adjusting P-values according to
|
A BSRInferenceComp object should be created by calling
"BSRInferenceComp"
A BSRInferenceComp object. The main application of this
method is to take a "universal" inference obtained by assigning
each gene to good logFC, P-values and expression levels whose role
is to find all the reachable targets per receptor/pathway, and
to update it by using actual logFC, P-values, and expression data.
The benefit is to save time when multiple sample comparisons are
performed, only one network exploration is necessary. Note that
if a restrictive logic such as positive.targets=TRUE
is used,
the result will be correct provided all the targets were in the
initial BSRInferenceComp object. If a restriction on the targets
was applied, then the update is likely to miss some targets, i.e.,
the statistical analysis will be wrong.
Note that correlations are set to 1 to avoid lengthy computations with scRNA-seq data and multiple cell populations.
The main function of this method is to support our SingleCellSignalR v2 package.
data(bsrdm.comp, package = "BulkSignalR") data(bsrinf.comp, package = "BulkSignalR") colA <- as.integer(1:2) colB <- as.integer(3:4) #bsrdm.comp <- as(bsrdm, "BSRDataModelComp") n <- nrow(ncounts(bsrdm.comp)) stats <- data.frame(pval = runif(n), logFC = rnorm(n, 0, 2), expr = runif(n, 0, 10)) rownames(stats) <- rownames(ncounts(bsrdm.comp)) # update stats$pval <- stats$pval / 100 stats$logFC <- stats$logFC + 0.5 bsrcc.2 <- BSRClusterComp(bsrdm.comp, colA, colB, stats) bsrinf.updated <- updateInference(bsrinf.comp, bsrcc.2, max.pval = 1, min.logFC = 0.1)
data(bsrdm.comp, package = "BulkSignalR") data(bsrinf.comp, package = "BulkSignalR") colA <- as.integer(1:2) colB <- as.integer(3:4) #bsrdm.comp <- as(bsrdm, "BSRDataModelComp") n <- nrow(ncounts(bsrdm.comp)) stats <- data.frame(pval = runif(n), logFC = rnorm(n, 0, 2), expr = runif(n, 0, 10)) rownames(stats) <- rownames(ncounts(bsrdm.comp)) # update stats$pval <- stats$pval / 100 stats$logFC <- stats$logFC + 0.5 bsrcc.2 <- BSRClusterComp(bsrdm.comp, colA, colB, stats) bsrinf.updated <- updateInference(bsrinf.comp, bsrcc.2, max.pval = 1, min.logFC = 0.1)