Skip to contents

Normalizes and transforms cell data in preparation for clustering (accepts dataframe, SingleCellExperiment and SpatialExperiment).

Usage

normalizeCells(
  cells,
  markers = NULL,
  assayIn = NULL,
  assayOut = "norm",
  imageID = "imageID",
  transformation = NULL,
  method = NULL,
  cores = 1
)

Arguments

cells

A Dataframe of SingleCellExperiment or SpatialExperiment containing cells and features to be normalized/transformed

markers

A list containing the names of cell markers which will be normalized and/or transformed.

assayIn

If input is a SingleCellExperiment or SpatialExperiment with multiple assays, specify the assay to be normalized and/or transformed.

assayOut

If input is a SingleCellExperiment or SpatialExperiment, the new of the normalized data.

imageID

If input is a SingleCellExperiment or SpatialExperiment, this is the name of the image ID variable in order to stratify. cells correctly

transformation

The transformation/s to be performed, default is NULL, accepted values: 'asinh' and 'sqrt'.

method

The normalization method/s to be performed, default is NULL, accepted values: 'mean', 'minMax', 'trim99', 'PC1'.

cores

The number or cores for parallel processing.

Value

returns a dataframe with individual cells as rows and features as columns.

Examples


library(cytomapper)
#> Loading required package: EBImage
#> Loading required package: SingleCellExperiment
#> Loading required package: SummarizedExperiment
#> 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
#> 
#> 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, intersect, is.unsorted,
#>     lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
#>     pmin.int, rank, rbind, rownames, sapply, setdiff, sort, table,
#>     tapply, union, 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
#> 
#> Attaching package: ‘IRanges’
#> The following objects are masked from ‘package:EBImage’:
#> 
#>     resize, tile
#> 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
#> The following object is masked from ‘package:EBImage’:
#> 
#>     channel
#> 
#> Attaching package: ‘cytomapper’
#> The following objects are masked from ‘package:Biobase’:
#> 
#>     channelNames, channelNames<-
data("pancreasSCE")
cells.normalized <- normalizeCells(
  cells = pancreasSCE,
  markers = c("CD99", "PIN", "CD8a", "CDH"),
  assayIn = "counts",
  assayOut = "normCounts",
  imageID = "ImageNb",
  transformation = "asinh",
  method = "trim99"
)