Class: Numo::GSL::SpMatrix

Inherits:
Object
  • Object
show all
Defined in:
ext/numo/gsl/spmatrix/gsl_spmatrix.c

Constant Summary

CCS =

This flag specifies compressed column storage.

INT2FIX(GSL_SPMATRIX_CCS)
CRS =

This flag specifies compressed row storage.

INT2FIX(GSL_SPMATRIX_CRS)
TRIPLET =

This flag specifies triplet storage.

INT2FIX(GSL_SPMATRIX_TRIPLET)

Class Method Summary collapse

Instance Method Summary collapse

Class Method Details

.d2sp(nary) ⇒ Numo::GSL::SpMatrix

This function converts the dense matrix A into sparse triplet format and stores the result in S.

Parameters:

  • nary (Numo::DFloat)

Returns:



536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 536

static VALUE
spmatrix_s_d2sp(VALUE klass, VALUE nary)
{
    gsl_spmatrix *s;
    gsl_matrix *m;
    narray_t *na;

    nary = cast_2d_contiguous(nary,cDF);
    GetNArray(nary,na);

    ALLOCA_GSL_MATRIX_FROM_NARRAY_R(nary,m);
    s = gsl_spmatrix_alloc(na->shape[0], na->shape[1]);

    gsl_spmatrix_d2sp(s, m);

    RB_GC_GUARD(nary);
    return TypedData_Wrap_Struct(cSpMatrix, &spmatrix_data_type, (void*)s);
}

.new(n1, n2, [nzmax,sptype]) ⇒ Object

This function allocates a sparse matrix of size n1-by-n2 and initializes it to all zeros. If the size of the matrix is not known at allocation time, both n1 and n2 may be set to 1, and they will automatically grow as elements are added to the matrix. The parameter nzmax specifies the maximum number of non-zero elements which will be added to the matrix. It does not need to be precisely known in advance, since storage space will automatically grow using gsl_spmatrix_realloc if nzmax is not large enough. Accurate knowledge of this parameter reduces the number of reallocation calls required. The parameter sptype specifies the storage format of the sparse matrix. Possible values are The allocated gsl_spmatrix structure is of size O(nzmax).

Parameters:

  • n1 (Integer)
  • n2 (Integer)
  • nzmax (Ingeger)
  • sptype (Ingeger)

    (default = GSL_SPMATRIX_TRIPLET)



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
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 112

static VALUE
spmatrix_s_new(int argc, VALUE *argv, VALUE klass)
{
    gsl_spmatrix *w;
    int narg;
    size_t sptype = GSL_SPMATRIX_TRIPLET;
    VALUE  n1, n2, v3, v4;

    narg = rb_scan_args(argc,argv,"22",&n1,&n2,&v3,&v4);
    switch(narg) {
    case 4:
        sptype = NUM2SIZET(v4);
    case 3:
        w = gsl_spmatrix_alloc_nzmax(NUM2SIZET(n1),NUM2SIZET(n2),NUM2SIZET(v3),sptype);
        break;
    case 2:
        w = gsl_spmatrix_alloc(NUM2SIZET(n1),NUM2SIZET(n2));
        break;
    default:
        rb_raise(rb_eArgError,"invalid number of argument: %d for 2..4",argc);
    }
    if (!w) {
        rb_raise(rb_eNoMemError,"fail to allocate struct");
    }
    return TypedData_Wrap_Struct(cSpMatrix, &spmatrix_data_type, (void*)w);
}

Instance Method Details

#add(other) ⇒ Numo::GSL::SpMatrix

This function computes the sum c = a + b. The three matrices must have the same dimensions and be stored in a compressed format.

Parameters:

Returns:



376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 376

static VALUE
spmatrix_add(VALUE self, VALUE other)
{
    gsl_spmatrix *a, *b, *c;
    VALUE result;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, a);
    TypedData_Get_Struct(other, gsl_spmatrix, &spmatrix_data_type, b);

    c = gsl_spmatrix_alloc_nzmax(a->size1, a->size2, a->nzmax+b->nzmax, a->sptype);
    result = TypedData_Wrap_Struct(cSpMatrix, &spmatrix_data_type, (void*)c);

    gsl_spmatrix_add(c,a,b);
    return result;
}

#ccsNumo::GSL::SpMatrix

This function creates a sparse matrix in compressed column format from the input sparse matrix T which must be in triplet format. A pointer to a newly allocated matrix is returned. The calling function should free the newly allocated matrix when it is no longer needed.

Returns:



486
487
488
489
490
491
492
493
494
495
496
497
498
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 486

static VALUE
spmatrix_ccs(VALUE self)
{
    gsl_spmatrix *x, *w;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, x);

    w = gsl_spmatrix_ccs(x);
    if (!w) {
        rb_raise(rb_eNoMemError,"fail to allocate struct");
    }
    return TypedData_Wrap_Struct(cSpMatrix, &spmatrix_data_type, (void*)w);
}

#crsNumo::GSL::SpMatrix

This function creates a sparse matrix in compressed row format from the input sparse matrix T which must be in triplet format. A pointer to a newly allocated matrix is returned. The calling function should free the newly allocated matrix when it is no longer needed.

Returns:



511
512
513
514
515
516
517
518
519
520
521
522
523
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 511

static VALUE
spmatrix_crs(VALUE self)
{
    gsl_spmatrix *x, *w;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, x);

    w = gsl_spmatrix_crs(x);
    if (!w) {
        rb_raise(rb_eNoMemError,"fail to allocate struct");
    }
    return TypedData_Wrap_Struct(cSpMatrix, &spmatrix_data_type, (void*)w);
}

#equal(b) ⇒ Numo::GSL::SpMatrix

This function returns 1 if the matrices a and b are equal (by comparison of element values) and 0 otherwise. The matrices a and b must be in the same sparse storage format for comparison.

Parameters:

Returns:



442
443
444
445
446
447
448
449
450
451
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 442

static VALUE
spmatrix_equal(VALUE self, VALUE other)
{
    gsl_spmatrix *w, *x;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);
    TypedData_Get_Struct(other, gsl_spmatrix, &spmatrix_data_type, x);
    gsl_spmatrix_equal(w,x);
    return self;
}

#get(i, j) ⇒ DFloat

This function returns element (i,j) of the matrix m. The matrix may be in triplet or compressed format.

Parameters:

  • i (UInt)
  • j (UInt)

Returns:

  • (DFloat)

    result



199
200
201
202
203
204
205
206
207
208
209
210
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 199

static VALUE
spmatrix_get(VALUE self, VALUE v1, VALUE v2)
{
    gsl_spmatrix *w;
    ndfunc_arg_in_t ain[2] = {{cSZ,0},{cSZ,0}};
    ndfunc_arg_out_t aout[1] = {{cDF,0}};
    ndfunc_t ndf = {iter_spmatrix_get, STRIDE_LOOP|NDF_EXTRACT, 2,1, ain,aout};

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);

    return na_ndloop3(&ndf, w, 2, v1, v2);
}

#memcpy(src) ⇒ Numo::GSL::SpMatrix

This function copies the elements of the sparse matrix src into dest. The two matrices must have the same dimensions and be in the same storage format.

Parameters:

Returns:



289
290
291
292
293
294
295
296
297
298
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 289

static VALUE
spmatrix_memcpy(VALUE self, VALUE other)
{
    gsl_spmatrix *w, *x;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);
    TypedData_Get_Struct(other, gsl_spmatrix, &spmatrix_data_type, x);
    gsl_spmatrix_memcpy(w,x);
    return self;
}

#minmaxArray

This function returns the minimum and maximum elements of the matrix m, storing them in min_out and max_out, and searching only the non-zero values.

Returns:

  • (Array)

    array of [min_out,max_out]



463
464
465
466
467
468
469
470
471
472
473
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 463

static VALUE
spmatrix_minmax(VALUE self)
{
    gsl_spmatrix *w;
    double d1, d2;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);

    gsl_spmatrix_minmax(w, &d1, &d2);
    return rb_assoc_new(DBL2NUM(d1),DBL2NUM(d2));
}

#nnzInteger

This function returns the number of non-zero elements in m.

Returns:

  • (Integer)


421
422
423
424
425
426
427
428
429
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 421

static VALUE
spmatrix_nnz(VALUE self)
{
    gsl_spmatrix *w;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);

    return SIZET2NUM(gsl_spmatrix_nnz(w));
}

#realloc(m) ⇒ Numo::GSL::SpMatrix

This function reallocates the storage space for m to accomodate nzmax non-zero elements. It is typically called internally by gsl_spmatrix_set if the user wants to add more elements to the sparse matrix than the previously specified nzmax.

Parameters:

  • m (Integer)

Returns:



151
152
153
154
155
156
157
158
159
160
161
162
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 151

static VALUE
spmatrix_realloc(VALUE self, VALUE v1)
{
    gsl_spmatrix *w;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);

    
    gsl_spmatrix_realloc(NUM2SIZET(v1),w);
    
    return self;
}

#scale(x) ⇒ Numo::GSL::SpMatrix

This function scales all elements of the matrix m by the constant factor x. The result m(i,j) \leftarrow x m(i,j) is stored in m.

Parameters:

  • x (Float)

Returns:



402
403
404
405
406
407
408
409
410
411
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 402

static VALUE
spmatrix_scale(VALUE self, VALUE v1)
{
    gsl_spmatrix *w;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);

    gsl_spmatrix_scale(w, NUM2DBL(v1));
    return self;
}

#set(i, j, j) ⇒ DFloat

This function sets element (i,j) of the matrix m to the value x. The matrix must be in triplet representation.

Parameters:

  • i (UInt)
  • j (UInt)
  • x (DFloat)

Returns:

  • (DFloat)

    self



248
249
250
251
252
253
254
255
256
257
258
259
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 248

static VALUE
spmatrix_set(VALUE self, VALUE v1, VALUE v2, VALUE v3)
{
    gsl_spmatrix *w;
    ndfunc_arg_in_t ain[3] = {{cSZ,0},{cSZ,0},{cDF,0}};
    ndfunc_t ndf = {iter_spmatrix_set, STRIDE_LOOP|NDF_EXTRACT, 3,0, ain,0};

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);

    na_ndloop3(&ndf, w, 3, v1, v2, v3);
    return self;
}

#set_zeroObject

This function sets (or resets) all the elements of the matrix m to zero.



268
269
270
271
272
273
274
275
276
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 268

static VALUE
spmatrix_set_zero(VALUE self)
{
    gsl_spmatrix *w;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);
    gsl_spmatrix_set_zero(w);
    return self;
}

#sp2dDFloat

This function converts the sparse matrix S into a dense matrix and stores the result in A. S must be in triplet format.

Returns:

  • (DFloat)

    result



564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 564

static VALUE
spmatrix_sp2d(VALUE self)
{
    gsl_spmatrix *s;
    gsl_matrix *m;
    size_t shape[2];
    VALUE result;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, s);

    shape[0] = s->size1;
    shape[1] = s->size2;
    result = rb_narray_new(cDF, 2, shape);

    ALLOCA_GSL_MATRIX_FROM_NARRAY_W(result,m);

    gsl_spmatrix_sp2d(m,s);
    return result;
}

#transposeObject

This function replaces the matrix m by its transpose, preserving the storage format of the input matrix. Currently, only triplet matrix inputs are supported.



332
333
334
335
336
337
338
339
340
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 332

static VALUE
spmatrix_transpose(VALUE self)
{
    gsl_spmatrix *w;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);
    gsl_spmatrix_transpose(w);
    return self;
}

#transpose2Object

This function replaces the matrix m by its transpose, but changes the storage format for compressed matrix inputs. Since compressed column storage is the transpose of compressed row storage, this function simply converts a CCS matrix to CRS and vice versa. This is the most efficient way to transpose a compressed storage matrix, but the user should note that the storage format of their compressed matrix will change on output. For triplet matrices, the output matrix is also in triplet storage.



356
357
358
359
360
361
362
363
364
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 356

static VALUE
spmatrix_transpose2(VALUE self)
{
    gsl_spmatrix *w;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);
    gsl_spmatrix_transpose2(w);
    return self;
}

#transpose_memcpy(src) ⇒ Numo::GSL::SpMatrix

This function copies the transpose of the sparse matrix src into dest. The dimensions of dest must match the transpose of the matrix src. Also, both matrices must use the same sparse storage format.

Parameters:

Returns:



312
313
314
315
316
317
318
319
320
321
# File 'ext/numo/gsl/spmatrix/gsl_spmatrix.c', line 312

static VALUE
spmatrix_transpose_memcpy(VALUE self, VALUE other)
{
    gsl_spmatrix *w, *x;

    TypedData_Get_Struct(self, gsl_spmatrix, &spmatrix_data_type, w);
    TypedData_Get_Struct(other, gsl_spmatrix, &spmatrix_data_type, x);
    gsl_spmatrix_transpose_memcpy(w,x);
    return self;
}