Title: | Advanced package of tensor decomposition based unsupervised feature extraction |
---|---|
Description: | This is an advanced version of TDbasedUFE, which is a comprehensive package to perform Tensor decomposition based unsupervised feature extraction. In contrast to TDbasedUFE which can perform simple the feature selection and the multiomics analyses, this package can perform more complicated and advanced features, but they are not so popularly required. Only users who require more specific features can make use of its functionality. |
Authors: | Y-h. Taguchi [aut, cre] |
Maintainer: | Y-h. Taguchi <[email protected]> |
License: | GPL-3 |
Version: | 1.7.0 |
Built: | 2024-11-19 04:32:37 UTC |
Source: | https://github.com/bioc/TDbasedUFEadv |
Title Perform SVD toward reduced matrix generated from a tensor with partial summation
computeSVD(matrix1, matrix2, dim = 10L, scale = TRUE)
computeSVD(matrix1, matrix2, dim = 10L, scale = TRUE)
matrix1 |
The first original matrix that generates a tensor |
matrix2 |
The second original matrix that generates a tensor |
dim |
The number of singular value vectors to be computed |
scale |
If matrix should be scaled or not |
Singular value vectors attributed to two sets of objects associated with singular value vectors attributed to features, by multiplying
matrix1 <- matrix(runif(200),20) matrix2 <- matrix(runif(400),20) SVD <- computeSVD(matrix1,matrix2)
matrix1 <- matrix(runif(200),20) matrix2 <- matrix(runif(400),20) SVD <- computeSVD(matrix1,matrix2)
Prepare condition matrix for expDrug
prepareCondDrugandDisease(expDrug)
prepareCondDrugandDisease(expDrug)
expDrug |
input gene expression profile |
Condition matrix for expDrug
library(RTCGA.rnaseq) Cancer_cell_lines <- list(ACC.rnaseq,BLCA.rnaseq,BRCA.rnaseq) Drug_and_Disease <- prepareexpDrugandDisease(Cancer_cell_lines) Cond <- prepareCondDrugandDisease(Drug_and_Disease$expDrug)
library(RTCGA.rnaseq) Cancer_cell_lines <- list(ACC.rnaseq,BLCA.rnaseq,BRCA.rnaseq) Drug_and_Disease <- prepareexpDrugandDisease(Cancer_cell_lines) Cond <- prepareCondDrugandDisease(Drug_and_Disease$expDrug)
Prepare Sample label for TCGA data
prepareCondTCGA( Multi_sample, Clinical, ID_column_of_Multi_sample, ID_column_of_Clinical )
prepareCondTCGA( Multi_sample, Clinical, ID_column_of_Multi_sample, ID_column_of_Clinical )
Multi_sample |
list of sample ids |
Clinical |
List of clinical data matrix from RTCGA.clinical |
ID_column_of_Multi_sample |
Column numbers used for conditions |
ID_column_of_Clinical |
Column numbers that include corresponding sample ids in clinical data |
list of sample labels
library(RTCGA.clinical) library(RTCGA.rnaseq) Clinical <- list(BLCA.clinical, BRCA.clinical, CESC.clinical, COAD.clinical) Multi_sample <- list( BLCA.rnaseq[seq_len(100), 1, drop = FALSE], BRCA.rnaseq[seq_len(100), 1, drop = FALSE], CESC.rnaseq[seq_len(100), 1, drop = FALSE], COAD.rnaseq[seq_len(100), 1, drop = FALSE] ) ID_column_of_Multi_sample <- c(770, 1482, 773, 791) ID_column_of_Clinical <- c(20, 20, 12, 14) cond <- prepareCondTCGA( Multi_sample, Clinical, ID_column_of_Multi_sample, ID_column_of_Clinical )
library(RTCGA.clinical) library(RTCGA.rnaseq) Clinical <- list(BLCA.clinical, BRCA.clinical, CESC.clinical, COAD.clinical) Multi_sample <- list( BLCA.rnaseq[seq_len(100), 1, drop = FALSE], BRCA.rnaseq[seq_len(100), 1, drop = FALSE], CESC.rnaseq[seq_len(100), 1, drop = FALSE], COAD.rnaseq[seq_len(100), 1, drop = FALSE] ) ID_column_of_Multi_sample <- c(770, 1482, 773, 791) ID_column_of_Clinical <- c(20, 20, 12, 14) cond <- prepareCondTCGA( Multi_sample, Clinical, ID_column_of_Multi_sample, ID_column_of_Clinical )
Generating gene expression of drug treated cell lines and a disease cell line
prepareexpDrugandDisease(Cancer_cell_lines)
prepareexpDrugandDisease(Cancer_cell_lines)
Cancer_cell_lines |
<- list(ACC.rnaseq,BLCA.rnaseq,BRCA.rnaseq) list that includes individual data set from RTCGA.rnaseq |
list of expDrug and expDisease
library(RTCGA.rnaseq) Cancer_cell_lines <- list(ACC.rnaseq,BLCA.rnaseq,BRCA.rnaseq) Drug_and_Disease <- prepareexpDrugandDisease(Cancer_cell_lines)
library(RTCGA.rnaseq) Cancer_cell_lines <- list(ACC.rnaseq,BLCA.rnaseq,BRCA.rnaseq) Drug_and_Disease <- prepareexpDrugandDisease(Cancer_cell_lines)
Prepare tensor from a list that includes multiple profiles
prepareTensorfromList(Multi, proj_dim)
prepareTensorfromList(Multi, proj_dim)
Multi |
a list that includes multiple profiles |
proj_dim |
the number of projection dimensions |
a tensor as a bundle of singular value vectors obtained by applying SVD to individual omics
library(MOFAdata) data("CLL_data") data("CLL_covariates") Z <- prepareTensorfromList(CLL_data,10L)
library(MOFAdata) data("CLL_data") data("CLL_covariates") Z <- prepareTensorfromList(CLL_data,10L)
Generate tensor from two matrices
prepareTensorfromMatrix(matrix1, matrix2)
prepareTensorfromMatrix(matrix1, matrix2)
matrix1 |
the first input matrix |
matrix2 |
the second input matrix |
A tensor generated from the first and second matrices
Z <- prepareTensorfromMatrix(matrix(runif(100),10),matrix(runif(100),10))
Z <- prepareTensorfromMatrix(matrix(runif(100),10),matrix(runif(100),10))
Prepare tensor generated from two matrices that share samples
prepareTensorRect( sample, feature, value, featureRange = GRanges(NULL), sampleData = list(NULL) )
prepareTensorRect( sample, feature, value, featureRange = GRanges(NULL), sampleData = list(NULL) )
sample |
Character vector of sample names |
feature |
list of features from two matrices |
value |
array, contents of |
featureRange |
Genomic Ranges to be associated with features |
sampleData |
List of conditional labeling associated with samples |
Tensor generated from two matrices that share samples
matrix1 <- matrix(runif(1000),200) #row features, column samples matrix2 <- matrix(runif(2000),400) #row features, column samples Z <- prepareTensorfromMatrix(t(matrix1),t(matrix2)) Z <- prepareTensorRect(sample=as.character(seq_len(50)), feature=list(as.character(seq_len(200)),as.character(seq_len(400))), sampleData=list(rep(seq_len(2),each=25)),value=Z)
matrix1 <- matrix(runif(1000),200) #row features, column samples matrix2 <- matrix(runif(2000),400) #row features, column samples Z <- prepareTensorfromMatrix(t(matrix1),t(matrix2)) Z <- prepareTensorRect(sample=as.character(seq_len(50)), feature=list(as.character(seq_len(200)),as.character(seq_len(400))), sampleData=list(rep(seq_len(2),each=25)),value=Z)
Select feature when projection strategy is employed for the case where features are shared with multiple omics profiles
selectFeatureProj( HOSVD, Multi, cond, de = 1e-04, p0 = 0.01, breaks = 100L, input_all = NULL )
selectFeatureProj( HOSVD, Multi, cond, de = 1e-04, p0 = 0.01, breaks = 100L, input_all = NULL )
HOSVD |
HOSVD |
Multi |
list of omics profiles, row: sample, column: feature |
cond |
list of conditions for individual omics profiles |
de |
initial value for optimization of standard deviation |
p0 |
Threshold P-value |
breaks |
The number of bins of histogram of P-values |
input_all |
The number of selected feature. if null, interactive mode is activated |
list composed of logical vector that represent which features are selected and p-values
library(TDbasedUFE) Multi <- list(matrix(runif(1000),10),matrix(runif(1000),10), matrix(runif(1000),10),matrix(runif(1000),10)) Z <- prepareTensorfromList(Multi,10L) Z <- aperm(Z,c(2,1,3)) Z <- PrepareSummarizedExperimentTensor(feature =as.character(1:10), sample=array("",1),value=Z) HOSVD <- computeHosvd(Z) cond <- rep(list(rep(1:2,each=5)),4) index <- selectFeatureProj(HOSVD,Multi,cond,de=0.1,input_all=2)
library(TDbasedUFE) Multi <- list(matrix(runif(1000),10),matrix(runif(1000),10), matrix(runif(1000),10),matrix(runif(1000),10)) Z <- prepareTensorfromList(Multi,10L) Z <- aperm(Z,c(2,1,3)) Z <- PrepareSummarizedExperimentTensor(feature =as.character(1:10), sample=array("",1),value=Z) HOSVD <- computeHosvd(Z) cond <- rep(list(rep(1:2,each=5)),4) index <- selectFeatureProj(HOSVD,Multi,cond,de=0.1,input_all=2)
Select features through the selection of singular value vectors
selectFeatureRect( SVD, cond, de = rep(1e-04, 2), p0 = 0.01, breaks = 100L, input_all = NULL )
selectFeatureRect( SVD, cond, de = rep(1e-04, 2), p0 = 0.01, breaks = 100L, input_all = NULL )
SVD |
SVD computed from matrix generated by partial summation of a tensor |
cond |
Condition to select singular value vectors |
de |
Initial values to be used for optimization of standard deviation |
p0 |
Threshold value for the significance |
breaks |
Number of bins of histogram of P-values |
input_all |
The ID of selected singular value vectors. If it is null, interactive mode is activated. |
List of lists that includes P-vales as well as if individual features selected.
set.seed(0) matrix1 <- matrix(runif(2000),200) matrix2 <- matrix(runif(4000),200) SVD <- computeSVD(matrix1,matrix2) index_all <- selectFeatureRect(SVD, list(NULL,rep(seq_len(2),each=5),rep(seq_len(2),each=10)),de=rep(0.5,2), input_all=1)
set.seed(0) matrix1 <- matrix(runif(2000),200) matrix2 <- matrix(runif(4000),200) SVD <- computeSVD(matrix1,matrix2) index_all <- selectFeatureRect(SVD, list(NULL,rep(seq_len(2),each=5),rep(seq_len(2),each=10)),de=rep(0.5,2), input_all=1)
Select features for a tensor generated from two matrices that share samples.
selectFeatureTransRect( HOSVD, cond, de = rep(1e-04, 2), p0 = 0.01, breaks = 100L, input_all = NULL )
selectFeatureTransRect( HOSVD, cond, de = rep(1e-04, 2), p0 = 0.01, breaks = 100L, input_all = NULL )
HOSVD |
HOSVD |
cond |
list of conditions |
de |
initial values for optimization of standard deviation |
p0 |
threshold value for the significance |
breaks |
number of bins of the histogram of P-values |
input_all |
The selected singular value vectors attributed to samples. if NULL, interactive mode |
list of logical vector that represent if the individual features are selected and P-values.
library(TDbasedUFE) set.seed(0) matrix1 <- matrix(runif(1000),20) #row features, column samples matrix2 <- matrix(runif(2000),40) #row features, column samples Z <- prepareTensorfromMatrix(t(matrix1),t(matrix2)) Z <- prepareTensorRect(sample=as.character(seq_len(50)), feature=list(as.character(seq_len(20)),as.character(seq_len(40))), sampleData=list(rep(seq_len(2),each=25)),value=Z) HOSVD <- computeHosvd(Z) cond <- list(attr(Z,"sampleData")[[1]],NULL,NULL) index_all <- selectFeatureTransRect(HOSVD,cond,de=c(0.1,0.1), input_all=2,p0=1e-10)
library(TDbasedUFE) set.seed(0) matrix1 <- matrix(runif(1000),20) #row features, column samples matrix2 <- matrix(runif(2000),40) #row features, column samples Z <- prepareTensorfromMatrix(t(matrix1),t(matrix2)) Z <- prepareTensorRect(sample=as.character(seq_len(50)), feature=list(as.character(seq_len(20)),as.character(seq_len(40))), sampleData=list(rep(seq_len(2),each=25)),value=Z) HOSVD <- computeHosvd(Z) cond <- list(attr(Z,"sampleData")[[1]],NULL,NULL) index_all <- selectFeatureTransRect(HOSVD,cond,de=c(0.1,0.1), input_all=2,p0=1e-10)
Class definitions
sample
character.
feature
list.
value
array.
featureRange
GRanges.
sampleData
list.
Convert SVD to that for the case where samples are shared between two matrices
transSVD(SVD)
transSVD(SVD)
SVD |
input SVD object generated from computeSVD function |
converted SVD objects
matrix1 <- matrix(runif(200),20) matrix2 <- matrix(runif(400),20) SVD <- computeSVD(matrix1,matrix2) SVD <- transSVD(SVD)
matrix1 <- matrix(runif(200),20) matrix2 <- matrix(runif(400),20) SVD <- computeSVD(matrix1,matrix2) SVD <- transSVD(SVD)