Install the released version from Bioconductor:
# Requires R >= 4.5.0
if (!require("BiocManager", quietly = TRUE)) install.packages("BiocManager")
BiocManager::install("BreastSubtypeR")
If you use BreastSubtypeR, please cite:
For BibTeX/LaTeX, run in R:
citation("BreastSubtypeR")
Breast cancer (BC) is a biologically heterogeneous disease with intrinsic molecular subtypes (e.g., Luminal A, Luminal B, HER2-enriched, Basal-like, Normal-like) that inform biological interpretation and clinical decision-making. While clinical assays such as Prosigna provide standardised subtyping in the clinic, research implementations have proliferated and diverge in pre-processing, gene mapping, and algorithmic assumptions—reducing reproducibility and complicating cross-cohort analyses.
BreastSubtypeR consolidates multiple published gene-expression signature classifiers into a unified, assumption-aware Bioconductor package with: - a unified multi-method API (run many classifiers in one call), - AUTO mode for cohort-aware method selection, - standardised, method-specific pre-processing for multiple input types (raw counts, FPKM, log2-processed arrays), - Entrez ID–based probe/gene mapping, - and a local Shiny app (iBreastSubtypeR) for non-programmers.
BS_Multi): execute several classifiers in a single call and compare results side by side.iBreastSubtypeR): point-and-click analysis; data stay on your machine.SummarizedExperiment compatibility.The package includes implementations of commonly used subtyping methods (NC-based and SSP-based):
| Method id | Short description | Group | Reference |
|---|---|---|---|
parker.original |
Original PAM50 by Parker et al., 2009 | NC-based | Parker et al., 2009 |
genefu.scale |
PAM50 implementation as in the genefu R package (scaled version) | NC-based | Gendoo et al., 2016 |
genefu.robust |
PAM50 implementation as in the genefu R package (robust version) | NC-based | Gendoo et al., 2016 |
cIHC |
Conventional ER-balancing using immunohistochemistry (IHC) | NC-based | Ciriello et al., 2015 |
cIHC.itr |
Iterative version of cIHC | NC-based | Curtis et al., 2012 |
PCAPAM50 |
Selects IHC-defined ER subsets, then uses Principal Component Analysis (PCA) to create ESR1 expression-based ER-balancing | NC-based | Raj-Kumar et al., 2019 |
ssBC |
Subgroup-specific gene-centering PAM50 | NC-based | Zhao et al., 2015 |
ssBC.v2 |
Updated subgroup-specific gene-centering PAM50 with refined quantiles | NC-based | Fernandez-Martinez et al., 2020 |
AIMS |
Absolute Intrinsic Molecular Subtyping (AIMS) method | SSP-based | Paquet & Hallett, 2015 |
sspbc |
Single-Sample Predictors for Breast Cancer (AIMS adaptation) | SSP-based | Staaf et al., 2022 |
The examples below use small example datasets shipped with the package. For your own data, provide a SummarizedExperiment with clinical metadata in colData (e.g., PatientID, ER/HER2; for ROR: TSIZE, NODE).
library(BreastSubtypeR)
# Example data
data("BreastSubtypeRobj")
data("OSLO2EMIT0obj")
1) Map & prepare (method-specific pre-processing + mapping)
# Pre-processing: automatically apply tailored normalisation, map probes/IDs to Entrez, and (optionally) impute missing values
data_input <- Mapping(
OSLO2EMIT0obj$se_obj,
method = "max", # mapping strategy (example)
RawCounts = FALSE,
impute = TRUE,
verbose = FALSE
)
Mapping() prepares expression inputs for downstream subtyping functions by:
2^x) for SSP-based methods.2^x) for SSP-based methods.method argument),BS_Multi or single-method callers.?Mapping for the full parameter list (e.g., RawCounts, method, impute, verbose) and Methods (Sections 2.3–2.4) in the paper for a complete description of the input/normalisation pipeline.2) Multi-method run (user-defined)
methods <- c("parker.original", "PCAPAM50", "sspbc")
res <- BS_Multi(
data_input = data_input,
methods = methods,
Subtype = FALSE,
hasClinical = FALSE
)
# Per-sample calls (methods × samples)
head(res$res_subtypes, 5)
3) AUTO mode (cohort-aware selection) + visualize
AUTO evaluates cohort diagnostics (for example, ER/HER2 distribution, subtype purity, and subgroup sizes) and selects methods compatible with the cohort. It disables classifiers whose distributional assumptions would likely be violated.
res_auto <- BS_Multi(
data_input = data_input,
methods = "AUTO",
Subtype = FALSE,
hasClinical = FALSE
)
# Visualise multi-method output and concordance
Vis_Multi(res_auto$res_subtypes)
4) Single-method run
PAM50 (NC-based)
res_pam <- BS_parker(
se_obj = data_input$se_NC, # object prepared for NC-based methods
calibration = "Internal",
internal = "medianCtr",
Subtype = FALSE,
hasClinical = FALSE
)
AIMS (SSP-based)
res_aims <- BS_AIMS(data_input$se_SSP)
BreastSubtypeR routes the supplied input to the appropriate, method-specific pre-processing pipeline automatically — see ?BS_Multi and Methods (Section 2.3) in the paper for details.AUTOmethods = "AUTO" (i.e. BS_Multi(methods = "AUTO", ...)) for exploratory datasets or cohorts of unknown / skewed composition.AUTO when you want the package to select only classifiers compatible with the cohort (it disables methods whose assumptions appear violated).BS_parker()).AUTO is designed to avoid misapplication of NC-based classifiers when cohort assumptions are violated; it does not produce a forced consensus label.For users new to R, we offer an intuitive Shiny app for interactive molecular subtyping.
BreastSubtypeR::iBreastSubtypeR() # interactive GUI (local)
If needed, install UI dependencies and re-run:
install.packages(c("shiny", "bslib"))
The app runs locally; no data leave your machine.
What you can do:
- Upload expression, clinical, and feature-annotation tables (clinical lives in colData). - Run single methods, or run multiple classifiers at once with BS_Multi and AUTO enabled for cohort-aware selection.
- Choose 5-class (incl. Normal-like) or 4-class (AIMS is 5-class only).
- Inspect per-sample concordance (entropy), heatmap and pie summaries.
- Export Calls-only or Full metrics. ROR is available for NC methods when TSIZE/NODE are present and numeric.
BreastSubtypeR harmonises many published, signature-based classifiers but has known limitations:
It is not a clinical-grade replacement for assays like Prosigna; clinical validation requires paired clinical assay data.
AUTO selects compatible methods; it does not perform consensus voting by default.
Yang Q., Hartman J., Sifakis E.G. (2025). BreastSubtypeR: A Unified R/Bioconductor Package for Intrinsic Molecular Subtyping in Breast Cancer Research. NAR Genomics and Bioinformatics. https://doi.org/10.1093/nargab/lqaf131
Parker JS, Mullins M, Cheang MCU, Leung S, Voduc D, Vickery T, et al. (2009). Supervised risk predictor of breast cancer based on intrinsic subtypes. J Clin Oncol, 27(8):1160–1167. https://doi.org/10.1200/JCO.2008.18.1370
Gendoo DMA, Ratanasirigulchai N, Schröder MS, Pare L, Parker JS, Prat A, Haibe-Kains B. (2016). Genefu: an R/Bioconductor package for computation of gene expression-based signatures in breast cancer. Bioinformatics, 32(7):1097–1099. https://doi.org/10.1093/bioinformatics/btv693
Ciriello G, Gatza ML, Beck AH, Wilkerson MD, Rhie SK, Pastore A, et al. (2015). Comprehensive molecular portraits of invasive lobular breast cancer. Cell, 163(2):506–519. https://doi.org/10.1016/j.cell.2015.09.033
Curtis C, Shah SP, Chin S-F, Turashvili G, Rueda OM, Dunning MJ, et al. (2012). The genomic and transcriptomic architecture of 2,000 breast tumours reveals novel subgroups. Nature, 486:346–352. https://doi.org/10.1038/nature10983
Raj-Kumar PK, Liu J, Hooke JA, Kovatich AJ, Kvecher L, Shriver CD, Hu H. (2019). PCA-PAM50 improves subtype assignment in ER-positive breast cancer. Sci Rep, 9:14386. https://doi.org/10.1038/s41598-019-44339-4
Zhao X, Rodland EA, Tibshirani R, Edvardsen H, Sauer T, Hovig E. (2015). Systematic evaluation of subtype prediction using gene expression profiles and intrinsic subtyping methods. Breast Cancer Res, 17:55. https://doi.org/10.1186/s13058-015-0520-4
Fernandez-Martinez A, Krop IE, Hillman DW, Polley M-YC, Parker JS, Huebner L, et al. (2020). Survival, pathologic response, and PAM50 subtype in stage II–III HER2-positive breast cancer treated with neoadjuvant chemotherapy and trastuzumab ± lapatinib. J Clin Oncol, 38(19):2140–2150. https://doi.org/10.1200/JCO.20.01276
Paquet ER, Hallett MT. (2015). Absolute assignment of breast cancer intrinsic molecular subtype. J Natl Cancer Inst, 107(1):357. https://doi.org/10.1093/jnci/dju357
Staaf J, Ringnér M, Vallon-Christersson J. (2022). Simple single-sample predictors for breast cancer subtype identification using gene expression data. npj Breast Cancer, 8:104. https://doi.org/10.1038/s41523-022-00465-3
sessionInfo()
#> R version 4.5.0 (2025-04-11 ucrt)
#> Platform: x86_64-w64-mingw32/x64
#> Running under: Windows 11 x64 (build 22631)
#>
#> Matrix products: default
#> LAPACK version 3.12.1
#>
#> locale:
#> [1] LC_COLLATE=C LC_CTYPE=English_Sweden.utf8
#> [3] LC_MONETARY=English_Sweden.utf8 LC_NUMERIC=C
#> [5] LC_TIME=English_Sweden.utf8
#>
#> time zone: Europe/Stockholm
#> tzcode source: internal
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] BiocStyle_2.36.0
#>
#> loaded via a namespace (and not attached):
#> [1] digest_0.6.37 R6_2.6.1 bookdown_0.44
#> [4] fastmap_1.2.0 xfun_0.53 cachem_1.1.0
#> [7] knitr_1.50 htmltools_0.5.8.1 rmarkdown_2.29
#> [10] lifecycle_1.0.4 cli_3.6.5 sass_0.4.10
#> [13] jquerylib_0.1.4 compiler_4.5.0 rstudioapi_0.17.1
#> [16] tools_4.5.0 evaluate_1.0.5 bslib_0.9.0
#> [19] yaml_2.3.10 BiocManager_1.30.26 jsonlite_2.0.0
#> [22] rlang_1.1.6