The decontam
functions isContaminant
and
isNotContaminant
are made available for
SummarizedExperiment
objects.
addContaminantQC(x, name = "isContaminant", ...)
addNotContaminantQC(x, name = "isNotContaminant", ...)
# S4 method for class 'SummarizedExperiment'
isContaminant(
seqtab,
assay.type = assay_name,
assay_name = "counts",
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",
control = NULL,
threshold = 0.5,
normalize = TRUE,
detailed = TRUE,
...
)
# S4 method for class 'SummarizedExperiment'
addContaminantQC(x, name = "contaminant", ...)
# S4 method for class 'SummarizedExperiment'
addNotContaminantQC(x, name = "not_contaminant", ...)
Character scalar
. A name for the column of the
colData
where results will be stored.
(Default: "contaminant"
or "not_contaminant"
)
for isContaminant
/ isNotContaminant
: arguments
passed on to decontam:isContaminant
or decontam:isNotContaminant
for addContaminantQC
/addNotContaminantQC
: arguments
passed on to isContaminant
/ isNotContaminant
Character scalar
. Specifies the name of assay
used in calculation. (Default: "counts"
)
Deprecated. Use assay.type
instead.
Character scalar
or NULL
. Defining
a column with numeric values from the colData
to use as
concentration information. (Default: NULL
)
Character scalar
or NULL
. Defining a
column with logical values from the colData
to define control and
non-control samples. (Default: NULL
)
Character scalar
or NULL
. Defining a
column with values interpretable as a factor from the colData
to use
as batch information. (Default: NULL
)
Numeric scalar
.. See
decontam:isContaminant
or
decontam:isNotContaminant
Logical scalar
. See
decontam:isContaminant
or
decontam:isNotContaminant
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
)
for isContaminant
/ isNotContaminant
a DataFrame
or for addContaminantQC
/addNotContaminantQC
a modified object
of class(x)
data(esophagus)
# setup of some mock data just for example
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.
rowData(esophagus)
#> 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
isNotContaminant(esophagus, control = "control")
#> Warning: Some batches have very few (<=4) samples.
#> DataFrame with 58 rows and 6 columns
#> freq prev p.freq p.prev p not.contaminant
#> <numeric> <integer> <numeric> <numeric> <numeric> <logical>
#> 59_8_22 0.11454876 3 NA 0.500000 0.500000 FALSE
#> 59_5_13 0.00261438 1 NA NA NA FALSE
#> 59_8_12 0.02027754 2 NA 0.166667 0.166667 TRUE
#> 65_3_22 0.00413645 2 NA 0.666667 0.666667 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 NA 0.500000 0.500000 FALSE
#> 65_9_9 0.00282926 2 NA 0.666667 0.666667 FALSE
#> 59_2_6 0.32771733 3 NA 0.500000 0.500000 FALSE
esophagus <- addNotContaminantQC(esophagus, control = "control")
#> Warning: Some batches have very few (<=4) samples.
#> Warning: The following values are already present in `rowData` and will be overwritten: 'freq', 'prev', 'p.freq', 'p.prev', 'p'. Consider using the 'name' argument to specify alternative names.
rowData(esophagus)
#> DataFrame with 58 rows and 7 columns
#> contaminant freq prev p.freq p.prev p
#> <logical> <numeric> <integer> <numeric> <numeric> <numeric>
#> 59_8_22 FALSE 0.11454876 3 NA 0.500000 0.500000
#> 59_5_13 FALSE 0.00261438 1 NA NA NA
#> 59_8_12 FALSE 0.02027754 2 NA 0.166667 0.166667
#> 65_3_22 FALSE 0.00413645 2 NA 0.666667 0.666667
#> 65_5_1 FALSE 0.00152207 1 NA NA NA
#> ... ... ... ... ... ... ...
#> 59_4_16 FALSE 0.00130719 1 NA NA NA
#> 59_8_3 FALSE 0.00522876 1 NA NA NA
#> 59_5_19 FALSE 0.09311117 3 NA 0.500000 0.500000
#> 65_9_9 FALSE 0.00282926 2 NA 0.666667 0.666667
#> 59_2_6 FALSE 0.32771733 3 NA 0.500000 0.500000
#> not.contaminant
#> <logical>
#> 59_8_22 FALSE
#> 59_5_13 FALSE
#> 59_8_12 TRUE
#> 65_3_22 FALSE
#> 65_5_1 FALSE
#> ... ...
#> 59_4_16 FALSE
#> 59_8_3 FALSE
#> 59_5_19 FALSE
#> 65_9_9 FALSE
#> 59_2_6 FALSE