R/preprocessingMinfiEwasWater_steps.R
assessSamplesMinfiEwasWater.RdCompute minfi QC metrics and detection P values, identify failed samples using the requested threshold, and return the assessment as a single object.
assessSamplesMinfiEwasWater(
rawData,
RGSet,
qcCutoff = 10.5,
detPtype = "m+u",
detPThreshold = 0.05,
verbose = FALSE,
logs = FALSE,
log_dir = NULL,
log_file = "log_assessSamplesMinfiEwasWater.txt"
)Object returned by buildRawMinfiEwasWater().
An RGChannelSet aligned with rawData.
Numeric. Cutoff passed to minfi::plotQC() when the QC plot
is drawn.
Character. Detection P-value mode passed to
minfi::detectionP(). Common values in minfi workflows are "m+u" and
"negative". The default used here is "m+u".
Numeric. Samples with mean detection P value above this threshold are marked as failed.
Logical. If TRUE, emit progress messages with message().
Logical. If TRUE, write the same messages to a log file.
Character or NULL. Directory used for the log file when
logs = TRUE.
Character. File name used when logs = TRUE.
A list with class "dnaEPICO_minfiEwasWater_assessment" containing
the QC object, detection P matrix, mean detection P values, and failed
sample identifiers.
ex <- dnaEPICO:::exampleMinfiBaseDataDnaEpico()
raw_data <- buildRawMinfiEwasWater(ex$RGSet, verbose = FALSE, logs = FALSE)
#> Loading required package: IlluminaHumanMethylation450kmanifest
#> Loading required package: minfi
#> 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, saveRDS, setdiff, table,
#> tapply, union, unique, unsplit, which.max, which.min
#> Loading required package: GenomicRanges
#> Loading required package: stats4
#> 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
#> Loading required package: GenomeInfoDb
#> 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: 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
#> Loading required package: Biostrings
#> Loading required package: XVector
#>
#> Attaching package: ‘Biostrings’
#> The following object is masked from ‘package:base’:
#>
#> strsplit
#> Loading required package: bumphunter
#> Loading required package: foreach
#> Loading required package: iterators
#> Loading required package: parallel
#> Loading required package: locfit
#> locfit 1.5-9.12 2025-03-05
#> Loading required package: IlluminaHumanMethylation450kanno.ilmn12.hg19
assessment <- assessSamplesMinfiEwasWater(
rawData = raw_data,
RGSet = ex$RGSet,
detPThreshold = 1,
verbose = FALSE,
logs = FALSE
)
names(assessment)
#> [1] "qc" "qcCutoff" "detP" "detPtype"
#> [5] "detPThreshold" "meanDetP" "failedSamples"