R/detectCellOutlier.R
detectCellOutlier.Rd
A wrapper function for isOutlier. Identify outliers from numeric vectors stored in the SingleCellExperiment object.
detectCellOutlier(
inSCE,
slotName,
itemName,
sample = NULL,
nmads = 3,
type = "both",
overwrite = TRUE
)
A SingleCellExperiment object.
Desired slot of SingleCellExperiment used for plotting. Possible options: "assays", "colData", "metadata", "reducedDims". Required.
Desired vector within the slot used for plotting. Required.
A single character specifying a name that can be found in
colData(inSCE)
to directly use the cell annotation; or a character
vector with as many elements as cells to indicates which sample each cell
belongs to. Default NULL. decontX will be run on cells from
each sample separately.
Integer. Number of median absolute deviation. Parameter may be adjusted for more lenient or stringent outlier cutoff. Default 3.
Character. Type/direction of outlier detection; whether the lower/higher outliers should be detected, or both. Options are "both", "lower", "higher".
Boolean. If TRUE, and this function has previously generated an outlier decision on the same itemName, the outlier decision will be overwritten. Default TRUE.
A SingleCellExperiment object with '' added to the colData slot. Additionally, the decontaminated counts will be added as an assay called 'decontXCounts'.
data(scExample, package = "singleCellTK")
sce <- subsetSCECols(sce, colData = "type != 'EmptyDroplet'")
sce <- runDecontX(sce[,sample(ncol(sce),20)])
#> Thu Apr 28 11:23:46 2022 ... Running 'DecontX'
#> --------------------------------------------------
#> Starting DecontX
#> --------------------------------------------------
#> Thu Apr 28 11:23:46 2022 .. Analyzing all cells
#> Thu Apr 28 11:23:46 2022 .... Generating UMAP and estimating cell types
#> Thu Apr 28 11:23:48 2022 .... Estimating contamination
#> Thu Apr 28 11:23:48 2022 ...... Completed iteration: 10 | converge: 0.01087
#> Thu Apr 28 11:23:48 2022 ...... Completed iteration: 20 | converge: 0.004081
#> Thu Apr 28 11:23:48 2022 ...... Completed iteration: 30 | converge: 0.002172
#> Thu Apr 28 11:23:48 2022 ...... Completed iteration: 40 | converge: 0.001166
#> Thu Apr 28 11:23:48 2022 ...... Completed iteration: 42 | converge: 0.0009294
#> Thu Apr 28 11:23:48 2022 .. Calculating final decontaminated matrix
#> --------------------------------------------------
#> Completed DecontX. Total time: 2.486478 secs
#> --------------------------------------------------
sce <- detectCellOutlier(sce, slotName = "colData", sample = sce$sample,
nmads = 4, itemName = "decontX_contamination", type = "both")