Progenetix is an open data resource that provides curated individual cancer copy number variation (CNV) profiles along with associated metadata sourced from published oncogenomic studies and various data repositories. This vignette provides a comprehensive guide on accessing genomic variant data within the Progenetix database.
If your focus lies in cancer cell lines, you can access data from cancercelllines.org by setting the domain
parameter to "https://cancercelllines.org"
in pgxLoader
function. This data repository originates from CNV profiling data of cell lines initially collected as part of Progenetix and currently includes additional types of genomic mutations.
library(pgxRpi)
library(SummarizedExperiment) # for pgxmatrix data
#> Loading required package: MatrixGenerics
#> Loading required package: matrixStats
#>
#> Attaching package: 'MatrixGenerics'
#> The following objects are masked from 'package:matrixStats':
#>
#> colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
#> colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
#> colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
#> colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
#> colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
#> colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
#> colWeightedMeans, colWeightedMedians, colWeightedSds,
#> colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
#> rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
#> rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
#> rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
#> rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
#> rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs, rowVars,
#> rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
#> rowWeightedSds, rowWeightedVars
#> Loading required package: GenomicRanges
#> Loading required package: stats4
#> Loading required package: BiocGenerics
#> Loading required package: generics
#>
#> Attaching package: 'generics'
#> The following objects are masked from 'package:base':
#>
#> as.difftime, as.factor, as.ordered, intersect, is.element, setdiff,
#> setequal, union
#>
#> Attaching package: 'BiocGenerics'
#> The following objects are masked from 'package:stats':
#>
#> IQR, mad, sd, var, xtabs
#> The following objects are masked from 'package:base':
#>
#> Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
#> as.data.frame, basename, cbind, colnames, dirname, do.call,
#> duplicated, eval, evalq, get, grep, grepl, is.unsorted, lapply,
#> mapply, match, mget, order, paste, pmax, pmax.int, pmin, pmin.int,
#> rank, rbind, rownames, sapply, saveRDS, table, tapply, unique,
#> unsplit, which.max, which.min
#> Loading required package: S4Vectors
#>
#> Attaching package: 'S4Vectors'
#> The following object is masked from 'package:utils':
#>
#> findMatches
#> The following objects are masked from 'package:base':
#>
#> I, expand.grid, unname
#> Loading required package: IRanges
#> Loading required package: GenomeInfoDb
#> Loading required package: Biobase
#> Welcome to Bioconductor
#>
#> Vignettes contain introductory material; view with
#> 'browseVignettes()'. To cite Bioconductor, see
#> 'citation("Biobase")', and for packages 'citation("pkgname")'.
#>
#> Attaching package: 'Biobase'
#> The following object is masked from 'package:MatrixGenerics':
#>
#> rowMedians
#> The following objects are masked from 'package:matrixStats':
#>
#> anyMissing, rowMedians
pgxLoader
functionThis function loads various data from Progenetix
database via the Beacon v2 API with some extensions (BeaconPlus).
The parameters of this function used in this tutorial:
type
: A string specifying output data type. "g_variants"
and "cnv_fraction"
are used in this tutorial.output
: A string specifying output file format. The available options depend on the type
parameter. When type
is "g_variants"
, available options are NULL
(default), "pgxseg"
, or "seg"
; when type
is "cnv_fraction"
, available options are NULL
(default) or "pgxmatrix"
.biosample_id
: Identifiers used in the query database for identifying biosamples.individual_id
: Identifiers used in the query database for identifying individuals.filters
: Identifiers used in public repositories, bio-ontology terms, or custom terms such as "NCIT:C7376"
. For more information about filters, see the documentation.codematches
: A logical value determining whether to exclude samples from child concepts
of specified filters in the ontology tree. If TRUE
, only samples exactly matching the specified filters will be included. Do not use this parameter when filters
include ontology-irrelevant filters such as pubmed or cohort identifiers. Default is FALSE
.limit
: Integer to specify the number of returned profiles. Default is 0
(return all).skip
: Integer to specify the number of skipped profiles. E.g. if skip = 2,limit=500
, the first 2*500=1000 profiles are skipped and the next 500 profiles are returned. Default is 0
(no skip).dataset
: A string specifying the dataset to query from the Beacon response. Default is NULL
, which includes results from all datasets.save_file
: A logical value determining whether to save variant data as a local file instead of direct return. Only used when the parameter type
is "g_variants"
. Default is FALSE
.filename
: A string specifying the path and name of the file to be saved. Only used if the parameter save_file
is TRUE
. Default is "variants.tsv"
, saved in the current working directory..num_cores
: An integer specifying the number of cores to use for parallel processing during Beacon v2 variant data queries from multiple biosamples. Default is 1
.domain
: A string specifying the domain of the query data resource. Default is "http://progenetix.org"
.entry_point
: A string specifying the entry point of the Beacon v2 API. Default is "beacon"
, resulting in the endpoint being "http://progenetix.org/beacon"
.Due to time-out constraints, variant data in Progenetix or cancercelllines can only be queried using biosample ids, rather than filters or individual ids. To speed up the process, you can use the num_cores
parameter to enable parallel processing. For more details about retrieving biosample ids, refer to the vignette Introduction_1_load_metadata.
# get 2 samples for demonstration
biosamples <- pgxLoader(type="biosamples", filters = "NCIT:C3512", limit=2)
biosample_id <- biosamples$biosample_id
biosample_id
#> [1] "pgxbs-m3io3tbz" "pgxbs-m3io6kha"
There are three output formats.
The default output format extracts variant data from the Beacon v2 response, containing variant id and associated analysis id, biosample id and individual id. The CNV data is represented as copy number change class following the GA4GH Variation Representation Specification (VRS).
variant_1 <- pgxLoader(type="g_variants", biosample_id = biosample_id)
head(variant_1)
#> variant_id analysis_id biosample_id individual_id
#> 1 pgxvar-673cbf8e725614730311c986 pgxcs-m3io3tbz pgxbs-m3io3tbz pgxind-m3io3tbz
#> 2 pgxvar-673cbf8e725614730311c987 pgxcs-m3io3tbz pgxbs-m3io3tbz pgxind-m3io3tbz
#> 3 pgxvar-673cbf8e725614730311c988 pgxcs-m3io3tbz pgxbs-m3io3tbz pgxind-m3io3tbz
#> 4 pgxvar-673cbf8e725614730311c989 pgxcs-m3io3tbz pgxbs-m3io3tbz pgxind-m3io3tbz
#> 5 pgxvar-673cbf8f725614730311c98a pgxcs-m3io3tbz pgxbs-m3io3tbz pgxind-m3io3tbz
#> 6 pgxvar-673cbf8f725614730311c98b pgxcs-m3io3tbz pgxbs-m3io3tbz pgxind-m3io3tbz
#> variant_internal_id variation_location_end
#> 1 1:120069370-149078383:EFO_0030068 149078383
#> 2 1:147810332-163328877:EFO_0030071 163328877
#> 3 1:163343802-247433974:EFO_0030071 247433974
#> 4 5:218387-181114186:EFO_0030071 181114186
#> 5 6:31980728-32013764:EFO_0030072 32013764
#> 6 6:117337256-117425595:EFO_0030068 117425595
#> variation_location_sequence_reference variation_location_start
#> 1 refseq:NC_000001.11 120069370
#> 2 refseq:NC_000001.11 147810332
#> 3 refseq:NC_000001.11 163343802
#> 4 refseq:NC_000005.10 218387
#> 5 refseq:NC_000006.12 31980728
#> 6 refseq:NC_000006.12 117337256
#> variation_location_type variant_copychange identifiers molecularAttributes
#> 1 SequenceLocation EFO:0030068 NA NA
#> 2 SequenceLocation EFO:0030071 NA NA
#> 3 SequenceLocation EFO:0030071 NA NA
#> 4 SequenceLocation EFO:0030071 NA NA
#> 5 SequenceLocation EFO:0030072 NA NA
#> 6 SequenceLocation EFO:0030068 NA NA
You can also query variant data from other Beacon v2 resources by setting the domain
and entry_point
parameters accordingly.
variant_2 <- pgxLoader(type="g_variants", biosample_id = "cellzbs-kftvksak", domain = "https://cancercelllines.org", entry_point="beacon")
head(variant_2)
#> variant_id analysis_id biosample_id
#> 1 cellzvar-63d118be42040c6f3561d4ac cellzcs-kftwtxov cellzbs-kftvksak
#> 2 cellzvar-63d118be42040c6f3561d4b4 cellzcs-kftwtxov cellzbs-kftvksak
#> 3 cellzvar-63d118be42040c6f3561d4b7 cellzcs-kftwtxov cellzbs-kftvksak
#> 4 cellzvar-63d118be42040c6f3561d4bd cellzcs-kftwtxov cellzbs-kftvksak
#> 5 cellzvar-63d118be42040c6f3561d4b3 cellzcs-kftwtxov cellzbs-kftvksak
#> 6 cellzvar-63d118be42040c6f3561d4ad cellzcs-kftwtxov cellzbs-kftvksak
#> individual_id variant_internal_id variation_location_end
#> 1 cellzind-B45e11EE 1:61735-16710048:DEL 16710048
#> 2 cellzind-B45e11EE 1:16710191-16788217:DUP 16788217
#> 3 cellzind-B45e11EE 1:16850550-16863509:DEL 16863509
#> 4 cellzind-B45e11EE 1:16864367-16935737:DEL 16935737
#> 5 cellzind-B45e11EE 1:16935752-29742865:DEL 29742865
#> 6 cellzind-B45e11EE 1:29745789-29748056:DEL 29748056
#> variation_location_sequence_reference variation_location_start
#> 1 refseq:NC_000001.11 61735
#> 2 refseq:NC_000001.11 16710191
#> 3 refseq:NC_000001.11 16850550
#> 4 refseq:NC_000001.11 16864367
#> 5 refseq:NC_000001.11 16935752
#> 6 refseq:NC_000001.11 29745789
#> variation_location_type variant_copychange identifiers molecularAttributes
#> 1 SequenceLocation EFO:0030067 NA NA
#> 2 SequenceLocation EFO:0030070 NA NA
#> 3 SequenceLocation EFO:0030067 NA NA
#> 4 SequenceLocation EFO:0030067 NA NA
#> 5 SequenceLocation EFO:0030067 NA NA
#> 6 SequenceLocation EFO:0030067 NA NA
output
= “pgxseg”)This format accesses data from Progenetix services API, which utilizes specialized resources within Progenetix. The ‘.pgxseg’ file format contains segment mean values (in log2
column), which are equal to log2(copy number of measured sample/copy number of control sample (usually 2)). A few variants are point mutations represented by columns reference_bases
and alternate_bases
.
variant_3 <- pgxLoader(type="g_variants", biosample_id = biosample_id,output = "pgxseg")
head(variant_3)
#> biosample_id reference_name start end log2 variant_type
#> 1 pgxbs-m3io3tbz 1 120069370 149078383 -0.4696 DEL
#> 2 pgxbs-m3io3tbz 1 147810332 163328877 0.4038 DUP
#> 3 pgxbs-m3io3tbz 1 163343802 247433974 0.3687 DUP
#> 4 pgxbs-m3io3tbz 5 218387 181114186 0.2860 DUP
#> 5 pgxbs-m3io3tbz 6 31980728 32013764 0.9103 DUP
#> 6 pgxbs-m3io3tbz 6 117337256 117425595 -0.4748 DEL
#> reference_sequence sequence variant_state_id variant_state_label
#> 1 . . EFO:0030068 low-level copy number loss
#> 2 . . EFO:0030071 low-level copy number gain
#> 3 . . EFO:0030071 low-level copy number gain
#> 4 . . EFO:0030071 low-level copy number gain
#> 5 . . EFO:0030072 high-level copy number gain
#> 6 . . EFO:0030068 low-level copy number loss
output
= “seg”)This format accesses data from Progenetix services API, which utilizes specialized resources within Progenetix. This format is compatible with IGV tool for visualization.
variant_4 <- pgxLoader(type="g_variants", biosample_id = biosample_id,output = "seg")
head(variant_4)
#> ID chrom loc.start loc.end num.mark seg.mean
#> 1 pgxbs-m3io3tbz 1 120069371 149078383 NA -0.4696
#> 2 pgxbs-m3io3tbz 1 147810333 163328877 NA 0.4038
#> 3 pgxbs-m3io3tbz 1 163343803 247433974 NA 0.3687
#> 4 pgxbs-m3io3tbz 5 218388 181114186 NA 0.2860
#> 5 pgxbs-m3io3tbz 6 31980729 32013764 NA 0.9103
#> 6 pgxbs-m3io3tbz 6 117337257 117425595 NA -0.4748
Setting save_file
to TRUE in the pgxLoader
function will save the retrieved variants data to a file
rather than returning it directly. By default, the data will be saved in the current working directory,
but you can specify a different path using the filename
parameter. This export functionality is
only available for variants data (when type='g_variants'
).
The following codes creates a ‘.pgxseg’ file with the name “variants.pgxseg” in “~/Downloads/” folder.
pgxLoader(type="g_variants", output="pgxseg", biosample_id=biosample_id, save_file=TRUE,
filename="~/Downloads/variants.pgxseg")
To visualize the ‘.pgxseg’ file, you can either upload it to this link or use the byconaut package for local visualization when dealing with a large number of samples.
The following codes creates a special ‘.seg’ file with the name “variants.seg” in “~/Downloads/” folder.
pgxLoader(type="g_variants", output="seg", biosample_id=biosample_id, save_file=TRUE,
filename="~/Downloads/variants.seg")
You can upload this ‘.seg’ file to IGV tool for visualization.
Because segment variants are not harmonized across samples, Progenetix provides processed CNV features, known as CNV fractions. These fractions represent the proportion of genomic regions overlapping one or more CNVs of a given type, facilitating sample-wise comparisons. The following query is based on filters, but biosample id and individual id are also available for sample-specific CNV fraction queries. For more information about filters, biosample id and individual id, as well as the use of parameters skip
, limit
, and codematches
, see the vignette Introduction_1_load_metadata.
CNV fraction data is retrieved through the Progenetix services API. When making a query, the domain of the data resource should be set to either "http://progenetix.org"
(by default) or "https://cancercelllines.org"
.
cnv_fraction <- pgxLoader(type="cnv_fraction", filters = "NCIT:C2948")
This includes CNV fraction across chromosome arms, whole chromosomes, or the whole genome.
names(cnv_fraction)
#> [1] "arm_cnv_frac" "chr_cnv_frac" "genome_cnv_frac"
The CNV fraction across chromosomal arms looks like this
head(cnv_fraction$arm_cnv_frac)[,c(1:4, 49:52)]
#> chr1p.dup chr1q.dup chr2p.dup chr2q.dup chr1p.del chr1q.del
#> pgxcs-kftvs0ri 0 0.000 0.000 0.000 0.000 0
#> pgxcs-kftvu9w2 0 0.000 0.000 0.000 0.000 0
#> pgxcs-kftvw1kw 0 0.000 0.000 0.000 0.000 0
#> pgxcs-kftvw1ld 0 0.000 0.000 0.000 0.000 0
#> pgxcs-kftvw1lu 0 0.000 0.000 0.000 0.225 0
#> pgxcs-kftvw1mb 0 0.979 0.989 0.003 0.000 0
#> chr2p.del chr2q.del
#> pgxcs-kftvs0ri 0 0
#> pgxcs-kftvu9w2 0 0
#> pgxcs-kftvw1kw 0 0
#> pgxcs-kftvw1ld 0 0
#> pgxcs-kftvw1lu 0 0
#> pgxcs-kftvw1mb 0 0
The row names are analyses ids from samples that belong to the input filter NCIT:C2948. There are 96 columns. The first 48 columns are duplication fraction across chromosomal arms, followed by deletion fraction. CNV fraction across whole chromosomes is similar, with the only difference in columns.
The CNV fraction across the genome (hg38) looks like this
head(cnv_fraction$genome_cnv_frac)
#> cnvfraction dupfraction delfraction
#> pgxcs-kftvs0ri 0.080 0.036 0.044
#> pgxcs-kftvu9w2 0.058 0.010 0.048
#> pgxcs-kftvw1kw 0.000 0.000 0.000
#> pgxcs-kftvw1ld 0.000 0.000 0.000
#> pgxcs-kftvw1lu 0.027 0.000 0.027
#> pgxcs-kftvw1mb 0.176 0.159 0.017
The first column is the total called fraction, followed by duplication fraction and deletion fraction.
cnvfrac_matrix <- pgxLoader(type="cnv_fraction", output="pgxmatrix", filters = "NCIT:C2948")
The returned data is stored in RangedSummarizedExperiment object, which is a matrix-like container where rows represent ranges of interest (as a GRanges object) and columns represent analyses derived from biosamples. The data looks like this
cnvfrac_matrix
#> class: RangedSummarizedExperiment
#> dim: 6212 53
#> metadata(0):
#> assays(1): cnv_matrix
#> rownames(6212): 1 2 ... 6211 6212
#> rowData names(1): type
#> colnames(53): pgxcs-kftvs0ri pgxcs-kftvu9w2 ... pgxcs-m3io7ksh
#> pgxcs-m3io7kt4
#> colData names(3): analysis_id biosample_id group_id
You could get the interval information by rowRanges
function from SummarizedExperiment package.
rowRanges(cnvfrac_matrix)
#> GRanges object with 6212 ranges and 1 metadata column:
#> seqnames ranges strand | type
#> <Rle> <IRanges> <Rle> | <character>
#> 1 chr1 0-400000 * | DUP
#> 2 chr1 400000-1400000 * | DUP
#> 3 chr1 1400000-2400000 * | DUP
#> 4 chr1 2400000-3400000 * | DUP
#> 5 chr1 3400000-4400000 * | DUP
#> ... ... ... ... . ...
#> 6208 chrY 52400000-53400000 * | DEL
#> 6209 chrY 53400000-54400000 * | DEL
#> 6210 chrY 54400000-55400000 * | DEL
#> 6211 chrY 55400000-56400000 * | DEL
#> 6212 chrY 56400000-57227415 * | DEL
#> -------
#> seqinfo: 24 sequences from an unspecified genome; no seqlengths
To access the CNV fraction matrix, use assay
accesssor from SummarizedExperiment package
assay(cnvfrac_matrix)[1:3,1:3]
#> pgxcs-kftvs0ri pgxcs-kftvu9w2 pgxcs-kftvw1kw
#> 1 0 0 0
#> 2 0 0 0
#> 3 0 0 0
The matrix contains 6212 rows, corresponding to genomic regions. These rows include 3,106 intervals with “gain status” and 3,106 intervals with “loss status”. The columns represent analysis profiles derived from samples that belong to the input filter.
The value is the fraction of the binned interval overlapping with one or more CNVs of the given type (DUP/DEL). For example, if the value in the second row, the first column is 0.2, it means that one or more duplication events overlapped with 20% of the genomic bin located in chromosome 1: 400000-1400000 in the first analysis profile.
To get associated biosample id and filters for analyses, use colData
function from SummarizedExperiment package:
colData(cnvfrac_matrix)
#> DataFrame with 53 rows and 3 columns
#> analysis_id biosample_id group_id
#> <character> <character> <character>
#> pgxcs-kftvs0ri pgxcs-kftvs0ri pgxbs-kftvh262 NCIT:C2948
#> pgxcs-kftvu9w2 pgxcs-kftvu9w2 pgxbs-kftvh9fp NCIT:C2948
#> pgxcs-kftvw1kw pgxcs-kftvw1kw pgxbs-kftvhf4h NCIT:C8893
#> pgxcs-kftvw1ld pgxcs-kftvw1ld pgxbs-kftvhf4i NCIT:C8893
#> pgxcs-kftvw1lu pgxcs-kftvw1lu pgxbs-kftvhf4k NCIT:C8893
#> ... ... ... ...
#> pgxcs-m3io7jp3 pgxcs-m3io7jp3 pgxbs-kftvkxqr NCIT:C2948
#> pgxcs-m3io7jq2 pgxcs-m3io7jq2 pgxbs-kftvl062 NCIT:C2948
#> pgxcs-m3io7kkt pgxcs-m3io7kkt pgxbs-kftvl35c NCIT:C2948
#> pgxcs-m3io7ksh pgxcs-m3io7ksh pgxbs-kftvl0ut NCIT:C2948
#> pgxcs-m3io7kt4 pgxcs-m3io7kt4 pgxbs-kftvl3oa NCIT:C2948
It is noted that the number of analysis profiles does not necessarily equal the number of samples. One biosample id may correspond to multiple analysis ids. group_id
corresponds to the meaning of filters
.
#> R version 4.5.0 (2025-04-11)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.2 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.22-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.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] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] SummarizedExperiment_1.39.0 Biobase_2.69.0
#> [3] GenomicRanges_1.61.0 GenomeInfoDb_1.45.3
#> [5] IRanges_2.43.0 S4Vectors_0.47.0
#> [7] BiocGenerics_0.55.0 generics_0.1.3
#> [9] MatrixGenerics_1.21.0 matrixStats_1.5.0
#> [11] future_1.40.0 pgxRpi_1.5.2
#> [13] BiocStyle_2.37.0
#>
#> loaded via a namespace (and not attached):
#> [1] tidyselect_1.2.1 dplyr_1.1.4 farver_2.1.2
#> [4] fastmap_1.2.0 digest_0.6.37 timechange_0.3.0
#> [7] lifecycle_1.0.4 survival_3.8-3 magrittr_2.0.3
#> [10] compiler_4.5.0 rlang_1.1.6 sass_0.4.10
#> [13] tools_4.5.0 yaml_2.3.10 data.table_1.17.0
#> [16] knitr_1.50 ggsignif_0.6.4 S4Arrays_1.9.0
#> [19] labeling_0.4.3 curl_6.2.2 DelayedArray_0.35.1
#> [22] RColorBrewer_1.1-3 abind_1.4-8 withr_3.0.2
#> [25] purrr_1.0.4 grid_4.5.0 ggpubr_0.6.0
#> [28] xtable_1.8-4 ggplot2_3.5.2 globals_0.17.0
#> [31] scales_1.4.0 dichromat_2.0-0.1 tinytex_0.57
#> [34] cli_3.6.5 crayon_1.5.3 rmarkdown_2.29
#> [37] future.apply_1.11.3 km.ci_0.5-6 httr_1.4.7
#> [40] survminer_0.5.0 cachem_1.1.0 splines_4.5.0
#> [43] parallel_4.5.0 BiocManager_1.30.25 XVector_0.49.0
#> [46] survMisc_0.5.6 vctrs_0.6.5 Matrix_1.7-3
#> [49] jsonlite_2.0.0 carData_3.0-5 bookdown_0.43
#> [52] car_3.1-3 rstatix_0.7.2 Formula_1.2-5
#> [55] listenv_0.9.1 magick_2.8.6 attempt_0.3.1
#> [58] tidyr_1.3.1 jquerylib_0.1.4 glue_1.8.0
#> [61] parallelly_1.43.0 codetools_0.2-20 lubridate_1.9.4
#> [64] gtable_0.3.6 UCSC.utils_1.5.0 tibble_3.2.1
#> [67] pillar_1.10.2 htmltools_0.5.8.1 R6_2.6.1
#> [70] KMsurv_0.1-5 evaluate_1.0.3 lattice_0.22-7
#> [73] backports_1.5.0 broom_1.0.8 bslib_0.9.0
#> [76] Rcpp_1.0.14 SparseArray_1.9.0 gridExtra_2.3
#> [79] xfun_0.52 zoo_1.8-14 pkgconfig_2.0.3