File: bind2.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 (429 lines) | stat: -rw-r--r-- 17,265 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
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
#### Containing all  cbind2() and rbind2() methods for all our Matrices

###-- General -----------------------------------------------------------

###-- Dense, incl Diagonal ----------------------------------------------

###-- Sparse ------------------------------------------------------------

setMethod("cbind2", signature(x = "sparseMatrix", y = "matrix"),
	  function(x, y, ...) cbind2(x, .Call(dense_to_Csparse, y)))
setMethod("cbind2", signature(x = "matrix", y = "sparseMatrix"),
	  function(x, y, ...) cbind2(.Call(dense_to_Csparse, x), y))
setMethod("rbind2", signature(x = "sparseMatrix", y = "matrix"),
	  function(x, y, ...) rbind2(x, .Call(dense_to_Csparse, y)))
setMethod("rbind2", signature(x = "matrix", y = "sparseMatrix"),
	  function(x, y, ...) rbind2(.Call(dense_to_Csparse, x), y))

## originally from ./Matrix.R : -------------------------------

## The trivial methods :
setMethod("cbind2", signature(x = "Matrix", y = "NULL"),
          function(x, y, ...) x)
setMethod("cbind2", signature(x = "Matrix", y = "missing"),
          function(x, y, ...) x)
setMethod("cbind2", signature(x = "NULL", y="Matrix"),
          function(x, y, ...) y)
## using "atomicVector" not just "numeric"
setMethod("cbind2", signature(x = "Matrix", y = "atomicVector"),
	  function(x, y, ...) cbind2(x, matrix(y, nrow = nrow(x))))
setMethod("cbind2", signature(x = "atomicVector", y = "Matrix"),
	  function(x, y, ...) cbind2(matrix(x, nrow = nrow(y)), y))
setMethod("cbind2", signature(x = "ANY", y = "Matrix"),
	  function(x, y, ...) .bail.out.2(.Generic, class(x), class(y)))
setMethod("cbind2", signature(x = "Matrix", y = "ANY"),
	  function(x, y, ...) .bail.out.2(.Generic, class(x), class(y)))

setMethod("rbind2", signature(x = "Matrix", y = "NULL"),
          function(x, y, ...) x)
setMethod("rbind2", signature(x = "Matrix", y = "missing"),
          function(x, y, ...) x)
setMethod("rbind2", signature(x = "NULL", y="Matrix"),
          function(x, y, ...) y)
setMethod("rbind2", signature(x = "Matrix", y = "atomicVector"),
	  function(x, y, ...) rbind2(x, matrix(y, ncol = ncol(x))))
setMethod("rbind2", signature(x = "atomicVector", y = "Matrix"),
	  function(x, y, ...) rbind2(matrix(x, ncol = ncol(y)), y))
setMethod("rbind2", signature(x = "ANY", y = "Matrix"),
	  function(x, y, ...) .bail.out.2(.Generic, class(x), class(y)))
setMethod("rbind2", signature(x = "Matrix", y = "ANY"),
	  function(x, y, ...) .bail.out.2(.Generic, class(x), class(y)))

## Makes sure one gets x decent error message for the unimplemented cases:
setMethod("cbind2", signature(x = "Matrix", y = "Matrix"),
	  function(x, y, ...) {
	      rowCheck(x,y)
	      .bail.out.2("cbind2", class(x), class(y))
          })

## Use a working fall back {particularly useful for sparse}:
## FIXME: implement rbind2 via "cholmod" for C* and Tsparse ones
setMethod("rbind2", signature(x = "Matrix", y = "Matrix"),
          function(x, y, ...) {
              colCheck(x,y)
              t(cbind2(t(x), t(y)))
          })

## originally from ./denseMatrix.R : -------------------------------

### cbind2
setMethod("cbind2", signature(x = "denseMatrix", y = "numeric"),
	  function(x, y, ...) {
	      d <- dim(x); nr <- d[1]; nc <- d[2]
	      y <- rep_len(y, nr) # 'silent procrustes'
	      ## beware of (packed) triangular, symmetric, ...
	      x <- as(x, geClass(x))
	      x@x <- c(x@x, as.double(y))
	      x@Dim[2] <- nc + 1L
	      if(is.character(dn <- x@Dimnames[[2]]))
		  x@Dimnames[[2]] <- c(dn, "")
	      x
	  })
## the same, (x,y) <-> (y,x):
setMethod("cbind2", signature(x = "numeric", y = "denseMatrix"),
	  function(x, y, ...) {
	      d <- dim(y); nr <- d[1]; nc <- d[2]
	      x <- rep_len(x, nr)
	      y <- as(y, geClass(y))
	      y@x <- c(as.double(x), y@x)
	      y@Dim[2] <- nc + 1L
	      if(is.character(dn <- y@Dimnames[[2]]))
		  y@Dimnames[[2]] <- c("", dn)
	      y
	  })


setMethod("cbind2", signature(x = "denseMatrix", y = "matrix"),
	  function(x, y, ...) cbind2(x, as_geSimpl(y)))
setMethod("cbind2", signature(x = "matrix", y = "denseMatrix"),
	  function(x, y, ...) cbind2(as_geSimpl(x), y))

cbind2DN <- function(dnx,dny, ncx,ncy) {
    ## R and S+ are different in which names they take
    ## if they differ -- but there's no warning in any case
    rn <- if(!is.null(dnx[[1]])) dnx[[1]] else dny[[1]]
    cx <- dnx[[2]] ; cy <- dny[[2]]
    cn <- if(is.null(cx) && is.null(cy)) NULL
	  else c(if(!is.null(cx)) cx else character(ncx),
		 if(!is.null(cy)) cy else character(ncy))
    list(rn, cn)
}

setMethod("cbind2", signature(x = "denseMatrix", y = "denseMatrix"),
	  function(x, y, ...) {
	      rowCheck(x,y)
	      ncx <- x@Dim[2]
	      ncy <- y@Dim[2]
	      ## beware of (packed) triangular, symmetric, ...
	      hasDN <- !is.null.DN(dnx <- dimnames(x)) | !is.null.DN(dny <- dimnames(y))
	      x <- as(x, geClass(x))
	      y <- as(y, geClass(y))
	      xx <- c(x@x, y@x)
	      ## be careful, e.g., if we have an 'n' and 'd'
	      if(identical((tr <- typeof(xx)), typeof(x@x))) {
		  x@x <- xx
                  x@Dim[2] <- ncx + ncy
		  if(hasDN) x@Dimnames <- cbind2DN(dnx,dny, ncx,ncy)
		  x
	      } else if(identical(tr, typeof(y@x))) {
		  y@x <- xx
		  y@Dim[2] <- ncx + ncy
		  if(hasDN) y@Dimnames <- cbind2DN(dnx,dny, ncx,ncy)
		  y
	      } else stop("resulting x-slot has different type than x's or y's")
	  })

### rbind2 -- analogous to cbind2 --- more to do for @x though:

setMethod("rbind2", signature(x = "denseMatrix", y = "numeric"),
	  function(x, y, ...) {
	      if(is.character(dn <- x@Dimnames[[1]])) dn <- c(dn, "")
	      y <- rbind2(as(x,"matrix"), y)
	      new(paste0(.M.kind(y), "geMatrix"), x = c(y),
                  Dim = x@Dim + 1:0, Dimnames = list(dn, x@Dimnames[[2]]))
	  })
## the same, (x,y) <-> (y,x):
setMethod("rbind2", signature(x = "numeric", y = "denseMatrix"),
	  function(x, y, ...) {
	      if(is.character(dn <- y@Dimnames[[1]])) dn <- c("", dn)
	      x <- rbind2(x, as(y,"matrix"))
	      new(paste0(.M.kind(x), "geMatrix"), x = c(x),
                  Dim = y@Dim + 1:0, Dimnames = list(dn, y@Dimnames[[2]]))
	  })

setMethod("rbind2", signature(x = "denseMatrix", y = "matrix"),
	  function(x, y, ...) rbind2(x, as_geSimpl(y)))
setMethod("rbind2", signature(x = "matrix", y = "denseMatrix"),
	  function(x, y, ...) rbind2(as_geSimpl(x), y))

rbind2DN <- function(dnx, dny, nrx,nry) {
    if(!is.null.DN(dnx) || !is.null.DN(dny)) {
	## R and S+ are different in which names they take
	## if they differ -- but there's no warning in any case
	list(if(is.null(rx <- dnx[[1]]) & is.null(ry <- dny[[1]]))
	     NULL else
	     c(if(!is.null(rx)) rx else character(nrx),
	       if(!is.null(ry)) ry else character(nry)),
	     if(!is.null(dnx[[2]])) dnx[[2]] else dny[[2]])
    } else list(NULL, NULL)
}

setMethod("rbind2", signature(x = "denseMatrix", y = "denseMatrix"),
	  function(x, y, ...) {
	      colCheck(x,y)
	      nrx <- x@Dim[1]
	      nry <- y@Dim[1]
	      ## beware of (packed) triangular, symmetric, ...
	      hasDN <- !is.null.DN(dnx <- dimnames(x)) | !is.null.DN(dny <- dimnames(y))
	      x <- as(x, geClass(x))
	      y <- as(y, geClass(y))
	      xx <- .Call(R_rbind2_vector, x, y)
	      ## be careful, e.g., if we have an 'n' and 'd'
	      if(identical((tr <- typeof(xx)), typeof(x@x))) {
		  x@x <- xx
		  x@Dim[1] <- nrx + nry
		  if(hasDN) x@Dimnames <- rbind2DN(dnx,dny, nrx,nry)
		  x
	      } else if(identical(tr, typeof(y@x))) {
		  y@x <- xx
		  y@Dim[1] <- nrx + nry
		  if(hasDN) y@Dimnames <- rbind2DN(dnx,dny, nrx,nry)
		  y
	      } else stop("resulting x-slot has different type than x's or y's")
	  })

## originally from ./diagMatrix.R : --------------------------------------

## For diagonalMatrix:  preserve sparseness {not always optimal, but "the law"}

## hack to suppress the obnoxious dispatch ambiguity warnings:
diag2Sp <- function(x) suppressWarnings(as(x, "CsparseMatrix"))

setMethod("cbind2", signature(x = "diagonalMatrix", y = "sparseMatrix"),
	  function(x, y, ...) cbind2(diag2Sp(x), as(y,"CsparseMatrix")))
setMethod("cbind2", signature(x = "sparseMatrix", y = "diagonalMatrix"),
	  function(x, y, ...) cbind2(as(x,"CsparseMatrix"), diag2Sp(y)))
setMethod("rbind2", signature(x = "diagonalMatrix", y = "sparseMatrix"),
	  function(x, y, ...) rbind2(diag2Sp(x), as(y,"CsparseMatrix")))
setMethod("rbind2", signature(x = "sparseMatrix", y = "diagonalMatrix"),
	  function(x, y, ...) rbind2(as(x,"CsparseMatrix"), diag2Sp(y)))

## in order to evade method dispatch ambiguity, but still remain "general"
## we use this hack instead of signature  x = "diagonalMatrix"
for(cls in names(getClass("diagonalMatrix")@subclasses)) {

 setMethod("cbind2", signature(x = cls, y = "matrix"),
	   function(x, y, ...) cbind2(diag2Sp(x), .Call(dense_to_Csparse, y)))
 setMethod("cbind2", signature(x = "matrix", y = cls),
	   function(x, y, ...) cbind2(.Call(dense_to_Csparse, x), diag2Sp(y)))
 setMethod("rbind2", signature(x = cls, y = "matrix"),
	   function(x, y, ...) rbind2(diag2Sp(x), .Call(dense_to_Csparse, y)))
 setMethod("rbind2", signature(x = "matrix", y = cls),
	   function(x, y, ...) rbind2(.Call(dense_to_Csparse, x), diag2Sp(y)))

 ## These are already defined for "Matrix"
 ## -- repeated here for method dispatch disambiguation	 {"design-FIXME" ?}
 setMethod("cbind2", signature(x = cls, y = "atomicVector"),
	   function(x, y, ...) cbind2(x, matrix(y, nrow = nrow(x))))
 setMethod("cbind2", signature(x = "atomicVector", y = cls),
	   function(x, y, ...) cbind2(matrix(x, nrow = nrow(y)), y))
 setMethod("rbind2", signature(x = cls, y = "atomicVector"),
	   function(x, y, ...) rbind2(x, matrix(y, ncol = ncol(x))))
 setMethod("rbind2", signature(x = "atomicVector", y = cls),
	   function(x, y, ...) rbind2(matrix(x, ncol = ncol(y)), y))
}


## originally from ./dsparseMatrix.R : --------------------------------

## FIXME: dimnames() handling should happen in C code
## ------> ../src/Csparse.c

## Fast - almost non-checking methods
.cbind2Csp <- function(x,y) .Call(Csparse_horzcat, as_Csp2(x), as_Csp2(y))
.rbind2Csp <- function(x,y) .Call(Csparse_vertcat, as_Csp2(x), as_Csp2(y))

cbind2sparse <- function(x,y) {
    ## beware of (packed) triangular, symmetric, ...
    if(identical(c(dnx <- dimnames(x),
		   dny <- dimnames(y)),
		 list(NULL,NULL,NULL,NULL)))
	## keep empty dimnames
	.cbind2Csp(x,y)
    else {
	## R and S+ are different in which names they take
	## if they differ -- but there's no warning in any case
	rn <- if(!is.null(dnx[[1]])) dnx[[1]] else dny[[1]]
	cx <- dnx[[2]] ; cy <- dny[[2]]
	cn <- if(is.null(cx) && is.null(cy)) NULL
	else c(if(!is.null(cx)) cx else character(ncol(x)),
	       if(!is.null(cy)) cy else character(ncol(y)))
	ans <- .cbind2Csp(x,y)
	ans@Dimnames <- list(rn, cn)
	ans
    }
}
setMethod("cbind2", signature(x = "sparseMatrix", y = "sparseMatrix"),
	  function(x, y, ...) {
	      rowCheck(x,y)
	      cbind2sparse(x,y)
	  })

rbind2sparse <- function(x,y) {
    ## beware of (packed) triangular, symmetric, ...
    if(identical(c(dnx <- dimnames(x),
		   dny <- dimnames(y)),
		 list(NULL,NULL,NULL,NULL)))
	## keep empty dimnames
	.rbind2Csp(x,y)
    else {
	## R and S+ are different in which names they take
	## if they differ -- but there's no warning in any case
	cn <- if(!is.null(dnx[[2]])) dnx[[2]] else dny[[2]]
	rx <- dnx[[1]] ; ry <- dny[[1]]
	rn <- if(is.null(rx) && is.null(ry)) NULL
	else c(if(!is.null(rx)) rx else character(nrow(x)),
	       if(!is.null(ry)) ry else character(nrow(y)))
	ans <- .rbind2Csp(x,y)
	ans@Dimnames <- list(rn, cn)
	ans
    }
}
setMethod("rbind2", signature(x = "sparseMatrix", y = "sparseMatrix"),
	  function(x, y, ...) {
	      colCheck(x,y)
	      rbind2sparse(x,y)
	  })

if(length(formals(cbind2)) >= 3) { ## newer R -- can use optional 'sparse = NA'

setMethod("cbind2", signature(x = "sparseMatrix", y = "denseMatrix"),
	  function(x, y, sparse = NA, ...) {
	      nr <- rowCheck(x,y)
	      if(is.na(sparse)) # result is sparse if "enough zeros" <==> sparseDefault() in Matrix()
		  sparse <- (nnzero(x,na.counted=TRUE)+nnzero(y,na.counted=TRUE)) * 2 <
		      as.double(nr) * (ncol(x)+ncol(y)) # as.double(): avoid integer overflow in '*'
	      if(sparse) cbind2sparse(x,y) else cbind2(as(x, "denseMatrix"), y)
	  })
setMethod("cbind2", signature(x = "denseMatrix", y = "sparseMatrix"),
	  function(x, y, sparse = NA, ...) {
	      nr <- rowCheck(x,y)
	      if(is.na(sparse)) # result is sparse if "enough zeros" <==> sparseDefault() in Matrix()
		  sparse <- (nnzero(x,na.counted=TRUE)+nnzero(y,na.counted=TRUE)) * 2 <
		      as.double(nr) * (ncol(x)+ncol(y))
	      if(sparse) cbind2sparse(x,y) else cbind2(x, as(y, "denseMatrix"))
	  })
setMethod("rbind2", signature(x = "sparseMatrix", y = "denseMatrix"),
	  function(x, y, sparse = NA, ...) {
	      nc <- colCheck(x,y)
	      if(is.na(sparse)) # result is sparse if "enough zeros" <==> sparseDefault() in Matrix()
		  sparse <- (nnzero(x,na.counted=TRUE)+nnzero(y,na.counted=TRUE)) * 2 <
		      (nrow(x)+nrow(y)) * as.double(nc)
	      if(sparse) rbind2sparse(x,y) else rbind2(as(x, "denseMatrix"), y)
	  })
setMethod("rbind2", signature(x = "denseMatrix", y = "sparseMatrix"),
	  function(x, y, sparse = NA, ...) {
	      nc <- colCheck(x,y)
	      if(is.na(sparse)) # result is sparse if "enough zeros" <==> sparseDefault() in Matrix()
		  sparse <- (nnzero(x,na.counted=TRUE)+nnzero(y,na.counted=TRUE)) * 2 <
		(nrow(x)+nrow(y)) * as.double(nc)
	      if(sparse) rbind2sparse(x,y) else rbind2(x, as(y, "denseMatrix"))
	  })

} else { ## older version of R -- cbind2() has no "..."

setMethod("cbind2", signature(x = "sparseMatrix", y = "denseMatrix"),
	  function(x, y, ...) {
	      nr <- rowCheck(x,y)
	      ## result is sparse if "enough zeros" <==> sparseDefault() in Matrix()
	      sparse <- (nnzero(x,na.counted=TRUE)+nnzero(y,na.counted=TRUE)) * 2 <
		  as.double(nr) * (ncol(x)+ncol(y))
	      if(sparse) cbind2sparse(x,y) else cbind2(as(x, "denseMatrix"), y)
	  })
setMethod("cbind2", signature(x = "denseMatrix", y = "sparseMatrix"),
	  function(x, y, ...) {
	      nr <- rowCheck(x,y)
	      ## result is sparse if "enough zeros" <==> sparseDefault() in Matrix()
	      sparse <- (nnzero(x,na.counted=TRUE)+nnzero(y,na.counted=TRUE)) * 2 <
		  as.double(nr) * (ncol(x)+ncol(y))
	      if(sparse) cbind2sparse(x,y) else cbind2(x, as(y, "denseMatrix"))
	  })
setMethod("rbind2", signature(x = "sparseMatrix", y = "denseMatrix"),
	  function(x, y, ...) {
	      nc <- colCheck(x,y)
	      ## result is sparse if "enough zeros" <==> sparseDefault() in Matrix()
	      sparse <- (nnzero(x,na.counted=TRUE)+nnzero(y,na.counted=TRUE)) * 2 <
		  (nrow(x)+nrow(y)) * as.double(nc)
	      if(sparse) rbind2sparse(x,y) else rbind2(as(x, "denseMatrix"), y)
	  })
setMethod("rbind2", signature(x = "denseMatrix", y = "sparseMatrix"),
	  function(x, y, ...) {
	      nc <- colCheck(x,y)
	      ## result is sparse if "enough zeros" <==> sparseDefault() in Matrix()
	      sparse <- (nnzero(x,na.counted=TRUE)+nnzero(y,na.counted=TRUE)) * 2 <
		  (nrow(x)+nrow(y)) * as.double(nc)
	      if(sparse) rbind2sparse(x,y) else rbind2(x, as(y, "denseMatrix"))
	  })
}# older R -- no "sparse = NA"




if(FALSE) {
    ## FIXME
    ##------------- maybe a bit faster --- but too much to maintain
    ## would have to be done for "rbind2" as well ...
setMethod("cbind2", signature(x = "sparseMatrix", y = "numeric"),
          function(x, y, ...) {
              d <- dim(x); nr <- d[1]; nc <- d[2]; cl <- class(x)
              x <- as(x, "CsparseMatrix")
              if(nr > 0) {
		  y <- rep_len(y, nr) # 'silent procrustes'
                  n0y <- y != 0
                  n.e <- length(x@i)
                  x@i <- c(x@i, (0:(nr-1))[n0y])
                  x@p <- c(x@p, n.e + sum(n0y))
                  x@x <- c(x@x, y[n0y])
              } else { ## nr == 0

              }
              x@Dim[2] <- nc + 1L
              if(is.character(dn <- x@Dimnames[[2]]))
                  x@Dimnames[[2]] <- c(dn, "")
              x
          })
## the same, (x,y) <-> (y,x):
setMethod("cbind2", signature(x = "numeric", y = "sparseMatrix"),
          function(x, y, ...) {
              d <- dim(y); nr <- d[1]; nc <- d[2]; cl <- class(y)
              y <-  as(y, "CsparseMatrix")
              if(nr > 0) {
		  x <- rep_len(x, nr) # 'silent procrustes'
                  n0x <- x != 0
                  y@i <- c((0:(nr-1))[n0x], y@i)
                  y@p <- c(0L, sum(n0x) + y@p)
                  y@x <- c(x[n0x], y@x)
              } else { ## nr == 0

              }
              y@Dim[2] <- nc + 1L
              if(is.character(dn <- y@Dimnames[[2]]))
                  y@Dimnames[[2]] <- c(dn, "")
              y
          })

}## -- no longer

## Can be made very efficient
setMethod("rbind2", signature(x = "indMatrix", y = "indMatrix"),
	  function(x, y, ...) {
	      dx <- x@Dim
	      dy <- y@Dim
	      if(dx[2] != dy[2])
		  stop(gettextf("Matrices must have same number of columns in %s",
				deparse(sys.call(sys.parent()))),
		       call. = FALSE, domain=NA)
	      new("indMatrix", Dim = c(dx[1]+dy[1], dx[2]),
		  perm = c(x@perm,y@perm),
		  Dimnames = rbind2DN(dimnames(x), dimnames(y), dx[1],dy[1]))
	  })