This function imputes assays values inside a MultiAssayExperiment
Source: R/imputeAssay.R
imputeAssay.Rd
These function allow the user to enter a
MultiAssayExperiment
and impute all the NA values inside assays.
Arguments
- multiassayexperiment
A
MultiAssayExperiment
with genes in the rows, samples in the columns- i
A numeric, logical, or character
vector
indicating the assays to perform imputation on (default 1L)- ...
Arguments passed on to
impute::impute.knn
data
An expression matrix with genes in the rows, samples in the columns
k
Number of neighbors to be used in the imputation (default=10)
rowmax
The maximum percent missing data allowed in any row (default 50%). For any rows with more than
rowmax
% missing are imputed using the overall mean per sample.colmax
The maximum percent missing data allowed in any column (default 80%). If any column has more than
colmax
% missing data, the program halts and reports an error.maxp
The largest block of genes imputed using the knn algorithm inside
impute.knn
(default 1500); larger blocks are divided by two-means clustering (recursively) prior to imputation. Ifmaxp=p
, only knn imputation is done.rng.seed
The seed used for the random number generator (default 362436069) for reproducibility.
Examples
example(getSubtypeMap)
#>
#> gtSbtM> library(curatedTCGAData)
#>
#> gtSbtM> gbm <- curatedTCGAData("GBM", c("RPPA*", "CNA*"), version = "2.0.1", FALSE)
#> Querying and downloading: GBM_CNACGH_CGH_hg_244a-20160128
#> see ?curatedTCGAData and browseVignettes('curatedTCGAData') for documentation
#> loading from cache
#> Querying and downloading: GBM_CNACGH_CGH_hg_415k_g4124a-20160128
#> see ?curatedTCGAData and browseVignettes('curatedTCGAData') for documentation
#> loading from cache
#> Querying and downloading: GBM_CNASNP-20160128
#> see ?curatedTCGAData and browseVignettes('curatedTCGAData') for documentation
#> loading from cache
#> Querying and downloading: GBM_RPPAArray-20160128
#> see ?curatedTCGAData and browseVignettes('curatedTCGAData') for documentation
#> loading from cache
#> Querying and downloading: GBM_colData-20160128
#> see ?curatedTCGAData and browseVignettes('curatedTCGAData') for documentation
#> loading from cache
#> Querying and downloading: GBM_metadata-20160128
#> see ?curatedTCGAData and browseVignettes('curatedTCGAData') for documentation
#> loading from cache
#> Querying and downloading: GBM_sampleMap-20160128
#> see ?curatedTCGAData and browseVignettes('curatedTCGAData') for documentation
#> loading from cache
#> harmonizing input:
#> removing 5922 sampleMap rows not in names(experiments)
#>
#> gtSbtM> getSubtypeMap(gbm)
#> GBM_annotations GBM_subtype
#> 1 Patient_ID Case
#> 2 methylation_subtypes MGMT promoter status
#> 3 mutation_subtypes IDH/codel subtype
#> 4 histological_subtypes Histology
#> 5 mrna_subtypes Original Subtype
#> 6 mrna_subtypes Transcriptome Subtype
#> 7 mrna_subtypes Pan-Glioma RNA Expression Cluster
#> 8 mrna_subtypes IDH-specific RNA Expression Cluster
#> 9 methylation_subtypes Pan-Glioma DNA Methylation Cluster
#> 10 methylation_subtypes IDH-specific DNA Methylation Cluster
#> 11 methylation_subtypes Supervised DNA Methylation Cluster
#> 12 methylation_subtypes Random Forest Sturm Cluster
#> 13 protein_subtypes RPPA cluster
#>
#> gtSbtM> sampleTables(gbm)
#> $`GBM_CNACGH_CGH_hg_244a-20160128`
#>
#> 01 10 11
#> 267 145 26
#>
#> $`GBM_CNACGH_CGH_hg_415k_g4124a-20160128`
#>
#> 01 10
#> 169 169
#>
#> $`GBM_CNASNP-20160128`
#>
#> 01 02 10 11
#> 577 13 488 26
#>
#> $`GBM_RPPAArray-20160128`
#>
#> 01 02
#> 233 11
#>
#>
#> gtSbtM> TCGAsplitAssays(gbm, c("01", "10"))
#> Warning: Some 'sampleCodes' not found in assays
#> Warning: Inconsistent barcode lengths: 28, 27
#> A MultiAssayExperiment object of 7 listed
#> experiments with user-defined names and respective classes.
#> Containing an ExperimentList class object of length 7:
#> [1] 01_GBM_CNACGH_CGH_hg_244a-20160128: RaggedExperiment with 81512 rows and 267 columns
#> [2] 10_GBM_CNACGH_CGH_hg_244a-20160128: RaggedExperiment with 81512 rows and 145 columns
#> [3] 01_GBM_CNACGH_CGH_hg_415k_g4124a-20160128: RaggedExperiment with 57975 rows and 169 columns
#> [4] 10_GBM_CNACGH_CGH_hg_415k_g4124a-20160128: RaggedExperiment with 57975 rows and 169 columns
#> [5] 01_GBM_CNASNP-20160128: RaggedExperiment with 602338 rows and 577 columns
#> [6] 10_GBM_CNASNP-20160128: RaggedExperiment with 602338 rows and 488 columns
#> [7] 01_GBM_RPPAArray-20160128: SummarizedExperiment with 208 rows and 233 columns
#> Functionality:
#> experiments() - obtain the ExperimentList instance
#> colData() - the primary/phenotype DataFrame
#> sampleMap() - the sample coordination DataFrame
#> `$`, `[`, `[[` - extract colData columns, subset, or experiment
#> *Format() - convert into a long or wide DataFrame
#> assays() - convert ExperimentList to a SimpleList of matrices
#> exportClass() - save data to flat files
#>
#> gtSbtM> getClinicalNames("COAD")
#> [1] "years_to_birth"
#> [2] "vital_status"
#> [3] "days_to_death"
#> [4] "days_to_last_followup"
#> [5] "tumor_tissue_site"
#> [6] "pathologic_stage"
#> [7] "pathology_T_stage"
#> [8] "pathology_N_stage"
#> [9] "pathology_M_stage"
#> [10] "gender"
#> [11] "date_of_initial_pathologic_diagnosis"
#> [12] "days_to_last_known_alive"
#> [13] "radiation_therapy"
#> [14] "histological_type"
#> [15] "residual_tumor"
#> [16] "number_of_lymph_nodes"
#> [17] "race"
#> [18] "ethnicity"
## convert data to matrix and add as experiment
gbm <-
c(gbm, RPPA_matrix = data.matrix(assay(gbm[["GBM_RPPAArray-20160128"]])))
imputeAssay(gbm, i = "RPPA_matrix")
#> Warning: 'experiments' dropped; see 'drops()'
#> harmonizing input:
#> removing 2124 sampleMap rows not in names(experiments)
#> removing 361 colData rownames not in sampleMap 'primary'
#> Warning: 5 rows with more than 50 % entries missing;
#> mean imputation used for these rows
#> A MultiAssayExperiment object of 5 listed
#> experiments with user-defined names and respective classes.
#> Containing an ExperimentList class object of length 5:
#> [1] GBM_CNACGH_CGH_hg_244a-20160128: RaggedExperiment with 81512 rows and 438 columns
#> [2] GBM_CNACGH_CGH_hg_415k_g4124a-20160128: RaggedExperiment with 57975 rows and 338 columns
#> [3] GBM_CNASNP-20160128: RaggedExperiment with 602338 rows and 1104 columns
#> [4] GBM_RPPAArray-20160128: SummarizedExperiment with 208 rows and 244 columns
#> [5] RPPA_matrix: matrix with 208 rows and 244 columns
#> Functionality:
#> experiments() - obtain the ExperimentList instance
#> colData() - the primary/phenotype DataFrame
#> sampleMap() - the sample coordination DataFrame
#> `$`, `[`, `[[` - extract colData columns, subset, or experiment
#> *Format() - convert into a long or wide DataFrame
#> assays() - convert ExperimentList to a SimpleList of matrices
#> exportClass() - save data to flat files