Skip to contents

NOTE: This man page is about DelayedArray internals and is provided for developers and advanced users only.

The DelayedAbind class provides a formal representation of a delayed abind() operation. It is a concrete subclass of the DelayedNaryOp virtual class, which itself is a subclass of the DelayedOp virtual class:


                          DelayedOp
                              ^
                              |
                        DelayedNaryOp
                              ^
                              |
                        DelayedAbind
  

DelayedAbind objects are used inside a DelayedArray object to represent the delayed abind() operations carried by the object. They're never exposed to the end user and are not intended to be manipulated directly.

Usage

# S4 method for class 'DelayedAbind'
is_noop(x)

# S4 method for class 'DelayedAbind'
summary(object, ...)

## ~ ~ ~ Seed contract ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~

# S4 method for class 'DelayedAbind'
dim(x)

# S4 method for class 'DelayedAbind'
dimnames(x)

# S4 method for class 'DelayedAbind'
extract_array(x, index)

## ~ ~ ~ Propagation of sparsity ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~

# S4 method for class 'DelayedAbind'
is_sparse(x)

# S4 method for class 'DelayedAbind'
extract_sparse_array(x, index)

Arguments

x, object

A DelayedAbind object.

index

See ?extract_array in the S4Arrays package for a description of the index argument.

...

Not used.

See also

Examples

## DelayedAbind extends DelayedNaryOp which extends DelayedOp:
extends("DelayedAbind")
#> [1] "DelayedAbind"  "DelayedNaryOp" "DelayedOp"     "Array"        

## ---------------------------------------------------------------------
## BASIC EXAMPLE
## ---------------------------------------------------------------------
m1 <- matrix(101:128, ncol=4)
m2 <- matrix(runif(16), ncol=4)
M1 <- DelayedArray(m1)
M2 <- DelayedArray(m2)
showtree(M1)
#> 7x4 integer: DelayedMatrix object
#> └─ 7x4 integer: [seed] matrix object
showtree(M2)
#> 4x4 double: DelayedMatrix object
#> └─ 4x4 double: [seed] matrix object

M3 <- rbind(M1, M2)
showtree(M3)
#> 11x4 double: DelayedMatrix object
#> └─ 11x4 double: Abind (along=1)
#>    ├─ 7x4 integer: [seed] matrix object
#>    └─ 4x4 double: [seed] matrix object
class(M3@seed)        # a DelayedAbind object
#> [1] "DelayedAbind"
#> attr(,"package")
#> [1] "DelayedArray"

M4 <- cbind(t(M1), M2)
showtree(M4)
#> 4x11 double: DelayedMatrix object
#> └─ 4x11 double: Abind (along=2)
#>    ├─ 4x7 integer: Aperm (perm=c(2,1))
#>    │  └─ 7x4 integer: [seed] matrix object
#>    └─ 4x4 double: [seed] matrix object
class(M4@seed)        # a DelayedAbind object
#> [1] "DelayedAbind"
#> attr(,"package")
#> [1] "DelayedArray"

## ---------------------------------------------------------------------
## PROPAGATION OF SPARSITY
## ---------------------------------------------------------------------
## DelayedAbind objects always propagate sparsity (granted that all the
## input arrays are sparse).

sm1 <- sparseMatrix(i=c(1, 1, 7, 7), j=c(1, 4, 1, 4),
                    x=c(11, 14, 71, 74), dims=c(7, 4))
SM1 <- DelayedArray(sm1)
sm2 <- sparseMatrix(i=c(1, 1, 4, 4), j=c(1, 4, 1, 4),
                    x=c(11, 14, 41, 44), dims=c(4, 4))
SM2 <- DelayedArray(sm2)
showtree(SM1)
#> 7x4 double, sparse: DelayedMatrix object
#> └─ 7x4 double, sparse: [seed] dgCMatrix object
showtree(SM2)
#> 4x4 double, sparse: DelayedMatrix object
#> └─ 4x4 double, sparse: [seed] dgCMatrix object
is_sparse(SM1)        # TRUE
#> [1] TRUE
is_sparse(SM2)        # TRUE
#> [1] TRUE

SM3 <- rbind(SM1, SM2)
showtree(SM3)
#> 11x4 double, sparse: DelayedMatrix object
#> └─ 11x4 double, sparse: Abind (along=1)
#>    ├─ 7x4 double, sparse: [seed] dgCMatrix object
#>    └─ 4x4 double, sparse: [seed] dgCMatrix object
class(SM3@seed)       # a DelayedAbind object
#> [1] "DelayedAbind"
#> attr(,"package")
#> [1] "DelayedArray"
is_sparse(SM3@seed)   # TRUE
#> [1] TRUE

SM4 <- cbind(SM2, t(SM1))
showtree(SM4)
#> 4x11 double, sparse: DelayedMatrix object
#> └─ 4x11 double, sparse: Abind (along=2)
#>    ├─ 4x4 double, sparse: [seed] dgCMatrix object
#>    └─ 4x7 double, sparse: Aperm (perm=c(2,1))
#>       └─ 7x4 double, sparse: [seed] dgCMatrix object
class(SM4@seed)       # a DelayedAbind object
#> [1] "DelayedAbind"
#> attr(,"package")
#> [1] "DelayedArray"
is_sparse(SM4@seed)   # TRUE
#> [1] TRUE

M5 <- rbind(SM2, M1)  # 2nd input array is not sparse!
showtree(M5)
#> 11x4 double: DelayedMatrix object
#> └─ 11x4 double: Abind (along=1)
#>    ├─ 4x4 double, sparse: [seed] dgCMatrix object
#>    └─ 7x4 integer: [seed] matrix object
class(M5@seed)        # a DelayedAbind object
#> [1] "DelayedAbind"
#> attr(,"package")
#> [1] "DelayedArray"
is_sparse(M5@seed)    # FALSE
#> [1] FALSE

## ---------------------------------------------------------------------
## SANITY CHECKS
## ---------------------------------------------------------------------
stopifnot(class(M3@seed) == "DelayedAbind")
stopifnot(class(M4@seed) == "DelayedAbind")
stopifnot(class(SM3@seed) == "DelayedAbind")
stopifnot(is_sparse(SM3@seed))
stopifnot(class(SM4@seed) == "DelayedAbind")
stopifnot(is_sparse(SM4@seed))
stopifnot(class(M5@seed) == "DelayedAbind")
stopifnot(!is_sparse(M5@seed))