ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
matrix.hpp
Go to the documentation of this file.
1/*
2 * $Id: matrix.hpp 313 2014-02-22 14:29:57Z jdl3 $
3 * Copyright (C) 2010, 2011, 2012, 2019, 2020, 2021, 2024 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_MATRIX_HPP
21#define CCGSL_MATRIX_HPP
22
23#include<gsl/gsl_matrix.h>
24#include<gsl/gsl_permute_matrix_double.h>
25#include<new>
26
27#include"exception.hpp"
28#include"vector.hpp"
29#include"permutation.hpp"
30
31// This file is used as a template
32
33namespace gsl {
34 namespace multiroot {
35 // declare function_fdf class
36 class function_fdf;
37 }
38 namespace multilarge {
39 namespace nlinear {
40 // declare class
41 class function_fdf;
42 }
43 }
44 namespace multiroot {
45 // declare function_fdf class
46 class function_fdf;
47 }
72 class matrix {
76 public:
81 owns_data = false;
82 ccgsl_pointer = 0;
83 count = 0; // initially nullptr will do
84 }
85 // Refines random access container
86 // Refines assignable
92 explicit matrix( size_t const n1, size_t const n2 ) : owns_data(true){
93 if( n1 > 0 and n2 > 0 ) ccgsl_pointer = gsl_matrix_alloc( n1, n2 );
94 else { ccgsl_pointer = new gsl_matrix; ccgsl_pointer->size1 = n1;
95 ccgsl_pointer->size2 = n2; ccgsl_pointer->data = 0; }
96 // just plausibly we could allocate matrix but not count
97 try { count = new size_t; } catch( std::bad_alloc& e ){
98 // try to tidy up before rethrowing
99 if( n1 > 0 and n2 > 0 ) gsl_matrix_free( ccgsl_pointer );
100 else delete ccgsl_pointer;
101 throw e;
102 }
103 *count = 1; // initially there is just one reference to ccgsl_pointer
104 }
110 explicit matrix( gsl_matrix* v ) : owns_data(true){
111 ccgsl_pointer = v;
112 // just plausibly we could fail to allocate count: no further action needed.
113 count = new size_t;
114 *count = 1; // initially there is just one reference to ccgsl_pointer
115 }
116#ifdef __GXX_EXPERIMENTAL_CXX0X__
121 matrix( std::initializer_list<std::initializer_list<double> > initializer_list ) : owns_data(true){
122 size_t const n1 = initializer_list.size();
123 size_t const n2 = initializer_list.begin()->size();
124 for( auto l : initializer_list ){
125 if( l.size() != n2 ){
126 gsl::exception e( "matrix rows have unequal sizes", __FILE__, __LINE__,
128 throw( e );
129 }
130 }
131 ccgsl_pointer = gsl_matrix_alloc( n1, n2 );
132 // just plausibly we could allocate matrix but not count
133 try { count = new size_t; } catch( std::bad_alloc& e ){
134 // try to tidy up before rethrowing
135 if( n1 > 0 and n2 > 0 ) gsl_matrix_free( ccgsl_pointer );
136 else delete ccgsl_pointer;
137 throw e;
138 }
139 *count = 1; // initially there is just one reference to ccgsl_pointer
140 // now copy
141 int r = 0;
142 for( auto row : initializer_list ){
143 size_t c = 0;
144 for( auto x : row ){ set( r, c, x ); ++c; }
145 ++r;
146 }
147 }
148#endif
149 // copy constructor
154 matrix( matrix const& v ) : owns_data(true),
156 if( count != 0 ) ++*count; // matrix is now shared.
157 }
158 // assignment operator
163 matrix& operator=( matrix const& v ){
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 and ccgsl_pointer->size2 > 0 ) gsl_matrix_free( ccgsl_pointer );
167 else delete ccgsl_pointer;
168 delete count;
169 }
170 // Then copy
173 count = v.count;
174 if( count != 0 ) ++*count; // block is now shared.
175 return *this;
176 }
177 // clone()
183 matrix clone() const {
184 matrix copy( get()->size1, get()->size2 );
185 // Now copy
186 gsl_matrix_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 and owns_data ){
198 if( ccgsl_pointer->size1 > 0 and ccgsl_pointer->size2 > 0 )
199 gsl_matrix_free( ccgsl_pointer );
200 else delete ccgsl_pointer; }
201 delete count;
202 }
203 }
204
205 // Allow possibility of assigning from gsl_matrix without sharing
215 if( count != 0 and --*count == 0 ){
216 // could have allocated null pointer
217 if( ccgsl_pointer != 0 ){
218 if( ccgsl_pointer->size1 != 0 and ccgsl_pointer->size1 != 0 )
219 gsl_matrix_free( ccgsl_pointer );
220 else delete ccgsl_pointer; }
221 }
222 ccgsl_pointer = v;
223 if(0 == count) count = new size_t;
224 *count = 1;
225 owns_data = false; // should never be able to delete ccgsl_pointer
226 }
230 void reset(){ matrix().swap( *this ); }
231#ifdef __GXX_EXPERIMENTAL_CXX0X__
237 ccgsl_pointer( v.ccgsl_pointer ), count( nullptr ){
238 std::swap( count, v.count );
239 v.ccgsl_pointer = nullptr;
240 }
247 matrix( std::move( v ) ).swap( *this );
248 return *this;
249 }
250#endif
251 private:
255 struct vector_ptr : public vector {
260 vector_ptr( vector const& v ) : vector( v ){}
265 vector& operator*(){ return *this; }
270 vector* operator->(){ return this; }
271 };
276 template<typename container,typename content,bool reverse_t> class iterator_base {
277 friend class vector;
278 public:
282 typedef std::bidirectional_iterator_tag iterator_category;
295 public:
296 // // operator*
302 // Always try to return something
303 static content something;
304 // First check that iterator is initialised.
305 if( v == 0 ){
306 gsl_error( "iterator not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
307 return something;
308 } else if( v->ccgsl_pointer == 0 ){
309 gsl_error( "matrix not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
310 return something;
311 }
312 // Check that position make sense
313 if( position >= v->size1() ){
314 gsl_error( "trying to dereference end()", __FILE__, __LINE__, exception::GSL_EFAILED );
315 return something;
316 }
317 // position and v are valid: return data
318 return v->row( position );
319 }
320 // // operator->
326 // Always try to return something
327 static content something_base;
328 static vector_ptr something( something_base );
329 // First check that iterator is initialised.
330 if( v == 0 ){
331 gsl_error( "iterator not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
332 return something;
333 } else if( v->ccgsl_pointer == 0 ){
334 gsl_error( "matrix not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
335 return something;
336 }
337 // Check that position make sense
338 if( position >= v->size1() ){
339 gsl_error( "trying to dereference end()", __FILE__, __LINE__, exception::GSL_EFAILED );
340 return something;
341 }
342 // position and v are valid: return data
343 vector_ptr ptr( v->row( position ) );
344#ifdef __GXX_EXPERIMENTAL_CXX0X__
345 return std::move( ptr );
346#else
347 return ptr;
348#endif
349 }
356 return this->v == i.v and this->position == i.position;
357 }
364 return not this->operator==( i );
365 }
366 protected:
371 void increment(){
372 // Only makes sense if v points to a matrix
373 if( v == 0 ){
374 gsl_error( "iterator not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
375 return;
376 } else if( v->ccgsl_pointer == 0 ){
377 gsl_error( "matrix not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
378 return;
379 }
380 // increment position and check against size
381 if( reverse_t ){ if( position >= 0 ) --position; }
382 else { if( position < v->size1() ) ++position; }
383 }
388 void decrement(){
389 // Only makes sense if v points to a matrix
390 if( v == 0 ){
391 gsl_error( "iterator not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
392 return;
393 } else if( v->ccgsl_pointer == 0 ){
394 gsl_error( "matrix not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
395 return;
396 }
397 // decrement position and check against size
398 if( reverse_t ){ if( position < v->size1() ) ++position; }
399 else { if( position >= 0 ) --position; }
400 }
404 iterator_base(){ v = 0; }
410 iterator_base( container* v, size_t position )
411 : v( v ), position( position ){}
415 container* v;
419 size_t position;
420 };
421 // Need to know about const_iterator_t
422 template<bool reverse_t> class const_iterator_t;
426 template<bool reverse_t> class iterator_t : public iterator_base<matrix,vector,reverse_t>{
427 public:
428 // // Refines output iterator
429 // // operator=
437 return *this;
438 }
439 // // Refines forward iterator
440 // // operator++ (both)
447 return *this;
448 }
454 // return value
457 return result;
458 }
459 // // Refines bidirectional iterator
460 // // operator-- (both)
467 return *this;
468 }
474 // return value
477 return result;
478 }
485 return this->v == i.v and this->position == i.position;
486 }
493 return not this->operator==( i );
494 }
499 protected:
500 friend class matrix;
501 // We need a constructor for matrix
508 : iterator_base<matrix,vector,reverse_t>( v, position ){}
509 };
513 template<bool reverse_t> class const_iterator_t
514 : public iterator_base<matrix const,vector const,reverse_t>{
515 public:
516 // // Refines output iterator
517 // // operator=
525 return *this;
526 }
527 // // Refines forward iterator
528 // // operator++ (both)
535 return *this;
536 }
542 // return value
545 return result;
546 }
547 // // Refines bidirectional iterator
548 // // operator-- (both)
555 return *this;
556 }
562 // return value
565 return result;
566 }
578 }
584 bool operator==( iterator_t<reverse_t> const& i ) const {
585 return this->v == i.v and this->position == i.position;
586 }
592 bool operator!=( iterator_t<reverse_t> const& i ) const {
593 return not this->operator==( i );
594 }
601 return this->v == i.v and this->position == i.position;
602 }
609 return not this->operator==( i );
610 }
611 protected:
612 // We need a constructor for matrix
613 friend class matrix;
620 : iterator_base<matrix const,vector const,reverse_t>( v, position ){}
621 };
622 public:
642 typedef size_t size_type;
643 // begin()
649 return iterator( this, 0 );
650 }
656 return const_iterator( this, 0 );
657 }
658 // end()
664 if( ccgsl_pointer == 0 ) return iterator( this, 0 );
665 return iterator( this, size1() );
666 }
672 if( ccgsl_pointer == 0 ) return const_iterator( this, 0 );
673 return const_iterator( this, size1() );
674 }
675 // rbegin()
681 if( ccgsl_pointer ==0 ) return reverse_iterator( this, 0 );
682 return reverse_iterator( this, size1() - 1 );
683 }
689 if( ccgsl_pointer == 0 ) return const_reverse_iterator( this, 0 );
690 return const_reverse_iterator( this, size1() - 1 );
691 }
692 // rend()
698 return reverse_iterator( this, -1 );
699 }
705 return const_reverse_iterator( this, -1 );
706 }
707 public:
708 // Sizes
713 size_t size1() const { return ccgsl_pointer == 0 ? 0 : ccgsl_pointer->size1; }
718 size_t size2() const { return ccgsl_pointer == 0 ? 0 : ccgsl_pointer->size2; }
726 double* data() {
727 if( ccgsl_pointer == 0 ) gsl_error( "null vector", __FILE__, __LINE__, GSL_EFAULT );
728#ifndef GSL_RANGE_CHECK_OFF
729 if( ccgsl_pointer->size2 != ccgsl_pointer->tda )
730 gsl_error( "matrix size2 and tda do not match", __FILE__, __LINE__, GSL_EBADLEN );
731#endif
732 return ccgsl_pointer->data; }
740 double const* data() const {
741 if( ccgsl_pointer == 0 ) gsl_error( "null vector", __FILE__, __LINE__, GSL_EFAULT );
742#ifndef GSL_RANGE_CHECK_OFF
743 if( ccgsl_pointer->size2 != ccgsl_pointer->tda )
744 gsl_error( "matrix size2 and tda do not match", __FILE__, __LINE__, GSL_EBADLEN );
745#endif
746 return ccgsl_pointer->data; }
752 void swap( matrix& m ){
753 std::swap( ccgsl_pointer, m.ccgsl_pointer );
754 std::swap( count, m.count );
755 }
762 void tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix const& src){
763 gsl_matrix_tricpy( Uplo, Diag, get(), src.get() );
764 }
771 void transpose_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix const& src){
772 gsl_matrix_transpose_tricpy( Uplo, Diag, get(), src.get() );
773 }
774 // view operations
783 matrix submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ){
784 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
785 *m = gsl_matrix_submatrix( get(), i, j, n1, n2 ).matrix;
786 return matrix( m );
787 }
793 vector row( size_t const i ){
794 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
795 *w = gsl_matrix_row( get(), i ).vector;
796 return vector( w );
797 }
803 vector column( size_t const j ){
804 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
805 *w = gsl_matrix_column( get(), j ).vector;
806 return vector( w );
807 }
812 vector
813 diagonal(){ gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
814 *w = gsl_matrix_diagonal( get() ).vector;
815 return vector( w );
816 }
822 vector subdiagonal( size_t const k ){
823 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
824 *w = gsl_matrix_subdiagonal( get(), k ).vector;
825 return vector( w );
826 }
832 vector superdiagonal( size_t const k ){
833 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
834 *w = gsl_matrix_superdiagonal( get(), k ).vector;
835 return vector( w );
836 }
844 vector subrow( size_t const i, size_t const offset, size_t const n ){
845 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
846 *w = gsl_matrix_subrow( get(), i, offset, n ).vector;
847 return vector( w );
848 }
856 vector subcolumn( size_t const j, size_t const offset, size_t const n ){
857 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
858 *w = gsl_matrix_subcolumn( get(), j, offset, n ).vector;
859 return vector( w );
860 }
868 static matrix view_array( double* base, size_t const n1, size_t const n2 ){
869 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
870 *m = gsl_matrix_view_array( base, n1, n2 ).matrix;
871 return matrix( m );
872 }
881 static matrix view_array_with_tda( double* base, size_t const n1, size_t const n2, size_t const tda ){
882 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
883 *m = gsl_matrix_view_array_with_tda( base, n1, n2, tda ).matrix;
884 return matrix( m );
885 }
893 static matrix view_vector( vector& v, size_t const n1, size_t const n2 ){
894 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
895 *m = gsl_matrix_view_vector( v.get(), n1, n2 ).matrix;
896 return matrix( m );
897 }
906 static matrix view_vector_with_tda( vector& v, size_t const n1, size_t const n2, size_t const tda ){
907 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
908 *m = gsl_matrix_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
909 return matrix( m );
910 }
911 // const versions ...
920 matrix const const_submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ) const {
921 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
922 *m = gsl_matrix_const_submatrix( get(), i, j, n1, n2 ).matrix;
923 return matrix( m );
924 }
930 vector const const_row( size_t const i ) const {
931 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
932 *w = gsl_matrix_const_row( get(), i ).vector;
933 return vector( w );
934 }
940 vector const const_column( size_t const j ) const {
941 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
942 *w = gsl_matrix_const_column( get(), j ).vector;
943 return vector( w );
944 }
949 vector const const_diagonal() const {
950 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
951 *w = gsl_matrix_const_diagonal( get() ).vector;
952 return vector( w );
953 }
959 vector const const_subdiagonal( size_t const k ) const {
960 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
961 *w = gsl_matrix_const_subdiagonal( get(), k ).vector;
962 return vector( w );
963 }
969 vector const const_superdiagonal( size_t const k ) const {
970 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
971 *w = gsl_matrix_const_superdiagonal( get(), k ).vector;
972 return vector( w );
973 }
981 vector const const_subrow( size_t const i, size_t const offset, size_t const n ) const {
982 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
983 *w = gsl_matrix_const_subrow( get(), i, offset, n ).vector;
984 return vector( w );
985 }
993 vector const const_subcolumn( size_t const j, size_t const offset, size_t const n ) const {
994 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
995 *w = gsl_matrix_const_subcolumn( get(), j, offset, n ).vector;
996 return vector( w );
997 }
1006 matrix const submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ) const {
1007 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1008 *m = gsl_matrix_const_submatrix( get(), i, j, n1, n2 ).matrix;
1009 return matrix( m );
1010 }
1016 vector const row( size_t const i ) const {
1017 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
1018 *w = gsl_matrix_const_row( get(), i ).vector;
1019 return vector( w );
1020 }
1026 vector const column( size_t const j ) const {
1027 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
1028 *w = gsl_matrix_const_column( get(), j ).vector;
1029 return vector( w );
1030 }
1035 vector const diagonal() const {
1036 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
1037 *w = gsl_matrix_const_diagonal( get() ).vector;
1038 return vector( w );
1039 }
1045 vector const subdiagonal( size_t const k ) const {
1046 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
1047 *w = gsl_matrix_const_subdiagonal( get(), k ).vector;
1048 return vector( w );
1049 }
1055 vector const superdiagonal( size_t const k ) const {
1056 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
1057 *w = gsl_matrix_const_superdiagonal( get(), k ).vector;
1058 return vector( w );
1059 }
1067 vector const subrow( size_t const i, size_t const offset, size_t const n ) const {
1068 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
1069 *w = gsl_matrix_const_subrow( get(), i, offset, n ).vector;
1070 return vector( w );
1071 }
1079 vector const subcolumn( size_t const j, size_t const offset, size_t const n ) const {
1080 gsl_vector* w = static_cast<gsl_vector*>( malloc( sizeof( gsl_vector ) ) );
1081 *w = gsl_matrix_const_subcolumn( get(), j, offset, n ).vector;
1082 return vector( w );
1083 }
1091 static matrix const const_view_array( double const* base, size_t const n1, size_t const n2 ){
1092 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1093 *m = gsl_matrix_const_view_array( base, n1, n2 ).matrix;
1094 return matrix( m );
1095 }
1104 static matrix const
1105 const_view_array_with_tda( double const* base, size_t const n1, size_t const n2, size_t const tda ){
1106 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1107 *m = gsl_matrix_const_view_array_with_tda( base, n1, n2, tda ).matrix;
1108 return matrix( m );
1109 }
1117 static matrix const const_view_vector( vector const& v, size_t const n1, size_t const n2 ){
1118 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1119 *m = gsl_matrix_const_view_vector( v.get(), n1, n2 ).matrix;
1120 return matrix( m );
1121 }
1130 static matrix const
1131 const_view_vector_with_tda( vector const& v, size_t const n1, size_t const n2, size_t const tda ){
1132 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1133 *m = gsl_matrix_const_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
1134 return matrix( m );
1135 }
1143 static matrix const view_array( double const* base, size_t const n1, size_t const n2 ){
1144 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1145 *m = gsl_matrix_const_view_array( base, n1, n2 ).matrix;
1146 return matrix( m );
1147 }
1156 static matrix const
1157 view_array_with_tda( double const* base, size_t const n1, size_t const n2, size_t const tda ){
1158 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1159 *m = gsl_matrix_const_view_array_with_tda( base, n1, n2, tda ).matrix;
1160 return matrix( m );
1161 }
1169 static matrix const view_vector( vector const& v, size_t const n1, size_t const n2 ){
1170 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1171 *m = gsl_matrix_const_view_vector( v.get(), n1, n2 ).matrix;
1172 return matrix( m );
1173 }
1182 static matrix const
1183 view_vector_with_tda( vector const& v, size_t const n1, size_t const n2, size_t const tda ){
1184 gsl_matrix* m = static_cast<gsl_matrix*>( malloc( sizeof( gsl_matrix ) ) );
1185 *m = gsl_matrix_const_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
1186 return matrix( m );
1187 }
1188 private:
1196 gsl_matrix* ccgsl_pointer;
1200 size_t* count;
1201 public:
1202 // shared reference functions
1207 gsl_matrix* get(){ return ccgsl_pointer; }
1212 gsl_matrix const* get() const { return ccgsl_pointer; }
1218 bool unique() const { return count != 0 and *count == 1; }
1223 size_t use_count() const { return count == 0 ? 0 : *count; }
1229#ifdef __GXX_EXPERIMENTAL_CXX0X__
1230 explicit
1231#endif
1232 operator bool() const { return ccgsl_pointer != 0; }
1233 // GSL functions
1241 static matrix calloc( size_t const n1, size_t const n2 ){ return matrix( gsl_matrix_calloc( n1, n2 ) ); }
1245 void set_zero(){ gsl_matrix_set_zero( get() ); }
1250 void set_all( double x ){ gsl_matrix_set_all( get(), x ); }
1256 int memcpy( matrix const& src ){ return gsl_matrix_memcpy( get(), src.get() ); }
1261 double max() const { return gsl_matrix_max( get() ); }
1266 double min() const { return gsl_matrix_min( get() ); }
1272 void minmax( double* min_out, double* max_out ) const {
1273 gsl_matrix_minmax( get(), min_out, max_out ); }
1279 void minmax( double& min_out, double& max_out ) const {
1280 gsl_matrix_minmax( get(), &min_out, &max_out ); }
1286 int add( matrix const& b ){ return gsl_matrix_add( get(), b.get() ); }
1292 int sub( matrix const& b ){ return gsl_matrix_sub( get(), b.get() ); }
1298 int scale( double const x ){ return gsl_matrix_scale( get(), x ); }
1304 int add_constant( double const x ){ return gsl_matrix_add_constant( get(), x ); }
1309 int isnull() const { return gsl_matrix_isnull( get() ); }
1314 int ispos() const { return gsl_matrix_ispos( get() ); }
1319 int isneg() const { return gsl_matrix_isneg( get() ); }
1324 int isnonneg() const { return gsl_matrix_isnonneg( get() ); }
1331 double get( size_t const i, size_t const j ) const { return gsl_matrix_get( get(), i, j ); }
1338 void set( size_t const i, size_t const j, double x ){ gsl_matrix_set( get(), i, j, x ); }
1345 double* ptr( size_t const i, size_t const j ){ return gsl_matrix_ptr( get(), i, j ); }
1352 double const* const_ptr( size_t const i, size_t const j ) const {
1353 return gsl_matrix_const_ptr( get(), i, j ); }
1359 int fread( FILE* stream ){ return gsl_matrix_fread( stream, get() ); }
1365 int fwrite( FILE* stream ) const { return gsl_matrix_fwrite( stream, get() ); }
1371 int fscanf( FILE* stream ){ return gsl_matrix_fscanf( stream, get() ); }
1378 int fprintf( FILE* stream, char const* format ) const {
1379 return gsl_matrix_fprintf( stream, get(), format ); }
1388 matrix( block& b, size_t const offset, size_t const n1, size_t const n2, size_t const d2 ){
1389 ccgsl_pointer = gsl_matrix_alloc_from_block( b.get(), offset, n1, n2, d2 );
1390 // just plausibly we could allocate vector but not count
1391 try { count = new size_t; } catch( std::bad_alloc& e ){
1392 // try to tidy up before rethrowing
1393 gsl_matrix_free( ccgsl_pointer );
1394 throw e;
1395 }
1396 *count = 1; // initially there is just one reference to ccgsl_pointer
1397 }
1406 matrix( matrix& m, size_t const k1, size_t const k2, size_t const n1, size_t const n2 ){
1407 ccgsl_pointer = gsl_matrix_alloc_from_matrix( m.get(), k1, k2, n1, n2 );
1408 // just plausibly we could allocate matrix but not count
1409 try { count = new size_t; } catch( std::bad_alloc& e ){
1410 // try to tidy up before rethrowing
1411 gsl_matrix_free( ccgsl_pointer );
1412 throw e;
1413 }
1414 *count = 1; // initially there is just one reference to ccgsl_pointer
1415 }
1416 // More functions
1420 void set_identity(){ gsl_matrix_set_identity( get() ); }
1427 int swap_rows( size_t const i, size_t const j ){ return gsl_matrix_swap_rows( get(), i, j ); }
1434 int swap_columns( size_t const i, size_t const j ){
1435 return gsl_matrix_swap_columns( get(), i, j ); }
1442 int swap_rowcol( size_t const i, size_t const j ){ return gsl_matrix_swap_rowcol( get(), i, j ); }
1447 int transpose(){ return gsl_matrix_transpose( get() ); }
1453 int transpose_memcpy( matrix const& src ){
1454 return gsl_matrix_transpose_memcpy( get(), src.get() ); }
1455#ifndef DOXYGEN_SKIP
1461 void max_index( size_t* imax, size_t* jmax ) const {
1462 gsl_matrix_max_index( get(), imax, jmax ); }
1468 void min_index( size_t* imin, size_t* jmin ) const {
1469 gsl_matrix_min_index( get(), imin, jmin ); }
1477 void minmax_index( size_t* imin, size_t* jmin, size_t* imax, size_t* jmax ) const {
1478 gsl_matrix_minmax_index( get(), imin, jmin, imax, jmax ); }
1479#endif
1485 void max_index( size_t& imax, size_t& jmax ) const {
1486 gsl_matrix_max_index( get(), &imax, &jmax ); }
1492 void min_index( size_t& imin, size_t& jmin ) const {
1493 gsl_matrix_min_index( get(), &imin, &jmin ); }
1501 void minmax_index( size_t& imin, size_t& jmin, size_t& imax, size_t& jmax ) const {
1502 gsl_matrix_minmax_index( get(), &imin, &jmin, &imax, &jmax ); }
1508 int
1510 return gsl_matrix_mul_elements( get(), b.get() ); }
1517 int div_elements( matrix const& b ){
1518 return gsl_matrix_div_elements( get(), b.get() ); }
1523 double norm1() const {
1524 return gsl_matrix_norm1( get() ); }
1530 int scale_rows( vector const& x ){
1531 return gsl_matrix_scale_rows( get(), x.get() ); }
1537 int scale_columns( vector const& x ){
1538 return gsl_matrix_scale_columns( get(), x.get() ); }
1544 int add_diagonal( double const x ){
1545 return gsl_matrix_add_diagonal( get(), x ); }
1552 int get_row( vector& v, size_t const i ) const {
1553 return gsl_matrix_get_row( v.get(), get(), i ); }
1560 int get_col( vector& v, size_t const j ) const {
1561 return gsl_matrix_get_col( v.get(), get(), j ); }
1568 int set_row( size_t const i, vector const& v ){
1569 return gsl_matrix_set_row( get(), i, v.get() ); }
1576 int set_col( size_t const j, vector const& v ){
1577 return gsl_matrix_set_col( get(), j, v.get() ); }
1578 // Extra functions for []
1586 vector operator[]( size_t const i ){
1587 // First check that iterator is initialised.
1588 if( ccgsl_pointer == 0 ){
1589 gsl_error( "matrix is null", __FILE__, __LINE__, exception::GSL_EFAULT );
1590 return vector();
1591 }
1592#ifndef GSL_RANGE_CHECK_OFF
1593#ifndef __GXX_EXPERIMENTAL_CXX0X__
1594 static vector something;
1595#endif
1596 //Check that position make sense
1597 if( i >= size1() ){
1598 gsl_error( "trying to read beyond last row of matrix",
1599 __FILE__, __LINE__, exception::GSL_EINVAL );
1600#ifdef __GXX_EXPERIMENTAL_CXX0X__
1601 return vector();
1602#else
1603 return something;
1604#endif
1605 }
1606 // n is a valid position
1607#endif
1608 vector v(0);
1609 gsl_vector_view w = gsl_matrix_row( ccgsl_pointer, i );
1610 v.wrap_gsl_vector_without_ownership(&(w.vector));
1611 return v;
1612 }
1620 vector const operator[]( size_t const i ) const {
1621 // First check that iterator is initialised.
1622 if( ccgsl_pointer == 0 ){
1623 gsl_error( "matrix is null", __FILE__, __LINE__, exception::GSL_EFAULT );
1624 return vector();
1625 }
1626 vector v(0);
1627 gsl_vector_view w = gsl_matrix_row( ccgsl_pointer, i );
1628 v.wrap_gsl_vector_without_ownership(&(w.vector));
1629 return v;
1630 }
1636 inline int
1638 return gsl_permute_matrix( p.get(), get() ); }
1639 };
1640
1641 // Extra functions for vector allocation from matrix objects
1643 inline vector vector::alloc_row_from_matrix( matrix& m, size_t const i ){
1644 return vector ( gsl_vector_alloc_row_from_matrix( m.get(), i ) ); }
1645 inline vector vector::alloc_col_from_matrix( matrix& m, size_t const i ){
1646 return vector ( gsl_vector_alloc_col_from_matrix( m.get(), i ) ); }
1648}
1649#endif
This class handles vectors as shared handles.
Definition: block.hpp:42
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_EFAILED
generic failure
Definition: exception.hpp:476
@ GSL_EINVAL
invalid argument supplied by user
Definition: exception.hpp:475
@ GSL_EBADLEN
matrix, vector lengths are not conformant
Definition: exception.hpp:490
@ GSL_EFAULT
invalid pointer
Definition: exception.hpp:474
Class that extends gsl_function_fdf so that it can be constructed from arbitrary function objects.
A class template for the const iterators.
Definition: matrix.hpp:514
bool operator==(iterator_t< reverse_t > const &i) const
Comparison with non-const iterator.
Definition: matrix.hpp:584
bool operator!=(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
Definition: matrix.hpp:608
const_iterator_t()
The default constructor.
Definition: matrix.hpp:570
bool operator!=(iterator_t< reverse_t > const &i) const
Comparison with non-const iterator.
Definition: matrix.hpp:592
const_iterator_t< reverse_t > & operator--()
The prefix – operator.
Definition: matrix.hpp:553
const_iterator_t< reverse_t > operator--(int)
The postfix – operator.
Definition: matrix.hpp:561
const_iterator_t(iterator_t< reverse_t > const &i)
A copy constructor.
Definition: matrix.hpp:575
const_iterator_t< reverse_t > & operator++()
The prefix ++ operator.
Definition: matrix.hpp:533
const_iterator_t< reverse_t > & operator=(const_iterator_t< reverse_t > const &i)
We can assign one output iterator from another.
Definition: matrix.hpp:522
const_iterator_t(matrix const *v, size_t position)
This constructor allows vector to create non-default iterators.
Definition: matrix.hpp:619
bool operator==(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
Definition: matrix.hpp:600
const_iterator_t< reverse_t > operator++(int)
The postfix ++ operator.
Definition: matrix.hpp:541
We create a suitable class for iterator types here.
Definition: matrix.hpp:276
void decrement()
Decrement the iterator.
Definition: matrix.hpp:388
reference operator*() const
Dereference the pointer.
Definition: matrix.hpp:301
iterator_base()
The iterator is default constructible.
Definition: matrix.hpp:404
bool operator==(iterator_base< container, content, reverse_t > const &i) const
The == operator.
Definition: matrix.hpp:355
vector_ptr pointer
An iterator must have a pointer typea.
Definition: matrix.hpp:290
bool operator!=(iterator_base< container, content, reverse_t > const &i) const
The != operator.
Definition: matrix.hpp:363
void increment()
Increment the iterator.
Definition: matrix.hpp:371
vector value_type
An iterator must have a value type.
Definition: matrix.hpp:286
container * v
Store a pointer to a matrix we can iterate over: 0 if no matrix.
Definition: matrix.hpp:415
std::bidirectional_iterator_tag iterator_category
An iterator must have an iterator category.
Definition: matrix.hpp:282
value_type reference
An iterator must have a reference type.
Definition: matrix.hpp:294
pointer operator->() const
Dereference the pointer.
Definition: matrix.hpp:325
iterator_base(container *v, size_t position)
This constructor allows vector to create non-default iterators.
Definition: matrix.hpp:410
size_t position
Mark position of iterator within matrix.
Definition: matrix.hpp:419
A class template for the two non-const iterators.
Definition: matrix.hpp:426
bool operator!=(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
Definition: matrix.hpp:492
iterator_t()
The default constructor.
Definition: matrix.hpp:498
iterator_t< reverse_t > & operator++()
The prefix ++ operator.
Definition: matrix.hpp:445
iterator_t< reverse_t > operator++(int)
The postfix ++ operator.
Definition: matrix.hpp:453
iterator_t< reverse_t > & operator--()
The prefix – operator.
Definition: matrix.hpp:465
iterator_t< reverse_t > operator--(int)
The postfix – operator.
Definition: matrix.hpp:473
bool operator==(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
Definition: matrix.hpp:484
iterator_t(matrix *v, size_t position)
This constructor allows vector to create non-default iterators.
Definition: matrix.hpp:507
iterator_t< reverse_t > & operator=(iterator_t< reverse_t > const &i)
We can assign one output iterator from another.
Definition: matrix.hpp:434
This class handles matrix objects as shared handles.
Definition: matrix.hpp:72
int transpose_memcpy(matrix const &src)
C++ version of gsl_matrix_transpose_memcpy().
Definition: matrix.hpp:1453
double * ptr(size_t const i, size_t const j)
C++ version of gsl_matrix_ptr().
Definition: matrix.hpp:1345
int sub(matrix const &b)
C++ version of gsl_matrix_sub().
Definition: matrix.hpp:1292
vector superdiagonal(size_t const k)
C++ version of gsl_matrix_superdiagonal().
Definition: matrix.hpp:832
int add(matrix const &b)
C++ version of gsl_matrix_add().
Definition: matrix.hpp:1286
vector subdiagonal(size_t const k)
C++ version of gsl_matrix_subdiagonal().
Definition: matrix.hpp:822
void reset()
Stop sharing ownership of the shared pointer.
Definition: matrix.hpp:230
static matrix const const_view_array(double const *base, size_t const n1, size_t const n2)
C++ version of gsl_matrix_const_view_array().
Definition: matrix.hpp:1091
~matrix()
The destructor only deletes the pointers if count reaches zero.
Definition: matrix.hpp:194
matrix(matrix &&v)
Move constructor.
Definition: matrix.hpp:236
iterator_t< true > reverse_iterator
The reverse_iterator type.
Definition: matrix.hpp:638
int isneg() const
C++ version of gsl_matrix_isneg().
Definition: matrix.hpp:1319
static matrix const const_view_array_with_tda(double const *base, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_const_view_array_with_tda().
Definition: matrix.hpp:1105
int add_diagonal(double const x)
C++ version of gsl_matrix_add_diagonal().
Definition: matrix.hpp:1544
const_iterator_t< false > const_iterator
The const_iterator type.
Definition: matrix.hpp:626
iterator end()
Get iterator pointing beyond last vector element.
Definition: matrix.hpp:663
const_iterator begin() const
Get iterator pointing to first vector element.
Definition: matrix.hpp:655
matrix & operator=(matrix &&v)
Move operator.
Definition: matrix.hpp:246
reverse_iterator rend()
Get iterator pointing beyond last vector element.
Definition: matrix.hpp:697
int transpose()
C++ version of gsl_matrix_transpose().
Definition: matrix.hpp:1447
static matrix view_array_with_tda(double *base, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_view_array_with_tda().
Definition: matrix.hpp:881
const_iterator end() const
Get iterator pointing beyond last vector element.
Definition: matrix.hpp:671
static matrix view_array(double *base, size_t const n1, size_t const n2)
C++ version of gsl_matrix_view_array().
Definition: matrix.hpp:868
double get(size_t const i, size_t const j) const
C++ version of gsl_matrix_get().
Definition: matrix.hpp:1331
double norm1() const
C++ version of gsl_matrix_norm1().
Definition: matrix.hpp:1523
vector const diagonal() const
Another C++ version of gsl_matrix_const_diagonal().
Definition: matrix.hpp:1035
int set_col(size_t const j, vector const &v)
C++ version of gsl_matrix_set_col().
Definition: matrix.hpp:1576
vector subrow(size_t const i, size_t const offset, size_t const n)
C++ version of gsl_matrix_subrow().
Definition: matrix.hpp:844
double * data()
Give access to the data block.
Definition: matrix.hpp:726
matrix(gsl_matrix *v)
Could construct from a gsl_matrix.
Definition: matrix.hpp:110
void minmax_index(size_t &imin, size_t &jmin, size_t &imax, size_t &jmax) const
C++ version of gsl_matrix_minmax_index().
Definition: matrix.hpp:1501
matrix(matrix const &v)
The copy constructor.
Definition: matrix.hpp:154
static matrix const view_array(double const *base, size_t const n1, size_t const n2)
Another C++ version of gsl_matrix_const_view_array().
Definition: matrix.hpp:1143
int swap_rowcol(size_t const i, size_t const j)
C++ version of gsl_matrix_swap_rowcol().
Definition: matrix.hpp:1442
void swap(matrix &m)
Swap two matrix objects.
Definition: matrix.hpp:752
vector const column(size_t const j) const
Another C++ version of gsl_matrix_const_column().
Definition: matrix.hpp:1026
int memcpy(matrix const &src)
C++ version of gsl_matrix_memcpy().
Definition: matrix.hpp:1256
static matrix const view_vector(vector const &v, size_t const n1, size_t const n2)
Another C++ version of gsl_matrix_const_view_vector().
Definition: matrix.hpp:1169
int fscanf(FILE *stream)
C++ version of gsl_matrix_fscanf().
Definition: matrix.hpp:1371
const_iterator_t< true > const_reverse_iterator
The const_reverse_t type.
Definition: matrix.hpp:634
void set_all(double x)
C++ version of gsl_matrix_set_all().
Definition: matrix.hpp:1250
int get_col(vector &v, size_t const j) const
C++ version of gsl_matrix_get_col().
Definition: matrix.hpp:1560
vector const const_subdiagonal(size_t const k) const
C++ version of gsl_matrix_const_subdiagonal().
Definition: matrix.hpp:959
static matrix view_vector(vector &v, size_t const n1, size_t const n2)
C++ version of gsl_matrix_view_vector().
Definition: matrix.hpp:893
size_t use_count() const
Find how many matrix objects share this pointer.
Definition: matrix.hpp:1223
void wrap_gsl_matrix_without_ownership(gsl_matrix *v)
This function is intended mainly for internal use.
Definition: matrix.hpp:214
int fwrite(FILE *stream) const
C++ version of gsl_matrix_fwrite().
Definition: matrix.hpp:1365
vector const row(size_t const i) const
Another C++ version of gsl_matrix_const_row().
Definition: matrix.hpp:1016
double max() const
C++ version of gsl_matrix_max().
Definition: matrix.hpp:1261
int scale(double const x)
C++ version of gsl_matrix_scale().
Definition: matrix.hpp:1298
int permute(permutation &p)
Permute the columns of this by permutation p.
Definition: matrix.hpp:1637
iterator begin()
Get iterator pointing to first vector element.
Definition: matrix.hpp:648
vector const superdiagonal(size_t const k) const
Another C++ version of gsl_matrix_const_superdiagonal().
Definition: matrix.hpp:1055
int set_row(size_t const i, vector const &v)
C++ version of gsl_matrix_set_row().
Definition: matrix.hpp:1568
static matrix const view_vector_with_tda(vector const &v, size_t const n1, size_t const n2, size_t const tda)
Another C++ version of gsl_matrix_const_view_vector_with_tda().
Definition: matrix.hpp:1183
int swap_columns(size_t const i, size_t const j)
C++ version of gsl_matrix_swap_columns().
Definition: matrix.hpp:1434
vector const operator[](size_t const i) const
This function allows us to use a matrix like an array.
Definition: matrix.hpp:1620
vector operator[](size_t const i)
This function allows us to use a matrix like an array.
Definition: matrix.hpp:1586
static matrix const const_view_vector(vector const &v, size_t const n1, size_t const n2)
C++ version of gsl_matrix_const_view_vector().
Definition: matrix.hpp:1117
int ispos() const
C++ version of gsl_matrix_ispos().
Definition: matrix.hpp:1314
matrix submatrix(size_t const i, size_t const j, size_t const n1, size_t const n2)
C++ version of gsl_matrix_submatrix().
Definition: matrix.hpp:783
vector subcolumn(size_t const j, size_t const offset, size_t const n)
C++ version of gsl_matrix_subcolumn().
Definition: matrix.hpp:856
vector const const_row(size_t const i) const
C++ version of gsl_matrix_const_row().
Definition: matrix.hpp:930
vector const subrow(size_t const i, size_t const offset, size_t const n) const
Another C++ version of gsl_matrix_const_subrow().
Definition: matrix.hpp:1067
vector row(size_t const i)
C++ version of gsl_matrix_row().
Definition: matrix.hpp:793
void set_identity()
C++ version of gsl_matrix_set_identity().
Definition: matrix.hpp:1420
gsl_matrix * ccgsl_pointer
The shared pointer.
Definition: matrix.hpp:1196
static matrix view_vector_with_tda(vector &v, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_view_vector_with_tda().
Definition: matrix.hpp:906
static matrix calloc(size_t const n1, size_t const n2)
C++ version of gsl_matrix_calloc().
Definition: matrix.hpp:1241
vector const subcolumn(size_t const j, size_t const offset, size_t const n) const
Another C++ version of gsl_matrix_const_subcolumn().
Definition: matrix.hpp:1079
const_reverse_iterator rend() const
Get iterator pointing beyond last vector element.
Definition: matrix.hpp:704
matrix(size_t const n1, size_t const n2)
This constructor creates a new matrix with n1 rows and n2 columns.
Definition: matrix.hpp:92
int div_elements(matrix const &b)
C++ version of gsl_matrix_div_elements().
Definition: matrix.hpp:1517
int scale_rows(vector const &x)
C++ version of gsl_matrix_scale_rows().
Definition: matrix.hpp:1530
double const * const_ptr(size_t const i, size_t const j) const
C++ version of gsl_matrix_const_ptr().
Definition: matrix.hpp:1352
static matrix const const_view_vector_with_tda(vector const &v, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_const_view_vector_with_tda().
Definition: matrix.hpp:1131
void tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix const &src)
Copy the upper or lower triangular part of matrix src to this.
Definition: matrix.hpp:762
void transpose_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix const &src)
Copy the upper or lower triangular part of matrix src to this.
Definition: matrix.hpp:771
matrix const submatrix(size_t const i, size_t const j, size_t const n1, size_t const n2) const
Another C++ version of gsl_matrix_const_submatrix().
Definition: matrix.hpp:1006
iterator_t< false > iterator
The iterator type.
Definition: matrix.hpp:630
size_t * count
The shared reference count.
Definition: matrix.hpp:1200
void minmax(double *min_out, double *max_out) const
C++ version of gsl_matrix_minmax().
Definition: matrix.hpp:1272
matrix & operator=(matrix const &v)
The assignment operator.
Definition: matrix.hpp:163
matrix(block &b, size_t const offset, size_t const n1, size_t const n2, size_t const d2)
C++ version of gsl_matrix_alloc_from_block().
Definition: matrix.hpp:1388
int fprintf(FILE *stream, char const *format) const
C++ version of gsl_matrix_fprintf().
Definition: matrix.hpp:1378
double min() const
C++ version of gsl_matrix_min().
Definition: matrix.hpp:1266
vector const subdiagonal(size_t const k) const
Another C++ version of gsl_matrix_const_subdiagonal().
Definition: matrix.hpp:1045
vector const const_superdiagonal(size_t const k) const
C++ version of gsl_matrix_const_superdiagonal().
Definition: matrix.hpp:969
bool unique() const
Find if this is the only object sharing the gsl_matrix.
Definition: matrix.hpp:1218
static matrix const view_array_with_tda(double const *base, size_t const n1, size_t const n2, size_t const tda)
Another C++ version of gsl_matrix_const_view_array_with_tda().
Definition: matrix.hpp:1157
reverse_iterator rbegin()
Get iterator pointing to first vector element.
Definition: matrix.hpp:680
int fread(FILE *stream)
C++ version of gsl_matrix_fread().
Definition: matrix.hpp:1359
bool owns_data
Used to allow a vector that does not own its data.
Definition: matrix.hpp:1192
void min_index(size_t &imin, size_t &jmin) const
C++ version of gsl_matrix_min_index().
Definition: matrix.hpp:1492
vector diagonal()
C++ version of gsl_matrix_diagonal().
Definition: matrix.hpp:813
vector const const_subrow(size_t const i, size_t const offset, size_t const n) const
C++ version of gsl_matrix_const_subrow().
Definition: matrix.hpp:981
int add_constant(double const x)
C++ version of gsl_matrix_add_constant().
Definition: matrix.hpp:1304
vector const const_diagonal() const
C++ version of gsl_matrix_const_diagonal().
Definition: matrix.hpp:949
const_reverse_iterator rbegin() const
Get iterator pointing to first vector element.
Definition: matrix.hpp:688
int isnonneg() const
C++ version of gsl_matrix_isnonneg().
Definition: matrix.hpp:1324
size_t size1() const
The number of rows of the matrix.
Definition: matrix.hpp:713
void set(size_t const i, size_t const j, double x)
C++ version of gsl_matrix_set().
Definition: matrix.hpp:1338
int swap_rows(size_t const i, size_t const j)
C++ version of gsl_matrix_swap_rows().
Definition: matrix.hpp:1427
int scale_columns(vector const &x)
C++ version of gsl_matrix_scale_columns().
Definition: matrix.hpp:1537
gsl_matrix * get()
Get the gsl_matrix.
Definition: matrix.hpp:1207
void minmax(double &min_out, double &max_out) const
C++ version of gsl_matrix_minmax().
Definition: matrix.hpp:1279
size_t size_type
A container must have a size_type.
Definition: matrix.hpp:642
int get_row(vector &v, size_t const i) const
C++ version of gsl_matrix_get_row().
Definition: matrix.hpp:1552
matrix(std::initializer_list< std::initializer_list< double > > initializer_list)
Could construct from a std::initializer_list in C++11.
Definition: matrix.hpp:121
int mul_elements(matrix const &b)
C++ version of gsl_matrix_mul_elements().
Definition: matrix.hpp:1509
gsl_matrix const * get() const
Get the gsl_matrix.
Definition: matrix.hpp:1212
size_t size2() const
The number of columns of the matrix.
Definition: matrix.hpp:718
vector const const_column(size_t const j) const
C++ version of gsl_matrix_const_column().
Definition: matrix.hpp:940
void max_index(size_t &imax, size_t &jmax) const
C++ version of gsl_matrix_max_index().
Definition: matrix.hpp:1485
void set_zero()
C++ version of gsl_matrix_set_zero().
Definition: matrix.hpp:1245
int isnull() const
C++ version of gsl_matrix_isnull().
Definition: matrix.hpp:1309
matrix(matrix &m, size_t const k1, size_t const k2, size_t const n1, size_t const n2)
C++ version of gsl_matrix_alloc_from_matrix().
Definition: matrix.hpp:1406
matrix()
The default constructor is only really useful for assigning to.
Definition: matrix.hpp:80
matrix clone() const
The clone function.
Definition: matrix.hpp:183
matrix const const_submatrix(size_t const i, size_t const j, size_t const n1, size_t const n2) const
C++ version of gsl_matrix_const_submatrix().
Definition: matrix.hpp:920
double const * data() const
Give access to the data block.
Definition: matrix.hpp:740
vector const const_subcolumn(size_t const j, size_t const offset, size_t const n) const
C++ version of gsl_matrix_const_subcolumn().
Definition: matrix.hpp:993
vector column(size_t const j)
C++ version of gsl_matrix_column().
Definition: matrix.hpp:803
Class that extends gsl_multifit_nlinear_fdf so that it can be constructed from arbitrary function obj...
Class that extends gsl_multilarge_nlinear_fdf so that it can be constructed from arbitrary function o...
Class that extends gsl_multiroot_function_fdf so that it can be constructed from arbitrary function o...
This class handles GSL permutation objects.
Definition: permutation.hpp:33
gsl_permutation * get()
Get the gsl_permutation.
This class handles vector objects as shared handles.
Definition: vector.hpp:74
gsl_vector * get()
Get the gsl_vector.
Definition: vector.hpp:1320
static vector alloc_row_from_matrix(matrix &m, size_t const i)
C++ version of gsl_vector_alloc_row_from_matrix().
void wrap_gsl_vector_without_ownership(gsl_vector *v)
This function is intended mainly for internal use.
Definition: vector.hpp:272
vector()
The default constructor is only really useful for assigning to.
Definition: vector.hpp:87
static vector alloc_col_from_matrix(matrix &m, size_t const j)
C++ version of gsl_vector_alloc_col_from_matrix().
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
Definition: rstat.hpp:299
double b(int order, double qq)
C++ version of gsl_sf_mathieu_b().
Definition: sf_mathieu.hpp:298
gsl_sf_result result
Typedef for gsl_sf_result.
Definition: sf_result.hpp:30
The gsl package creates an interface to the GNU Scientific Library for C++.
Definition: blas.hpp:34
This is a pointer-like type for iterator return values.
Definition: matrix.hpp:255
vector_ptr(vector const &v)
Typically we have to construct from a vector.
Definition: matrix.hpp:260
vector * operator->()
Dereference operator.
Definition: matrix.hpp:270
vector & operator*()
Dereference operator.
Definition: matrix.hpp:265