| fastMisc {Matrix} | R Documentation |
“Low Level” Coercions and Methods
Description
“Semi-API” functions used internally by Matrix,
often to bypass S4 dispatch and avoid the associated overhead.
These are exported to provide this capability to expert users.
Typical users should continue to rely on S4 generic functions
to dispatch suitable methods, by calling,
e.g., as(., <class>) for coercions.
Usage
.M2kind(from, kind = ".", sparse = NA)
.M2gen(from, kind = ".")
.M2sym(from, ...)
.M2tri(from, ...)
.M2diag(from)
.M2v(from)
.M2m(from)
.M2unpacked(from)
.M2packed(from)
.M2C(from)
.M2R(from)
.M2T(from)
.M2V(from)
.m2V(from, kind = ".")
.sparse2dense(from, packed = FALSE)
.diag2dense(from, kind = ".", shape = "t", packed = FALSE, uplo = "U")
.ind2dense(from, kind = "n")
.m2dense(from, class = ".ge", uplo = "U", diag = "N", trans = FALSE)
.dense2sparse(from, repr = "C")
.diag2sparse(from, kind = ".", shape = "t", repr = "C", uplo = "U")
.ind2sparse(from, kind = "n", repr = ".")
.m2sparse(from, class = ".gC", uplo = "U", diag = "N", trans = FALSE)
.tCRT(x, lazy = TRUE)
.diag.dsC(x, Chx = Cholesky(x, LDL = TRUE), res.kind = "diag")
.solve.dgC.lu (a, b, tol = .Machine$double.eps, check = TRUE)
.solve.dgC.qr (a, b, order = 3L, check = TRUE)
.solve.dgC.chol(a, b, check = TRUE)
.updateCHMfactor(object, parent, mult = 0)
Arguments
from, x, a, b |
a |
kind |
a string ( |
shape |
a string ( |
repr |
a string ( |
packed |
a logical indicating if the result should
inherit from |
sparse |
a logical indicating if the result should inherit
from |
uplo |
a string ( |
diag |
a string ( |
trans |
a logical indicating if the result should be a 1-row
matrix rather than a 1-column matrix where |
class |
a string whose first three characters specify the class
of the result. It should match the pattern
|
... |
optional arguments passed to |
lazy |
a logical indicating if the transpose should be constructed with minimal allocation, but possibly without preserving representation. |
Chx |
optionally, the |
res.kind |
a string in |
tol |
see |
order |
see |
check |
a logical indicating if the first argument should be
tested for inheritance from |
object |
a Cholesky factorization inheriting from virtual class
|
parent |
|
mult |
a numeric vector of postive length. Only the first element is used, and that must be finite. |
Details
Functions with names of the form .<A>2<B> implement coercions
from virtual class A to the “nearest” non-virtual subclass of
virtual class B, where the virtual classes are abbreviated as follows:
MVmmatrix
vvector
denseunpackedpackedsparseCRTgensymtridiagind
Abbreviations should be seen as a guide, rather than as an
exact description of behaviour. Notably, .m2dense,
.m2sparse, and .m2V accept vectors that are
not matrices.
.tCRT(x)
If lazy = TRUE, then .tCRT constructs the transpose
of x using the most efficient representation,
which for ‘CRT’ is ‘RCT’. If lazy = FALSE,
then .tCRT preserves the representation of x,
behaving as the corresponding methods for generic function t.
.diag.dsC(x)
.diag.dsC computes (or uses if Chx is supplied)
the Cholesky factorization of x as L D L' in order
to calculate one of several possible statistics from the diagonal
entries of D. See res.kind under ‘Arguments’.
.solve.dgC.*(a, b)
.solve.dgC.lu(a, b) needs a square matrix a.
.solve.dgC.qr(a, b) needs a “long” matrix a,
with nrow(a) >= ncol(a).
.solve.dgC.chol(a, b) needs a “wide” matrix a,
with nrow(a) <= ncol(a).
All three may be used to solve sparse linear systems directly.
Only .solve.dgC.qr and .solve.dgC.chol be used
to solve sparse least squares problems.
.updateCHMfactor(object, parent, mult)
.updateCHMfactor updates object with the result
of Cholesky factorizing
F(parent) + mult[1] * diag(nrow(parent)),
i.e., F(parent) plus mult[1] times the identity matrix,
where F = identity if parent is a dsCMatrix
and F = tcrossprod if parent is a dgCMatrix.
The nonzero pattern of F(parent) must match
that of S if object = Cholesky(S, ...).
Examples
D. <- diag(x = c(1, 1, 2, 3, 5, 8))
D.0 <- Diagonal(x = c(0, 0, 0, 3, 5, 8))
S. <- toeplitz(as.double(1:6))
C. <- new("dgCMatrix", Dim = c(3L, 4L),
p = c(0L, 1L, 1L, 1L, 3L), i = c(1L, 0L, 2L), x = c(-8, 2, 3))
stopifnot(exprs = {
identical(.M2tri (D.), as(D., "triangularMatrix"))
identical(.M2sym (D.), as(D., "symmetricMatrix"))
identical(.M2diag(D.), as(D., "diagonalMatrix"))
identical(.M2kind(C., "l"),
as(C., "lMatrix"))
identical(.M2kind(.sparse2dense(C.), "l"),
as(as(C., "denseMatrix"), "lMatrix"))
identical(.diag2sparse(D.0, ".", "t", "C"),
.dense2sparse(.diag2dense(D.0, ".", "t", TRUE), "C"))
identical(.M2gen(.diag2dense(D.0, ".", "s", FALSE)),
.sparse2dense(.M2gen(.diag2sparse(D.0, ".", "s", "T"))))
identical(S.,
.M2m(.m2sparse(S., ".sR")))
identical(S. * lower.tri(S.) + diag(1, 6L),
.M2m(.m2dense (S., ".tr", "L", "U")))
identical(.M2R(C.), .M2R(.M2T(C.)))
identical(.tCRT(C.), .M2R(t(C.)))
})
A <- tcrossprod(C.)/6 + Diagonal(3, 1/3); A[1,2] <- 3; A
stopifnot(exprs = {
is.numeric( x. <- c(2.2, 0, -1.2) )
all.equal(x., .solve.dgC.lu(A, c(1,0,0), check=FALSE))
all.equal(x., .solve.dgC.qr(A, c(1,0,0), check=FALSE))
})
## Solving sparse least squares:
X <- rbind(A, Diagonal(3)) # design matrix X (for L.S.)
Xt <- t(X) # *transposed* X (for L.S.)
(y <- drop(crossprod(Xt, 1:3)) + c(-1,1)/1000) # small rand.err.
str(solveCh <- .solve.dgC.chol(Xt, y, check=FALSE)) # Xt *is* dgC..
stopifnot(exprs = {
all.equal(solveCh$coef, 1:3, tol = 1e-3)# rel.err ~ 1e-4
all.equal(solveCh$coef, drop(solve(tcrossprod(Xt), Xt %*% y)))
all.equal(solveCh$coef, .solve.dgC.qr(X, y, check=FALSE))
})