tidyFlowCore
R
is an open-source statistical environment which can be
easily modified to enhance its functionality via packages. tidyFlowCore
is an R
package available via Bioconductor, an open-source
repository for R packages related to biostatistics and biocomputing.
R
can be installed on any operating system from CRAN, after which you can install
tidyFlowCore
by using the following commands in your R
session:
tidyFlowCore
adopts the so-called “tidy” functional programming paradigm developed by
Wickham et al. in the tidyverse
ecosystem of R packages
(Wickham, François, Henry, Müller, and Vaughan, 2023). For information
about the tidyverse
ecosytem broadly, feel free to
reference the (free) R for Data
Science book, the tidyverse
website, or this
paper describing the larger tidyomics
project.
tidyFlowCore
integrates the flowCore
Bioconductor package’s data
analysis capabilities with those of the tidyverse
. If
you’re relatively unfamiliar with the Bioconductor project, you might be
interested in this
blog post.
Learning to use R
and Bioconductor
can be
challenging, so it’s important to know where to get help. The main place
to ask questions about tidyFlowCore
is the Bioconductor support site.
Use the tidyFlowCore
tag there and look at previous
posts.
You can also ask questions on GitHub or Twitter. But remember, if you’re asking for help, follow the posting guidelines. Make sure to include a simple example that reproduces your issue (a “reprex”) and your session info to help developers understand and solve your problem.
tidyFlowCore
If you use tidyFlowCore for your research, please use the following citation.
citation("tidyFlowCore")
#> Warning in person1(given = given[[i]], family = family[[i]], middle =
#> middle[[i]], : It is recommended to use 'given' instead of 'middle'.
#> To cite package 'tidyFlowCore' in publications use:
#>
#> Keyes TJ (2024). _tidyFlowCore: Bringing flowCore to the tidyverse_.
#> doi:10.18129/B9.bioc.tidyFlowCore
#> <https://doi.org/10.18129/B9.bioc.tidyFlowCore>,
#> https://github.com/keyes-timothy/tidyflowCore/tidyFlowCore - R
#> package version 1.1.0,
#> <http://www.bioconductor.org/packages/tidyFlowCore>.
#>
#> A BibTeX entry for LaTeX users is
#>
#> @Manual{,
#> title = {tidyFlowCore: Bringing flowCore to the tidyverse},
#> author = {Timothy J Keyes},
#> year = {2024},
#> url = {http://www.bioconductor.org/packages/tidyFlowCore},
#> note = {https://github.com/keyes-timothy/tidyflowCore/tidyFlowCore - R package version 1.1.0},
#> doi = {10.18129/B9.bioc.tidyFlowCore},
#> }
tidyFlowCore
quick starttidyFlowCore
allows you to treat flowCore
data structures like tidy data.frame
s or
tibbles
. It does so by implementing dplyr, tidyr, and
ggplot2 verbs that can be deployed directly on the
flowFrame
and flowSet
S4 classes.
In this section, we give a brief example of how
tidyFlowCore
can enable a data analysis pipeline to use all
the useful functions of the flowCore
package and many of
the functions of the dplyr
, tidyr
, and
ggplot2
packages.
For our example here, we download some publicly available mass
cytometry (CyTOF) data downloadable through the (Weber, M, Soneson, and
Charlotte, 2019) package. These data are made available as a
flowCore::flowSet
S4 object, Bioconductor
’s
standard data structure for cytometry data.
# read data from the HDCytoData package
bcr_flowset <- HDCytoData::Bodenmiller_BCR_XL_flowSet()
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'SummarizedExperiment'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'IRanges'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'GenomeInfoDb'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'GenomicRanges'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'XVector'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'S4Arrays'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'DelayedArray'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'SparseArray'
#> Warning: replacing previous import 'S4Arrays::read_block' by
#> 'DelayedArray::read_block' when loading 'SummarizedExperiment'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'ExperimentHub'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'AnnotationHub'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'AnnotationDbi'
#> Warning: replacing previous import 'BiocGenerics::setequal' by
#> 'S4Vectors::setequal' when loading 'Biostrings'
#> Warning: multiple methods tables found for 'setequal'
#> see ?HDCytoData and browseVignettes('HDCytoData') for documentation
#> downloading 1 resources
#> retrieving 1 resource
#> loading from cache
#> Warning in updateObjectFromSlots(object, ..., verbose = verbose): dropping
#> slot(s) 'colnames' from object = 'flowSet'
To read more about this dataset, run the following command:
The flowCore
package natively supports multiple types of
data preprocessing and transformations for cytometry data through the
use of its tranform
class.
For example, if we want to apply the standard arcsinh transformation often used for CyTOF data to our current dataset, we could use the following code:
asinh_transformation <- flowCore::arcsinhTransform(a = 0, b = 1/5, c = 0)
transformation_list <-
flowCore::transformList(
colnames(bcr_flowset),
asinh_transformation
)
transformed_bcr_flowset <- flowCore::transform(bcr_flowset, transformation_list)
Alternatively, we can also use the tidyverse
’s
functional programming paradigm to perform the same transformation. For
this, we use the mutate
-across
framework via
tidyFlowCore
:
Suppose we’re interested in counting the number of cells that belong
to each cell type (encoded in the population_id
column of
bcr_flowset
) in our dataset. Using standard
flowCore
functions, we could perform this calculation in a
few steps:
# extract all expression matrices from our flowSet
combined_matrix <- flowCore::fsApply(bcr_flowset, exprs)
# take out the concatenated population_id column
combined_population_id <- combined_matrix[, 'population_id']
# perform the calculation
table(combined_population_id)
#> combined_population_id
#> 1 2 3 4 5 6 7 8
#> 3265 6651 62890 51150 1980 18436 24518 3901
tidyFlowCore
allows us to perform the same operation
simply using the dplyr
package’s count
function, with output provided in the convenient form of a
tibble
:
bcr_flowset |>
dplyr::count(population_id)
#> # A tibble: 8 × 2
#> population_id n
#> <dbl> <int>
#> 1 1 3265
#> 2 2 6651
#> 3 3 62890
#> 4 4 51150
#> 5 5 1980
#> 6 6 18436
#> 7 7 24518
#> 8 8 3901
tidyFlowCore
also makes it easy to perform the counting
after grouping by other variables in our metadata. For example, we can
see that bcr_flowset
contains data from multiple .FCS
files, each of which is associated with a file name.
flowCore::pData(object = bcr_flowset)
#> name
#> PBMC8_30min_patient1_BCR-XL.fcs PBMC8_30min_patient1_BCR-XL.fcs
#> PBMC8_30min_patient1_Reference.fcs PBMC8_30min_patient1_Reference.fcs
#> PBMC8_30min_patient2_BCR-XL.fcs PBMC8_30min_patient2_BCR-XL.fcs
#> PBMC8_30min_patient2_Reference.fcs PBMC8_30min_patient2_Reference.fcs
#> PBMC8_30min_patient3_BCR-XL.fcs PBMC8_30min_patient3_BCR-XL.fcs
#> PBMC8_30min_patient3_Reference.fcs PBMC8_30min_patient3_Reference.fcs
#> PBMC8_30min_patient4_BCR-XL.fcs PBMC8_30min_patient4_BCR-XL.fcs
#> PBMC8_30min_patient4_Reference.fcs PBMC8_30min_patient4_Reference.fcs
#> PBMC8_30min_patient5_BCR-XL.fcs PBMC8_30min_patient5_BCR-XL.fcs
#> PBMC8_30min_patient5_Reference.fcs PBMC8_30min_patient5_Reference.fcs
#> PBMC8_30min_patient6_BCR-XL.fcs PBMC8_30min_patient6_BCR-XL.fcs
#> PBMC8_30min_patient6_Reference.fcs PBMC8_30min_patient6_Reference.fcs
#> PBMC8_30min_patient7_BCR-XL.fcs PBMC8_30min_patient7_BCR-XL.fcs
#> PBMC8_30min_patient7_Reference.fcs PBMC8_30min_patient7_Reference.fcs
#> PBMC8_30min_patient8_BCR-XL.fcs PBMC8_30min_patient8_BCR-XL.fcs
#> PBMC8_30min_patient8_Reference.fcs PBMC8_30min_patient8_Reference.fcs
tidyFlowCore
makes it easy to perform grouped tidy
operations, like counting, using information in our
flowSet
’s metadata:
bcr_flowset |>
# use the .tidyFlowCore_identifier pronoun to access the name of
# each experiment in the flowSet
dplyr::count(.tidyFlowCore_identifier, population_id)
#> # A tibble: 128 × 3
#> .tidyFlowCore_identifier population_id n
#> <chr> <dbl> <int>
#> 1 PBMC8_30min_patient1_BCR-XL.fcs 1 31
#> 2 PBMC8_30min_patient1_BCR-XL.fcs 2 112
#> 3 PBMC8_30min_patient1_BCR-XL.fcs 3 761
#> 4 PBMC8_30min_patient1_BCR-XL.fcs 4 1307
#> 5 PBMC8_30min_patient1_BCR-XL.fcs 5 5
#> 6 PBMC8_30min_patient1_BCR-XL.fcs 6 127
#> 7 PBMC8_30min_patient1_BCR-XL.fcs 7 444
#> 8 PBMC8_30min_patient1_BCR-XL.fcs 8 51
#> 9 PBMC8_30min_patient1_Reference.fcs 1 52
#> 10 PBMC8_30min_patient1_Reference.fcs 2 132
#> # ℹ 118 more rows
flowFrame
and flowSet
data objects have a
clear relationship with one another in the flowCore
application programming interface (API). Essentially,
flowSet
s are nested flowFrame
s - or, in other
words, flowSet
s are made up of multiple
flowFrame
s!
tidyFlowCore
provides a useful API for converting
between flowSet
and flowFrame
data structures
at various degrees of nesting using the
group
/nest
and
ungroup
/unnest
verbs. Note that in the
dplyr
and tidyr
APIs,
group
/nest
and
ungroup
/unnest
are not
synonyms (grouped data.frames
are different from nested
data.frames
). However, because of how
flowFrame
s and flowSet
s are structured,
tidyFlowCore
’s group
/nest
and
ungroup
/unnest
functions have identical
behavior, respectively.
# unnesting a flowSet results in a flowFrame with an additional column,
# 'tidyFlowCore_name` that identifies cells based on which experiment in the
# original flowSet they come from
bcr_flowset |>
dplyr::ungroup()
#> flowFrame object 'file199e9927af3'
#> with 172791 cells and 40 observables:
#> name desc range minRange maxRange
#> $P1 Time Time 2399633 0.0000 2399632
#> $P2 Cell_length Cell_length 69 0.0000 68
#> $P3 CD3(110:114)Dd CD3(110:114)Dd 9383 -61.6796 9382
#> $P4 CD45(In115)Dd CD45(In115)Dd 5035 0.0000 5034
#> $P5 BC1(La139)Dd BC1(La139)Dd 14306 -100.8797 14305
#> ... ... ... ... ... ...
#> $P36 group_id group_id 3 0 2
#> $P37 patient_id patient_id 9 0 8
#> $P38 sample_id sample_id 17 0 16
#> $P39 population_id population_id 9 0 8
#> $P40 .tidyFlowCore_name .tidyFlowCore_name 17 0 16
#> 297 keywords are stored in the 'description' slot
# flowSets can be unnested and re-nested for various analyses
bcr_flowset |>
dplyr::ungroup() |>
# group_by cell type
dplyr::group_by(population_id) |>
# calculate the mean HLA-DR expression of each cell population
dplyr::summarize(mean_hladr_expression = mean(`HLA-DR(Yb174)Dd`)) |>
dplyr::select(population_id, mean_hladr_expression)
#> # A tibble: 8 × 2
#> population_id mean_hladr_expression
#> <dbl> <dbl>
#> 1 3 3.67
#> 2 7 3.33
#> 3 4 4.33
#> 4 2 87.1
#> 5 6 88.2
#> 6 8 3.12
#> 7 1 51.4
#> 8 5 18.0
tidyFlowCore
also provides a direct interface between
ggplot2
and flowFrame
or flowSet
data objects. For example…
# cell population names, from the HDCytoData documentation
population_names <-
c(
"B-cells IgM-",
"B-cells IgM+",
"CD4 T-cells",
"CD8 T-cells",
"DC",
"monocytes",
"NK cells",
"surface-"
)
# calculate mean CD20 expression across all cells
mean_cd20_expression <-
bcr_flowset |>
dplyr::ungroup() |>
dplyr::summarize(mean_expression = mean(asinh(`CD20(Sm147)Dd` / 5))) |>
dplyr::pull(mean_expression)
# calculate mean CD4 expression across all cells
mean_cd4_expression <-
bcr_flowset |>
dplyr::ungroup() |>
dplyr::summarize(mean_expression = mean(asinh(`CD4(Nd145)Dd` / 5))) |>
dplyr::pull(mean_expression)
bcr_flowset |>
# preprocess all columns that represent protein measurements
dplyr::mutate(dplyr::across(-ends_with("_id"), \(.x) asinh(.x / 5))) |>
# plot a CD4 vs. CD45 scatterplot
ggplot2::ggplot(ggplot2::aes(x = `CD20(Sm147)Dd`, y = `CD4(Nd145)Dd`)) +
# add some reference lines
ggplot2::geom_hline(
yintercept = mean_cd4_expression,
color = "red",
linetype = "dashed"
) +
ggplot2::geom_vline(
xintercept = mean_cd20_expression,
color = "red",
linetype = "dashed"
) +
ggplot2::geom_point(size = 0.1, alpha = 0.1) +
# facet by cell population
ggplot2::facet_wrap(
facets = ggplot2::vars(population_id),
labeller =
ggplot2::as_labeller(
\(population_id) population_names[as.numeric(population_id)]
)
) +
# axis labels
ggplot2::labs(
x = "CD20 expression (arcsinh)",
y = "CD4 expression (arcsinh)"
)
Using some standard functions from the ggplot2
library,
we can create a scatterplot of CD4 vs. CD20 expression in the different
cell populations included in the bcr_flowset
flowSet
. We can see, unsurprisingly, that both B-cell
populations are highest for CD20 expression, whereas CD4+ T-helper cells
are highest for CD4 expression.
The tidyFlowCore package (Keyes, 2024) was made possible thanks to the following:
This package was developed using biocthis.
Code for creating the vignette
## Create the vignette
library("rmarkdown")
system.time(render("tidyFlowCore.Rmd", "BiocStyle::html_document"))
## Extract the R code
library("knitr")
knit("tidyFlowCore.Rmd", tangle = TRUE)
Date the vignette was generated.
#> [1] "2024-11-01 06:27:32 UTC"
Wallclock time spent generating the vignette.
#> Time difference of 25.73 secs
R
session information.
#> ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
#> setting value
#> version R version 4.4.1 (2024-06-14)
#> os Ubuntu 24.04.1 LTS
#> system x86_64, linux-gnu
#> ui X11
#> language (EN)
#> collate C
#> ctype en_US.UTF-8
#> tz Etc/UTC
#> date 2024-11-01
#> pandoc 3.2.1 @ /usr/local/bin/ (via rmarkdown)
#>
#> ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
#> package * version date (UTC) lib source
#> abind 1.4-8 2024-09-12 [2] RSPM (R 4.4.0)
#> AnnotationDbi 1.69.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> AnnotationHub * 3.15.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> backports 1.5.0 2024-05-23 [2] RSPM (R 4.4.0)
#> bibtex 0.5.1 2023-01-26 [2] RSPM (R 4.4.0)
#> Biobase * 2.67.0 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> BiocFileCache * 2.15.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> BiocGenerics * 0.53.1 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> BiocManager 1.30.25 2024-08-28 [2] RSPM (R 4.4.0)
#> BiocStyle * 2.35.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> BiocVersion 3.21.1 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> Biostrings 2.75.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> bit 4.5.0 2024-09-20 [2] RSPM (R 4.4.0)
#> bit64 4.5.2 2024-09-22 [2] RSPM (R 4.4.0)
#> blob 1.2.4 2023-03-17 [2] RSPM (R 4.4.0)
#> bslib 0.8.0 2024-07-29 [2] RSPM (R 4.4.0)
#> buildtools 1.0.0 2024-10-31 [3] local (/pkg)
#> cachem 1.1.0 2024-05-16 [2] RSPM (R 4.4.0)
#> cli 3.6.3 2024-06-21 [2] RSPM (R 4.4.0)
#> colorspace 2.1-1 2024-07-26 [2] RSPM (R 4.4.0)
#> crayon 1.5.3 2024-06-20 [2] RSPM (R 4.4.0)
#> curl 5.2.3 2024-09-20 [2] RSPM (R 4.4.0)
#> cytolib 2.19.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> DBI 1.2.3 2024-06-02 [2] RSPM (R 4.4.0)
#> dbplyr * 2.5.0 2024-03-19 [2] RSPM (R 4.4.0)
#> DelayedArray 0.33.1 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> digest 0.6.37 2024-08-19 [2] RSPM (R 4.4.0)
#> dplyr 1.1.4 2023-11-17 [2] RSPM (R 4.4.0)
#> evaluate 1.0.1 2024-10-10 [2] RSPM (R 4.4.0)
#> ExperimentHub * 2.15.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> fansi 1.0.6 2023-12-08 [2] RSPM (R 4.4.0)
#> farver 2.1.2 2024-05-13 [2] RSPM (R 4.4.0)
#> fastmap 1.2.0 2024-05-15 [2] RSPM (R 4.4.0)
#> filelock 1.0.3 2023-12-11 [2] RSPM (R 4.4.0)
#> flowCore * 2.19.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> generics * 0.1.3 2022-07-05 [2] RSPM (R 4.4.0)
#> GenomeInfoDb * 1.43.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> GenomeInfoDbData 1.2.13 2024-11-01 [2] Bioconductor
#> GenomicRanges * 1.59.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> ggplot2 3.5.1 2024-04-23 [2] RSPM (R 4.4.0)
#> glue 1.8.0 2024-09-30 [2] RSPM (R 4.4.0)
#> gtable 0.3.6 2024-10-25 [2] RSPM (R 4.4.0)
#> HDCytoData * 1.25.0 2024-05-02 [2] Bioconductor 3.20 (R 4.4.1)
#> highr 0.11 2024-05-26 [2] RSPM (R 4.4.0)
#> htmltools 0.5.8.1 2024-04-04 [2] RSPM (R 4.4.0)
#> httr 1.4.7 2023-08-15 [2] RSPM (R 4.4.0)
#> IRanges * 2.41.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> jquerylib 0.1.4 2021-04-26 [2] RSPM (R 4.4.0)
#> jsonlite 1.8.9 2024-09-20 [2] RSPM (R 4.4.0)
#> KEGGREST 1.47.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> knitr 1.48 2024-07-07 [2] RSPM (R 4.4.0)
#> labeling 0.4.3 2023-08-29 [2] RSPM (R 4.4.0)
#> lattice 0.22-6 2024-03-20 [2] RSPM (R 4.4.0)
#> lifecycle 1.0.4 2023-11-07 [2] RSPM (R 4.4.0)
#> lubridate 1.9.3 2023-09-27 [2] RSPM (R 4.4.0)
#> magrittr 2.0.3 2022-03-30 [2] RSPM (R 4.4.0)
#> maketools 1.3.1 2024-10-31 [3] Github (jeroen/maketools@d46f92c)
#> Matrix 1.7-1 2024-10-18 [2] RSPM (R 4.4.0)
#> MatrixGenerics * 1.19.0 2024-10-30 [2] https://bioc.r-universe.dev (R 4.4.1)
#> matrixStats * 1.4.1 2024-09-08 [2] RSPM (R 4.4.0)
#> memoise 2.0.1 2021-11-26 [2] RSPM (R 4.4.0)
#> mime 0.12 2021-09-28 [2] RSPM (R 4.4.0)
#> munsell 0.5.1 2024-04-01 [2] RSPM (R 4.4.0)
#> pillar 1.9.0 2023-03-22 [2] RSPM (R 4.4.0)
#> pkgconfig 2.0.3 2019-09-22 [2] RSPM (R 4.4.0)
#> plyr 1.8.9 2023-10-02 [2] RSPM (R 4.4.0)
#> png 0.1-8 2022-11-29 [2] RSPM (R 4.4.0)
#> purrr 1.0.2 2023-08-10 [2] RSPM (R 4.4.0)
#> R6 2.5.1 2021-08-19 [2] RSPM (R 4.4.0)
#> rappdirs 0.3.3 2021-01-31 [2] RSPM (R 4.4.0)
#> Rcpp 1.0.13 2024-07-17 [2] RSPM (R 4.4.0)
#> RefManageR * 1.4.0 2022-09-30 [2] RSPM (R 4.4.0)
#> rlang 1.1.4 2024-06-04 [2] RSPM (R 4.4.0)
#> rmarkdown 2.28 2024-08-17 [2] RSPM (R 4.4.0)
#> RProtoBufLib 2.19.0 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> RSQLite 2.3.7 2024-05-27 [2] RSPM (R 4.4.0)
#> S4Arrays 1.7.1 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> S4Vectors * 0.45.0 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> sass 0.4.9 2024-03-15 [2] RSPM (R 4.4.0)
#> scales 1.3.0 2023-11-28 [2] RSPM (R 4.4.0)
#> sessioninfo * 1.2.2 2021-12-06 [2] RSPM (R 4.4.0)
#> SparseArray 1.7.0 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> stringi 1.8.4 2024-05-06 [2] RSPM (R 4.4.0)
#> stringr 1.5.1 2023-11-14 [2] RSPM (R 4.4.0)
#> SummarizedExperiment * 1.37.0 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> sys 3.4.3 2024-10-04 [2] RSPM (R 4.4.0)
#> tibble 3.2.1 2023-03-20 [2] RSPM (R 4.4.0)
#> tidyFlowCore * 1.1.0 2024-11-01 [1] https://bioc.r-universe.dev (R 4.4.1)
#> tidyr 1.3.1 2024-01-24 [2] RSPM (R 4.4.0)
#> tidyselect 1.2.1 2024-03-11 [2] RSPM (R 4.4.0)
#> timechange 0.3.0 2024-01-18 [2] RSPM (R 4.4.0)
#> UCSC.utils 1.3.0 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> utf8 1.2.4 2023-10-22 [2] RSPM (R 4.4.0)
#> vctrs 0.6.5 2023-12-01 [2] RSPM (R 4.4.0)
#> withr 3.0.2 2024-10-28 [2] RSPM (R 4.4.0)
#> xfun 0.49 2024-10-31 [2] CRAN (R 4.4.1)
#> xml2 1.3.6 2023-12-04 [2] RSPM (R 4.4.0)
#> XVector 0.47.0 2024-10-31 [2] https://bioc.r-universe.dev (R 4.4.1)
#> yaml 2.3.10 2024-07-26 [2] RSPM (R 4.4.0)
#> zlibbioc 1.51.2 2024-10-21 [2] Bioconductor 3.20 (R 4.4.1)
#>
#> [1] /tmp/RtmpD59kDQ/Rinst18e84c583340
#> [2] /github/workspace/pkglib
#> [3] /usr/local/lib/R/site-library
#> [4] /usr/lib/R/site-library
#> [5] /usr/lib/R/library
#>
#> ──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
This vignette was generated using BiocStyle (Oleś, 2024) with knitr (Xie, 2024) and rmarkdown (Allaire, Xie, Dervieux et al., 2024) running behind the scenes.
Citations made with RefManageR (McLean, 2017).
[1] J. Allaire, Y. Xie, C. Dervieux, et al. rmarkdown: Dynamic Documents for R. R package version 2.28. 2024. URL: https://github.com/rstudio/rmarkdown.
[2] T. J. Keyes. tidyFlowCore: Bringing flowCore to the tidyverse. https://github.com/keyes-timothy/tidyflowCore/tidyFlowCore - R package version 1.1.0. 2024. DOI: 10.18129/B9.bioc.tidyFlowCore. URL: http://www.bioconductor.org/packages/tidyFlowCore.
[3] M. W. McLean. “RefManageR: Import and Manage BibTeX and BibLaTeX References in R”. In: The Journal of Open Source Software (2017). DOI: 10.21105/joss.00338.
[4] A. Oleś. BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.35.0. 2024. URL: https://github.com/Bioconductor/BiocStyle.
[5] R Core Team. R: A Language and Environment for Statistical Computing. R Foundation for Statistical Computing. Vienna, Austria, 2024. URL: https://www.R-project.org/.
[6] Weber, L. M, Soneson, et al. “HDCytoData: Collection of high-dimensional cytometry benchmark datasets in Bioconductor object formats”. In: F1000Research 8.v2 (2019), p. 1459.
[7] H. Wickham. “testthat: Get Started with Testing”. In: The R Journal 3 (2011), pp. 5–10. URL: https://journal.r-project.org/archive/2011-1/RJournal_2011-1_Wickham.pdf.
[8] H. Wickham, W. Chang, R. Flight, et al. sessioninfo: R Session Information. R package version 1.2.2, https://r-lib.github.io/sessioninfo/. 2021. URL: https://github.com/r-lib/sessioninfo#readme.
[9] H. Wickham, R. François, L. Henry, et al. dplyr: A Grammar of Data Manipulation. R package version 1.1.4, https://github.com/tidyverse/dplyr. 2023. URL: https://dplyr.tidyverse.org.
[10] Y. Xie. knitr: A General-Purpose Package for Dynamic Report Generation in R. R package version 1.48. 2024. URL: https://yihui.org/knitr/.