If you are new to this package, then please consult the vignette Introduction: microarray quality assessment with arrayQualityMetrics.
This vignette addresses advanced topics. It explains how to customize the report by selecting specific modules and sections, or by adding your own ones. Furthermore, we will see how to (programmatically) postprocess the results of the outlier detection, or how to adapt the detection criteria to your needs.
Terminology: In the documentation of this package, we refer
to a module as a self-contained element of a report that
investigates one particular quality metric. A module consists of a
figure and an explanatory text. It may also contain a data structure (an
object of class outlierDetection) that marks a subset of
the arrays in the dataset as outliers according to the particular metric
investigated in the module. We refer to a section as a
collection of one or more modules that are thematically related.
For the following examples, let us load the needed packages and some data.
Some of the computations that are needed for the modules are common
between several modules, and thus we perform them once, beforehand.
These functions are called prepdata() and
prepaffy(), and we refer to their manual page for details.
For example,
The arguments intgroup and do.logtransform
are the same as for arrayQualityMetrics(), but in
prepdata() they have no defaults, so we need to set them
explicitely.
The package contains a variety of functions that compute modules, and they are listed on a manual page which you can access by typing:
Here, let us create a report with only two quality metrics modules: boxplots and density plots.
bo = aqm.boxplot(preparedData)
de = aqm.density(preparedData)
qm = list("Boxplot" = bo, "Density" = de)The objects bo and de are of class
aqmReportModule; please consult the class’ manual page for
more information.
If you want to create your own modules, please have a look at the
code for the various existing functions for this purpose, and adapt it.
The function aqm.pca() is a good place to start.
Some of the modules perform outlier detection. For instance, in the
default report produced by arrayQualityMetrics(), the
module headed Boxplots is followed by one headed Outlier
detection for Boxplots. In the corresponding
aqmReportModule object, this is reflected by a non-trivial
value for the slot named outliers:
## An object of class "outlierDetection"
## Slot "statistic":
## JD-ALD009-v5-U133A.CEL JD-ALD051-v5-U133A.CEL JD-ALD052-v5-U133A.CEL
## 0.0685725 0.1006875 0.1625925
## JD-ALD057-v5-U133A.CEL JD-ALD078-v5-U133A.CEL JD-ALD180-v5-U133A.CEL
## 0.0569800 0.1190800 0.2310725
## JD-ALD226-NA-v5-U133A.CEL JD-ALD232-v5-U133A.CEL JD-ALD237-v5-U133A.CEL
## 0.2181700 0.0923500 0.1526125
## JD-ALD244-v5-U133A.CEL JD-ALD294-v5-U133A.CEL JD-ALD380-v5-U133A.CEL
## 0.1327825 0.3755475 0.0524125
## JD-ALD381-v5-U133A.CEL JD-ALD384-v5-U133A.CEL JD-ALD385-v5-U133A.CEL
## 0.0934375 0.4407825 0.1645400
## JD-ALD420-v5-U133A.CEL JD-ALD421-v5-U133A.CEL JD-ALD431-v5-U133A.CEL
## 0.1054900 0.1232900 0.1928375
## JD-ALD433-v5-U133A.CEL JD-ALD520-v5-U133A.CEL
## 0.1501825 0.0579050
##
## Slot "threshold":
## JD-ALD385-v5-U133A.CEL
## 0.3073812
##
## Slot "which":
## JD-ALD294-v5-U133A.CEL JD-ALD384-v5-U133A.CEL
## 11 14
##
## Slot "description":
## [1] "Kolmogorov-Smirnov statistic <i>K<sub>a</sub></i>"
## [2] "data-driven"
The slot named statistic contains, for each array, a
single number based on which outlier detection is performed. For
instance, in the case of bo, the slot
bo@outliers@statistic is the Kolmogorov-Smirnov statistic
for the comparison between each array’s intensity distribution and the
distribution of the pooled data. The slot threshold
contains the threshold against which the valuess of
statistic were compared. Arrays with a value of
statistic greater than threshold are called
outliers. Their indices are listed in the vector which.
Finally, the slot description contains a textual
description of the definition of statistic and indicates
how the threshold was chosen. The actual details of outlier
detection are different for each report module, and are documented in
the figure caption of the report module. For more information, please
look at the code of the report module generating function of interest –
for instance, at the first few lines of the boxplot
function. The code there uses the helper functions outliers
and boxplotOutliers, which are documented in their manual
pages.
A report is rendered by calling the function
aqm.writereport() with a list of
aqmReportModule objects.
outdir = tempdir()
aqm.writereport(modules = qm, reporttitle = "My example", outdir = outdir,
arrayTable = pData(MLL.A))
outdir## [1] "/tmp/RtmpxV1N9j"
Point your browser to the index.html file in that
directory.
## R version 4.5.3 (2026-03-11)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=en_US.UTF-8
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## time zone: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] ALLMLL_1.51.0 affy_1.89.0
## [3] Biobase_2.71.0 BiocGenerics_0.57.1
## [5] generics_0.1.4 arrayQualityMetrics_3.67.4
## [7] BiocStyle_2.39.0
##
## loaded via a namespace (and not attached):
## [1] affyPLM_1.87.0 farver_2.1.2 blob_1.3.0
## [4] Biostrings_2.79.5 S7_0.2.1-1 fastmap_1.2.0
## [7] XML_3.99-0.23 digest_0.6.39 rpart_4.1.27
## [10] lifecycle_1.0.5 cluster_2.1.8.2 statmod_1.5.1
## [13] survival_3.8-6 KEGGREST_1.51.1 RSQLite_2.4.6
## [16] magrittr_2.0.5 genefilter_1.93.0 compiler_4.5.3
## [19] setRNG_2024.2-1 rlang_1.2.0 Hmisc_5.2-5
## [22] sass_0.4.10 tools_4.5.3 yaml_2.3.12
## [25] beadarray_2.61.3 data.table_1.18.2.1 knitr_1.51
## [28] htmlwidgets_1.6.4 interp_1.1-6 bit_4.6.0
## [31] plyr_1.8.9 RColorBrewer_1.1-3 foreign_0.8-91
## [34] hwriter_1.3.2.1 sys_3.4.3 nnet_7.3-20
## [37] grid_4.5.3 stats4_4.5.3 latticeExtra_0.6-31
## [40] preprocessCore_1.73.0 xtable_1.8-8 colorspace_2.1-2
## [43] ggplot2_4.0.2 scales_1.4.0 cli_3.6.6
## [46] rmarkdown_2.31 crayon_1.5.3 rstudioapi_0.18.0
## [49] reshape2_1.4.5 httr_1.4.8 gridSVG_1.7-7
## [52] DBI_1.3.0 cachem_1.1.0 stringr_1.6.0
## [55] splines_4.5.3 AnnotationDbi_1.73.1 vsn_3.79.6
## [58] BiocManager_1.30.27 XVector_0.51.0 matrixStats_1.5.0
## [61] base64enc_0.1-6 vctrs_0.7.3 Matrix_1.7-5
## [64] jsonlite_2.0.0 IRanges_2.45.0 S4Vectors_0.49.2
## [67] bit64_4.8.0 Formula_1.2-5 htmlTable_2.4.3
## [70] jpeg_0.1-11 systemfonts_1.3.2 maketools_1.3.2
## [73] limma_3.67.3 hexbin_1.28.5 jquerylib_0.1.4
## [76] affyio_1.81.0 annotate_1.89.0 glue_1.8.1
## [79] gcrma_2.83.0 codetools_0.2-20 stringi_1.8.7
## [82] gtable_0.3.6 deldir_2.0-4 GenomicRanges_1.63.2
## [85] htmltools_0.5.9 Seqinfo_1.1.0 BeadDataPackR_1.63.0
## [88] R6_2.6.1 textshaping_1.0.5 evaluate_1.0.5
## [91] lattice_0.22-9 png_0.1-9 backports_1.5.1
## [94] memoise_2.0.1 bslib_0.10.0 Rcpp_1.1.1-1
## [97] svglite_2.2.2 gridExtra_2.3 checkmate_2.3.4
## [100] xfun_0.57 MatrixGenerics_1.23.0 buildtools_1.0.0