Title: | Tensor Decomposition Based Unsupervised Feature Extraction |
---|---|
Description: | This is a comprehensive package to perform Tensor decomposition based unsupervised feature extraction. It can perform unsupervised feature extraction. It uses tensor decomposition. It is applicable to gene expression, DNA methylation, and histone modification etc. It can perform multiomics analysis. It is also potentially applicable to single cell omics data sets. |
Authors: | Y-h. Taguchi [aut, cre] |
Maintainer: | Y-h. Taguchi <[email protected]> |
License: | GPL-3 |
Version: | 1.7.0 |
Built: | 2024-12-19 04:08:56 UTC |
Source: | https://github.com/bioc/TDbasedUFE |
Title Compute higher order singular value decomposition
computeHosvd(Z, dims = c(10, dim(attr(Z, "value"))[-1]), scale = TRUE)
computeHosvd(Z, dims = c(10, dim(attr(Z, "value"))[-1]), scale = TRUE)
Z |
array that includes omics data |
dims |
dimensions to be computed by HOSVD |
scale |
If value is scaled |
List that includes output from HOSVD
Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(6)),c(3,2)), feature=as.character(seq_len(10)), value=array(runif(10*3*2),c(10,3,2))) HOSVD <- computeHosvd(Z)
Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(6)),c(3,2)), feature=as.character(seq_len(10)), value=array(runif(10*3*2),c(10,3,2))) HOSVD <- computeHosvd(Z)
Title Compute higher order singular value decomposition from the tensor generated from squared matrix
computeHosvdSqure( Z, dims = unlist(lapply(dim(attr(Z, "value")), function(x) { min(10, x) })), scale = TRUE )
computeHosvdSqure( Z, dims = unlist(lapply(dim(attr(Z, "value")), function(x) { min(10, x) })), scale = TRUE )
Z |
A tensor including sample names, feature values, associated with featureRange and sample properties |
dims |
dimensions to be computed by HOSVD |
scale |
If value is scaled |
List that includes output from HOSVD
omics1 <- matrix(runif(100),10) dimnames(omics1) <- list(seq_len(10),seq_len(10)) omics2 <- matrix(runif(100),10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- PrepareSummarizedExperimentTensorSquare( sample=matrix(colnames(omics1),1), feature=list(omics1=rownames(omics1), omics2=rownames(omics2)), value=convertSquare(Multi), sampleData=list(NA)) HOSVD <- computeHosvdSqure(Z)
omics1 <- matrix(runif(100),10) dimnames(omics1) <- list(seq_len(10),seq_len(10)) omics2 <- matrix(runif(100),10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- PrepareSummarizedExperimentTensorSquare( sample=matrix(colnames(omics1),1), feature=list(omics1=rownames(omics1), omics2=rownames(omics2)), value=convertSquare(Multi), sampleData=list(NA)) HOSVD <- computeHosvdSqure(Z)
Generate squared tensor from multiomics data
convertSquare(Multi)
convertSquare(Multi)
Multi |
A list that include multiomics data |
A tensor computed from multiomics data
omics1 <- matrix(runif(100),10) dimnames(omics1) <- list(seq_len(10),seq_len(10)) omics2 <- matrix(runif(100),10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- convertSquare(Multi)
omics1 <- matrix(runif(100),10) dimnames(omics1) <- list(seq_len(10),seq_len(10)) omics2 <- matrix(runif(100),10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- convertSquare(Multi)
Title Generate feature values formatted as a tensor format
PrepareSummarizedExperimentTensor( sample, feature, value, featureRange = GRanges(NULL), sampleData = list(NULL) )
PrepareSummarizedExperimentTensor( sample, feature, value, featureRange = GRanges(NULL), sampleData = list(NULL) )
sample |
Sample names |
feature |
Feature id names |
value |
Feature values |
featureRange |
Genomic coordinate attributed to feature id (if any) |
sampleData |
Sample property (labels etc) |
A tensor including sample names, feature id, feature values, associated with featureRange and sample properties
require(GenomicRanges) Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(6)),c(3,2)), feature=as.character(seq_len(10)), value=array(runif(10*3*2),c(10,3,2)))
require(GenomicRanges) Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(6)),c(3,2)), feature=as.character(seq_len(10)), value=array(runif(10*3*2),c(10,3,2)))
Title Generate feature values formatted as a tensor format from Squared matrix
PrepareSummarizedExperimentTensorSquare( sample = list(NULL), feature, value, featureRange = GRanges(NULL), sampleData = list(NULL) )
PrepareSummarizedExperimentTensorSquare( sample = list(NULL), feature, value, featureRange = GRanges(NULL), sampleData = list(NULL) )
sample |
Sample names |
feature |
Feature id names |
value |
Squared Feature values |
featureRange |
Genomic coordinate attributed to feature id (if any) |
sampleData |
Sample property (labels etc) |
A tensor including sample names, feature values, associated with featureRange and sample properties
omics1 <- matrix(runif(100),10) dimnames(omics1) <- list(seq_len(10),seq_len(10)) omics2 <- matrix(runif(100),10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- PrepareSummarizedExperimentTensorSquare( sample=matrix(colnames(omics1),1), feature=list(omics1=rownames(omics1), omics2=rownames(omics2)), value=convertSquare(Multi), sampleData=list(NA))
omics1 <- matrix(runif(100),10) dimnames(omics1) <- list(seq_len(10),seq_len(10)) omics2 <- matrix(runif(100),10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- PrepareSummarizedExperimentTensorSquare( sample=matrix(colnames(omics1),1), feature=list(omics1=rownames(omics1), omics2=rownames(omics2)), value=convertSquare(Multi), sampleData=list(NA))
Title Select features
selectFeature(HOSVD, input_all, de = 1e-04, p0 = 0.01, breaks = 100)
selectFeature(HOSVD, input_all, de = 1e-04, p0 = 0.01, breaks = 100)
HOSVD |
output from HOSVD |
input_all |
Selected singular value IDs |
de |
Initial value for optimization of standard deviation |
p0 |
Threshold P-value |
breaks |
The number of bins |
List that includes selected features and computed P-value
set.seed(2) require(rTensor) HOSVD <- hosvd(as.tensor(array(runif(10000*3*3),c(10000,3,3))),c(10,3,3)) input_all <- c(2,2) index <- selectFeature(HOSVD,input_all,de=0.01,p0=0.01)
set.seed(2) require(rTensor) HOSVD <- hosvd(as.tensor(array(runif(10000*3*3),c(10000,3,3))),c(10,3,3)) input_all <- c(2,2) index <- selectFeature(HOSVD,input_all,de=0.01,p0=0.01)
Title Select features (for tensor generated from squared matrix)
selectFeatureSquare( HOSVD, input_all, Multi, de = rep(1e-04, dim(HOSVD$U[[3]])[2]), p0 = 0.01, breaks = 100, interact = TRUE )
selectFeatureSquare( HOSVD, input_all, Multi, de = rep(1e-04, dim(HOSVD$U[[3]])[2]), p0 = 0.01, breaks = 100, interact = TRUE )
HOSVD |
output from HOSVD applied to tensor generated from squared matrix |
input_all |
Selected singular value vector IDs |
Multi |
Multiomics data |
de |
Initial value for optimization of standard deviation |
p0 |
Threshold P-value |
breaks |
The number of bins |
interact |
if interact mode or not |
List that includes selected features and computed P-value
omics1 <- matrix(runif(100000),ncol=10) dimnames(omics1) <- list(seq_len(10000),seq_len(10)) omics2 <- matrix(runif(100000),ncol=10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- PrepareSummarizedExperimentTensorSquare( sample=matrix(colnames(omics1),1), feature=list(omics1=rownames(omics1), omics2=rownames(omics2)), value=convertSquare(Multi), sampleData=list(NA)) HOSVD <- computeHosvdSqure(Z) cond <- list(0,rep(seq_len(2),each=5),c("A","B")) input_all <- selectSingularValueVectorLarge(HOSVD,cond,input_all=c(1,1)) index <- selectFeatureSquare(HOSVD,input_all,Multi,de=c(0.1,0.1), interact=FALSE)
omics1 <- matrix(runif(100000),ncol=10) dimnames(omics1) <- list(seq_len(10000),seq_len(10)) omics2 <- matrix(runif(100000),ncol=10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- PrepareSummarizedExperimentTensorSquare( sample=matrix(colnames(omics1),1), feature=list(omics1=rownames(omics1), omics2=rownames(omics2)), value=convertSquare(Multi), sampleData=list(NA)) HOSVD <- computeHosvdSqure(Z) cond <- list(0,rep(seq_len(2),each=5),c("A","B")) input_all <- selectSingularValueVectorLarge(HOSVD,cond,input_all=c(1,1)) index <- selectFeatureSquare(HOSVD,input_all,Multi,de=c(0.1,0.1), interact=FALSE)
Title Select singular value vectors from HOSVD (boxplot version)
selectSingularValueVectorLarge(HOSVD, cond, input_all = NULL)
selectSingularValueVectorLarge(HOSVD, cond, input_all = NULL)
HOSVD |
output from HOSVD |
cond |
Labels to select singular value vector number |
input_all |
if list is not null, no interactive mode is activated but provided values are used. |
Selected singular value vector IDs
Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(6)),c(3,2)), feature=as.character(seq_len(10)), value=array(runif(10*3*2),c(10,3,2))) HOSVD <- computeHosvd(Z) cond <- list(0,c("A","B","C"),c("A","B")) input_all <- selectSingularValueVectorLarge(HOSVD,cond,input_all=c(1,1))
Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(6)),c(3,2)), feature=as.character(seq_len(10)), value=array(runif(10*3*2),c(10,3,2))) HOSVD <- computeHosvd(Z) cond <- list(0,c("A","B","C"),c("A","B")) input_all <- selectSingularValueVectorLarge(HOSVD,cond,input_all=c(1,1))
Title Select singular value vectors from HOSVD
selectSingularValueVectorSmall(HOSVD, input_all = NULL)
selectSingularValueVectorSmall(HOSVD, input_all = NULL)
HOSVD |
output from HOSVD |
input_all |
if ist is no null, no interactive mode is activated but provided values are used. |
Selected singular value vector IDs
Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(6)),c(3,2)), feature=as.character(seq_len(10)), value=array(runif(10*3*2),c(10,3,2))) HOSVD <- computeHosvd(Z) input_all <- selectSingularValueVectorSmall(HOSVD,input_all=c(1,1))
Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(6)),c(3,2)), feature=as.character(seq_len(10)), value=array(runif(10*3*2),c(10,3,2))) HOSVD <- computeHosvd(Z) input_all <- selectSingularValueVectorSmall(HOSVD,input_all=c(1,1))
Title Show selected features as Table
tableFeatures(Z, index)
tableFeatures(Z, index)
Z |
Tensor of features |
index |
List that includes selected features and P-values |
Table list of selected features
set.seed(2) require(rTensor) HOSVD <- hosvd(as.tensor(array(runif(10000*3*3),c(10000,3,3))),c(10,3,3)) input_all <- c(2,2) index <- selectFeature(HOSVD,input_all,de=0.01,p0=0.01) index$index[seq_len(100)] <- TRUE Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(9)),c(3,3)), feature=as.character(seq_len(10000)), value=array(runif(10000*3*3),c(10,3,3))) head(tableFeatures(Z,index))
set.seed(2) require(rTensor) HOSVD <- hosvd(as.tensor(array(runif(10000*3*3),c(10000,3,3))),c(10,3,3)) input_all <- c(2,2) index <- selectFeature(HOSVD,input_all,de=0.01,p0=0.01) index$index[seq_len(100)] <- TRUE Z <- PrepareSummarizedExperimentTensor( sample=matrix(as.character(seq_len(9)),c(3,3)), feature=as.character(seq_len(10000)), value=array(runif(10000*3*3),c(10,3,3))) head(tableFeatures(Z,index))
Title Show selected features as Table (for Squared one)
tableFeaturesSquare(Z, index, id)
tableFeaturesSquare(Z, index, id)
Z |
Tensor of features |
index |
List that includes selected features and P-values |
id |
feature to be shown |
Table list of selected features
omics1 <- matrix(runif(100000),ncol=10) dimnames(omics1) <- list(seq_len(10000),seq_len(10)) omics2 <- matrix(runif(100000),ncol=10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- PrepareSummarizedExperimentTensorSquare( sample=matrix(colnames(omics1),1), feature=list(omics1=rownames(omics1), omics2=rownames(omics2)), value=convertSquare(Multi), sampleData=list(NA)) HOSVD <- computeHosvdSqure(Z) cond <- list(0,rep(seq_len(2),each=5),c("A","B")) input_all <- selectSingularValueVectorLarge(HOSVD,cond,input_all=c(1,1)) index <- selectFeatureSquare(HOSVD,input_all,Multi,de=c(0.1,0.1), interact=FALSE) index[[1]]$index[1:100]<-TRUE index[[1]]$p.value[1:100] <- 1e-3 tableFeaturesSquare(Z,index,1)
omics1 <- matrix(runif(100000),ncol=10) dimnames(omics1) <- list(seq_len(10000),seq_len(10)) omics2 <- matrix(runif(100000),ncol=10) dimnames(omics2) <- dimnames(omics1) Multi <- list(omics1,omics2) Z <- PrepareSummarizedExperimentTensorSquare( sample=matrix(colnames(omics1),1), feature=list(omics1=rownames(omics1), omics2=rownames(omics2)), value=convertSquare(Multi), sampleData=list(NA)) HOSVD <- computeHosvdSqure(Z) cond <- list(0,rep(seq_len(2),each=5),c("A","B")) input_all <- selectSingularValueVectorLarge(HOSVD,cond,input_all=c(1,1)) index <- selectFeatureSquare(HOSVD,input_all,Multi,de=c(0.1,0.1), interact=FALSE) index[[1]]$index[1:100]<-TRUE index[[1]]$p.value[1:100] <- 1e-3 tableFeaturesSquare(Z,index,1)