DelayedUnaryIsoOpStack-class | R Documentation |
NOTE: This man page is about DelayedArray internals and is provided for developers and advanced users only.
The DelayedUnaryIsoOpStack class provides a formal representation of a stack of delayed unary isometric operations, that is, of a group of delayed unary isometric operations stacked (a.k.a. piped) together. It is a concrete subclass of the DelayedUnaryIsoOp virtual class, which itself is a subclass of the DelayedUnaryOp virtual class, which itself is a subclass of the DelayedOp virtual class:
DelayedOp ^ | DelayedUnaryOp ^ | DelayedUnaryIsoOp ^ | DelayedUnaryIsoOpStack
DelayedUnaryIsoOpStack objects are used inside a DelayedArray object to represent groups of delayed unary isometric operations carried by the object. They're never exposed to the end user and are not intended to be manipulated directly.
## S4 method for signature 'DelayedUnaryIsoOpStack'
summary(object, ...)
## ~ ~ ~ Seed contract ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
## DelayedUnaryIsoOpStack objects inherit the default dim()
## and dimnames() methods defined for DelayedUnaryIsoOp
## derivatives, but overwite their extract_array() method.
## S4 method for signature 'DelayedUnaryIsoOpStack'
extract_array(x, index)
## ~ ~ ~ Propagation of sparsity ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~
## S4 method for signature 'DelayedUnaryIsoOpStack'
is_sparse(x)
## S4 method for signature 'DelayedUnaryIsoOpStack'
extract_sparse_array(x, index)
x , object |
A DelayedUnaryIsoOpStack object. |
index |
See |
... |
Not used. |
A DelayedUnaryIsoOpStack object is used to represent the delayed version of an operation of the form:
out <- a |> OP1 |> OP2 |> ... |> OPk
where:
OP1
, OP2
, ..., OPk
are isometric array
transformations i.e. operations that return an array with the
same dimensions as the input array.
a
is the input array.
The output (out
) is an array of same dimensions as a
.
In addition, each operation (OP
) in the pipe must satisfy the
property that each value in the output array must be determined **solely**
by the corresponding value in the input array. In other words:
a |> OP |> `[`(i_1, i_2, ..., i_n) # i.e. OP(a)[i_1, i_2, ..., i_n]
must be equal to:
a |> `[`(i_1, i_2, ..., i_n) |> OP # i.e. OP(a[i_1, i_2, ..., i_n])
for any valid multidimensional index (i_1, i_2, ..., i_n).
We refer to this property as the locality principle.
Concrete examples:
Things like is.na()
, is.finite()
, logical negation
(!
), nchar()
, tolower()
.
Most functions in the Math and Math2
groups e.g. log()
, sqrt()
, abs()
,
ceiling()
, round()
, etc...
Notable exceptions are the cum*()
functions (cummin()
,
cummax()
, cumsum()
, and cumprod()
): they don't
satisfy the locality principle.
Operations in the Ops group when one operand is
an array and the other a scalar e.g. a + 10
, 2 ^ a
,
a <= 0.5
, etc...
DelayedOp objects.
showtree
to visualize the nodes and access the
leaves in the tree of delayed operations carried by a
DelayedArray object.
extract_array in the S4Arrays package.
extract_sparse_array
in the
SparseArray package.
## DelayedUnaryIsoOpStack extends DelayedUnaryIsoOp, which extends
## DelayedUnaryOp, which extends DelayedOp:
extends("DelayedUnaryIsoOpStack")
## ---------------------------------------------------------------------
## BASIC EXAMPLE
## ---------------------------------------------------------------------
m0 <- matrix(runif(12), ncol=3)
M0 <- DelayedArray(m0)
showtree(M0)
M <- log(1 + M0) / 10
showtree(M)
class(M@seed) # a DelayedUnaryIsoOpStack object
## ---------------------------------------------------------------------
## PROPAGATION OF SPARSITY
## ---------------------------------------------------------------------
sm0 <- sparseMatrix(i=c(1, 4), j=c(1, 3), x=c(11, 43), dims=4:3)
SM0 <- DelayedArray(sm0)
showtree(SM0)
is_sparse(SM0) # TRUE
M1 <- SM0 - 11
showtree(M1)
class(M1@seed) # a DelayedUnaryIsoOpStack object
is_sparse(M1@seed) # FALSE
SM2 <- 10 * SM0
showtree(SM2)
class(SM2@seed) # a DelayedUnaryIsoOpStack object
is_sparse(SM2@seed) # TRUE
M3 <- SM0 / 0
showtree(M3)
class(M3@seed) # a DelayedUnaryIsoOpStack object
is_sparse(M3@seed) # FALSE
SM4 <- log(1 + SM0) / 10
showtree(SM4)
class(SM4@seed) # a DelayedUnaryIsoOpStack object
is_sparse(SM4@seed) # TRUE
SM5 <- 2 ^ SM0 - 1
showtree(SM5)
class(SM5@seed) # a DelayedUnaryIsoOpStack object
is_sparse(SM5@seed) # TRUE
## ---------------------------------------------------------------------
## SANITY CHECKS
## ---------------------------------------------------------------------
stopifnot(class(M@seed) == "DelayedUnaryIsoOpStack")
stopifnot(class(M1@seed) == "DelayedUnaryIsoOpStack")
stopifnot(!is_sparse(M1@seed))
stopifnot(class(SM2@seed) == "DelayedUnaryIsoOpStack")
stopifnot(is_sparse(SM2@seed))
stopifnot(class(M3@seed) == "DelayedUnaryIsoOpStack")
stopifnot(!is_sparse(M3@seed))
stopifnot(class(SM4@seed) == "DelayedUnaryIsoOpStack")
stopifnot(is_sparse(SM4@seed))
stopifnot(class(SM5@seed) == "DelayedUnaryIsoOpStack")
stopifnot(is_sparse(SM5@seed))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.