The decontam functions isContaminant and isNotContaminant are made available for SummarizedExperiment objects.

# S4 method for class 'SummarizedExperiment'
isContaminant(
  seqtab,
  assay.type = assay_name,
  assay_name = "counts",
  name = "isContaminant",
  concentration = NULL,
  control = NULL,
  batch = NULL,
  threshold = 0.1,
  normalize = TRUE,
  detailed = TRUE,
  ...
)

# S4 method for class 'SummarizedExperiment'
isNotContaminant(
  seqtab,
  assay.type = assay_name,
  assay_name = "counts",
  name = "isNotContaminant",
  control = NULL,
  threshold = 0.5,
  normalize = TRUE,
  detailed = FALSE,
  ...
)

addContaminantQC(x, name = "isContaminant", ...)

# S4 method for class 'SummarizedExperiment'
addContaminantQC(x, name = "isContaminant", ...)

addNotContaminantQC(x, name = "isNotContaminant", ...)

# S4 method for class 'SummarizedExperiment'
addNotContaminantQC(x, name = "isNotContaminant", ...)

Arguments

seqtab, x

a SummarizedExperiment

assay.type

Character scalar. Specifies which assay to use for calculation. (Default: "counts")

assay_name

Deprecated. Use assay.type instead.

name

Character scalar. A name for the column of the colData where results will be stored. (Default: "isContaminant")

concentration

Character scalar or NULL. Defining a column with numeric values from the colData to use as concentration information. (Default: NULL)

control

Character scalar or NULL. Defining a column with logical values from the colData to define control and non-control samples. (Default: NULL)

batch

Character scalar or NULL. Defining a column with values interpretable as a factor from the colData to use as batch information. (Default: NULL)

threshold

Numeric scalar.. See decontam:isContaminant or decontam:isNotContaminant

normalize

Logical scalar. See decontam:isContaminant or decontam:isNotContaminant

detailed

Logical scalar. If TRUE, the return value is a data.frame containing diagnostic information on the contaminant decision. If FALSE, the return value is a logical vector containing the binary contaminant classifications. (Default: TRUE)

...

Value

for isContaminant/ isNotContaminant a DataFrame or for addContaminantQC/addNotContaminantQC a modified object of class(x)

Examples

data(esophagus)
# setup of some mock data
colData(esophagus)$concentration <- c(1,2,3)
colData(esophagus)$control <- c(FALSE,FALSE,TRUE)

isContaminant(esophagus,
              method = "frequency",
              concentration = "concentration")
#> Warning: Some batches have very few (<=4) samples.
#> DataFrame with 58 rows and 6 columns
#>               freq      prev    p.freq    p.prev         p contaminant
#>          <numeric> <integer> <numeric> <logical> <numeric>   <logical>
#> 59_8_22 0.11454876         3  0.227201        NA  0.227201       FALSE
#> 59_5_13 0.00261438         1        NA        NA        NA       FALSE
#> 59_8_12 0.02027754         2  0.609754        NA  0.609754       FALSE
#> 65_3_22 0.00413645         2  0.156237        NA  0.156237       FALSE
#> 65_5_1  0.00152207         1        NA        NA        NA       FALSE
#> ...            ...       ...       ...       ...       ...         ...
#> 59_4_16 0.00130719         1        NA        NA        NA       FALSE
#> 59_8_3  0.00522876         1        NA        NA        NA       FALSE
#> 59_5_19 0.09311117         3  0.588847        NA  0.588847       FALSE
#> 65_9_9  0.00282926         2  0.830389        NA  0.830389       FALSE
#> 59_2_6  0.32771733         3  0.724753        NA  0.724753       FALSE
esophagus <- addContaminantQC(esophagus,
                              method = "frequency",
                              concentration = "concentration")
#> Warning: Some batches have very few (<=4) samples.
colData(esophagus)
#> DataFrame with 3 rows and 2 columns
#>   concentration   control
#>       <numeric> <logical>
#> B             1     FALSE
#> C             2     FALSE
#> D             3      TRUE

isNotContaminant(esophagus, control = "control")
#> Warning: Some batches have very few (<=4) samples.
#>  [1] FALSE FALSE  TRUE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE
#> [13] FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE
#> [25] FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE
#> [37] FALSE FALSE  TRUE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE FALSE
#> [49] FALSE FALSE  TRUE FALSE FALSE FALSE FALSE FALSE FALSE FALSE
#> attr(,"metadata")
#> attr(,"metadata")$neg
#> [1] FALSE FALSE  TRUE
#> 
#> attr(,"metadata")$threshold
#> [1] 0.5
#> 
#> attr(,"metadata")$normalize
#> [1] TRUE
#> 
#> attr(,"metadata")$detailed
#> [1] FALSE
#> 
esophagus <- addNotContaminantQC(esophagus, control = "control")
#> Warning: Some batches have very few (<=4) samples.
colData(esophagus)
#> DataFrame with 3 rows and 2 columns
#>   concentration   control
#>       <numeric> <logical>
#> B             1     FALSE
#> C             2     FALSE
#> D             3      TRUE