Title: | xcore expression regulators inference |
---|---|
Description: | xcore is an R package for transcription factor activity modeling based on known molecular signatures and user's gene expression data. Accompanying xcoredata package provides a collection of molecular signatures, constructed from publicly available ChiP-seq experiments. xcore use ridge regression to model changes in expression as a linear combination of molecular signatures and find their unknown activities. Obtained, estimates can be further tested for significance to select molecular signatures with the highest predicted effect on the observed expression changes. |
Authors: | Maciej Migdał [aut, cre] , Bogumił Kaczkowski [aut] |
Maintainer: | Maciej Migdał <[email protected]> |
License: | GPL-2 |
Version: | 1.11.0 |
Built: | 2024-10-31 06:36:40 UTC |
Source: | https://github.com/bioc/xcore |
addSignatures
extends mae
by adding to it new experiments.
Rows consistency is ensured by taking an intersection of rows after new
experiments are added.
addSignatures(mae, ..., intersect_rows = TRUE)
addSignatures(mae, ..., intersect_rows = TRUE)
mae |
MultiAssayExperiment object. |
... |
named experiments to be added to |
intersect_rows |
logical flag indicating if only common rows across
experiments should be included. Only set to |
MultiAssayExperiment object with new experiments added.
data("rinderpest_mini", "remap_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- prepareCountsForRegression( counts = rinderpest_mini, design = design, base_lvl = base_lvl) mae <- addSignatures(mae, remap = remap_mini)
data("rinderpest_mini", "remap_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- prepareCountsForRegression( counts = rinderpest_mini, design = design, base_lvl = base_lvl) mae <- addSignatures(mae, remap = remap_mini)
Returns a array obtained by applying a function to rows of submatrices of the input matrix, where the submatrices are divided into specified groups of columns.
applyOverColumnGroups(mat, groups, f, ...)
applyOverColumnGroups(mat, groups, f, ...)
mat |
a matrix. |
groups |
a vector giving columns grouping. |
f |
function to be applied. |
... |
optional arguments to |
a matrix of dimensions nrow(mat) x nlevels(groups)
.
applyOverDFList
operates on a list of data frames where all data frames
has the same size and columns. Column of interest is extracted from each data
frame and column binded in groups
, next fun
is applied over
rows. Final result is a matrix with result for each group on a separate column.
Function is parallelized over groups.
applyOverDFList(list_of_df, col_name, fun, groups)
applyOverDFList(list_of_df, col_name, fun, groups)
list_of_df |
list of |
col_name |
string specifying column in |
fun |
function to apply, should take a single vector as a argument. |
groups |
factor defining how elements of |
matrix with nrow(list_of_df[[1]])
rows and
nlevels(groups)
columns.
Transform design matrix to factor
design2factor(design)
design2factor(design)
design |
design matrix |
factor
## Not run: design <- matrix(data = c(1, 1, 0, 0, 0, 0, 1, 1), nrow = 4, ncol = 2, dimnames = list(c(paste("sample", 1:4)), c("gr1", "gr2"))) design2factor(design) ## End(Not run)
## Not run: design <- matrix(data = c(1, 1, 0, 0, 0, 0, 1, 1), nrow = 4, ncol = 2, dimnames = list(c(paste("sample", 1:4)), c("gr1", "gr2"))) design2factor(design) ## End(Not run)
Estimate goodness of fit statistic of penalized linear regression models. Works with different goodness of fit statistic functions.
estimateStat(x, y, u, s, method = "cv", nfold = 10, statistic = rsq, alpha = 0)
estimateStat(x, y, u, s, method = "cv", nfold = 10, statistic = rsq, alpha = 0)
x |
input matrix, of dimension nobs x nvars; each row is an observation
vector. Can be in sparse matrix format (inherit from class
|
y |
response variable. Quantitative for |
u |
offset vector as in |
s |
user supplied lambda. |
method |
currently only cross-validation is implemented. |
nfold |
number of fold to use in cross-validation. |
statistic |
function computing goodness of fit statistic. Should accept
|
alpha |
The elasticnet mixing parameter, with
|
numeric vector of statistic
estimates.
Filter signatures overlapping low or high number of promoters. Useful to get rid of signatures that have very low variance.
filterSignatures( mae, min = 0.05, max = 0.95, ref_experiment = "Y", omit_experiments = c("Y", "U") )
filterSignatures( mae, min = 0.05, max = 0.95, ref_experiment = "Y", omit_experiments = c("Y", "U") )
mae |
MultiAssayExperiment object. |
min |
length one numeric between 0 and 1 defining minimum promoter coverage for the signature to pass filtering. |
max |
length one numeric between 0 and 1 defining maximum promoter coverage for the signature to pass filtering. |
ref_experiment |
string giving name of experiment to use for inferring total number of promoters. |
omit_experiments |
character giving names of experiments to exclude from filtering. |
MultiAssayExperiment object with selected experiments filtered.
data("rinderpest_mini", "remap_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- prepareCountsForRegression( counts = rinderpest_mini, design = design, base_lvl = base_lvl) mae <- addSignatures(mae, remap = remap_mini) mae <- filterSignatures(mae)
data("rinderpest_mini", "remap_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- prepareCountsForRegression( counts = rinderpest_mini, design = design, base_lvl = base_lvl) mae <- addSignatures(mae, remap = remap_mini) mae <- filterSignatures(mae)
Fisher's method is a meta-analysis technique used to combine the results from independent statistical tests with the same hypothesis (Wikipedia article).
fisherMethod(p.value, lower.tail = FALSE, log.p = TRUE)
fisherMethod(p.value, lower.tail = FALSE, log.p = TRUE)
p.value |
a numeric vector of p-values to combine. |
lower.tail |
logical; if TRUE (default), probabilities are
|
log.p |
logical; if TRUE, probabilities p are given as log(p). |
a number giving combined p-value.
getCoverage
calculates coverage of regions (rows in interaction
matrix) by features (columns). It is possible to specify features grouping
variable gr
then coverage tells how many distinct groups the region
overlap with.
getCoverage(mat, gr)
getCoverage(mat, gr)
mat |
dgCMatrix interaction matrix such as produced by
|
gr |
factor specifying features groups. Must have length equal
to number of columns in |
Numeric vector.
data("remap_mini") y <- colnames(remap_mini) # simple coverage gr <- seq_along(y) %>% as.factor() getCoverage(remap_mini, gr) # per cell type coverage gr <- sub(".*\\.", "", y) %>% as.factor() getCoverage(remap_mini, gr)
data("remap_mini") y <- colnames(remap_mini) # simple coverage gr <- seq_along(y) %>% as.factor() getCoverage(remap_mini, gr) # per cell type coverage gr <- sub(".*\\.", "", y) %>% as.factor() getCoverage(remap_mini, gr)
getInteractionMatrix
construct interaction matrix between two Granges
objects. Names of object a
became row names and names of b
column names.
getInteractionMatrix(a, b, ext = 500, count = FALSE)
getInteractionMatrix(a, b, ext = 500, count = FALSE)
a |
GRanges object. |
b |
GRanges object. |
ext |
Integer specifying number of base pairs the |
count |
Logical indicating if matrix should hold number of overlaps
between |
Sparse matrix of class dgCMatrix, with rows corresponding to
a
and columns to b
. Each cell holds a number indicating
how many times a
and b
overlapped.
a <- GenomicRanges::GRanges( seqnames = c("chr20", "chr4"), ranges = IRanges::IRanges( start = c(62475984L, 173530220L), end = c(62476001L, 173530236L)), strand = c("-", "-"), name = c("hg19::chr20:61051039..61051057,-;hg_188273.1", "hg19::chr4:174451370..174451387,-;hg_54881.1")) b <- GenomicRanges::GRanges( seqnames = c("chr4", "chr20"), ranges = IRanges::IRanges( start = c(173530229L, 63864270L), end = c(173530236L, 63864273L)), strand = c("-", "-"), name = c("HAND2", "GATA5")) getInteractionMatrix(a, b)
a <- GenomicRanges::GRanges( seqnames = c("chr20", "chr4"), ranges = IRanges::IRanges( start = c(62475984L, 173530220L), end = c(62476001L, 173530236L)), strand = c("-", "-"), name = c("hg19::chr20:61051039..61051057,-;hg_188273.1", "hg19::chr4:174451370..174451387,-;hg_54881.1")) b <- GenomicRanges::GRanges( seqnames = c("chr4", "chr20"), ranges = IRanges::IRanges( start = c(173530229L, 63864270L), end = c(173530236L, 63864273L)), strand = c("-", "-"), name = c("HAND2", "GATA5")) getInteractionMatrix(a, b)
Calculate variance weighted average coefficients matrix
getVarianceWeightedAvgCoeff(pvalues, groups)
getVarianceWeightedAvgCoeff(pvalues, groups)
pvalues |
list of |
groups |
factor giving the grouping. |
variance weighted average coefficients matrix
Check if argument is a binary flag
isTRUEorFALSE(x)
isTRUEorFALSE(x)
x |
object to test |
binary flag
Calculate Mean Absolute Error
mae(y, yhat, ...)
mae(y, yhat, ...)
y |
numeric vector of observed expression values. |
yhat |
numeric vector of predicted expression values. |
... |
not used. |
numeric vector
Helper summarizing MAE object
maeSummary(mae)
maeSummary(mae)
mae |
MultiAssayExperiment object. |
named list giving number of rows and columns, overall mean and
standard deviation in mae
's experiments.
modelGeneExpression
uses parallelization if parallel backend is
registered. For that reason we advise against passing parallel
argument
to internally called cv.glmnet
routine.
modelGeneExpression( mae, yname = "Y", uname = "U", xnames, design = NULL, standardize = TRUE, parallel = FALSE, pvalues = TRUE, precalcmodels = NULL, ... )
modelGeneExpression( mae, yname = "Y", uname = "U", xnames, design = NULL, standardize = TRUE, parallel = FALSE, pvalues = TRUE, precalcmodels = NULL, ... )
mae |
MultiAssayExperiment object such as produced by
|
yname |
string indicating experiment in |
uname |
string indicating experiment in |
xnames |
character indicating experiments in |
design |
matrix giving the design matrix for the samples. Default
( |
standardize |
logical flag indicating if the molecular signatures should
be scaled. Advised to be set to |
parallel |
parallel argument to internally used
|
pvalues |
logical flag indicating if significance testing for the estimated molecular signatures activities should be performed. |
precalcmodels |
optional list of precomputed |
... |
arguments passed to glmnet::cv.glmnet. |
For speeding up the calculations consider lowering number of folds used in
internally run cv.glmnet
by specifying nfolds
argument. By default 10 fold cross validation is used.
The relationship between the expression (Y) and molecular signatures (X) is described using linear model formulation. The pipeline attempts to model the change in expression between basal expression level (u) and each sample, with the goal of finding the unknown molecular signatures activities. Linear models are fit using popular ridge regression implementation glmnet (Friedman, Hastie, and Tibshirani 2010).
If pvalues
is set to TRUE
the significance of the estimated
molecular signatures activities is tested using methodology introduced by
(Cule, Vineis, and De Iorio 2011) which original implementation can be found
in ridge-package.
If replicates are available the signatures activities estimates and their standard error estimates can be combined. This is done by averaging signatures activities estimates and pooling their significance estimates using Stouffer's method for the Z-scores and Fisher's method for the p-values.
For detailed pipeline description we refer interested user to paper accompanying this package.
Nested list with following elements
Named list with elements corresponding to
signatures specified in xnames
. Each of these is a list holding
'cv.glmnet'
objects corresponding to each sample.
Named list with elements corresponding to
signatures specified in xnames
. Each of these is a list holding
data.frame
of signature's p-values and test statistics
estimated for each sample.
Named list with elements corresponding to
signatures specified in xnames
. Each of these is a matrix
holding replicate average Z-scores with columns corresponding to groups
in the design.
Named list with elements corresponding to
signatures specified in xnames
. Each of these is a matrix
holding replicate averaged signatures activities with columns
corresponding to groups in the design.
Named list of a data.frame
s holding replicate
average molecular signatures, overall molecular signatures Z-score and
p-values calculated over groups using Stouffer's and Fisher's methods.
data("rinderpest_mini", "remap_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- prepareCountsForRegression( counts = rinderpest_mini, design = design, base_lvl = base_lvl) mae <- addSignatures(mae, remap = remap_mini) mae <- filterSignatures(mae) res <- modelGeneExpression( mae = mae, xnames = "remap", nfolds = 5)
data("rinderpest_mini", "remap_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- prepareCountsForRegression( counts = rinderpest_mini, design = design, base_lvl = base_lvl) mae <- addSignatures(mae, remap = remap_mini) mae <- filterSignatures(mae) res <- modelGeneExpression( mae = mae, xnames = "remap", nfolds = 5)
Internal function used in modelGeneExpression
. It runs ridge regression
parallelly across signatures and samples as specified by experiment design.
modelGeneExpression_ridge_regression_wraper( mae, yname, uname, xnames, groups, standardize, parallel, precalcmodels, ... )
modelGeneExpression_ridge_regression_wraper( mae, yname, uname, xnames, groups, standardize, parallel, precalcmodels, ... )
mae |
MultiAssayExperiment object such as produced by
|
yname |
string indicating experiment in |
uname |
string indicating experiment in |
xnames |
character indicating experiments in |
groups |
factor representation of design matrix. |
standardize |
logical flag indicating if the molecular signatures should
be scaled. Advised to be set to |
parallel |
parallel argument to internally used
|
precalcmodels |
optional list of precomputed |
... |
arguments passed to glmnet::cv.glmnet. |
Named list with elements corresponding to signatures specified in
xnames
. Each of these is a list holding 'cv.glmnet'
objects
corresponding to each sample.
Internal function used in modelGeneExpression
. It runs ridgePvals
parallelly across signatures and samples as specified by experiment design.
modelGeneExpression_significance_testing_wraper( mae, yname, uname, xnames, groups, standardize, regression_models )
modelGeneExpression_significance_testing_wraper( mae, yname, uname, xnames, groups, standardize, regression_models )
mae |
MultiAssayExperiment object such as produced by
|
yname |
string indicating experiment in |
uname |
string indicating experiment in |
xnames |
character indicating experiments in |
groups |
factor representation of design matrix. |
standardize |
logical flag indicating if the molecular signatures should
be scaled. Advised to be set to |
regression_models |
Named list with elements corresponding to signatures
specified in |
Named list with elements corresponding to signatures specified in
xnames
. Each of these is a list holding data.frame
of
signature's p-values and test statistics estimated for each sample.
Calculate Mean Squared Error
mse(y, yhat, ...)
mse(y, yhat, ...)
y |
numeric vector of observed expression values. |
yhat |
numeric vector of predicted expression values. |
... |
not used. |
numeric vector
Expression counts are processed using edgeR following
User's Guide.
Shortly, counts for each sample are filtered for lowly expressed promoters,
normalized for the library size and transformed into counts per million (CPM).
Optionally, CPM are log2 transformed with addition of pseudo count. Basal
level expression is calculated by averaging base_lvl
samples
expression values.
prepareCountsForRegression( counts, design, base_lvl, log2 = TRUE, pseudo_count = 1L, drop_base_lvl = TRUE )
prepareCountsForRegression( counts, design, base_lvl, log2 = TRUE, pseudo_count = 1L, drop_base_lvl = TRUE )
counts |
matrix of read counts. |
design |
matrix giving the design matrix for the samples. Columns corresponds to samples groups and rows to samples names. |
base_lvl |
string indicating group in |
log2 |
logical flag indicating if counts should be log2(counts per million) should be returned. |
pseudo_count |
integer count to be added before taking log2. |
drop_base_lvl |
logical flag indicating if |
MultiAssayExperiment object with two experiments:
matrix giving expression values averaged over basal level samples
matrix of expression values
design with base_lvl
dropped is stored in metadata and directly
available for modelGeneExpression
.
data("rinderpest_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- prepareCountsForRegression( counts = rinderpest_mini, design = design, base_lvl = base_lvl)
data("rinderpest_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- prepareCountsForRegression( counts = rinderpest_mini, design = design, base_lvl = base_lvl)
regressionData
orgnize expression data and experiment design into
MultiAssayExperiment object that can be further used in xcore
framework. Additionally, function calculate basal expression level, for
latter use in expression modeling, by averaging base_lvl
samples
expression values.
regressionData(expr_mat, design, base_lvl, drop_base_lvl = TRUE)
regressionData(expr_mat, design, base_lvl, drop_base_lvl = TRUE)
expr_mat |
matrix of expression values. |
design |
matrix giving the design matrix for the samples. Columns corresponds to samples groups and rows to samples names. |
base_lvl |
string indicating group in |
drop_base_lvl |
logical flag indicating if |
Note that regressionData
does not apply any normalization or
transformation to the input data! Use prepareCountsForRegression
if you want to start with raw expression counts.
MultiAssayExperiment object with two experiments:
matrix giving expression values averaged over basal level samples
matrix of expression values
design with base_lvl
dropped is stored in metadata and directly
available for modelGeneExpression
.
data("rinderpest_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- regressionData( expr_mat = rinderpest_mini, design = design, base_lvl = base_lvl)
data("rinderpest_mini") base_lvl <- "00hr" design <- matrix( data = c(1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1), ncol = 3, nrow = 9, byrow = TRUE, dimnames = list(colnames(rinderpest_mini), c("00hr", "12hr", "24hr"))) mae <- regressionData( expr_mat = rinderpest_mini, design = design, base_lvl = base_lvl)
Molecular signatures data intended for use in xcore vignette and examples.
It is build ReMap2020 molecular signatures constructed against FANTOM5
annotation, which can be found in xcoredata
package.
Here the data is only a subset limited to core promoters
(promoters_f5_core
) and randomly selected 600 signatures.
data(remap_mini)
data(remap_mini)
A dgCMatrix
with 14191 rows and 600 columns holding
interaction matrix for subset of ReMap2020 molecular signatures against
FANTOM5 annotation. Rows corresponds to FANTOM5 promoters and columns to
signatures.
Replicate averaged Z-scores is calculated by dividing replicate average coefficient by replicate pooled standard error.
repVarianceWeightedAvgZscore(pvalues, groups)
repVarianceWeightedAvgZscore(pvalues, groups)
pvalues |
Data frame with |
groups |
Factor giving group membership for samples in |
Numeric matrix of averaged Z-scores. Columns correspond to groups and rows to predictors.
Standard error estimation and significance testing for coefficients
estimated in linear ridge regression. ridgePvals
re-implement
original method by (Cule et al. BMC Bioinformatics 2011.) found in
ridge-package. This function is intended to use with
cv.glmnet
output.
ridgePvals(x, y, beta, lambda, standardizex = TRUE, svdX = NULL)
ridgePvals(x, y, beta, lambda, standardizex = TRUE, svdX = NULL)
x |
input matrix, same as used in |
y |
response variable, same as used in |
beta |
matrix of coefficients, estimated using
|
lambda |
lambda value for which |
standardizex |
logical flag for x variable standardization, should be
set to same value as |
svdX |
optional singular-value decomposition of |
a data.frame with columns
beta
's names
beta
's standard errors
beta
's test statistic
beta
's p-values
Expression data intended for use in xcore vignette and examples. It is build
from FANTOM5's 293SLAM rinderpest infection time course dataset. Here the
data is only a subset limited to core promoters (promoters_f5_core
).
data(rinderpest_mini)
data(rinderpest_mini)
A matrix
with 14191 rows and 6 columns holding expression
counts from CAGE-seq experiment. Rows corresponds to FANTOM5 promoters and
columns to time points at which expression was measured 0 and 24 hours post
infection.
Calculate $R^2$
rsq(y, yhat, offset)
rsq(y, yhat, offset)
y |
numeric vector of observed expression values. |
yhat |
numeric vector of predicted expression values. |
offset |
numeric vector giving basal expression level. |
numeric vector
Simplify Interaction Matrix
simplifyInteractionMatrix(mat, alpha = 0.5, colname = NA)
simplifyInteractionMatrix(mat, alpha = 0.5, colname = NA)
mat |
dgCMatrix interaction matrix such as produced by
|
alpha |
Number between 0 and 1 specifying voting threshold. Eg. for 3 column matrix alpha 0.5 will give voting criteria >= 2. |
colname |
character giving new column name. |
dgCMatrix
Stouffer's Z-score method is a meta-analysis technique used to combine the results from independent statistical tests with the same hypothesis. It is closely related to Fisher's method, but operates on Z-scores instead of p-values (Wikipedia article).
stoufferZMethod(z)
stoufferZMethod(z)
z |
a numeric vector of Z-score to combine. |
a number giving combined Z-score.
Subset matrix keeping unmatched rows as NA.
subsetWithMissing(mat, rows)
subsetWithMissing(mat, rows)
mat |
matrix |
rows |
character |
a matrix
translateCounts
renames counts matrix rownames according to supplied
dict
ionary. Function can handle many to one assignments by taking a
sum or an average over counts
rows. Other types of ambiguous
assignments are not supported.
translateCounts(counts, dict)
translateCounts(counts, dict)
counts |
matrix of expression values. |
dict |
named character vector mapping |
matrix of expression values with new rownames.
counts <- matrix( data = c(5, 4, 3, 2), nrow = 2, dimnames = list( c("ENSG00000130700", "ENSG00000089225"), c("treatment", "control") ) ) dict <- c(ENSG00000130700 = "GATA5", ENSG00000089225 = "TBX5") translateCounts(counts, dict)
counts <- matrix( data = c(5, 4, 3, 2), nrow = 2, dimnames = list( c("ENSG00000130700", "ENSG00000089225"), c("treatment", "control") ) ) dict <- c(ENSG00000130700 = "GATA5", ENSG00000089225 = "TBX5") translateCounts(counts, dict)