File: denseMatrix.R

package info (click to toggle)
rmatrix 1.3-2-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 7,024 kB
  • sloc: ansic: 42,435; makefile: 330; sh: 180
file content (272 lines) | stat: -rw-r--r-- 9,981 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
### Simple fallback methods for all dense matrices
### These are "cheap" to program, but potentially far from efficient;
### Methods for specific subclasses will overwrite these:

setAs("ANY", "denseMatrix", function(from) Matrix(from, sparse=FALSE, doDiag=FALSE))
## Conceivably, could write
## setAs("matrix", "denseMatrix", ....) which was slightly more efficient than
##  Matrix(.)  but would have many things in common

setAs("denseMatrix", "generalMatrix", as_geSimpl)

## dense to sparse:
## : if we do this, do it "right", i.e. preserve symmetric/triangular!
## setAs("denseMatrix", "dsparseMatrix",
## ## MM thought that  as() will take the ``closest'' match; but that fails!
## ##      function(from) as(as(from, "dgeMatrix"), "dsparseMatrix"))
##       function(from) as(as(from, "dgeMatrix"), "dgCMatrix"))

.dense2C <- function(from, kind = NA, uplo = "U", symDimnames = FALSE) {
    useK <- is.character(kind) && length(kind) == 1 &&
        kind %in% c("gen", "sym", "tri")
    if(!useK) {
        cl <- class(from)
        cld <- getClassDef(cl) ## get it once (speedup)
    }
    r <- .Call(dense_to_Csparse, from)# goes via "generalMatrix"
    ## FIXME: for symmetric / triangular matrices, this is a waste, notably if packed
    if (useK && kind == "gen"  ||  !useK && extends(cld, "generalMatrix"))
	r
    else if(useK && kind == "sym" || !useK && extends(cld, "symmetricMatrix"))
	forceCspSymmetric(r, uplo, isTri = FALSE, symDimnames=symDimnames)
    else if(!useK && extends(cld, "diagonalMatrix"))
	stop("diagonalMatrix in .dense2C() -- should never happen, please report!")
    else { ## we have "triangular" :
        if(useK) {
            cl <- class(from)
            cld <- getClassDef(cl) ## get it once (speedup)
        }
	if	(extends(cld,"dMatrix")) as(r, "dtCMatrix")
        else if (extends(cld,"lMatrix")) as(r, "ltCMatrix")
        else if (extends(cld,"nMatrix")) as(r, "ntCMatrix")
        else if (extends(cld,"zMatrix")) as(r, "ztCMatrix")
	else stop(gettextf("undefined method for class %s", dQuote(cl)), domain=NA)
    }
}

setAs("denseMatrix", "CsparseMatrix", function(from) .dense2C(from))

## This sometimes fails (eg. for "lsyMatrix"), and we really want to
## use the generic ``go via Csparse'' (top of ./sparseMatrix.R) instead
## setAs("denseMatrix",  "sparseMatrix",
##       function(from) {
## 	  cl <- class(from)
## 	  cld <- getClassDef(cl)
## 	  if (extends(cld, "generalMatrix"))
## 	      .Call(dense_to_Csparse, from)
## 	  else ## i.e. triangular | symmetric
## 	      as_Csparse(from, cld)
##       })

setAs("denseMatrix", "TsparseMatrix",
      function(from) as(.dense2C(from), "TsparseMatrix"))


setMethod("show", signature(object = "denseMatrix"),
          function(object) prMatrix(object))
##- ## FIXME: The following is only for the "dMatrix" objects that are not
##- ##	      "dense" nor "sparse" -- i.e. "packed" ones :
##- ## But these could be printed better -- "." for structural zeros.
##- setMethod("show", signature(object = "dMatrix"), prMatrix)
##- ## and improve this as well:
##- setMethod("show", signature(object = "pMatrix"), prMatrix)
##- ## this should now be superfluous [keep for safety for the moment]:

setMethod("dim<-", signature(x = "denseMatrix", value = "ANY"),
	  function(x, value) {
	      if(!is.numeric(value) || length(value) != 2)
		  stop("dim(.) value must be numeric of length 2")
	      if(prod(dim(x)) != prod(value <- as.integer(value)))
		  stop("dimensions don't match the number of cells")
	      clx <- as.character(MatrixClass(class(x))) # as.*(): drop attr
	      if(substring(clx,2) == "geMatrix") {
		  x@Dim <- value
		  if(length(x@factors) > 0)
		      x@factors <- list()
		  x
	      } else { ## other "denseMatrix"
		  x <- as_geSimpl2(x, clx)
		  dim(x) <- value
                  x
	      }
          })



## Using "index" for indices should allow
## integer (numeric), logical, or character (names!) indices :

## use geClass() when 'i' or 'j' are missing:
## since  symmetric, triangular, .. will not be preserved anyway:
setMethod("[", signature(x = "denseMatrix", i = "index", j = "missing",
			 drop = "logical"),
	  function (x, i, j, ..., drop) {
	      if((na <- nargs()) == 3)
		  r <- as(x, "matrix")[i, drop=drop]
	      else if(na == 4)
		  r <- as(x, "matrix")[i, , drop=drop]
	      else stop(gettextf("invalid nargs()= %d", na), domain=NA)
	      if(is.null(dim(r))) r else as(r, geClass(x))
	  })

setMethod("[", signature(x = "denseMatrix", i = "missing", j = "index",
			 drop = "logical"),
	  function (x, i, j, ..., drop) {
	      r <- as(x, "matrix")[, j, drop=drop]
	      if(is.null(dim(r))) r else as(r, geClass(x))
	  })

setMethod("[", signature(x = "denseMatrix", i = "index", j = "index",
			 drop = "logical"),
	  function (x, i, j, ..., drop) {
	      r <- callGeneric(x = as(x, "matrix"), i=i, j=j, drop=drop)
	      if(is.null(dim(r)))
		  r
	      else {
		  cld <- getClassDef(cl <- class(x))
		  if(extends(cld, "symmetricMatrix") &&
		     length(i) == length(j) && isTRUE(all(i == j)))
                      ## keep original symmetric class (but not "dpo")
                      as(r, class2(cl, .M.kindC(cld)))

		  else as_smartClass(r, cl)
	      }
	  })

.dense.sub.i.2col <- function(x, i, j, ..., drop) {
    r <- as(x, "matrix")[ i ]
    if(is.null(dim(r))) r else as(r, geClass(x))
}
setMethod("[", signature(x = "denseMatrix", i = "matrix", j = "missing"),#drop="ANY"
	  .dense.sub.i.2col)
setMethod("[", signature(x = "denseMatrix", i = "matrix", j = "missing", drop="missing"),
	  .dense.sub.i.2col)


## Now the "[<-" ones --- see also those in ./Matrix.R
## It's recommended to use setReplaceMethod() rather than setMethod("[<-",.)
## even though the former is currently just a wrapper for the latter

## x[] <- value :
setReplaceMethod("[", signature(x = "denseMatrix", i = "missing", j = "missing",
				value = "ANY"),## double/logical/...
	  function (x, value) {
	      x <- as(x, "generalMatrix")
	      x@x[] <- value
	      validObject(x)# check if type and lengths above match
	      x
	  })

## FIXME: 1) These are far from efficient
## -----
setReplaceMethod("[", signature(x = "denseMatrix", i = "index", j = "missing",
				value = "replValue"),
		 function (x, i, j, ..., value) {
		     r <- as(x, "matrix")
## 		     message("`[<-` with nargs()= ",nargs())
		     if((na <- nargs()) == 3)
			 r[i] <- value
		     else if(na == 4)
			 r[i, ] <- value
		     else stop(gettextf("invalid nargs()= %d", na), domain=NA)
		     as(r, geClass(x))
		 })

setReplaceMethod("[", signature(x = "denseMatrix", i = "missing", j = "index",
				value = "replValue"),
		 function (x, i, j, ..., value) {
		     r <- as(x, "matrix")
		     r[, j] <- value
		     as(r, geClass(x))
		 })

setReplaceMethod("[", signature(x = "denseMatrix", i = "index", j = "index",
				value = "replValue"),
		 function (x, i, j, ..., value) {
		     r <- as(x, "matrix")
		     r[i, j] <- value
		     as_smartClass(r, class(x)) ## was as(r, class(x))
		 })

setReplaceMethod("[", signature(x = "denseMatrix", i = "matrix",  # 2-col.matrix
				j = "missing", value = "replValue"),
		 function(x, i, j, ..., value) {
		     r <- as(x, "matrix")
		     r[ i ] <- value
		     as(r, geClass(x))
		 })


setMethod("isSymmetric", signature(object = "denseMatrix"),
	  function(object, tol = 100*.Machine$double.eps, tol1 = 8*tol, ...) {
	      ## pretest: is it square?
	      d <- dim(object)
              if((n <- d[1L]) != d[2L]) return(FALSE)
	      if(n <= 1L) return(TRUE)
	      ## else: square (n x n) matrix, n >= 2 :
              is.z <- is(object, "zMatrix")
	      ## initial tests, fast for large non-symmetric:
	      if(length(tol1)) {
		  ## initial pre-tests, fast for large non-symmetric:
		  Cj <- if(is.z) Conj else identity
		  for(i in unique(c(1L, 2L, n-1L, n)))
		      if(is.character(all.equal(object[i, ], Cj(object[, i]),
						tolerance = tol1, ...))) return(FALSE)
	      }
	      ## else slower test
	      if (is(object,"dMatrix"))
		  isTRUE(all.equal(as(  object,  "dgeMatrix"),
				   as(t(object), "dgeMatrix"), tolerance = tol, ...))
	      else if (is(object, "nMatrix"))
		  identical(as(  object,  "ngeMatrix"),
			    as(t(object), "ngeMatrix"))
	      else if (is(object, "lMatrix"))# not possible currently
		  ## test for exact equality; FIXME(?): identical() too strict?
		  identical(as(  object,  "lgeMatrix"),
			    as(t(object), "lgeMatrix"))
	      else if (is.z) ## will error out here
		  isTRUE(all.equal(as(       object,   "zgeMatrix"),
				   as(Conj(t(object)), "zgeMatrix"),
				   tolerance = tol, ...))
	      else if (is(object, "iMatrix")) ## will error out here
		  identical(as(object, "igeMatrix"),
			    as(t(object), "igeMatrix"))
	  })

## rather methods in ./triangularMatrix.R
## setMethod("isTriangular", signature(object = "triangularMatrix"),
## 	  function(object, ...) TRUE)

setMethod("isTriangular", signature(object = "denseMatrix"), isTriMat)

setMethod("isDiagonal", signature(object = "denseMatrix"), .is.diagonal)

setMethod("rcond", signature(x = "denseMatrix", norm = "character"),
	  function(x, norm, ...)
	  rcond(as(as(x, "dMatrix"), "dgeMatrix"), norm=norm, ...))

setMethod("symmpart", signature(x = "denseMatrix"),
	  function(x) symmpart(as(x, "dMatrix")))
setMethod("skewpart", signature(x = "denseMatrix"),
	  function(x) skewpart(as(x, "dMatrix")))

setMethod("is.na", signature(x = "denseMatrix"),
	  function(x) {
	      if(any((inax <- is.na(x@x)))) {
		  r <- as(x, "lMatrix")#-> logical x-slot
		  r@x <- inax
		  as(r, "nMatrix")
	      } else {
		  d <- x@Dim
		  new("ngCMatrix", Dim = d, Dimnames = dimnames(x),
		      i = integer(0), p = rep.int(0L, d[2]+1L))
	      }
	  })

if(.Matrix.avoiding.as.matrix) {
setMethod("qr", signature(x = "ddenseMatrix"),
	  function(x, ...) qr.default(ge2mat(..2dge(x)), ...))
setMethod("qr", signature(x = "denseMatrix"),
	  function(x, ...) qr(as(x, "ddenseMatrix"), ...))
}