View source: R/pipeline_functions.R
update_eset.feature | R Documentation |
update_eset.feature
reassigns the featureData slot of ExpressionSet object based on user's demand. It is mainly used for gene ID conversion.
update_eset.feature(
use_eset = NULL,
use_feature_info = NULL,
from_feature = NULL,
to_feature = NULL,
merge_method = "median",
distribute_method = "equal"
)
use_eset |
ExpressionSet class object. |
use_feature_info |
data.frame, a data frame contains feature information, it can be obtained by calling |
from_feature |
character, orginal ID. Must be one of the column names in |
to_feature |
character, target ID. Must be one of the column names in |
merge_method |
character, the agglomeration method to be used for merging gene expression value. This should be one of, "median", "mean", "max" or "min". Default is "median". |
distribute_method |
character, the agglomeration method to be used for distributing the gene expression value. This should be one of, "mean" or "equal". Default is "equal". |
User can pass a conversion table to use_eset
for the ID conversion. A conversion table can be obtained from the original featureDta slot
(if one called the load.exp.GEO
function and set getGPL==TRUE) or by running the get_IDtransfer
function.
The mapping between original ID and target ID can be summerised into 4 categories.
1) One-to-one, simply replaces the original ID with target ID;
2) Many-to-one, the expression value for the target ID will be merged from its original ID;
3) One-to-many, the expression value for the original ID will be distributed to the matched target IDs;
4) Many-to-many, apply part 3) first, then part 2).
Return an ExressionSet object with updated feature information.
mat1 <- matrix(rnorm(10000),nrow=1000,ncol=10)
colnames(mat1) <- paste0('Sample',1:ncol(mat1))
rownames(mat1) <- paste0('Gene',1:nrow(mat1))
eset <- generate.eset(exp_mat=mat1)
test_transfer_table <- data.frame(
'Gene'=c('Gene1','Gene1','Gene2','Gene3','Gene4'),
'Transcript'=c('T11','T12','T2','T3','T3'))
new_eset <- update_eset.feature(use_eset=eset,
use_feature_info=test_transfer_table,
from_feature='Gene',
to_feature='Transcript',
merge_method='median',
distribute_method='equal'
)
print(Biobase::exprs(eset)[test_transfer_table$Gene,])
print(Biobase::exprs(new_eset))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.