SparseArray-misc-methods | R Documentation |
This man page documents various base array operations that are supported by SparseArray derivatives, and that didn't belong to any of the groups of operations documented in the other man pages of the SparseArray package.
# --- unary isometric array transformations ---
## S4 method for signature 'COO_SparseArray'
is.na(x)
## S4 method for signature 'SVT_SparseArray'
is.na(x)
## S4 method for signature 'COO_SparseArray'
is.nan(x)
## S4 method for signature 'SVT_SparseArray'
is.nan(x)
## S4 method for signature 'COO_SparseArray'
is.infinite(x)
## S4 method for signature 'SVT_SparseArray'
is.infinite(x)
## S4 method for signature 'COO_SparseArray'
tolower(x)
## S4 method for signature 'COO_SparseArray'
toupper(x)
## S4 method for signature 'COO_SparseArray'
nchar(x, type="chars", allowNA=FALSE, keepNA=NA)
# --- N-ary isometric array transformations ---
## S4 method for signature 'SparseArray'
pmin(..., na.rm=FALSE)
## S4 method for signature 'SparseArray'
pmax(..., na.rm=FALSE)
x |
A SparseArray derivative. |
type , allowNA , keepNA |
See |
... |
SparseArray derivatives. |
na.rm |
See |
More operations will be added in the future.
See man pages of the corresponding base functions (e.g.
?base::is.na
, ?base::nchar
,
?base::pmin
, etc...) for the value returned by
these methods.
Note that, like the base functions, the methods documented in this man page are endomorphisms i.e. they return an array-like object of the same class as the input.
base::is.na
and
base::is.infinite
in base R.
base::tolower
in base R.
base::nchar
in base R.
base::pmin
in base R.
SparseArray objects.
Ordinary array objects in base R.
a <- array(c(0, 2.77, NA, 0, NaN, -Inf), dim=5:3)
svt <- SparseArray(a) # SVT_SparseArray object
class(svt)
is.na(svt) # SVT_SparseArray object of type "logical"
is.nan(svt) # SVT_SparseArray object of type "logical"
is.infinite(svt) # SVT_SparseArray object of type "logical"
svt1 <- poissonSparseMatrix(500, 20, density=0.2)
svt2 <- poissonSparseMatrix(500, 20, density=0.25) * 0.77
pmin(svt1, svt2)
pmax(svt1, svt2)
## Sanity checks:
res <- is.na(svt)
stopifnot(is(res, "SVT_SparseArray"), type(res) == "logical",
identical(as.array(res), is.na(a)))
res <- is.nan(svt)
stopifnot(is(res, "SVT_SparseArray"), type(res) == "logical",
identical(as.array(res), is.nan(a)))
res <- is.infinite(svt)
stopifnot(is(res, "SVT_SparseArray"), type(res) == "logical",
identical(as.array(res), is.infinite(a)))
res <- pmin(svt1, svt2)
stopifnot(is(res, "SVT_SparseArray"),
identical(as.array(res), pmin(as.array(svt1), as.array(svt2))))
res <- pmax(svt1, svt2)
stopifnot(is(res, "SVT_SparseArray"),
identical(as.array(res), pmax(as.array(svt1), as.array(svt2))))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.