Last data update: 2014.03.03

R: Apply a function on subsets of data frames
bpaggregateR Documentation

Apply a function on subsets of data frames

Description

This is a parallel version of aggregate.

Usage


## S4 method for signature 'formula,BiocParallelParam'
bpaggregate(x, data, FUN, ..., 
    BPREDO=list(), BPPARAM=bpparam())

## S4 method for signature 'data.frame,BiocParallelParam'
bpaggregate(x, by, FUN, ..., 
    simplify=TRUE, BPREDO=list(), BPPARAM=bpparam())

## S4 method for signature 'matrix,BiocParallelParam'
bpaggregate(x, by, FUN, ..., 
    simplify=TRUE, BPREDO=list(), BPPARAM=bpparam())

## S4 method for signature 'ANY,missing'
bpaggregate(x, ..., BPREDO=list(), BPPARAM=bpparam())

Arguments

x

A data.frame, matrix or a formula.

by

A list of factors by which x is split; applicable when x is data.frame or matrix.

data

A data.frame; applicable when x is a formula.

FUN

Function to apply.

...

Additional arguments for FUN.

simplify

If set to TRUE, the return values of FUN will be simplified using simplify2array.

BPPARAM

An optional BiocParallelParam instance determining the parallel back-end to be used during evaluation.

BPREDO

A list of output from bpaggregate with one or more failed elements. When a list is given in BPREDO, bpok is used to identify errors, tasks are rerun and inserted into the original results.

Details

bpaggregate is a generic with methods for data.frame matrix and formula objects. x is divided into subsets according to factors in by. Data chunks are sent to the workers, FUN is applied and results are returned as a data.frame.

The function is similar in spirit to aggregate from the stats package but aggregate is not explicitly called. The bpaggregate formula method reformulates the call and dispatches to the data.frame method which in turn distributes data chunks to workers with bplapply.

Value

See aggregate.

Author(s)

Martin Morgan mailto:mtmorgan@fhcrc.org.

Examples


if (require(Rsamtools) && require(GenomicAlignments)) {

  fl <- system.file("extdata", "ex1.bam", package="Rsamtools")
  param <- ScanBamParam(what = c("flag", "mapq"))
  gal <- readGAlignments(fl, param=param) 

  ## Report the mean map quality by range cutoff:
  cutoff <- rep(0, length(gal))
  cutoff[start(gal) > 1000 & start(gal) < 1500] <- 1
  cutoff[start(gal) > 1500] <- 2 
  bpaggregate(as.data.frame(mcols(gal)$mapq), list(cutoff = cutoff), mean)

}

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(BiocParallel)
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/BiocParallel/bpaggregate.Rd_%03d_medium.png", width=480, height=480)
> ### Name: bpaggregate
> ### Title: Apply a function on subsets of data frames
> ### Aliases: bpaggregate bpaggregate,formula,BiocParallelParam-method
> ###   bpaggregate,matrix,BiocParallelParam-method
> ###   bpaggregate,data.frame,BiocParallelParam-method
> ###   bpaggregate,ANY,missing-method
> 
> ### ** Examples
> 
> 
> if (require(Rsamtools) && require(GenomicAlignments)) {
+ 
+   fl <- system.file("extdata", "ex1.bam", package="Rsamtools")
+   param <- ScanBamParam(what = c("flag", "mapq"))
+   gal <- readGAlignments(fl, param=param) 
+ 
+   ## Report the mean map quality by range cutoff:
+   cutoff <- rep(0, length(gal))
+   cutoff[start(gal) > 1000 & start(gal) < 1500] <- 1
+   cutoff[start(gal) > 1500] <- 2 
+   bpaggregate(as.data.frame(mcols(gal)$mapq), list(cutoff = cutoff), mean)
+ 
+ }
Loading required package: Rsamtools
Loading required package: GenomeInfoDb
Loading required package: stats4
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Loading required package: S4Vectors

Attaching package: 'S4Vectors'

The following objects are masked from 'package:base':

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: IRanges
Loading required package: GenomicRanges
Loading required package: Biostrings
Loading required package: XVector
Loading required package: GenomicAlignments
Loading required package: SummarizedExperiment
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")'.

  cutoff mcols(gal)$mapq
1      0        93.98635
2      1        94.30386
3      2        59.49398
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>