gs_pop_get_stats | R Documentation |
Extract stats from populations(or nodes)
gs_pop_get_stats(x, ...)
gh_pop_get_stats(
x,
nodes = NULL,
type = "count",
xml = FALSE,
inverse.transform = FALSE,
stats.fun.arg = list(),
...
)
x |
a GatingSet or GatingHierarchy |
... |
arguments passed to gs_get_pop_paths method. |
nodes |
the character vector specifies the populations of interest. default is all available nodes |
type |
the character vector specifies the type of pop stats or a function used to compute population stats. when character, it is expected to be either "count" or "percent". Default is "count" (total number of events in the populations). when a function, it takes a flowFrame object through 'fr' argument and return the stats as a named vector. |
xml |
whether to extract xml stats or openCyto stats |
inverse.transform |
logical flag . Whether inverse transform the data before computing the stats. |
stats.fun.arg |
a list of arguments passed to ‘type' when ’type' is a function. |
a data.table that contains stats values (if MFI, for each marker per column) along with 'pop' column and 'sample' column (when used on a 'GatingSet')
## Not run:
dataDir <- system.file("extdata",package="flowWorkspaceData")
suppressMessages(gs <- load_gs(list.files(dataDir, pattern = "gs_manual",full = TRUE)))
# get stats all nodes
dt <- gs_pop_get_stats(gs) #default is "count"
nodes <- c("CD4", "CD8")
gs_pop_get_stats(gs, nodes, "percent")
# pass a build-in function
gs_pop_get_stats(gs, nodes, type = pop.MFI)
# compute the stats based on the raw data scale
gs_pop_get_stats(gs, nodes, type = pop.MFI, inverse.transform = TRUE)
# supply user-defined stats fun
pop.quantiles <- function(fr){
chnls <- colnames(fr)
res <- matrixStats::colQuantiles(exprs(fr), probs = 0.75)
names(res) <- chnls
res
}
gs_pop_get_stats(gs, nodes, type = pop.quantiles)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.