Title: | Delayed Arrays of Random Values |
---|---|
Description: | Implements a DelayedArray of random values where the realization of the sampled values is delayed until they are needed. Reproducible sampling within any subarray is achieved by chunking where each chunk is initialized with a different random seed and stream. The usual distributions in the stats package are supported, along with scalar, vector and arrays for the parameters. |
Authors: | Aaron Lun [aut, cre] |
Maintainer: | Aaron Lun <[email protected]> |
License: | GPL-3 |
Version: | 1.15.0 |
Built: | 2024-10-30 05:26:13 UTC |
Source: | https://github.com/bioc/DelayedRandomArray |
The RandomArraySeed is a DelayedArray seed that performs reproducible, on-demand sampling of randomly distributed values. Note that this is a virtual class; the intention is to define concrete subclasses corresponding to specific parameterized distributions.
The array is conceptually partitioned into contiguous chunks of the same shape. The random values in each chunk are initialized with a different seed and stream via the PCG32 pseudo-random number generator (see the dqrng package for details). This design allows us to rapidly access any given subarray without having to do jump-aheads from the start of the stream.
The default chunking dimensions are set to the square root of the array dimensions - or 100, whichever is larger.
This scheme provides decent though suboptimal performance along any dimension.
If the access pattern is known beforehand, a better chunking scheme can often be chosen and passed to the chunkdim
argument.
Note that changing the chunking dimensions will change the ordering of array values, even if the seeds are unchanged. This may be unexpected, given that chunking in real datasets will never change the data, only the performance of access operations. However, it is largely unavoidable in this context as the random number stream is rearranged within the array.
The chunkdim(x)
method will return the chunk dimensions of a RandomArraySeed instance x
.
This will be used by the DelayedArray machinery to optimize block processing by extracting whole chunks where possible.
To sample from a specific distribution, we can implement a concrete subclass of the RandomArraySeed.
This is done by implementing methods for sampleDistrFun
and sampleDistrParam
.
In the code chunks below, x
is an instance of a RandomArraySeed subclass:
sampleDistrFun(x)
returns a quantile function that accepts a vector of cumulative probabilities p
and returns a numeric vector of quantiles.
A typical example is qnorm
, though similar functions from the stats package can also be used.
The output vector should be the same length as p
; any other distributional parameters should be recycled to the length of p
.
sampleDistrParam(x)
returns a character vector specifying the names of the distributional parameters as slots of x
.
For example, for a subclass that samples from a normal distribution, this might be "mean"
and "sd"
.
Each distributional parameter is expected to be numeric.
The extract_array
method for the RandomArraySeed will automatically use both of the above methods to sample from the specified distribution.
This is achieved by randomly sampling from a standard uniform distribution, treating the values as probabilities and converting them into quantiles.
Distributional parameters are passed to the relevant quantile function to obtain a random value from the desired distribution. Each parameter can be:
A numeric scalar, which is used throughout the array.
A numeric vector, which is recycled along the length of the array. This traverses the array along the first dimension, then the second, then the third, and so on; for matrices, this is equivalent to column-major ordering.
A numeric array-like object of the same dimensions as dim
,
where each entry contains the parameter value for the corresponding entry of the output array.
This can be another DelayedArray object.
For certain distributions, we may expect a large number of zeroes in the random output.
We provide the option to treat the sampled values as being sparse, by setting sparse=TRUE
in the constructors of the relevant subclasses.
This is optional as most distributions will not yield sparse arrays for most of their possible parameter space.
When sparse=TRUE
, the block processing machinery in DelayedArray will return a sparse array.
This gives downstream applications the opportunity to use more efficient sparse algorithms when relevant.
However, this option does not affect the sampling itself; the result is always the same as a dense array, just that the output is coerced into a SVT_SparseArray.
We can determine whether a RandomArraySeed x
has a sparse interpretation with is_sparse(x)
.
Aaron Lun
The RandomUnifArraySeed class, which implements sampling from a uniform distribution.
The RandomPoisArraySeed class, which implements sampling from a Poisson distribution.
A DelayedArray subclass that performs on-the-fly sampling of beta-distributed values.
RandomBetaArraySeed(dim, shape1, shape2, ncp = 0, chunkdim = NULL) ## S4 method for signature 'RandomBetaArraySeed' DelayedArray(seed) RandomBetaArray(dim, shape1, shape2, ncp = 0, chunkdim = NULL)
RandomBetaArraySeed(dim, shape1, shape2, ncp = 0, chunkdim = NULL) ## S4 method for signature 'RandomBetaArraySeed' DelayedArray(seed) RandomBetaArray(dim, shape1, shape2, ncp = 0, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
shape1 , shape2 , ncp
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
seed |
A RandomBetaArraySeed object. |
All constructors return an instance of a RandomBetaArray object, containing random draws from a beta distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomBetaArraySeed(c(1e5, 1e5), shape1=1, shape2=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomBetaArraySeed(c(1e5, 1e5), shape1=runif(1e5), shape2=2) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) s1 <- rsparsematrix(1e5, 1e5, density=0.00001) s1 <- abs(DelayedArray(s1)) + 1 X3 <- RandomBetaArraySeed(c(1e5, 1e5), shape1=s1, shape2=s1+1) Y3 <- DelayedArray(X3) Y3
X <- RandomBetaArraySeed(c(1e5, 1e5), shape1=1, shape2=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomBetaArraySeed(c(1e5, 1e5), shape1=runif(1e5), shape2=2) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) s1 <- rsparsematrix(1e5, 1e5, density=0.00001) s1 <- abs(DelayedArray(s1)) + 1 X3 <- RandomBetaArraySeed(c(1e5, 1e5), shape1=s1, shape2=s1+1) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of binomial-distributed values.
RandomBinomArraySeed(dim, size, prob, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomBinomArraySeed' DelayedArray(seed) RandomBinomArray(dim, size, prob, chunkdim = NULL, sparse = FALSE)
RandomBinomArraySeed(dim, size, prob, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomBinomArraySeed' DelayedArray(seed) RandomBinomArray(dim, size, prob, chunkdim = NULL, sparse = FALSE)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
size , prob
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
sparse |
Logical scalar indicating whether the sampled array should be treated as sparse. |
seed |
A RandomBinomArraySeed object. |
All constructors return an instance of a RandomBinomArray object, containing random draws from a binomial distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomBinomArraySeed(c(1e5, 1e5), size=10, prob=0.5) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomBinomArraySeed(c(1e5, 1e5), size=10, prob=runif(1e5, 0, 0.1), sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) size <- rsparsematrix(1e5, 1e5, density=0.00001) size <- round(abs(DelayedArray(size)) * 10) X3 <- RandomBinomArraySeed(c(1e5, 1e5), size=size, prob=0.5) Y3 <- DelayedArray(X3) Y3
X <- RandomBinomArraySeed(c(1e5, 1e5), size=10, prob=0.5) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomBinomArraySeed(c(1e5, 1e5), size=10, prob=runif(1e5, 0, 0.1), sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) size <- rsparsematrix(1e5, 1e5, density=0.00001) size <- round(abs(DelayedArray(size)) * 10) X3 <- RandomBinomArraySeed(c(1e5, 1e5), size=size, prob=0.5) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of Cauchy-distributed values.
RandomCauchyArraySeed(dim, location = 0, scale = 1, chunkdim = NULL) ## S4 method for signature 'RandomCauchyArraySeed' DelayedArray(seed) RandomCauchyArray(dim, location = 0, scale = 1, chunkdim = NULL)
RandomCauchyArraySeed(dim, location = 0, scale = 1, chunkdim = NULL) ## S4 method for signature 'RandomCauchyArraySeed' DelayedArray(seed) RandomCauchyArray(dim, location = 0, scale = 1, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
location , scale
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
seed |
A RandomCauchyArraySeed object. |
All constructors return an instance of a RandomCauchyArray object, containing random draws from a Cauchy distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomCauchyArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomCauchyArraySeed(c(1e5, 1e5), location=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) loc <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomCauchyArraySeed(c(1e5, 1e5), location=loc) Y3 <- DelayedArray(X3) Y3
X <- RandomCauchyArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomCauchyArraySeed(c(1e5, 1e5), location=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) loc <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomCauchyArraySeed(c(1e5, 1e5), location=loc) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of chi-squared-distributed values.
RandomChisqArraySeed(dim, df, ncp = 0, chunkdim = NULL) ## S4 method for signature 'RandomChisqArraySeed' DelayedArray(seed) RandomChisqArray(dim, df, ncp = 0, chunkdim = NULL)
RandomChisqArraySeed(dim, df, ncp = 0, chunkdim = NULL) ## S4 method for signature 'RandomChisqArraySeed' DelayedArray(seed) RandomChisqArray(dim, df, ncp = 0, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
df , ncp
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
seed |
A RandomChisqArraySeed object. |
All constructors return an instance of a RandomChisqArray object, containing random draws from a chi-squared distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomChisqArraySeed(c(1e5, 1e5), df=5) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomChisqArraySeed(c(1e5, 1e5), df=runif(1e5)*20) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) df <- rsparsematrix(1e5, 1e5, density=0.00001) df <- abs(DelayedArray(df) + 1) * 10 X3 <- RandomChisqArraySeed(c(1e5, 1e5), df=df) Y3 <- DelayedArray(X3) Y3
X <- RandomChisqArraySeed(c(1e5, 1e5), df=5) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomChisqArraySeed(c(1e5, 1e5), df=runif(1e5)*20) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) df <- rsparsematrix(1e5, 1e5, density=0.00001) df <- abs(DelayedArray(df) + 1) * 10 X3 <- RandomChisqArraySeed(c(1e5, 1e5), df=df) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of exponentially distributed values.
RandomExpArraySeed(dim, rate = 1, chunkdim = NULL) ## S4 method for signature 'RandomExpArraySeed' DelayedArray(seed) RandomExpArray(dim, rate = 1, chunkdim = NULL)
RandomExpArraySeed(dim, rate = 1, chunkdim = NULL) ## S4 method for signature 'RandomExpArraySeed' DelayedArray(seed) RandomExpArray(dim, rate = 1, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
rate |
Numeric vector used as |
chunkdim |
Integer vector of length equal to |
seed |
A RandomExpArraySeed object. |
All constructors return an instance of a RandomExpArray object, containing random draws from a exponential distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomExpArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomExpArraySeed(c(1e5, 1e5), rate=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) rate <- rsparsematrix(1e5, 1e5, density=0.00001) rate <- abs(DelayedArray(rate)) + 1 X3 <- RandomExpArraySeed(c(1e5, 1e5), rate=rate) Y3 <- DelayedArray(X3) Y3
X <- RandomExpArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomExpArraySeed(c(1e5, 1e5), rate=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) rate <- rsparsematrix(1e5, 1e5, density=0.00001) rate <- abs(DelayedArray(rate)) + 1 X3 <- RandomExpArraySeed(c(1e5, 1e5), rate=rate) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of F-distributed values.
RandomFArraySeed(dim, df1, df2, ncp, chunkdim = NULL) ## S4 method for signature 'RandomFArraySeed' DelayedArray(seed) RandomFArray(dim, df1, df2, ncp, chunkdim = NULL)
RandomFArraySeed(dim, df1, df2, ncp, chunkdim = NULL) ## S4 method for signature 'RandomFArraySeed' DelayedArray(seed) RandomFArray(dim, df1, df2, ncp, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
df1 , df2 , ncp
|
Numeric vector used as the argument of the same name in If |
chunkdim |
Integer vector of length equal to |
seed |
A RandomFArraySeed object. |
All constructors return an instance of a RandomFArray object, containing random draws from a exponential distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomFArraySeed(c(1e5, 1e5), df1=1, df2=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomFArraySeed(c(1e5, 1e5), df1=runif(1e5)*10, df2=10) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) ncp <- rsparsematrix(1e5, 1e5, density=0.00001) ncp <- abs(DelayedArray(ncp)) + 1 X3 <- RandomFArraySeed(c(1e5, 1e5), df1=1, df2=10, ncp=ncp) Y3 <- DelayedArray(X3) Y3
X <- RandomFArraySeed(c(1e5, 1e5), df1=1, df2=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomFArraySeed(c(1e5, 1e5), df1=runif(1e5)*10, df2=10) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) ncp <- rsparsematrix(1e5, 1e5, density=0.00001) ncp <- abs(DelayedArray(ncp)) + 1 X3 <- RandomFArraySeed(c(1e5, 1e5), df1=1, df2=10, ncp=ncp) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of gamma-distributed values.
RandomGammaArraySeed(dim, shape, rate = 1, scale = 1/rate, chunkdim = NULL) ## S4 method for signature 'RandomGammaArraySeed' DelayedArray(seed) RandomGammaArray(dim, shape, rate = 1, scale = 1/rate, chunkdim = NULL)
RandomGammaArraySeed(dim, shape, rate = 1, scale = 1/rate, chunkdim = NULL) ## S4 method for signature 'RandomGammaArraySeed' DelayedArray(seed) RandomGammaArray(dim, shape, rate = 1, scale = 1/rate, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
shape , rate , scale
|
Numeric vector used as the argument of the same name in If |
chunkdim |
Integer vector of length equal to |
seed |
A RandomGammaArraySeed object. |
All constructors return an instance of a RandomGammaArray object, containing random draws from a gamma distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomGammaArraySeed(c(1e5, 1e5), shape=1, rate=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomGammaArraySeed(c(1e5, 1e5), shape=runif(1e5), rate=2) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) s1 <- rsparsematrix(1e5, 1e5, density=0.00001) s1 <- abs(DelayedArray(s1)) + 1 X3 <- RandomGammaArraySeed(c(1e5, 1e5), shape=s1, rate=s1+1) Y3 <- DelayedArray(X3) Y3
X <- RandomGammaArraySeed(c(1e5, 1e5), shape=1, rate=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomGammaArraySeed(c(1e5, 1e5), shape=runif(1e5), rate=2) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) s1 <- rsparsematrix(1e5, 1e5, density=0.00001) s1 <- abs(DelayedArray(s1)) + 1 X3 <- RandomGammaArraySeed(c(1e5, 1e5), shape=s1, rate=s1+1) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of geometric-distributed values.
RandomGeomArraySeed(dim, prob, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomGeomArraySeed' DelayedArray(seed) RandomGeomArray(dim, prob, chunkdim = NULL, sparse = FALSE)
RandomGeomArraySeed(dim, prob, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomGeomArraySeed' DelayedArray(seed) RandomGeomArray(dim, prob, chunkdim = NULL, sparse = FALSE)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
prob |
Numeric vector used as |
chunkdim |
Integer vector of length equal to |
sparse |
Logical scalar indicating whether the sampled array should be treated as sparse. |
seed |
A RandomGeomArraySeed object. |
All constructors return an instance of a RandomGeomArray object, containing random draws from a geometric distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomGeomArraySeed(c(1e5, 1e5), prob=0.5) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomGeomArraySeed(c(1e5, 1e5), prob=runif(1e5, 0, 0.1), sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) prob <- RandomUnifArray(c(1e5, 1e5)) X3 <- RandomGeomArraySeed(c(1e5, 1e5), prob=prob) Y3 <- DelayedArray(X3) Y3
X <- RandomGeomArraySeed(c(1e5, 1e5), prob=0.5) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomGeomArraySeed(c(1e5, 1e5), prob=runif(1e5, 0, 0.1), sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) prob <- RandomUnifArray(c(1e5, 1e5)) X3 <- RandomGeomArraySeed(c(1e5, 1e5), prob=prob) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of hypergeometric-distributed values.
RandomHyperArraySeed(dim, m, n, k, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomHyperArraySeed' DelayedArray(seed) RandomHyperArray(dim, m, n, k, chunkdim = NULL, sparse = FALSE)
RandomHyperArraySeed(dim, m, n, k, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomHyperArraySeed' DelayedArray(seed) RandomHyperArray(dim, m, n, k, chunkdim = NULL, sparse = FALSE)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
m , n , k
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
sparse |
Logical scalar indicating whether the sampled array should be treated as sparse. |
seed |
A RandomHyperArraySeed object. |
All constructors return an instance of a RandomHyperArray object, containing random draws from a hypergeometric distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomHyperArraySeed(c(1e5, 1e5), m=10, n=20, k=15) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomHyperArraySeed(c(1e5, 1e5), m=round(runif(1e5, 10, 20)), n=20, k=15, sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) m <- round(RandomUnifArray(c(1e5, 1e5), 10, 20)) X3 <- RandomHyperArraySeed(c(1e5, 1e5), m=m, n=50, k=20) Y3 <- DelayedArray(X3) Y3
X <- RandomHyperArraySeed(c(1e5, 1e5), m=10, n=20, k=15) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomHyperArraySeed(c(1e5, 1e5), m=round(runif(1e5, 10, 20)), n=20, k=15, sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) m <- round(RandomUnifArray(c(1e5, 1e5), 10, 20)) X3 <- RandomHyperArraySeed(c(1e5, 1e5), m=m, n=50, k=20) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of log-normally distributed values.
RandomLnormArraySeed(dim, meanlog = 0, sdlog = 1, chunkdim = NULL) ## S4 method for signature 'RandomLnormArraySeed' DelayedArray(seed) RandomLnormArray(dim, meanlog = 0, sdlog = 1, chunkdim = NULL)
RandomLnormArraySeed(dim, meanlog = 0, sdlog = 1, chunkdim = NULL) ## S4 method for signature 'RandomLnormArraySeed' DelayedArray(seed) RandomLnormArray(dim, meanlog = 0, sdlog = 1, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
meanlog , sdlog
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
seed |
A RandomLnormArraySeed object. |
All constructors return an instance of a RandomLnormArray object, containing random draws from a log-normal distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomLnormArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomLnormArraySeed(c(1e5, 1e5), meanlog=runif(1e5), sdlog=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) meanlog <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomLnormArraySeed(c(1e5, 1e5), meanlog=meanlog) Y3 <- DelayedArray(X3) Y3
X <- RandomLnormArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomLnormArraySeed(c(1e5, 1e5), meanlog=runif(1e5), sdlog=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) meanlog <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomLnormArraySeed(c(1e5, 1e5), meanlog=meanlog) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of log-normally distributed values.
RandomLogisArraySeed(dim, location = 0, scale = 1, chunkdim = NULL) ## S4 method for signature 'RandomLogisArraySeed' DelayedArray(seed) RandomLogisArray(dim, location = 0, scale = 1, chunkdim = NULL)
RandomLogisArraySeed(dim, location = 0, scale = 1, chunkdim = NULL) ## S4 method for signature 'RandomLogisArraySeed' DelayedArray(seed) RandomLogisArray(dim, location = 0, scale = 1, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
location , scale
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
seed |
A RandomLogisArraySeed object. |
All constructors return an instance of a RandomLogisArray object, containing random draws from a log-normal distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomLogisArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomLogisArraySeed(c(1e5, 1e5), location=runif(1e5), scale=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) location <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomLogisArraySeed(c(1e5, 1e5), location=location) Y3 <- DelayedArray(X3) Y3
X <- RandomLogisArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomLogisArraySeed(c(1e5, 1e5), location=runif(1e5), scale=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) location <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomLogisArraySeed(c(1e5, 1e5), location=location) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of negative binomial-distributed values.
RandomNbinomArraySeed( dim, prob = prob, size = size, mu = mu, chunkdim = NULL, sparse = FALSE ) ## S4 method for signature 'RandomNbinomArraySeed' DelayedArray(seed) RandomNbinomArray(dim, prob, size, mu, chunkdim = NULL, sparse = FALSE)
RandomNbinomArraySeed( dim, prob = prob, size = size, mu = mu, chunkdim = NULL, sparse = FALSE ) ## S4 method for signature 'RandomNbinomArraySeed' DelayedArray(seed) RandomNbinomArray(dim, prob, size, mu, chunkdim = NULL, sparse = FALSE)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
prob , size , mu
|
Numeric vector used as the argument of the same name in Exactly one of |
chunkdim |
Integer vector of length equal to |
sparse |
Logical scalar indicating whether the sampled array should be treated as sparse. |
seed |
A RandomNbinomArraySeed object. |
All constructors return an instance of a RandomNbinomArray object, containing random draws from a negative binomial distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomNbinomArraySeed(c(1e5, 1e5), size=10, mu=20) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomNbinomArraySeed(c(1e5, 1e5), size=10, mu=runif(1e5), sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) lambda <- rsparsematrix(1e5, 1e5, density=0.00001) lambda <- abs(DelayedArray(lambda)) + 0.1 X3 <- RandomNbinomArraySeed(c(1e5, 1e5), size=1, mu=lambda) Y3 <- DelayedArray(X3) Y3
X <- RandomNbinomArraySeed(c(1e5, 1e5), size=10, mu=20) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomNbinomArraySeed(c(1e5, 1e5), size=10, mu=runif(1e5), sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) lambda <- rsparsematrix(1e5, 1e5, density=0.00001) lambda <- abs(DelayedArray(lambda)) + 0.1 X3 <- RandomNbinomArraySeed(c(1e5, 1e5), size=1, mu=lambda) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of normally distributed values.
RandomNormArraySeed(dim, mean = 0, sd = 1, chunkdim = NULL) ## S4 method for signature 'RandomNormArraySeed' DelayedArray(seed) RandomNormArray(dim, mean = 0, sd = 1, chunkdim = NULL)
RandomNormArraySeed(dim, mean = 0, sd = 1, chunkdim = NULL) ## S4 method for signature 'RandomNormArraySeed' DelayedArray(seed) RandomNormArray(dim, mean = 0, sd = 1, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
mean , sd
|
Numeric vector used as |
chunkdim |
Integer vector of length equal to |
seed |
A RandomNormArraySeed object. |
All constructors return an instance of a RandomNormArray object, containing random draws from a normal distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomNormArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomNormArraySeed(c(1e5, 1e5), mean=runif(1e5), sd=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) mean <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomNormArraySeed(c(1e5, 1e5), mean=mean) Y3 <- DelayedArray(X3) Y3
X <- RandomNormArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomNormArraySeed(c(1e5, 1e5), mean=runif(1e5), sd=runif(1e5)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) mean <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomNormArraySeed(c(1e5, 1e5), mean=mean) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of Poisson-distributed values.
RandomPoisArraySeed(dim, lambda, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomPoisArraySeed' DelayedArray(seed) RandomPoisArray(dim, lambda, chunkdim = NULL, sparse = FALSE)
RandomPoisArraySeed(dim, lambda, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomPoisArraySeed' DelayedArray(seed) RandomPoisArray(dim, lambda, chunkdim = NULL, sparse = FALSE)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
lambda |
Numeric vector used as |
chunkdim |
Integer vector of length equal to |
sparse |
Logical scalar indicating whether the sampled array should be treated as sparse. |
seed |
A RandomPoisArraySeed object. |
All constructors return an instance of a RandomPoisArray object, containing random draws from a Poisson distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomPoisArraySeed(c(1e5, 1e5), lambda=2) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomPoisArraySeed(c(1e5, 1e5), lambda=runif(1e5), sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) lambda <- rsparsematrix(1e5, 1e5, density=0.00001) lambda <- abs(DelayedArray(lambda)) + 0.1 X3 <- RandomPoisArraySeed(c(1e5, 1e5), lambda=lambda) Y3 <- DelayedArray(X3) Y3
X <- RandomPoisArraySeed(c(1e5, 1e5), lambda=2) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomPoisArraySeed(c(1e5, 1e5), lambda=runif(1e5), sparse=TRUE) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) lambda <- rsparsematrix(1e5, 1e5, density=0.00001) lambda <- abs(DelayedArray(lambda)) + 0.1 X3 <- RandomPoisArraySeed(c(1e5, 1e5), lambda=lambda) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of F-distributed values.
RandomTArraySeed(dim, df, ncp, chunkdim = NULL) ## S4 method for signature 'RandomTArraySeed' DelayedArray(seed) RandomTArray(dim, df, ncp, chunkdim = NULL)
RandomTArraySeed(dim, df, ncp, chunkdim = NULL) ## S4 method for signature 'RandomTArraySeed' DelayedArray(seed) RandomTArray(dim, df, ncp, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
df , ncp
|
Numeric vector used as the argument of the same name in If |
chunkdim |
Integer vector of length equal to |
seed |
A RandomTArraySeed object. |
All constructors return an instance of a RandomTArray object, containing random draws from a exponential distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomTArraySeed(c(1e5, 1e5), df=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomTArraySeed(c(1e5, 1e5), df=sample(20, 1e5, replace=TRUE)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) ncp <- rsparsematrix(1e5, 1e5, density=0.00001) ncp <- abs(DelayedArray(ncp)) + 1 X3 <- RandomTArraySeed(c(1e5, 1e5), df=10, ncp=ncp) Y3 <- DelayedArray(X3) Y3
X <- RandomTArraySeed(c(1e5, 1e5), df=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomTArraySeed(c(1e5, 1e5), df=sample(20, 1e5, replace=TRUE)) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) ncp <- rsparsematrix(1e5, 1e5, density=0.00001) ncp <- abs(DelayedArray(ncp)) + 1 X3 <- RandomTArraySeed(c(1e5, 1e5), df=10, ncp=ncp) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of uniformly distributed values.
RandomUnifArraySeed(dim, min = 0, max = 1, chunkdim = NULL) ## S4 method for signature 'RandomUnifArraySeed' DelayedArray(seed) RandomUnifArray(dim, min = 0, max = 1, chunkdim = NULL)
RandomUnifArraySeed(dim, min = 0, max = 1, chunkdim = NULL) ## S4 method for signature 'RandomUnifArraySeed' DelayedArray(seed) RandomUnifArray(dim, min = 0, max = 1, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
min , max
|
Numeric vector used as |
chunkdim |
Integer vector of length equal to |
seed |
A RandomUnifArraySeed object. |
All constructors return an instance of a RandomUnifArray object, containing random draws from a uniform distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomUnifArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomUnifArraySeed(c(1e5, 1e5), min=1:1e5, max=1:1e5*2) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) min <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomUnifArraySeed(c(1e5, 1e5), min=min, max=DelayedArray(min)+1) Y3 <- DelayedArray(X3) Y3
X <- RandomUnifArraySeed(c(1e5, 1e5)) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomUnifArraySeed(c(1e5, 1e5), min=1:1e5, max=1:1e5*2) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) min <- rsparsematrix(1e5, 1e5, density=0.00001) X3 <- RandomUnifArraySeed(c(1e5, 1e5), min=min, max=DelayedArray(min)+1) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of Weibull-distributed values.
RandomWeibullArraySeed(dim, shape, scale = 1, chunkdim = NULL) ## S4 method for signature 'RandomWeibullArraySeed' DelayedArray(seed) RandomWeibullArray(dim, shape, scale = 1, chunkdim = NULL)
RandomWeibullArraySeed(dim, shape, scale = 1, chunkdim = NULL) ## S4 method for signature 'RandomWeibullArraySeed' DelayedArray(seed) RandomWeibullArray(dim, shape, scale = 1, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
shape , scale
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
seed |
A RandomWeibullArraySeed object. |
All constructors return an instance of a RandomWeibullArray object, containing random draws from a Weibull distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomWeibullArraySeed(c(1e5, 1e5), shape=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomWeibullArraySeed(c(1e5, 1e5), shape=round(runif(1e5, 10, 20))) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) shape <- round(RandomUnifArray(c(1e5, 1e5), 10, 20)) X3 <- RandomWeibullArraySeed(c(1e5, 1e5), shape=shape) Y3 <- DelayedArray(X3) Y3
X <- RandomWeibullArraySeed(c(1e5, 1e5), shape=10) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomWeibullArraySeed(c(1e5, 1e5), shape=round(runif(1e5, 10, 20))) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) shape <- round(RandomUnifArray(c(1e5, 1e5), 10, 20)) X3 <- RandomWeibullArraySeed(c(1e5, 1e5), shape=shape) Y3 <- DelayedArray(X3) Y3
A DelayedArray subclass that performs on-the-fly sampling of Wilcoxon-distributed values.
RandomWilcoxArraySeed(dim, m, n, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomWilcoxArraySeed' DelayedArray(seed) RandomWilcoxArray(dim, m, n, chunkdim = NULL)
RandomWilcoxArraySeed(dim, m, n, chunkdim = NULL, sparse = FALSE) ## S4 method for signature 'RandomWilcoxArraySeed' DelayedArray(seed) RandomWilcoxArray(dim, m, n, chunkdim = NULL)
dim |
Integer vector of positive length, specifying the dimensions of the array. |
m , n
|
Numeric vector used as the argument of the same name in |
chunkdim |
Integer vector of length equal to |
sparse |
Logical scalar indicating whether the sampled array should be treated as sparse. |
seed |
A RandomWilcoxArraySeed object. |
All constructors return an instance of a RandomWilcoxArray object, containing random draws from a Wilcox distribution with the specified parameters.
Aaron Lun
The RandomArraySeed class, for details on chunking and the distributional parameters.
X <- RandomWilcoxArraySeed(c(1e5, 1e5), m=10, n=20) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomWilcoxArraySeed(c(1e5, 1e5), m=round(runif(1e5, 10, 20)), n=20) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) m <- round(RandomUnifArray(c(1e5, 1e5), 10, 20)) X3 <- RandomWilcoxArraySeed(c(1e5, 1e5), m=m, n=50) Y3 <- DelayedArray(X3) Y3
X <- RandomWilcoxArraySeed(c(1e5, 1e5), m=10, n=20) Y <- DelayedArray(X) Y # Fiddling with the distribution parameters: X2 <- RandomWilcoxArraySeed(c(1e5, 1e5), m=round(runif(1e5, 10, 20)), n=20) Y2 <- DelayedArray(X2) Y2 # Using another array as input: library(Matrix) m <- round(RandomUnifArray(c(1e5, 1e5), 10, 20)) X3 <- RandomWilcoxArraySeed(c(1e5, 1e5), m=m, n=50) Y3 <- DelayedArray(X3) Y3