Method to run the FlowSOM clustering algorithm. This function runs FlowSOM on a data.table with cells (rows) vs markers (columns) with new columns for FlowSOM clusters and metaclusters. Output data will be "flowsom.res.original" (for clusters) and "flowsom.res.meta" (for metaclusters). Uses the R packages "FlowSOM" for clustering, "flowCore" for handling .fcs files, "Biobase" for creating a flow frame, "data.table" for handling data.table format.
Usage
run.flowsom(dat, use.cols, xdim=14, ydim=14, meta.k='auto',
max.meta=20, clust.seed=42, meta.seed=42, clust.name="FlowSOM_cluster",
meta.clust.name="FlowSOM_metacluster", mem.ctrl=TRUE, verbose=TRUE)
Arguments
- dat
NO DEFAULT. data.frame. Input sample.
- use.cols
NO DEFAULT. Vector of column names to use for clustering.
- xdim
DEFAULT = 14. Numeric. Number of first level clusters across the x-axis. xdim x ydim = total number of first level FlowSOM clusters.
- ydim
DEFAULT = 14. Numeric. Number of first level clusters across the y-axis. xdim x ydim = total number of first level FlowSOM clusters.
- meta.k
DEFAULT = 'auto'. If set to 'auto', then number of metaclusters will be determined automatically. Alternatively, can specify the desired number of metaclusters to create. If set to zero (0), no metaclusters will be created.
- max.meta
DEFAULT = 20. Only used if meta.k is set to 'auto'. This parameter indicates the maximum number of metaclusters FlowSOM will try out when determining the optimal number of metaclusters for the dataset.
- clust.seed
DEFAULT = 42 Numeric. Clustering seed for reproducibility.
- meta.seed
DEFAULT = 42 Numeric. Metaclustering seed for reproducibility.
- clust.name
DEFAULT = "FlowSOM_cluster". Character. Name of the resulting 'cluster' parameter.
- meta.clust.name
DEFAULT = "FlowSOM_metacluster". Character. Name of the resulting 'metacluster' parameter.
- mem.ctrl
DEFAULT = TRUE. Runs gc() (garbage collection) after a number of steps to free up memory that hasn't been released quickly enough.
- verbose
DEFAULT = TRUE. Logical. Whether to print progress messages.
Author
Thomas Ashhurst, thomas.ashhurst@sydney.edu.au Felix Marsh-Wakefield, felix.marsh-wakefield@sydney.edu.au
Examples
# Run FlowSOM on demonstration dataset
res <- Spectre::run.flowsom(Spectre::demo.clustered,
use.cols = c("NK11_asinh", "CD3_asinh",
"CD45_asinh", "Ly6G_asinh", "CD11b_asinh",
"B220_asinh", "CD8a_asinh", "Ly6C_asinh",
"CD4_asinh"))
#> Preparing data
#> Starting FlowSOM
#> Freeing up memory
#> Building SOM
#> Mapping data to SOM
#> Building MST
#> Adding cluster labels to input dataset
#> Warning: Column 'FlowSOM_cluster' already exists and will be replaced.
#> Starting metaclustering
#> Adding metacluster labels to input dataset
#> Warning: Column 'FlowSOM_metacluster' already exists and will be replaced.