Best practise to download and prepare TCGA-BRCA dataset for BOBaFIT’s analysis
BOBaFIT 1.6.0
The data preparation is an important step before the BOBaFIT analysis. In this vignete we will explain how to download the TCGA-BRCA dataset(Tomczak, Czerwińska, and Wiznerowicz 2015) the package TCGAbiolink (Colaprico et al. 2015) and how to add information like chromosomal arm and CN value of each segments, which operating principle of the package. Further, here we show the column names of the input file for all the BOBaFIT function.
To download the TCGCA-BRCA(Tomczak, Czerwińska, and Wiznerowicz 2015), we used the R package TCGAbiolinks (Colaprico et al. 2015)and we construct the query. The query includs Breast Cancer samples analyzed by SNParray method (GenomeWide_SNP6), obtaining their Copy Number (CN) profile.
BiocManager::install("TCGAbiolinks")
library(TCGAbiolinks)
query <- GDCquery(project = "TCGA-BRCA",
data.category = "Copy Number Variation",
data.type = "Copy Number Segment",
sample.type = "Primary Tumor"
)
#Selecting first 100 samples using the TCGA barcode
subset <- query[[1]][[1]]
barcode <- subset$cases[1:100]
TCGA_BRCA_CN_segments <- GDCquery(project = "TCGA-BRCA",
data.category = "Copy Number Variation",
data.type = "Copy Number Segment",
sample.type = "Primary Tumor",
barcode = barcode
)
GDCdownload(TCGA_BRCA_CN_segments, method = "api", files.per.chunk = 50)
#prepare a data.frame where working
data <- GDCprepare(TCGA_BRCA_CN_segments, save = TRUE,
save.filename= "TCGA_BRCA_CN_segments.txt")
In the last step, a dataframe with the segments of all samples is prepared. However some information are missing, so the dataset is not ready as BOBaFIT input.
Further, here we show the column names of the input file for all the BOBaFIT function.
names(data)
BOBaFIT_names <- c("ID", "chr", "start", "end", "Num_Probes",
"Segment_Mean","Sample")
names(data)<- BOBaFIT_names
names(data)
The arm column is an very important information that support the diploid region check ofDRrefit
and the chromosome list computation of ComputeNormalChromosome
. As it is lacking in the TCGA-BRCAdataset, the function Popeye
has been specially designed to calculate which chromosomal arm the segment belongs to. Thanks to this algorithm, not only the TCGA-BRCA dataset, but any database you want to analyze can be analyzed by any function of BOBaFIT.
library(BOBaFIT)
segments <- Popeye(data)
chr | start | end | width | strand | ID | Num_Probes | Segment_Mean | Sample | arm | chrarm |
---|---|---|---|---|---|---|---|---|---|---|
1 | 62920 | 21996664 | 21933745 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 12088 | -0.4756 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p |
1 | 22001786 | 22002025 | 240 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 2 | -7.4436 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p |
1 | 22004046 | 22010750 | 6705 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 2 | -2.1226 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p |
1 | 22011632 | 25256850 | 3245219 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 1914 | -0.4808 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p |
1 | 25266637 | 25320198 | 53562 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 22 | -2.1144 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p |
1 | 25320253 | 30316360 | 4996108 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 2434 | -0.4905 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p |
The last step is the computation of the copy number value from the “Segment_Mean
” column (logR), with the following formula. At this point the data is ready to be analyzed by the whole package.
#When data coming from SNParray platform are used, the user have to apply the
#compression factor in the formula (0.55). In case of WGS/WES data, the
#correction factor is equal to 1.
compression_factor <- 0.55
segments$CN <- 2^(segments$Segment_Mean/compression_factor + 1)
chr | start | end | width | strand | ID | Num_Probes | Segment_Mean | Sample | arm | chrarm | CN |
---|---|---|---|---|---|---|---|---|---|---|---|
1 | 62920 | 21996664 | 21933745 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 12088 | -0.4756 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p | 1.0983004 |
1 | 22001786 | 22002025 | 240 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 2 | -7.4436 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p | 0.0001686 |
1 | 22004046 | 22010750 | 6705 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 2 | -2.1226 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p | 0.1378076 |
1 | 22011632 | 25256850 | 3245219 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 1914 | -0.4808 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p | 1.0911264 |
1 | 25266637 | 25320198 | 53562 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 22 | -2.1144 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p | 0.1392391 |
1 | 25320253 | 30316360 | 4996108 | * | 01428281-1653-4839-b5cf-167bc62eb147 | 2434 | -0.4905 | TCGA-BH-A18R-01A-11D-A12A-01 | p | 1p | 1.0778690 |
## R version 4.3.1 (2023-06-16)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## time zone: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] GenomicRanges_1.54.0 IRanges_2.36.0 S4Vectors_0.40.0
## [4] dplyr_1.1.3 BOBaFIT_1.6.0 BiocStyle_2.30.0
##
## loaded via a namespace (and not attached):
## [1] RColorBrewer_1.1-3 rstudioapi_0.15.0
## [3] jsonlite_1.8.7 magrittr_2.0.3
## [5] magick_2.8.1 GenomicFeatures_1.54.0
## [7] farver_2.1.1 rmarkdown_2.25
## [9] BiocIO_1.12.0 zlibbioc_1.48.0
## [11] vctrs_0.6.4 memoise_2.0.1
## [13] Rsamtools_2.18.0 RCurl_1.98-1.12
## [15] base64enc_0.1-3 htmltools_0.5.6.1
## [17] S4Arrays_1.2.0 progress_1.2.2
## [19] curl_5.1.0 SparseArray_1.2.0
## [21] Formula_1.2-5 sass_0.4.7
## [23] bslib_0.5.1 htmlwidgets_1.6.2
## [25] plyr_1.8.9 cachem_1.0.8
## [27] GenomicAlignments_1.38.0 lifecycle_1.0.3
## [29] pkgconfig_2.0.3 Matrix_1.6-1.1
## [31] R6_2.5.1 fastmap_1.1.1
## [33] GenomeInfoDbData_1.2.11 MatrixGenerics_1.14.0
## [35] digest_0.6.33 colorspace_2.1-0
## [37] GGally_2.1.2 reshape_0.8.9
## [39] OrganismDbi_1.44.0 AnnotationDbi_1.64.0
## [41] Hmisc_5.1-1 RSQLite_2.3.1
## [43] labeling_0.4.3 filelock_1.0.2
## [45] fansi_1.0.5 polyclip_1.10-6
## [47] httr_1.4.7 abind_1.4-5
## [49] compiler_4.3.1 withr_2.5.1
## [51] bit64_4.0.5 htmlTable_2.4.1
## [53] backports_1.4.1 BiocParallel_1.36.0
## [55] DBI_1.1.3 ggforce_0.4.1
## [57] biomaRt_2.58.0 MASS_7.3-60
## [59] rappdirs_0.3.3 DelayedArray_0.28.0
## [61] rjson_0.2.21 tools_4.3.1
## [63] foreign_0.8-85 nnet_7.3-19
## [65] glue_1.6.2 restfulr_0.0.15
## [67] grid_4.3.1 checkmate_2.2.0
## [69] cluster_2.1.4 reshape2_1.4.4
## [71] generics_0.1.3 gtable_0.3.4
## [73] BSgenome_1.70.0 tidyr_1.3.0
## [75] ensembldb_2.26.0 data.table_1.14.8
## [77] hms_1.1.3 xml2_1.3.5
## [79] utf8_1.2.4 XVector_0.42.0
## [81] BiocGenerics_0.48.0 pillar_1.9.0
## [83] stringr_1.5.0 tweenr_2.0.2
## [85] BiocFileCache_2.10.0 lattice_0.22-5
## [87] rtracklayer_1.62.0 bit_4.0.5
## [89] biovizBase_1.50.0 RBGL_1.78.0
## [91] tidyselect_1.2.0 Biostrings_2.70.0
## [93] knitr_1.44 gridExtra_2.3
## [95] bookdown_0.36 ggbio_1.50.0
## [97] ProtGenerics_1.34.0 SummarizedExperiment_1.32.0
## [99] stats4_4.3.1 xfun_0.40
## [101] Biobase_2.62.0 matrixStats_1.0.0
## [103] stringi_1.7.12 lazyeval_0.2.2
## [105] yaml_2.3.7 evaluate_0.22
## [107] codetools_0.2-19 NbClust_3.0.1
## [109] tibble_3.2.1 graph_1.80.0
## [111] BiocManager_1.30.22 cli_3.6.1
## [113] rpart_4.1.21 munsell_0.5.0
## [115] jquerylib_0.1.4 dichromat_2.0-0.1
## [117] Rcpp_1.0.11 GenomeInfoDb_1.38.0
## [119] dbplyr_2.3.4 png_0.1-8
## [121] XML_3.99-0.14 parallel_4.3.1
## [123] ggplot2_3.4.4 blob_1.2.4
## [125] prettyunits_1.2.0 AnnotationFilter_1.26.0
## [127] plyranges_1.22.0 bitops_1.0-7
## [129] VariantAnnotation_1.48.0 scales_1.2.1
## [131] purrr_1.0.2 crayon_1.5.2
## [133] rlang_1.1.1 KEGGREST_1.42.0
Colaprico, Antonio, Tiago C. Silva, Catharina Olsen, Luciano Garofano, Claudia Cava, Davide Garolini, Thais S. Sabedot, et al. 2015. “TCGAbiolinks: An R/Bioconductor Package for Integrative Analysis of Tcga Data.” Nucleic Acids Research 44 (8): e71–e71. https://doi.org/10.1093/nar/gkv1507.
Tomczak, Katarzyna, Patrycja Czerwińska, and Maciej Wiznerowicz. 2015. “Review the Cancer Genome Atlas (Tcga): An Immeasurable Source of Knowledge.” Współczesna Onkologia 1A: 68–77. https://doi.org/10.5114/wo.2014.47136.