SCM

SCM Repository

[matrix] Diff of /pkg/Matrix/src/Csparse.c
ViewVC logotype

Diff of /pkg/Matrix/src/Csparse.c

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 2646, Mon Feb 21 10:57:49 2011 UTC revision 2677, Sat Jun 25 19:18:12 2011 UTC
# Line 1  Line 1 
1                          /* Sparse matrices in compressed column-oriented form */                          /* Sparse matrices in compressed column-oriented form */
2    
3    #include <stdint.h> // C99 for int64_t
4  #include "Csparse.h"  #include "Csparse.h"
5  #include "Tsparse.h"  #include "Tsparse.h"
6  #include "chm_common.h"  #include "chm_common.h"
# Line 136  Line 138 
138              }              }
139      }      }
140      if (!sorted)      if (!sorted)
141          /* cannot easily use cholmod_l_sort(.) ... -> "error out" :*/          /* cannot easily use cholmod_sort(.) ... -> "error out" :*/
142          return mkString(_("slot j is not increasing inside a column"));          return mkString(_("slot j is not increasing inside a column"));
143      else if(!strictly) /* sorted, but not strictly */      else if(!strictly) /* sorted, but not strictly */
144          return mkString(_("slot j is not *strictly* increasing inside a column"));          return mkString(_("slot j is not *strictly* increasing inside a column"));
# Line 154  Line 156 
156      /* This loses the symmetry property, since cholmod_dense has none,      /* This loses the symmetry property, since cholmod_dense has none,
157       * BUT, much worse (FIXME!), it also transforms CHOLMOD_PATTERN ("n") matrices       * BUT, much worse (FIXME!), it also transforms CHOLMOD_PATTERN ("n") matrices
158       * to numeric (CHOLMOD_REAL) ones : */       * to numeric (CHOLMOD_REAL) ones : */
159      CHM_DN chxd = cholmod_l_sparse_to_dense(chxs, &c);      CHM_DN chxd = cholmod_sparse_to_dense(chxs, &c);
160      int Rkind = (chxs->xtype == CHOLMOD_PATTERN)? -1 : Real_kind(x);      int Rkind = (chxs->xtype == CHOLMOD_PATTERN)? -1 : Real_kind(x);
161      R_CheckStack();      R_CheckStack();
162    
# Line 165  Line 167 
167  SEXP Csparse_to_nz_pattern(SEXP x, SEXP tri)  SEXP Csparse_to_nz_pattern(SEXP x, SEXP tri)
168  {  {
169      CHM_SP chxs = AS_CHM_SP__(x);      CHM_SP chxs = AS_CHM_SP__(x);
170      CHM_SP chxcp = cholmod_l_copy(chxs, chxs->stype, CHOLMOD_PATTERN, &c);      CHM_SP chxcp = cholmod_copy(chxs, chxs->stype, CHOLMOD_PATTERN, &c);
171      int tr = asLogical(tri);      int tr = asLogical(tri);
172      R_CheckStack();      R_CheckStack();
173    
# Line 230  Line 232 
232    
233  SEXP Csparse_to_matrix(SEXP x)  SEXP Csparse_to_matrix(SEXP x)
234  {  {
235      return chm_dense_to_matrix(cholmod_l_sparse_to_dense(AS_CHM_SP__(x), &c),      return chm_dense_to_matrix(cholmod_sparse_to_dense(AS_CHM_SP__(x), &c),
236                                 1 /*do_free*/, GET_SLOT(x, Matrix_DimNamesSym));                                 1 /*do_free*/, GET_SLOT(x, Matrix_DimNamesSym));
237  }  }
238    
239  SEXP Csparse_to_Tsparse(SEXP x, SEXP tri)  SEXP Csparse_to_Tsparse(SEXP x, SEXP tri)
240  {  {
241      CHM_SP chxs = AS_CHM_SP__(x);      CHM_SP chxs = AS_CHM_SP__(x);
242      CHM_TR chxt = cholmod_l_sparse_to_triplet(chxs, &c);      CHM_TR chxt = cholmod_sparse_to_triplet(chxs, &c);
243      int tr = asLogical(tri);      int tr = asLogical(tri);
244      int Rkind = (chxs->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;      int Rkind = (chxs->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
245      R_CheckStack();      R_CheckStack();
# Line 257  Line 259 
259    
260      if (!(chx->stype))      if (!(chx->stype))
261          error(_("Nonsymmetric matrix in Csparse_symmetric_to_general"));          error(_("Nonsymmetric matrix in Csparse_symmetric_to_general"));
262      chgx = cholmod_l_copy(chx, /* stype: */ 0, chx->xtype, &c);      chgx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c);
263      /* xtype: pattern, "real", complex or .. */      /* xtype: pattern, "real", complex or .. */
264      return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",      return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",
265                                GET_SLOT(x, Matrix_DimNamesSym));                                GET_SLOT(x, Matrix_DimNamesSym));
# Line 270  Line 272 
272      int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;      int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
273      R_CheckStack();      R_CheckStack();
274    
275      chgx = cholmod_l_copy(chx, /* stype: */ uploT, chx->xtype, &c);      chgx = cholmod_copy(chx, /* stype: */ uploT, chx->xtype, &c);
276      /* xtype: pattern, "real", complex or .. */      /* xtype: pattern, "real", complex or .. */
277      return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",      return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",
278                                GET_SLOT(x, Matrix_DimNamesSym));                                GET_SLOT(x, Matrix_DimNamesSym));
# Line 282  Line 284 
284       *       since cholmod (& cs) lacks sparse 'int' matrices */       *       since cholmod (& cs) lacks sparse 'int' matrices */
285      CHM_SP chx = AS_CHM_SP__(x);      CHM_SP chx = AS_CHM_SP__(x);
286      int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;      int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
287      CHM_SP chxt = cholmod_l_transpose(chx, chx->xtype, &c);      CHM_SP chxt = cholmod_transpose(chx, chx->xtype, &c);
288      SEXP dn = PROTECT(duplicate(GET_SLOT(x, Matrix_DimNamesSym))), tmp;      SEXP dn = PROTECT(duplicate(GET_SLOT(x, Matrix_DimNamesSym))), tmp;
289      int tr = asLogical(tri);      int tr = asLogical(tri);
290      R_CheckStack();      R_CheckStack();
# Line 301  Line 303 
303      CHM_SP      CHM_SP
304          cha = AS_CHM_SP(a),          cha = AS_CHM_SP(a),
305          chb = AS_CHM_SP(b),          chb = AS_CHM_SP(b),
306          chc = cholmod_l_ssmult(cha, chb, /*out_stype:*/ 0,          chc = cholmod_ssmult(cha, chb, /*out_stype:*/ 0,
307                                 /* values:= is_numeric (T/F) */ cha->xtype > 0,                                 /* values:= is_numeric (T/F) */ cha->xtype > 0,
308                                 /*out sorted:*/ 1, &c);                                 /*out sorted:*/ 1, &c);
309      const char *cl_a = class_P(a), *cl_b = class_P(b);      const char *cl_a = class_P(a), *cl_b = class_P(b);
# Line 352  Line 354 
354      SEXP dn = PROTECT(allocVector(VECSXP, 2));      SEXP dn = PROTECT(allocVector(VECSXP, 2));
355      R_CheckStack();      R_CheckStack();
356    
357      chTr = cholmod_l_transpose((tr) ? chb : cha, chb->xtype, &c);      chTr = cholmod_transpose((tr) ? chb : cha, chb->xtype, &c);
358      chc = cholmod_l_ssmult((tr) ? cha : chTr, (tr) ? chTr : chb,      chc = cholmod_ssmult((tr) ? cha : chTr, (tr) ? chTr : chb,
359                           /*out_stype:*/ 0, cha->xtype, /*out sorted:*/ 1, &c);                           /*out_stype:*/ 0, cha->xtype, /*out sorted:*/ 1, &c);
360      cholmod_l_free_sparse(&chTr, &c);      cholmod_free_sparse(&chTr, &c);
361    
362      /* Preserve triangularity and unit-triangularity if appropriate;      /* Preserve triangularity and unit-triangularity if appropriate;
363       * see Csparse_Csparse_prod() for comments */       * see Csparse_Csparse_prod() for comments */
# Line 381  Line 383 
383      CHM_SP cha = AS_CHM_SP(a);      CHM_SP cha = AS_CHM_SP(a);
384      SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));      SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));
385      CHM_DN chb = AS_CHM_DN(b_M);      CHM_DN chb = AS_CHM_DN(b_M);
386      CHM_DN chc = cholmod_l_allocate_dense(cha->nrow, chb->ncol, cha->nrow,      CHM_DN chc = cholmod_allocate_dense(cha->nrow, chb->ncol, cha->nrow,
387                                          chb->xtype, &c);                                          chb->xtype, &c);
388      SEXP dn = PROTECT(allocVector(VECSXP, 2));      SEXP dn = PROTECT(allocVector(VECSXP, 2));
389      double one[] = {1,0}, zero[] = {0,0};      double one[] = {1,0}, zero[] = {0,0};
# Line 398  Line 400 
400          SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;          SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;
401          cha = AS_CHM_SP(da);          cha = AS_CHM_SP(da);
402      }      }
403      cholmod_l_sdmult(cha, 0, one, zero, chb, chc, &c);      cholmod_sdmult(cha, 0, one, zero, chb, chc, &c);
404      SET_VECTOR_ELT(dn, 0,       /* establish dimnames */      SET_VECTOR_ELT(dn, 0,       /* establish dimnames */
405                     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0)));                     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0)));
406      SET_VECTOR_ELT(dn, 1,      SET_VECTOR_ELT(dn, 1,
# Line 412  Line 414 
414      CHM_SP cha = AS_CHM_SP(a);      CHM_SP cha = AS_CHM_SP(a);
415      SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));      SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));
416      CHM_DN chb = AS_CHM_DN(b_M);      CHM_DN chb = AS_CHM_DN(b_M);
417      CHM_DN chc = cholmod_l_allocate_dense(cha->ncol, chb->ncol, cha->ncol,      CHM_DN chc = cholmod_allocate_dense(cha->ncol, chb->ncol, cha->ncol,
418                                          chb->xtype, &c);                                          chb->xtype, &c);
419      SEXP dn = PROTECT(allocVector(VECSXP, 2)); int nprot = 2;      SEXP dn = PROTECT(allocVector(VECSXP, 2)); int nprot = 2;
420      double one[] = {1,0}, zero[] = {0,0};      double one[] = {1,0}, zero[] = {0,0};
# Line 422  Line 424 
424          SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;          SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;
425          cha = AS_CHM_SP(da);          cha = AS_CHM_SP(da);
426      }      }
427      cholmod_l_sdmult(cha, 1, one, zero, chb, chc, &c);      cholmod_sdmult(cha, 1, one, zero, chb, chc, &c);
428      SET_VECTOR_ELT(dn, 0,       /* establish dimnames */      SET_VECTOR_ELT(dn, 0,       /* establish dimnames */
429                     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 1)));                     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 1)));
430      SET_VECTOR_ELT(dn, 1,      SET_VECTOR_ELT(dn, 1,
# Line 445  Line 447 
447  #endif  #endif
448      CHM_SP chcp, chxt,      CHM_SP chcp, chxt,
449          chx = (trip ?          chx = (trip ?
450                 cholmod_l_triplet_to_sparse(cht, cht->nnz, &c) :                 cholmod_triplet_to_sparse(cht, cht->nnz, &c) :
451                 AS_CHM_SP(x));                 AS_CHM_SP(x));
452      SEXP dn = PROTECT(allocVector(VECSXP, 2));      SEXP dn = PROTECT(allocVector(VECSXP, 2));
453      R_CheckStack();      R_CheckStack();
454    
455      if (!tr) chxt = cholmod_l_transpose(chx, chx->xtype, &c);      if (!tr) chxt = cholmod_transpose(chx, chx->xtype, &c);
456      chcp = cholmod_l_aat((!tr) ? chxt : chx, (int *) NULL, 0, chx->xtype, &c);      chcp = cholmod_aat((!tr) ? chxt : chx, (int *) NULL, 0, chx->xtype, &c);
457      if(!chcp) {      if(!chcp) {
458          UNPROTECT(1);          UNPROTECT(1);
459          error(_("Csparse_crossprod(): error return from cholmod_l_aat()"));          error(_("Csparse_crossprod(): error return from cholmod_aat()"));
460      }      }
461      cholmod_l_band_inplace(0, chcp->ncol, chcp->xtype, chcp, &c);      cholmod_band_inplace(0, chcp->ncol, chcp->xtype, chcp, &c);
462      chcp->stype = 1;      chcp->stype = 1;
463      if (trip) cholmod_l_free_sparse(&chx, &c);      if (trip) cholmod_free_sparse(&chx, &c);
464      if (!tr) cholmod_l_free_sparse(&chxt, &c);      if (!tr) cholmod_free_sparse(&chxt, &c);
465      SET_VECTOR_ELT(dn, 0,       /* establish dimnames */      SET_VECTOR_ELT(dn, 0,       /* establish dimnames */
466                     duplicate(VECTOR_ELT(GET_SLOT(x, Matrix_DimNamesSym),                     duplicate(VECTOR_ELT(GET_SLOT(x, Matrix_DimNamesSym),
467                                          (tr) ? 0 : 1)));                                          (tr) ? 0 : 1)));
# Line 480  Line 482 
482      /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */      /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
483      int tr = (cl[1] == 't');      int tr = (cl[1] == 't');
484      CHM_SP chx = AS_CHM_SP__(x);      CHM_SP chx = AS_CHM_SP__(x);
485      CHM_SP ans = cholmod_l_copy(chx, chx->stype, chx->xtype, &c);      CHM_SP ans = cholmod_copy(chx, chx->stype, chx->xtype, &c);
486      double dtol = asReal(tol);      double dtol = asReal(tol);
487      int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;      int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
488      R_CheckStack();      R_CheckStack();
489    
490      if(!cholmod_l_drop(dtol, ans, &c))      if(!cholmod_drop(dtol, ans, &c))
491          error(_("cholmod_l_drop() failed"));          error(_("cholmod_drop() failed"));
492      return chm_sparse_to_SEXP(ans, 1,      return chm_sparse_to_SEXP(ans, 1,
493                                tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,                                tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
494                                Rkind, tr ? diag_P(x) : "",                                Rkind, tr ? diag_P(x) : "",
# Line 502  Line 504 
504      R_CheckStack();      R_CheckStack();
505    
506      /* TODO: currently drops dimnames - and we fix at R level */      /* TODO: currently drops dimnames - and we fix at R level */
507      return chm_sparse_to_SEXP(cholmod_l_horzcat(chx, chy, 1, &c),      return chm_sparse_to_SEXP(cholmod_horzcat(chx, chy, 1, &c),
508                                1, 0, Rkind, "", R_NilValue);                                1, 0, Rkind, "", R_NilValue);
509  }  }
510    
# Line 515  Line 517 
517      R_CheckStack();      R_CheckStack();
518    
519      /* TODO: currently drops dimnames - and we fix at R level */      /* TODO: currently drops dimnames - and we fix at R level */
520      return chm_sparse_to_SEXP(cholmod_l_vertcat(chx, chy, 1, &c),      return chm_sparse_to_SEXP(cholmod_vertcat(chx, chy, 1, &c),
521                                1, 0, Rkind, "", R_NilValue);                                1, 0, Rkind, "", R_NilValue);
522  }  }
523    
# Line 523  Line 525 
525  {  {
526      CHM_SP chx = AS_CHM_SP__(x);      CHM_SP chx = AS_CHM_SP__(x);
527      int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;      int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
528      CHM_SP ans = cholmod_l_band(chx, asInteger(k1), asInteger(k2), chx->xtype, &c);      CHM_SP ans = cholmod_band(chx, asInteger(k1), asInteger(k2), chx->xtype, &c);
529      R_CheckStack();      R_CheckStack();
530    
531      return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "",      return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "",
# Line 541  Line 543 
543      }      }
544      else { /* unit triangular (diag='U'): "fill the diagonal" & diag:= "N" */      else { /* unit triangular (diag='U'): "fill the diagonal" & diag:= "N" */
545          CHM_SP chx = AS_CHM_SP__(x);          CHM_SP chx = AS_CHM_SP__(x);
546          CHM_SP eye = cholmod_l_speye(chx->nrow, chx->ncol, chx->xtype, &c);          CHM_SP eye = cholmod_speye(chx->nrow, chx->ncol, chx->xtype, &c);
547          double one[] = {1, 0};          double one[] = {1, 0};
548          CHM_SP ans = cholmod_l_add(chx, eye, one, one, TRUE, TRUE, &c);          CHM_SP ans = cholmod_add(chx, eye, one, one, TRUE, TRUE, &c);
549          int uploT = (*uplo_P(x) == 'U') ? 1 : -1;          int uploT = (*uplo_P(x) == 'U') ? 1 : -1;
550          int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;          int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
551    
552          R_CheckStack();          R_CheckStack();
553          cholmod_l_free_sparse(&eye, &c);          cholmod_free_sparse(&eye, &c);
554          return chm_sparse_to_SEXP(ans, 1, uploT, Rkind, "N",          return chm_sparse_to_SEXP(ans, 1, uploT, Rkind, "N",
555                                    GET_SLOT(x, Matrix_DimNamesSym));                                    GET_SLOT(x, Matrix_DimNamesSym));
556      }      }
# Line 602  Line 604 
604    
605      if (chx->stype) /* symmetricMatrix */      if (chx->stype) /* symmetricMatrix */
606          /* for now, cholmod_submatrix() only accepts "generalMatrix" */          /* for now, cholmod_submatrix() only accepts "generalMatrix" */
607          chx = cholmod_l_copy(chx, /* stype: */ 0, chx->xtype, &c);          chx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c);
608    
609      return chm_sparse_to_SEXP(cholmod_l_submatrix(chx,      return chm_sparse_to_SEXP(cholmod_submatrix(chx,
610                                  (rsize < 0) ? NULL : INTEGER(i), rsize,                                  (rsize < 0) ? NULL : INTEGER(i), rsize,
611                                  (csize < 0) ? NULL : INTEGER(j), csize,                                  (csize < 0) ? NULL : INTEGER(j), csize,
612                                                    TRUE, TRUE, &c),                                                    TRUE, TRUE, &c),
# Line 612  Line 614 
614                                /* FIXME: drops dimnames */ R_NilValue);                                /* FIXME: drops dimnames */ R_NilValue);
615  }  }
616    
617    /**
618     * Subassignment:  x[i,j]  <- value
619     *
620     * @param x
621     * @param i_ integer row    index 0-origin vector (as returned from R .ind.prep2())
622     * @param j_ integer column index 0-origin vector
623     * @param value currently must be a dsparseVector {which is recycled if needed}
624     *
625     * @return a Csparse matrix like x, but with the values replaced
626     */
627    SEXP Csparse_subassign(SEXP x, SEXP i_, SEXP j_, SEXP value)
628    {
629        // TODO: for other classes consider using a trick as  RallocedReal() in ./chm_common.c
630        static const char
631            *valid_cM [] = {"dgCMatrix",// the only one, for "the moment", .. TODO
632                            ""},
633            *valid_spv[] = {"dsparseVector",
634                            ""};
635    
636        int ctype = Matrix_check_class_etc(x, valid_cM);
637        if (ctype < 0)
638            error(_("invalid class of 'x' in Csparse_subassign()"));
639        // value: assume a  "dsparseVector" for now -- slots: (i, length, x)
640        ctype = Matrix_check_class_etc(value, valid_spv);
641        if (ctype < 0)
642            error(_("invalid class of 'value' in Csparse_subassign()"));
643    
644        SEXP
645            islot   = GET_SLOT(x, Matrix_iSym),
646            dimslot = GET_SLOT(x, Matrix_DimSym),
647            i_cp = PROTECT(coerceVector(i_, INTSXP)),
648            j_cp = PROTECT(coerceVector(j_, INTSXP));
649            // for d.CMatrix and l.CMatrix  but not n.CMatrix:
650    
651        int *dims = INTEGER(dimslot),
652            ncol = dims[1], /* nrow = dims[0], */
653            *i = INTEGER(i_cp), len_i = LENGTH(i_cp),
654            *j = INTEGER(j_cp), len_j = LENGTH(j_cp),
655            k,
656            nnz_x = LENGTH(islot);
657        int nnz = nnz_x;
658    
659    #define MATRIX_SUBASSIGN_VERBOSE
660    // Temporary hack for debugging --- remove eventually -- FIXME
661    #ifdef MATRIX_SUBASSIGN_VERBOSE
662        Rboolean verbose = i[0] < 0;
663        if(verbose) i[0] = -i[0];
664    #endif
665    
666        SEXP val_i_slot;
667        PROTECT(val_i_slot = coerceVector(GET_SLOT(value, Matrix_iSym), REALSXP));
668        double *val_i = REAL(val_i_slot);
669        int nnz_val =  LENGTH(GET_SLOT(value, Matrix_iSym));
670        // for dsparseVector only:
671        double *val_x =   REAL (GET_SLOT(value, Matrix_xSym));
672        int64_t len_val = (int64_t) asReal(GET_SLOT(value, Matrix_lengthSym));
673        /* llen_i = (int64_t) len_i; */
674    
675        double z_ans = 0.;
676    
677        SEXP ans;
678        /* Instead of simple "duplicate": PROTECT(ans = duplicate(x)) , build up: */
679        ans = PROTECT(NEW_OBJECT(MAKE_CLASS("dgCMatrix" /* <- TODO*/)));
680        SET_SLOT(ans, Matrix_DimSym,      duplicate(dimslot));
681        SET_SLOT(ans, Matrix_DimNamesSym, duplicate(GET_SLOT(x, Matrix_DimNamesSym)));
682        SET_SLOT(ans, Matrix_pSym,        duplicate(GET_SLOT(x, Matrix_pSym)));
683        SEXP r_pslot = GET_SLOT(ans, Matrix_pSym);
684        // and assign the i- and x- slots at the end, as they are potentially modified
685        // not just in content, but also in their *length*
686        int *rp = INTEGER(r_pslot),
687            *ri = Calloc(nnz_x, int);       // to contain the final i - slot
688        // for d.CMatrix only:
689        double *rx = Calloc(nnz_x, double); // to contain the final x - slot
690        Memcpy(ri, INTEGER(islot), nnz_x);
691        Memcpy(rx, REAL(GET_SLOT(x, Matrix_xSym)), nnz_x);
692        // NB:  nnz_x : will always be the "current allocated length" of (i, x) slots
693        // --   nnz   : the current *used* length; always   nnz <= nnz_x
694    
695        int jj, j_val = 0; // in "running" conceptionally through all value[i+ jj*len_i]
696        // values, we are "below"/"before" the (j_val)-th non-zero one.
697        // e.g. if value = (0,0,...,0), have nnz_val == 0, j_val must remain == 0
698        int64_t ii_val;// == "running" index (i + jj*len_i) % len_val for value[]
699        for(jj = 0, ii_val=0; jj < len_j; jj++) {
700            int j__ = j[jj];
701            /* int64_t j_l = jj * llen_i; */
702            R_CheckUserInterrupt();
703            for(int ii = 0; ii < len_i; ii++, ii_val++) {
704                int i__ = i[ii], p1, p2;
705                if(nnz_val && ii_val >= len_val) { // "recycle" indexing into value[]
706                    ii_val -= len_val; // = (ii + jj*len_i) % len_val
707                    j_val = 0;
708                }
709                int64_t ii_v1;//= ii_val + 1;
710                double v, /* := value[(ii + j_l) % len_val]
711                             = dsparseVector_sub((ii + j_l) % len_val,
712                                                 nnz_val, val_i, val_x, len_val)
713                          */
714                    M_ij;
715                int ind;
716                Rboolean have_entry = FALSE;
717    
718                // note that rp[]'s may have *changed* even when 'j' remained!
719                // "FIXME": do this only *when* rp[] has changed
720                p1 = rp[j__], p2 = rp[j__ + 1];
721    
722                // v :=  value[(ii + j_l) % len_val] = value[ii_val]
723                v = z_ans;
724                if(j_val < nnz_val) { // maybe find v := non-zero value[ii_val]
725                    ii_v1 = ii_val + 1;
726                    if(ii_v1 < val_i[j_val]) { // typical case: are still in zero-stretch
727                        v = z_ans; // v = 0
728                    } else if(ii_v1 == val_i[j_val]) { // have a match
729                        v = val_x[j_val];
730                        j_val++;// from now on, look at the next non-zero entry
731                    } else { //  ii_v1 > val_i[j_val]
732                        REprintf("programming thinko in Csparse_subassign(*, i=%d,j=%d): ii_v=%d, v@i[j_val=%ld]=%g\n",
733                                 i__,j__, ii_v1, j_val, val_i[j_val]);
734                        j_val++;// from now on, look at the next non-zero entry
735                    }
736                }
737                // --------------- M_ij := getM(i., j.) --------------------------------
738                M_ij = z_ans; // as in  ./t_sparseVector.c
739                for(ind = p1; ind < p2; ind++) {
740                    if(ri[ind] >= i__) {
741                        if(ri[ind] == i__) {
742                            M_ij = rx[ind];
743    #ifdef MATRIX_SUBASSIGN_VERBOSE
744                            if(verbose) REprintf("have entry x[%d, %d] = %g\n",
745                                                 i__, j__, M_ij);
746    #endif
747                            have_entry = TRUE;
748                        } else { // ri[ind] > i__
749    #ifdef MATRIX_SUBASSIGN_VERBOSE
750                            if(verbose)
751                                REprintf("@i > i__ = %d --> ind-- = %d\n", i__, ind);
752    #endif
753                        }
754                        break;
755                    }
756                }
757    
758                //-- R:  if(getM(i., j.) != (v <- getV(ii, jj)))
759    
760                if(M_ij != v) { // contents differ ==> value needs to be changed
761    #ifdef MATRIX_SUBASSIGN_VERBOSE
762                    if(verbose)
763                        REprintf("setting x[%d, %d] <- %g", i__,j__, v);
764    #endif
765                    // (otherwise: nothing to do):
766                    // setM(i__, j__, v)
767                    // ----------------------------------------------------------
768    
769                    // Case I --------------------------------------------
770    /*              if(v == z_ans) { // remove x[i, j] = M_ij  which we know is *non*-zero */
771    /*                  // we know : have_entry = TRUE ; */
772    /*                  //  ri[ind] == i__; M_ij = rx[ind]; */
773    /* #ifdef MATRIX_SUBASSIGN_VERBOSE */
774    /*                  if(verbose) */
775    /*                      REprintf(" rm ind=%d\n", ind); */
776    /* #endif */
777    /*                  // remove the 'ind'-th element from x@i and x@x : */
778    /*                  nnz-- ; */
779    /*                  for(k=ind; k < nnz; k++) { */
780    /*                      ri[k] = ri[k+1]; */
781    /*                      rx[k] = rx[k+1]; */
782    /*                  } */
783    /*                  for(k=j__ + 1; k <= ncol; k++) { */
784    /*                      rp[k] = rp[k] - 1; */
785    /*                  } */
786    /*              } */
787    /*              else  */
788                    if(have_entry) {
789                        // Case II ----- replace (non-empty) x[i,j] by v -------
790    #ifdef MATRIX_SUBASSIGN_VERBOSE
791                        if(verbose)
792                            REprintf(" repl.  ind=%d\n", ind);
793    #endif
794                        rx[ind] = v;
795                    } else {
796                        // Case III ---- v != 0 : insert v into "empty" x[i,j] ----
797    
798                        // extend the  i  and  x  slot by one entry : ---------------------
799    
800                        if(nnz+1 > nnz_x) { // need to reallocate:
801    #ifdef MATRIX_SUBASSIGN_VERBOSE
802                            if(verbose) REprintf(" Realloc()ing: nnz_x=%d", nnz_x);
803    #endif
804                            // do it "only" 1x,..4x at the very most increasing by the
805                            // nnz-length of "value":
806                            nnz_x += (1 + nnz_val / 4);
807    #ifdef MATRIX_SUBASSIGN_VERBOSE
808                            if(verbose) REprintf("(nnz_v=%d) --> %d ", nnz_val, nnz_x);
809    #endif
810                            // C doc on realloc() says that the old content is *preserve*d
811                            ri = Realloc(ri, nnz_x, int);
812                            rx = Realloc(rx, nnz_x, double);
813                        }
814    
815                        // 3) fill them ...
816    
817                        int i1 = ind;
818    #ifdef MATRIX_SUBASSIGN_VERBOSE
819                        if(verbose)
820                            REprintf(" INSERT p12=(%d,%d) -> ind=%d -> i1 = %d\n",
821                                     p1,p2, ind, i1);
822    #endif
823    
824                        // shift the "upper values" *before* the insertion:
825                        for(int l = nnz-1; l >= i1; l--) {
826                            ri[l+1] = ri[l];
827                            rx[l+1] = rx[l];
828                        }
829                        ri[i1] = i__;
830                        rx[i1] = v;
831                        nnz++;
832    
833                        // the columns j "right" of the current one :
834                        for(k=j__ + 1; k <= ncol; k++)
835                            rp[k]++;
836                    }
837                }
838    #ifdef MATRIX_SUBASSIGN_VERBOSE
839                else if(verbose) REprintf("M_ij == v = %g\n", v);
840    #endif
841            }// for( ii )
842        }// for( jj )
843    
844        // now assign the i- and x- slots,  free memory and return :
845        Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_iSym,  INTSXP, nnz)), ri, nnz);
846        Memcpy(   REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), rx, nnz);
847        Free(ri);
848        Free(rx);
849        UNPROTECT(4);
850        return ans;
851    }
852    
853  SEXP Csparse_MatrixMarket(SEXP x, SEXP fname)  SEXP Csparse_MatrixMarket(SEXP x, SEXP fname)
854  {  {
855      FILE *f = fopen(CHAR(asChar(fname)), "w");      FILE *f = fopen(CHAR(asChar(fname)), "w");
# Line 619  Line 857 
857      if (!f)      if (!f)
858          error(_("failure to open file \"%s\" for writing"),          error(_("failure to open file \"%s\" for writing"),
859                CHAR(asChar(fname)));                CHAR(asChar(fname)));
860      if (!cholmod_l_write_sparse(f, AS_CHM_SP(x),      if (!cholmod_write_sparse(f, AS_CHM_SP(x),
861                                (CHM_SP)NULL, (char*) NULL, &c))                                (CHM_SP)NULL, (char*) NULL, &c))
862          error(_("cholmod_l_write_sparse returned error code"));          error(_("cholmod_write_sparse returned error code"));
863      fclose(f);      fclose(f);
864      return R_NilValue;      return R_NilValue;
865  }  }
# Line 835  Line 1073 
1073      if (cls[1] != 'g')      if (cls[1] != 'g')
1074          error(_("Only 'g'eneral sparse matrix types allowed"));          error(_("Only 'g'eneral sparse matrix types allowed"));
1075                                  /* allocate and populate the triplet */                                  /* allocate and populate the triplet */
1076      T = cholmod_l_allocate_triplet((size_t)nrow, (size_t)ncol, (size_t)nnz, 0,      T = cholmod_allocate_triplet((size_t)nrow, (size_t)ncol, (size_t)nnz, 0,
1077                                      xtype, &c);                                      xtype, &c);
1078      T->x = x;      T->x = x;
1079      tri = (int*)T->i;      tri = (int*)T->i;
# Line 845  Line 1083 
1083          trj[ii] = j[ii] - ((!mj && index1) ? 1 : 0);          trj[ii] = j[ii] - ((!mj && index1) ? 1 : 0);
1084      }      }
1085                                  /* create the cholmod_sparse structure */                                  /* create the cholmod_sparse structure */
1086      A = cholmod_l_triplet_to_sparse(T, nnz, &c);      A = cholmod_triplet_to_sparse(T, nnz, &c);
1087      cholmod_l_free_triplet(&T, &c);      cholmod_free_triplet(&T, &c);
1088                                  /* copy the information to the SEXP */                                  /* copy the information to the SEXP */
1089      ans = PROTECT(NEW_OBJECT(MAKE_CLASS(cls)));      ans = PROTECT(NEW_OBJECT(MAKE_CLASS(cls)));
1090  /* FIXME: This has been copied from chm_sparse_to_SEXP in chm_common.c */  /* FIXME: This has been copied from chm_sparse_to_SEXP in chm_common.c */
1091                                  /* allocate and copy common slots */                                  /* allocate and copy common slots */
1092      nnz = cholmod_l_nnz(A, &c);      nnz = cholmod_nnz(A, &c);
1093      dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2));      dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2));
1094      dims[0] = A->nrow; dims[1] = A->ncol;      dims[0] = A->nrow; dims[1] = A->ncol;
1095      Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, INTSXP, A->ncol + 1)), (int*)A->p, A->ncol + 1);      Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, INTSXP, A->ncol + 1)), (int*)A->p, A->ncol + 1);
# Line 864  Line 1102 
1102          error(_("code not yet written for cls = \"lgCMatrix\""));          error(_("code not yet written for cls = \"lgCMatrix\""));
1103      }      }
1104  /* FIXME: dimnames are *NOT* put there yet (if non-NULL) */  /* FIXME: dimnames are *NOT* put there yet (if non-NULL) */
1105      cholmod_l_free_sparse(&A, &c);      cholmod_free_sparse(&A, &c);
1106      UNPROTECT(1);      UNPROTECT(1);
1107      return ans;      return ans;
1108  }  }

Legend:
Removed from v.2646  
changed lines
  Added in v.2677

root@r-forge.r-project.org
ViewVC Help
Powered by ViewVC 1.0.0  
Thanks to:
Vienna University of Economics and Business University of Wisconsin - Madison Powered By FusionForge