ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
spmatrix.hpp
Go to the documentation of this file.
1/*
2 * $Id: spmatrix.hpp 313 2014-02-22 14:29:57Z jdl3 $
3 * Copyright (C) 2010, 2011, 2012, 2019, 2020, 2021 John D Lamb
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or (at
8 * your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
18 */
19
20#ifndef CCGSL_SPMATRIX_HPP
21#define CCGSL_SPMATRIX_HPP
22
23#include<cstddef>
24#include<gsl/gsl_spmatrix.h>
25#include<new>
26
27#include"exception.hpp"
28#include"vector.hpp"
29#include"matrix.hpp"
30
31namespace gsl {
38 class spmatrix {
39 public:
43 enum {
45 COO = GSL_SPMATRIX_COO,
47 CSC = GSL_SPMATRIX_CSC,
49 CSR = GSL_SPMATRIX_CSR,
51 TRIPLET = GSL_SPMATRIX_TRIPLET,
53 CCS = GSL_SPMATRIX_CCS,
55 CRS = GSL_SPMATRIX_CRS
56 };
61 ccgsl_pointer = 0;
62 count = 0; // initially nullptr will do
63 }
64 // Refines random access container
65 // Refines assignable
71 explicit spmatrix( size_t const n1, size_t const n2 ){
72 if( n1 > 0 and n2 > 0 ) ccgsl_pointer = gsl_spmatrix_alloc( n1, n2 );
73 else { ccgsl_pointer = new gsl_spmatrix; ccgsl_pointer->size1 = n1;
74 ccgsl_pointer->size2 = n2; ccgsl_pointer->data = 0; }
75 // just plausibly we could allocate spmatrix but not count
76 try { count = new size_t; } catch( std::bad_alloc& e ){
77 // try to tidy up before rethrowing
78 if( n1 > 0 and n2 > 0 ) gsl_spmatrix_free( ccgsl_pointer );
79 else delete ccgsl_pointer;
80 throw e;
81 }
82 *count = 1; // initially there is just one reference to ccgsl_pointer
83 }
92 explicit spmatrix( size_t const n1, size_t const n2, size_t nzmax, size_t sptype ){
93 if( n1 > 0 and n2 > 0 )
94 ccgsl_pointer = gsl_spmatrix_alloc_nzmax( n1, n2, nzmax, sptype );
95 else { ccgsl_pointer = new gsl_spmatrix; ccgsl_pointer->size1 = n1;
96 ccgsl_pointer->size2 = n2; ccgsl_pointer->data = 0; }
97 // just plausibly we could allocate spmatrix but not count
98 try { count = new size_t; } catch( std::bad_alloc& e ){
99 // try to tidy up before rethrowing
100 if( n1 > 0 and n2 > 0 ) gsl_spmatrix_free( ccgsl_pointer );
101 else delete ccgsl_pointer;
102 throw e;
103 }
104 *count = 1; // initially there is just one reference to ccgsl_pointer
105 }
111 explicit spmatrix( gsl_spmatrix* v ){
112 ccgsl_pointer = v;
113 // just plausibly we could fail to allocate count: no further action needed.
114 count = new size_t;
115 *count = 1; // initially there is just one reference to ccgsl_pointer
116 }
117#ifdef __GXX_EXPERIMENTAL_CXX0X__
122 spmatrix( std::initializer_list<std::initializer_list<double> > initializer_list ){
123 size_t const n1 = initializer_list.size();
124 size_t const n2 = initializer_list.begin()->size();
125 for( auto l : initializer_list ){
126 if( l.size() != n2 ){
127 gsl::exception e( "matrix rows have unequal sizes", __FILE__, __LINE__,
129 throw( e );
130 }
131 }
132 ccgsl_pointer = gsl_spmatrix_alloc( n1, n2 );
133 // just plausibly we could allocate spmatrix but not count
134 try { count = new size_t; } catch( std::bad_alloc& e ){
135 // try to tidy up before rethrowing
136 if( n1 > 0 and n2 > 0 ) gsl_spmatrix_free( ccgsl_pointer );
137 else delete ccgsl_pointer;
138 throw e;
139 }
140 *count = 1; // initially there is just one reference to ccgsl_pointer
141 // now copy
142 int r = 0;
143 for( auto row : initializer_list ){
144 size_t c = 0;
145 for( auto x : row ){ set( r, c, x ); ++c; }
146 ++r;
147 }
148 }
149#endif
150 // copy constructor
156 if( count != 0 ) ++*count; // spmatrix is now shared.
157 }
158 // assignment operator
164 // first, possibly delete anything pointed to by this
165 if( count == 0 or --*count == 0 ){
166 if( ccgsl_pointer != 0 and ccgsl_pointer->size1 > 0
167 and ccgsl_pointer->size2 > 0 ) gsl_spmatrix_free( ccgsl_pointer );
168 else delete ccgsl_pointer;
169 delete count;
170 }
171 // Then copy
173 count = v.count;
174 if( count != 0 ) ++*count; // block is now shared.
175 return *this;
176 }
177 // clone()
183 spmatrix clone() const {
184 spmatrix copy( get()->size1, get()->size2 );
185 // Now copy
186 gsl_spmatrix_memcpy( copy.get(), get() );
187 // return new object
188 return copy;
189 }
190 // destructor
195 if( count != 0 and --*count == 0 ){
196 // could have allocated null pointer
197 if( ccgsl_pointer != 0 ){
198 if( ccgsl_pointer->size1 > 0 and ccgsl_pointer->size2 > 0 )
199 gsl_spmatrix_free( ccgsl_pointer );
200 else delete ccgsl_pointer; }
201 delete count;
202 }
203 }
204 // Allow possibility of assigning from gsl_spmatrix without sharing
213 void wrap_gsl_spmatrix_without_ownership( gsl_spmatrix* v ){
214 if( count != 0 and --*count == 0 ){
215 // could have allocated null pointer
216 if( ccgsl_pointer != 0 ){
217 if( ccgsl_pointer->size1 != 0 and ccgsl_pointer->size1 != 0 )
218 gsl_spmatrix_free( ccgsl_pointer );
219 else delete ccgsl_pointer; }
220 }
221 ccgsl_pointer = v;
222 if(0 == count) count = new size_t;
223 *count = 2; // should never be able to delete ccgsl_pointer
224 }
228 void reset(){ spmatrix().swap( *this ); }
229#ifdef __GXX_EXPERIMENTAL_CXX0X__
235 std::swap( count, v.count );
236 v.ccgsl_pointer = nullptr;
237 }
244 spmatrix( std::move( v ) ).swap( *this );
245 return *this;
246 }
247#endif
248 public:
249 // Sizes
254 size_t size1() const { return ccgsl_pointer == 0 ? 0 : ccgsl_pointer->size1; }
259 size_t size2() const { return ccgsl_pointer == 0 ? 0 : ccgsl_pointer->size2; }
264 size_t nzmax() const { return ccgsl_pointer == 0 ? 0 : ccgsl_pointer->nzmax; }
269 int sptype() const { return ccgsl_pointer == 0 ? 0 : ccgsl_pointer->sptype; }
274 double* data() {
275 if( ccgsl_pointer == 0 ) gsl_error( "null vector", __FILE__, __LINE__, GSL_EFAULT );
276 return ccgsl_pointer->data; }
281 double const* data() const {
282 if( ccgsl_pointer == 0 ) gsl_error( "null vector", __FILE__, __LINE__, GSL_EFAULT );
283 return ccgsl_pointer->data; }
289 void swap( spmatrix& m ){
290 std::swap( ccgsl_pointer, m.ccgsl_pointer );
291 std::swap( count, m.count );
292 }
298 int realloc( size_t const nzmax ){ return gsl_spmatrix_realloc( nzmax, get() ); }
303 size_t nnz() const { return gsl_spmatrix_nnz( get() ); }
308 char const* type() const { return gsl_spmatrix_type( get() ); }
313 int set_zero(){ return gsl_spmatrix_set_zero( get() ); }
318 int tree_rebuild(){ return gsl_spmatrix_tree_rebuild( get() ); }
324 spmatrix csc() const {
325 spmatrix dest { size1(), size2() };
326 gsl_spmatrix_csc( dest.get(), get() );
327 return dest;
328 }
334 spmatrix csr() const {
335 spmatrix dest { size1(), size2() };
336 gsl_spmatrix_csr( dest.get(), get() );
337 return dest;
338 }
344 spmatrix compress( int const sptype ) const {
345 spmatrix dest { gsl_spmatrix_compress( get(), sptype ) };
346 return dest; }
352 spmatrix dest { gsl_spmatrix_compcol( get() ) };
353 return dest; }
358 spmatrix ccs() const {
359 spmatrix dest { gsl_spmatrix_ccs( get() ) };
360 return dest; }
365 spmatrix crs() const {
366 spmatrix dest { gsl_spmatrix_crs( get() ) };
367 return dest; }
373 int memcpy( spmatrix& dest ) const { return gsl_spmatrix_memcpy( dest.get(), get() ); }
380 int fprintf( FILE* stream, char const* format ) const {
381 return gsl_spmatrix_fprintf( stream, get(), format ); }
387 int fwrite( FILE* stream ) const {
388 return gsl_spmatrix_fwrite( stream, get() ); }
394 int fread( FILE* stream ){
395 return gsl_spmatrix_fread( stream, get() ); }
402 double get( size_t const i, size_t const j ) const {
403 return gsl_spmatrix_get( get(), i, j ); }
411 int set( size_t const i, size_t const j, double const x ){
412 return gsl_spmatrix_set( get(), i, j, x ); }
419 double* ptr( size_t const i, size_t const j ) const {
420 return gsl_spmatrix_ptr( get(), i, j ); }
421 // BEGIN REAL ONLY
428 int minmax( double& min_out, double& max_out ) const {
429 return gsl_spmatrix_minmax( get(), &min_out, &max_out ); }
436 int min_index( size_t& imin_out, size_t& jmin_out ) const {
437 return gsl_spmatrix_min_index( get(), &imin_out, &jmin_out ); }
438 // END REAL ONLY
444 int scale( double const x ){ return gsl_spmatrix_scale( get(), x ); }
450 int scale_columns( gsl::vector const& x ){
451 return gsl_spmatrix_scale_columns( get(), x.get() ); }
457 int scale_rows( gsl::vector const& x ){
458 return gsl_spmatrix_scale_rows( get(), x.get() ); }
464 int equal( spmatrix const& b ) const { return gsl_spmatrix_equal( get(), b.get() ); }
469 int transpose(){ return gsl_spmatrix_transpose( get() ); }
474 int transpose2(){ return gsl_spmatrix_transpose2( get() ); }
475 // Now the static functions
483 static int realloc( size_t const nzmax, gsl::spmatrix& spMatrix ){
484 return gsl_spmatrix_realloc( nzmax, spMatrix.get() ); }
490 static size_t nnz( spmatrix const& spMatrix ){
491 return gsl_spmatrix_nnz( spMatrix.get() ); }
497 static char const* type( spmatrix const& spMatrix ){
498 return gsl_spmatrix_type( spMatrix.get() ); }
504 static int set_zero( spmatrix& spMatrix ){
505 return gsl_spmatrix_set_zero( spMatrix.get() ); }
511 static int tree_rebuild( spmatrix& spMatrix ){
512 return gsl_spmatrix_tree_rebuild( spMatrix.get() ); }
519 static spmatrix csc( spmatrix const& spMatrix ){
520 spmatrix dest { spMatrix.size1(), spMatrix.size2() };
521 gsl_spmatrix_csc( dest.get(), spMatrix.get() );
522 return dest; }
529 static spmatrix csr( spmatrix const& spMatrix ){
530 spmatrix dest { spMatrix.size1(), spMatrix.size2() };
531 gsl_spmatrix_csr( dest.get(), spMatrix.get() );
532 return dest; }
539 static spmatrix compress( spmatrix const& spMatrix, int const sptype ){
540 spmatrix dest { gsl_spmatrix_compress( spMatrix.get(), sptype ) };
541 return dest; }
547 static spmatrix compcol( spmatrix const& spMatrix ){
548 spmatrix dest { gsl_spmatrix_compcol( spMatrix.get() ) };
549 return dest; }
555 static spmatrix ccs( spmatrix const& spMatrix ){
556 spmatrix dest { gsl_spmatrix_ccs( spMatrix.get() ) };
557 return dest; }
563 static spmatrix crs( spmatrix const& spMatrix ){
564 spmatrix dest { gsl_spmatrix_crs( spMatrix.get() ) };
565 return dest; }
572 static int memcpy( spmatrix& dest, spmatrix const& src ){
573 return gsl_spmatrix_memcpy( dest.get(), src.get() ); }
581 static int fprintf( FILE* stream, spmatrix const& spMatrix, char const* format ){
582 return gsl_spmatrix_fprintf( stream, spMatrix.get(), format ); }
588 static spmatrix fscanf( FILE* stream ){
589 return spmatrix { gsl_spmatrix_fscanf( stream ) }; }
596 static int fwrite( FILE* stream, spmatrix& spMatrix ){
597 return gsl_spmatrix_fwrite( stream, spMatrix.get() ); }
604 static int fread( FILE* stream, spmatrix& spMatrix ){
605 return gsl_spmatrix_fread( stream, spMatrix.get() ); }
613 static int add( spmatrix& c, spmatrix const& a, spmatrix const& b ){
614 return gsl_spmatrix_add( c.get(), a.get(), b.get() ); }
621 static spmatrix add( spmatrix const& a, spmatrix const& b ){
622 spmatrix c { a.size1(), a.size2() };
623 gsl_spmatrix_add( c.get(), a.get(), b.get() );
624 return c;
625 }
626#ifndef GSL_DISABLE_DEPRECATED
633 static int add_to_dense( gsl::matrix& a, spmatrix const&b ){
634 return gsl_spmatrix_add_to_dense( a.get(), b.get() ); }
635#endif
642 static int dense_add( gsl::matrix& a, spmatrix const&b ){
643 return gsl_spmatrix_dense_add( a.get(), b.get() ); }
650 static int dense_sub( gsl::matrix& a, spmatrix const&b ){
651 return gsl_spmatrix_dense_sub( a.get(), b.get() ); }
657 static double norm1( gsl::spmatrix const& a ){
658 return gsl_spmatrix_norm1( a.get() ); }
665 static int d2sp( gsl::spmatrix& S, gsl::matrix const& A ){
666 return gsl_spmatrix_d2sp( S.get(), A.get() ); }
673 static int sp2d( gsl::matrix& A, spmatrix const& S ){
674 return gsl_spmatrix_sp2d( A.get(), S.get() ); }
681 static int equal( spmatrix const& a, spmatrix const& b ){
682 return gsl_spmatrix_equal( a.get(), b.get() ); }
688 static int transpose( spmatrix& a ){ return gsl_spmatrix_transpose( a.get() ); }
694 static int transpose2( spmatrix& a ){ return gsl_spmatrix_transpose2( a.get() ); }
701 static int transpose_memcpy( spmatrix& dest, spmatrix const& src ){
702 return gsl_spmatrix_transpose_memcpy( dest.get(), src.get() ); }
703 private:
707 gsl_spmatrix* ccgsl_pointer;
711 size_t* count;
712 public:
713 // shared reference functions
718 gsl_spmatrix* get(){ return ccgsl_pointer; }
723 gsl_spmatrix const* get() const { return ccgsl_pointer; }
729 bool unique() const { return count != 0 and *count == 1; }
734 size_t use_count() const { return count == 0 ? 0 : *count; }
740#ifdef __GXX_EXPERIMENTAL_CXX0X__
741 explicit
742#endif
743 operator bool() const { return ccgsl_pointer != 0; }
744 };
745}
746#endif
This class is used to handle gsl exceptions so that gsl can use these rather than the GSL error handl...
Definition: exception.hpp:387
@ GSL_EBADLEN
matrix, vector lengths are not conformant
Definition: exception.hpp:490
This class handles matrix objects as shared handles.
Definition: matrix.hpp:72
gsl_matrix * get()
Get the gsl_matrix.
Definition: matrix.hpp:1207
This class handles sparse matrix objects as shared handles.
Definition: spmatrix.hpp:38
spmatrix ccs() const
C++ version of gsl_spmatrix_ccs().
Definition: spmatrix.hpp:358
spmatrix compcol() const
C++ version of gsl_spmatrix_compcol().
Definition: spmatrix.hpp:351
static double norm1(gsl::spmatrix const &a)
C++ version of gsl_spmatrix_norm1().
Definition: spmatrix.hpp:657
double const * data() const
Give access to the data block.
Definition: spmatrix.hpp:281
int memcpy(spmatrix &dest) const
C++ version of gsl_spmatrix_memcpy().
Definition: spmatrix.hpp:373
char const * type() const
C++ version of gsl_spmatrix_type().
Definition: spmatrix.hpp:308
spmatrix(size_t const n1, size_t const n2, size_t nzmax, size_t sptype)
This constructor creates a new matrix with n1 rows and n2 columns.
Definition: spmatrix.hpp:92
static int sp2d(gsl::matrix &A, spmatrix const &S)
C++ version of gsl_spmatrix_sp2d().
Definition: spmatrix.hpp:673
int scale_rows(gsl::vector const &x)
C++ version of gsl_spmatrix_scale_rows().
Definition: spmatrix.hpp:457
static spmatrix fscanf(FILE *stream)
C++ version of gsl_spmatrix_fscanf().
Definition: spmatrix.hpp:588
int scale(double const x)
C++ version of gsl_spmatrix_scale().
Definition: spmatrix.hpp:444
int scale_columns(gsl::vector const &x)
C++ version of gsl_spmatrix_scale_columns().
Definition: spmatrix.hpp:450
int sptype() const
The type of the spmatrix.
Definition: spmatrix.hpp:269
static int tree_rebuild(spmatrix &spMatrix)
C++ version of gsl_spmatrix_tree_rebuild().
Definition: spmatrix.hpp:511
static int realloc(size_t const nzmax, gsl::spmatrix &spMatrix)
C++ version of gsl_spmatrix_realloc(): matches a function in gsl that is intended for internal use.
Definition: spmatrix.hpp:483
static int fread(FILE *stream, spmatrix &spMatrix)
C++ version of gsl_spmatrix_fread().
Definition: spmatrix.hpp:604
int realloc(size_t const nzmax)
C++ version of gsl_spmatrix_realloc().
Definition: spmatrix.hpp:298
gsl_spmatrix const * get() const
Get the gsl_spmatrix.
Definition: spmatrix.hpp:723
static int transpose_memcpy(spmatrix &dest, spmatrix const &src)
C++ version of gsl_spmatrix_transpose_memcpy().
Definition: spmatrix.hpp:701
static spmatrix crs(spmatrix const &spMatrix)
C++ version of gsl_spmatrix_crs().
Definition: spmatrix.hpp:563
double * ptr(size_t const i, size_t const j) const
C++ version of gsl_spmatrix_ptr().
Definition: spmatrix.hpp:419
static spmatrix compcol(spmatrix const &spMatrix)
C++ version of gsl_spmatrix_compcol().
Definition: spmatrix.hpp:547
size_t nzmax() const
The maximum number of nonzero elements.
Definition: spmatrix.hpp:264
int tree_rebuild()
C++ version of gsl_spmatrix_tree_rebuild().
Definition: spmatrix.hpp:318
int set_zero()
C++ version of gsl_spmatrix_set_zero().
Definition: spmatrix.hpp:313
gsl_spmatrix * ccgsl_pointer
The shared pointer.
Definition: spmatrix.hpp:707
spmatrix(size_t const n1, size_t const n2)
This constructor creates a new matrix with n1 rows and n2 columns.
Definition: spmatrix.hpp:71
double * data()
Give access to the data block.
Definition: spmatrix.hpp:274
int transpose2()
C++ version of gsl_spmatrix_transpose2().
Definition: spmatrix.hpp:474
static int dense_sub(gsl::matrix &a, spmatrix const &b)
C++ version of gsl_spmatrix_dense_sub().
Definition: spmatrix.hpp:650
size_t size2() const
The number of columns of the matrix.
Definition: spmatrix.hpp:259
spmatrix(gsl_spmatrix *v)
Could construct from a gsl_spmatrix.
Definition: spmatrix.hpp:111
int set(size_t const i, size_t const j, double const x)
C++ version of gsl_spmatrix_set().
Definition: spmatrix.hpp:411
static int memcpy(spmatrix &dest, spmatrix const &src)
C++ version of gsl_spmatrix_memcpy().
Definition: spmatrix.hpp:572
static size_t nnz(spmatrix const &spMatrix)
C++ version of gsl_spmatrix_nnz().
Definition: spmatrix.hpp:490
gsl_spmatrix * get()
Get the gsl_spmatrix.
Definition: spmatrix.hpp:718
int minmax(double &min_out, double &max_out) const
C++ version of gsl_spmatrix_minmax().
Definition: spmatrix.hpp:428
spmatrix & operator=(spmatrix &&v)
Move operator.
Definition: spmatrix.hpp:243
~spmatrix()
The destructor only deletes the pointers if count reaches zero.
Definition: spmatrix.hpp:194
int fprintf(FILE *stream, char const *format) const
C++ version of gsl_spmatrix_fprintf().
Definition: spmatrix.hpp:380
static int transpose2(spmatrix &a)
C++ version of gsl_spmatrix_transpose2().
Definition: spmatrix.hpp:694
void wrap_gsl_spmatrix_without_ownership(gsl_spmatrix *v)
This function is intended mainly for internal use.
Definition: spmatrix.hpp:213
spmatrix(spmatrix const &v)
The copy constructor.
Definition: spmatrix.hpp:155
size_t use_count() const
Find how many spmatrix objects share this pointer.
Definition: spmatrix.hpp:734
spmatrix csc() const
C++ version of gsl_spmatrix_csc().
Definition: spmatrix.hpp:324
spmatrix(std::initializer_list< std::initializer_list< double > > initializer_list)
Could construct from a std::initializer_list in C++11 and later.
Definition: spmatrix.hpp:122
spmatrix compress(int const sptype) const
C++ version of gsl_spmatrix_compress().
Definition: spmatrix.hpp:344
static spmatrix compress(spmatrix const &spMatrix, int const sptype)
C++ version of gsl_spmatrix_compress().
Definition: spmatrix.hpp:539
size_t size1() const
The number of rows of the matrix.
Definition: spmatrix.hpp:254
static spmatrix csc(spmatrix const &spMatrix)
C++ version of gsl_spmatrix_csc().
Definition: spmatrix.hpp:519
int fread(FILE *stream)
C++ version of gsl_spmatrix_fread().
Definition: spmatrix.hpp:394
static int equal(spmatrix const &a, spmatrix const &b)
C++ version of gsl_spmatrix_equal().
Definition: spmatrix.hpp:681
spmatrix crs() const
C++ version of gsl_spmatrix_crs().
Definition: spmatrix.hpp:365
spmatrix & operator=(spmatrix const &v)
The assignment operator.
Definition: spmatrix.hpp:163
static int d2sp(gsl::spmatrix &S, gsl::matrix const &A)
C++ version of gsl_spmatrix_d2sp().
Definition: spmatrix.hpp:665
double get(size_t const i, size_t const j) const
C++ version of gsl_spmatrix_get().
Definition: spmatrix.hpp:402
static spmatrix add(spmatrix const &a, spmatrix const &b)
C++ version of gsl_spmatrix_add().
Definition: spmatrix.hpp:621
int fwrite(FILE *stream) const
C++ version of gsl_spmatrix_fwrite().
Definition: spmatrix.hpp:387
static int fwrite(FILE *stream, spmatrix &spMatrix)
C++ version of gsl_spmatrix_fwrite().
Definition: spmatrix.hpp:596
int equal(spmatrix const &b) const
C++ version of gsl_spmatrix_equal().
Definition: spmatrix.hpp:464
spmatrix clone() const
The clone function.
Definition: spmatrix.hpp:183
bool unique() const
Find if this is the only object sharing the gsl_spmatrix.
Definition: spmatrix.hpp:729
static int set_zero(spmatrix &spMatrix)
C++ version of gsl_spmatrix_set_zero().
Definition: spmatrix.hpp:504
spmatrix()
The default constructor is only really useful for assigning to.
Definition: spmatrix.hpp:60
spmatrix csr() const
C++ version of gsl_spmatrix_csr().
Definition: spmatrix.hpp:334
static int transpose(spmatrix &a)
C++ version of gsl_spmatrix_transpose().
Definition: spmatrix.hpp:688
spmatrix(spmatrix &&v)
Move constructor.
Definition: spmatrix.hpp:234
static char const * type(spmatrix const &spMatrix)
C++ version of gsl_spmatrix_type().
Definition: spmatrix.hpp:497
int min_index(size_t &imin_out, size_t &jmin_out) const
C++ version of gsl_spmatrix_min_index().
Definition: spmatrix.hpp:436
static spmatrix csr(spmatrix const &spMatrix)
C++ version of gsl_spmatrix_csr().
Definition: spmatrix.hpp:529
static int dense_add(gsl::matrix &a, spmatrix const &b)
C++ version of gsl_spmatrix_dense_add().
Definition: spmatrix.hpp:642
static int fprintf(FILE *stream, spmatrix const &spMatrix, char const *format)
C++ version of gsl_spmatrix_fprintf().
Definition: spmatrix.hpp:581
static int add_to_dense(gsl::matrix &a, spmatrix const &b)
C++ version of gsl_spmatrix_add_to_dense().
Definition: spmatrix.hpp:633
size_t nnz() const
C++ version of gsl_spmatrix_nnz().
Definition: spmatrix.hpp:303
int transpose()
C++ version of gsl_spmatrix_transpose().
Definition: spmatrix.hpp:469
static int add(spmatrix &c, spmatrix const &a, spmatrix const &b)
C++ version of gsl_spmatrix_add().
Definition: spmatrix.hpp:613
void swap(spmatrix &m)
Swap two spmatrix objects.
Definition: spmatrix.hpp:289
void reset()
Stop sharing ownership of the shared pointer.
Definition: spmatrix.hpp:228
@ CCS
GSL_SPMATRIX_CCS.
Definition: spmatrix.hpp:53
@ TRIPLET
GSL_SPMATRIX_TRIPLET.
Definition: spmatrix.hpp:51
@ COO
GSL_SPMATRIX_COO.
Definition: spmatrix.hpp:45
@ CRS
GSL_SPMATRIX_CRS.
Definition: spmatrix.hpp:55
@ CSR
GSL_SPMATRIX_CSR.
Definition: spmatrix.hpp:49
@ CSC
GSL_SPMATRIX_CSC.
Definition: spmatrix.hpp:47
static spmatrix ccs(spmatrix const &spMatrix)
C++ version of gsl_spmatrix_ccs().
Definition: spmatrix.hpp:555
size_t * count
The shared reference count.
Definition: spmatrix.hpp:711
This class handles vector objects as shared handles.
Definition: vector.hpp:74
gsl_vector * get()
Get the gsl_vector.
Definition: vector.hpp:1320
double b(int order, double qq)
C++ version of gsl_sf_mathieu_b().
Definition: sf_mathieu.hpp:298
double a(int order, double qq)
C++ version of gsl_sf_mathieu_a().
Definition: sf_mathieu.hpp:272
The gsl package creates an interface to the GNU Scientific Library for C++.
Definition: blas.hpp:34