ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
matrix_complex.hpp
Go to the documentation of this file.
1/*
2 * $Id: matrix_complex.hpp 313 2014-02-22 14:29:57Z jdl3 $
3 * Copyright (C) 2010, 2011, 2012, 2019, 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_COMPLEX_HPP
21#define CCGSL_MATRIX_COMPLEX_HPP
22
23#include<gsl/gsl_matrix.h>
24#include<new>
25#ifdef __GXX_EXPERIMENTAL_CXX0X__
26#include<complex>
27#endif
28#include<gsl/gsl_permute_matrix_complex_double.h>
29
30#include"exception.hpp"
31#include"vector_complex.hpp"
32#include"permutation.hpp"
33
34// This file is a template
35#define CCGSL_MTY 2
36
37namespace gsl {
63 public:
68 ccgsl_pointer = 0;
69 count = 0; // initially nullptr will do
70 }
71 // Refines random access container
72 // Refines assignable
78 explicit matrix_complex( size_t const n1, size_t const n2 ){
79 if( n1 > 0 and n2 > 0 ) ccgsl_pointer = gsl_matrix_complex_alloc( n1, n2 );
80 else { ccgsl_pointer = new gsl_matrix_complex; ccgsl_pointer->size1 = n1;
81 ccgsl_pointer->size2 = n2; ccgsl_pointer->data = 0; }
82 // just plausibly we could allocate matrix but not count
83 try { count = new size_t; } catch( std::bad_alloc& e ){
84 // try to tidy up before rethrowing
85 if( n1 > 0 and n2 > 0 ) gsl_matrix_complex_free( ccgsl_pointer );
86 else delete ccgsl_pointer;
87 throw e;
88 }
89 *count = 1; // initially there is just one reference to ccgsl_pointer
90 }
96 explicit matrix_complex( gsl_matrix_complex* v ){
97 ccgsl_pointer = v;
98 // just plausibly we could fail to allocate count: no further action needed.
99 count = new size_t;
100 *count = 1; // initially there is just one reference to ccgsl_pointer
101 }
102#ifdef __GXX_EXPERIMENTAL_CXX0X__
107 matrix_complex( std::initializer_list<std::initializer_list<std::complex<double> > > initializer_list ){
108 size_t const n1 = initializer_list.size();
109 size_t const n2 = initializer_list.begin()->size();
110 for( auto l : initializer_list ){
111 if( l.size() != n2 ){
112 gsl::exception e( "matrix rows have unequal sizes", __FILE__, __LINE__,
114 throw( e );
115 }
116 }
117 ccgsl_pointer = gsl_matrix_complex_alloc( n1, n2 );
118 // just plausibly we could allocate matrix but not count
119 try { count = new size_t; } catch( std::bad_alloc& e ){
120 // try to tidy up before rethrowing
121 if( n1 > 0 and n2 > 0 ) gsl_matrix_complex_free( ccgsl_pointer );
122 else delete ccgsl_pointer;
123 throw e;
124 }
125 *count = 1; // initially there is just one reference to ccgsl_pointer
126 // now copy
127 int r = 0;
128 for( auto row : initializer_list ){
129 size_t c = 0;
130 for( auto x : row ){ set( r, c, x ); ++c; }
131 ++r;
132 }
133 }
134#endif
135 // copy constructor
141 if( count != 0 ) ++*count; // matrix_complex is now shared.
142 }
143 // assignment operator
149 // first, possibly delete anything pointed to by this
150 if( count == 0 or --*count == 0 ){
151 if( ccgsl_pointer != 0 and ccgsl_pointer->size1 > 0 and ccgsl_pointer->size2 > 0 ) gsl_matrix_complex_free( ccgsl_pointer );
152 else delete ccgsl_pointer;
153 delete count;
154 }
155 // Then copy
157 count = v.count;
158 if( count != 0 ) ++*count; // block_complex is now shared.
159 return *this;
160 }
161 // clone()
168 matrix_complex copy( get()->size1, get()->size2 );
169 // Now copy
170 gsl_matrix_complex_memcpy( copy.get(), get() );
171 // return new object
172 return copy;
173 }
174 // destructor
179 if( count != 0 and --*count == 0 ){
180 // could have allocated null pointer
181 if( ccgsl_pointer != 0 ){
182 if( ccgsl_pointer->size1 > 0 and ccgsl_pointer->size2 > 0 ) gsl_matrix_complex_free( ccgsl_pointer );
183 else delete ccgsl_pointer; }
184 delete count;
185 }
186 }
190 void reset(){ matrix_complex().swap( *this ); }
191#ifdef __GXX_EXPERIMENTAL_CXX0X__
197 std::swap( count, v.count );
198 v.ccgsl_pointer = nullptr;
199 }
206 matrix_complex( std::move( v ) ).swap( *this );
207 return *this;
208 }
209#endif
210 private:
224 vector_complex& operator*(){ return *this; }
229 vector_complex* operator->(){ return this; }
230 };
235 template<typename container,typename content,bool reverse_t> class iterator_base {
236 friend class vector_complex;
237 public:
241 typedef std::bidirectional_iterator_tag iterator_category;
254 public:
255 // // operator*
261 // Always try to return something
262 static content something;
263 // First check that iterator is initialised.
264 if( v == 0 ){
265 gsl_error( "iterator not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
266 return something;
267 } else if( v->ccgsl_pointer == 0 ){
268 gsl_error( "matrix not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
269 return something;
270 }
271 // Check that position make sense
272 if( position >= v->size1() ){
273 gsl_error( "trying to dereference end()", __FILE__, __LINE__, exception::GSL_EFAILED );
274 return something;
275 }
276 // position and v are valid: return data
277 return v->row( position );
278 }
279 // // operator->
285 // Always try to return something
286 static content something_base;
287 static vector_complex_ptr something( something_base );
288 // First check that iterator is initialised.
289 if( v == 0 ){
290 gsl_error( "iterator not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
291 return something;
292 } else if( v->ccgsl_pointer == 0 ){
293 gsl_error( "matrix not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
294 return something;
295 }
296 // Check that position make sense
297 if( position >= v->size1() ){
298 gsl_error( "trying to dereference end()", __FILE__, __LINE__, exception::GSL_EFAILED );
299 return something;
300 }
301 // position and v are valid: return data
302 vector_complex_ptr ptr( v->row( position ) );
303#ifdef __GXX_EXPERIMENTAL_CXX0X__
304 return std::move( ptr );
305#else
306 return ptr;
307#endif
308 }
315 return this->v == i.v and this->position == i.position;
316 }
323 return not this->operator==( i );
324 }
325 protected:
330 void increment(){
331 // Only makes sense if v points to a matrix
332 if( v == 0 ){
333 gsl_error( "iterator not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
334 return;
335 } else if( v->ccgsl_pointer == 0 ){
336 gsl_error( "matrix not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
337 return;
338 }
339 // increment position and check against size
340 if( reverse_t ){ if( position >= 0 ) --position; }
341 else { if( position < v->size1() ) ++position; }
342 }
347 void decrement(){
348 // Only makes sense if v points to a matrix
349 if( v == 0 ){
350 gsl_error( "iterator not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
351 return;
352 } else if( v->ccgsl_pointer == 0 ){
353 gsl_error( "matrix not initialised", __FILE__, __LINE__, exception::GSL_EFAILED );
354 return;
355 }
356 // decrement position and check against size
357 if( reverse_t ){ if( position < v->size1() ) ++position; }
358 else { if( position >= 0 ) --position; }
359 }
363 iterator_base(){ v = 0; }
369 iterator_base( container* v, size_t position )
370 : v( v ), position( position ){}
374 container* v;
378 size_t position;
379 };
380 // Need to know about const_iterator_t
381 template<bool reverse_t> class const_iterator_t;
385 template<bool reverse_t> class iterator_t : public iterator_base<matrix_complex,vector_complex,reverse_t>{
386 public:
387 // // Refines output iterator
388 // // operator=
396 return *this;
397 }
398 // // Refines forward iterator
399 // // operator++ (both)
406 return *this;
407 }
413 // return value
416 return result;
417 }
418 // // Refines bidirectional iterator
419 // // operator-- (both)
426 return *this;
427 }
433 // return value
436 return result;
437 }
444 return this->v == i.v and this->position == i.position;
445 }
452 return not this->operator==( i );
453 }
458 protected:
459 friend class matrix_complex;
460 // We need a constructor for matrix_complex
468 };
472 template<bool reverse_t> class const_iterator_t
473 : public iterator_base<matrix_complex const,vector_complex const,reverse_t>{
474 public:
475 // // Refines output iterator
476 // // operator=
484 return *this;
485 }
486 // // Refines forward iterator
487 // // operator++ (both)
494 return *this;
495 }
501 // return value
504 return result;
505 }
506 // // Refines bidirectional iterator
507 // // operator-- (both)
514 return *this;
515 }
521 // return value
524 return result;
525 }
537 }
543 bool operator==( iterator_t<reverse_t> const& i ) const {
544 return this->v == i.v and this->position == i.position;
545 }
551 bool operator!=( iterator_t<reverse_t> const& i ) const {
552 return not this->operator==( i );
553 }
560 return this->v == i.v and this->position == i.position;
561 }
568 return not this->operator==( i );
569 }
570 protected:
571 // We need a constructor for matrix_complex
572 friend class matrix_complex;
579 : iterator_base<matrix_complex const,vector_complex const,reverse_t>( v, position ){}
580 };
581 public:
601 typedef size_t size_type;
602 // begin()
608 return iterator( this, 0 );
609 }
615 return const_iterator( this, 0 );
616 }
617 // end()
623 if( ccgsl_pointer == 0 ) return iterator( this, 0 );
624 return iterator( this, size1() );
625 }
631 if( ccgsl_pointer == 0 ) return const_iterator( this, 0 );
632 return const_iterator( this, size1() );
633 }
634 // rbegin()
640 if( ccgsl_pointer ==0 ) return reverse_iterator( this, 0 );
641 return reverse_iterator( this, size1() - 1 );
642 }
648 if( ccgsl_pointer == 0 ) return const_reverse_iterator( this, 0 );
649 return const_reverse_iterator( this, size1() - 1 );
650 }
651 // rend()
657 return reverse_iterator( this, -1 );
658 }
664 return const_reverse_iterator( this, -1 );
665 }
666 public:
667 // Sizes
672 size_t size1() const { return ccgsl_pointer == 0 ? 0 : ccgsl_pointer->size1; }
677 size_t size2() const { return ccgsl_pointer == 0 ? 0 : ccgsl_pointer->size2; }
684 gsl_matrix_complex* tmp = ccgsl_pointer; ccgsl_pointer = m.ccgsl_pointer; m.ccgsl_pointer = tmp;
685 size_t* tmp2 = count; count = m.count; m.count = tmp2;
686 }
693 void tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex const& src){
694 gsl_matrix_complex_tricpy( Uplo, Diag, get(), src.get() );
695 }
702 void transpose_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex const& src){
703 gsl_matrix_complex_transpose_tricpy( Uplo, Diag, get(), src.get() );
704 }
705 // view operations
714 matrix_complex submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ){
715 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
716 *m = gsl_matrix_complex_submatrix( get(), i, j, n1, n2 ).matrix;
717 return matrix_complex( m );
718 }
724 vector_complex row( size_t const i ){
725 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
726 *w = gsl_matrix_complex_row( get(), i ).vector;
727 return vector_complex( w );
728 }
734 vector_complex column( size_t const j ){
735 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
736 *w = gsl_matrix_complex_column( get(), j ).vector;
737 return vector_complex( w );
738 }
744 diagonal(){ gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
745 *w = gsl_matrix_complex_diagonal( get() ).vector;
746 return vector_complex( w );
747 }
753 vector_complex subdiagonal( size_t const k ){
754 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
755 *w = gsl_matrix_complex_subdiagonal( get(), k ).vector;
756 return vector_complex( w );
757 }
763 vector_complex superdiagonal( size_t const k ){
764 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
765 *w = gsl_matrix_complex_superdiagonal( get(), k ).vector;
766 return vector_complex( w );
767 }
775 vector_complex subrow( size_t const i, size_t const offset, size_t const n ){
776 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
777 *w = gsl_matrix_complex_subrow( get(), i, offset, n ).vector;
778 return vector_complex( w );
779 }
787 vector_complex subcolumn( size_t const j, size_t const offset, size_t const n ){
788 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
789 *w = gsl_matrix_complex_subcolumn( get(), j, offset, n ).vector;
790 return vector_complex( w );
791 }
799 static matrix_complex view_array( double* base, size_t const n1, size_t const n2 ){
800 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
801 *m = gsl_matrix_complex_view_array( base, n1, n2 ).matrix;
802 return matrix_complex( m );
803 }
812 static matrix_complex view_array_with_tda( double* base, size_t const n1, size_t const n2, size_t const tda ){
813 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
814 *m = gsl_matrix_complex_view_array_with_tda( base, n1, n2, tda ).matrix;
815 return matrix_complex( m );
816 }
824 static matrix_complex view_vector( vector_complex& v, size_t const n1, size_t const n2 ){
825 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
826 *m = gsl_matrix_complex_view_vector( v.get(), n1, n2 ).matrix;
827 return matrix_complex( m );
828 }
837 static matrix_complex view_vector_with_tda( vector_complex& v, size_t const n1, size_t const n2, size_t const tda ){
838 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
839 *m = gsl_matrix_complex_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
840 return matrix_complex( m );
841 }
842 // const versions ...
851 matrix_complex const const_submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ) const {
852 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
853 *m = gsl_matrix_complex_const_submatrix( get(), i, j, n1, n2 ).matrix;
854 return matrix_complex( m );
855 }
861 vector_complex const const_row( size_t const i ) const {
862 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
863 *w = gsl_matrix_complex_const_row( get(), i ).vector;
864 return vector_complex( w );
865 }
871 vector_complex const const_column( size_t const j ) const {
872 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
873 *w = gsl_matrix_complex_const_column( get(), j ).vector;
874 return vector_complex( w );
875 }
881 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
882 *w = gsl_matrix_complex_const_diagonal( get() ).vector;
883 return vector_complex( w );
884 }
890 vector_complex const const_subdiagonal( size_t const k ) const {
891 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
892 *w = gsl_matrix_complex_const_subdiagonal( get(), k ).vector;
893 return vector_complex( w );
894 }
900 vector_complex const const_superdiagonal( size_t const k ) const {
901 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
902 *w = gsl_matrix_complex_const_superdiagonal( get(), k ).vector;
903 return vector_complex( w );
904 }
912 vector_complex const const_subrow( size_t const i, size_t const offset, size_t const n ) const {
913 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
914 *w = gsl_matrix_complex_const_subrow( get(), i, offset, n ).vector;
915 return vector_complex( w );
916 }
924 vector_complex const const_subcolumn( size_t const j, size_t const offset, size_t const n ) const {
925 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
926 *w = gsl_matrix_complex_const_subcolumn( get(), j, offset, n ).vector;
927 return vector_complex( w );
928 }
936 static matrix_complex const const_view_array( double const* base, size_t const n1, size_t const n2 ){
937 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
938 *m = gsl_matrix_complex_const_view_array( base, n1, n2 ).matrix;
939 return matrix_complex( m );
940 }
949 static matrix_complex const
950 const_view_array_with_tda( double const* base, size_t const n1, size_t const n2, size_t const tda ){
951 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
952 *m = gsl_matrix_complex_const_view_array_with_tda( base, n1, n2, tda ).matrix;
953 return matrix_complex( m );
954 }
962 static matrix_complex const const_view_vector( vector_complex const& v, size_t const n1, size_t const n2 ){
963 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
964 *m = gsl_matrix_complex_const_view_vector( v.get(), n1, n2 ).matrix;
965 return matrix_complex( m );
966 }
975 static matrix_complex const
976 const_view_vector_with_tda( vector_complex const& v, size_t const n1, size_t const n2, size_t const tda ){
977 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
978 *m = gsl_matrix_complex_const_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
979 return matrix_complex( m );
980 }
989 matrix_complex const submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ) const {
990 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
991 *m = gsl_matrix_complex_const_submatrix( get(), i, j, n1, n2 ).matrix;
992 return matrix_complex( m );
993 }
999 vector_complex const row( size_t const i ) const {
1000 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
1001 *w = gsl_matrix_complex_const_row( get(), i ).vector;
1002 return vector_complex( w );
1003 }
1009 vector_complex const column( size_t const j ) const {
1010 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
1011 *w = gsl_matrix_complex_const_column( get(), j ).vector;
1012 return vector_complex( w );
1013 }
1018 vector_complex const diagonal() const {
1019 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
1020 *w = gsl_matrix_complex_const_diagonal( get() ).vector;
1021 return vector_complex( w );
1022 }
1028 vector_complex const subdiagonal( size_t const k ) const {
1029 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
1030 *w = gsl_matrix_complex_const_subdiagonal( get(), k ).vector;
1031 return vector_complex( w );
1032 }
1038 vector_complex const superdiagonal( size_t const k ) const {
1039 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
1040 *w = gsl_matrix_complex_const_superdiagonal( get(), k ).vector;
1041 return vector_complex( w );
1042 }
1050 vector_complex const subrow( size_t const i, size_t const offset, size_t const n ) const {
1051 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
1052 *w = gsl_matrix_complex_const_subrow( get(), i, offset, n ).vector;
1053 return vector_complex( w );
1054 }
1062 vector_complex const subcolumn( size_t const j, size_t const offset, size_t const n ) const {
1063 gsl_vector_complex* w = static_cast<gsl_vector_complex*>( malloc( sizeof( gsl_vector_complex ) ) );
1064 *w = gsl_matrix_complex_const_subcolumn( get(), j, offset, n ).vector;
1065 return vector_complex( w );
1066 }
1074 static matrix_complex const view_array( double const* base, size_t const n1, size_t const n2 ){
1075 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
1076 *m = gsl_matrix_complex_const_view_array( base, n1, n2 ).matrix;
1077 return matrix_complex( m );
1078 }
1087 static matrix_complex const
1088 view_array_with_tda( double const* base, size_t const n1, size_t const n2, size_t const tda ){
1089 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
1090 *m = gsl_matrix_complex_const_view_array_with_tda( base, n1, n2, tda ).matrix;
1091 return matrix_complex( m );
1092 }
1100 static matrix_complex const view_vector( vector_complex const& v, size_t const n1, size_t const n2 ){
1101 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
1102 *m = gsl_matrix_complex_const_view_vector( v.get(), n1, n2 ).matrix;
1103 return matrix_complex( m );
1104 }
1113 static matrix_complex const
1114 view_vector_with_tda( vector_complex const& v, size_t const n1, size_t const n2, size_t const tda ){
1115 gsl_matrix_complex* m = static_cast<gsl_matrix_complex*>( malloc( sizeof( gsl_matrix_complex ) ) );
1116 *m = gsl_matrix_complex_const_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
1117 return matrix_complex( m );
1118 }
1119 private:
1123 gsl_matrix_complex* ccgsl_pointer;
1127 size_t* count;
1128 public:
1129 // shared reference functions
1134 gsl_matrix_complex* get(){ return ccgsl_pointer; }
1139 gsl_matrix_complex const* get() const { return ccgsl_pointer; }
1145 bool unique() const { return count != 0 and *count == 1; }
1150 size_t use_count() const { return count == 0 ? 0 : *count; }
1156#ifdef __GXX_EXPERIMENTAL_CXX0X__
1157 explicit
1158#endif
1159 operator bool() const { return ccgsl_pointer != 0; }
1160 // GSL functions
1168 static matrix_complex calloc( size_t const n1, size_t const n2 ){ return matrix_complex( gsl_matrix_complex_calloc( n1, n2 ) ); }
1172 void set_zero(){ gsl_matrix_complex_set_zero( get() ); }
1177 void set_all( complex x ){ gsl_matrix_complex_set_all( get(), x ); }
1183 int memcpy( matrix_complex const& src ){ return gsl_matrix_complex_memcpy( get(), src.get() ); }
1189 int add( matrix_complex const& b ){ return gsl_matrix_complex_add( get(), b.get() ); }
1195 int sub( matrix_complex const& b ){ return gsl_matrix_complex_sub( get(), b.get() ); }
1201 int scale( complex const x ){ return gsl_matrix_complex_scale( get(), x ); }
1207 int add_constant( complex const x ){ return gsl_matrix_complex_add_constant( get(), x ); }
1212 int isnull() const { return gsl_matrix_complex_isnull( get() ); }
1217 int ispos() const { return gsl_matrix_complex_ispos( get() ); }
1222 int isneg() const { return gsl_matrix_complex_isneg( get() ); }
1227 int isnonneg() const { return gsl_matrix_complex_isnonneg( get() ); }
1234 complex get( size_t const i, size_t const j ) const { return gsl_matrix_complex_get( get(), i, j ); }
1241 void set( size_t const i, size_t const j, complex x ){ gsl_matrix_complex_set( get(), i, j, x ); }
1248 complex_ptr ptr( size_t const i, size_t const j ){
1249 if( i >= ccgsl_pointer->size1 )
1250 gsl_error( "Index out of range", __FILE__, __LINE__, exception::GSL_EINVAL );
1251 if( j >= ccgsl_pointer->size2 )
1252 gsl_error( "Index out of range", __FILE__, __LINE__, exception::GSL_EINVAL );
1253 return complex_ptr( get()->data + CCGSL_MTY * (i * get()->tda + j ) ); }
1260 complex_ptr const const_ptr( size_t const i, size_t const j ) const {
1261 if( i >= ccgsl_pointer->size1 )
1262 gsl_error( "Index out of range", __FILE__, __LINE__, exception::GSL_EINVAL );
1263 if( j >= ccgsl_pointer->size2 )
1264 gsl_error( "Index out of range", __FILE__, __LINE__, exception::GSL_EINVAL );
1265 return complex_ptr( get()->data + CCGSL_MTY * (i * get()->tda + j ) ); }
1271 int fread( FILE* stream ){ return gsl_matrix_complex_fread( stream, get() ); }
1277 int fwrite( FILE* stream ) const { return gsl_matrix_complex_fwrite( stream, get() ); }
1283 int fscanf( FILE* stream ){ return gsl_matrix_complex_fscanf( stream, get() ); }
1290 int fprintf( FILE* stream, char const* format ) const {
1291 return gsl_matrix_complex_fprintf( stream, get(), format ); }
1300 matrix_complex( block_complex& b, size_t const offset, size_t const n1, size_t const n2, size_t const d2 ){
1301 ccgsl_pointer = gsl_matrix_complex_alloc_from_block( b.get(), offset, n1, n2, d2 );
1302 // just plausibly we could allocate vector_complex but not count
1303 try { count = new size_t; } catch( std::bad_alloc& e ){
1304 // try to tidy up before rethrowing
1305 gsl_matrix_complex_free( ccgsl_pointer );
1306 throw e;
1307 }
1308 *count = 1; // initially there is just one reference to ccgsl_pointer
1309 }
1318 matrix_complex( matrix_complex& m, size_t const k1, size_t const k2, size_t const n1, size_t const n2 ){
1319 ccgsl_pointer = gsl_matrix_complex_alloc_from_matrix( m.get(), k1, k2, n1, n2 );
1320 // just plausibly we could allocate matrix_complex but not count
1321 try { count = new size_t; } catch( std::bad_alloc& e ){
1322 // try to tidy up before rethrowing
1323 gsl_matrix_complex_free( ccgsl_pointer );
1324 throw e;
1325 }
1326 *count = 1; // initially there is just one reference to ccgsl_pointer
1327 }
1328 // More functions
1332 void set_identity(){ gsl_matrix_complex_set_identity( get() ); }
1339 int swap_rows( size_t const i, size_t const j ){ return gsl_matrix_complex_swap_rows( get(), i, j ); }
1346 int swap_columns( size_t const i, size_t const j ){
1347 return gsl_matrix_complex_swap_columns( get(), i, j ); }
1355 int swap_rowcol( size_t const i, size_t const j ){ return gsl_matrix_complex_swap_rowcol( get(), i, j ); }
1360 int transpose(){ return gsl_matrix_complex_transpose( get() ); }
1367 return gsl_matrix_complex_transpose_memcpy( get(), src.get() ); }
1373 int
1375 return gsl_matrix_complex_mul_elements( get(), b.get() ); }
1383 return gsl_matrix_complex_div_elements( get(), b.get() ); }
1390 return gsl_matrix_complex_conjtrans_memcpy( get(), src.get() ); }
1397 return gsl_matrix_complex_scale_rows( get(), x.get() ); }
1404 return gsl_matrix_complex_scale_columns( get(), x.get() ); }
1410 int add_diagonal( complex const x ){
1411 return gsl_matrix_complex_add_diagonal( get(), x ); }
1418 int get_row( vector_complex& v, size_t const i ) const {
1419 return gsl_matrix_complex_get_row( v.get(), get(), i ); }
1426 int get_col( vector_complex& v, size_t const j ) const {
1427 return gsl_matrix_complex_get_col( v.get(), get(), j ); }
1434 int set_row( size_t const i, vector_complex const& v ){
1435 return gsl_matrix_complex_set_row( get(), i, v.get() ); }
1442 int set_col( size_t const j, vector_complex const& v ){
1443 return gsl_matrix_complex_set_col( get(), j, v.get() ); }
1444 // Extra functions for []
1452 vector_complex operator[]( size_t const i ){
1453 // First check that iterator is initialised.
1454 if( ccgsl_pointer == 0 ){
1455 gsl_error( "matrix_complex is null", __FILE__, __LINE__, exception::GSL_EFAILED );
1456 return vector_complex();
1457 }
1458 vector_complex v(0);
1459 gsl_vector_complex_view w = gsl_matrix_complex_row( ccgsl_pointer, i );
1461 return v;
1462 }
1470 vector_complex const operator[]( size_t const i ) const {
1471 // First check that iterator is initialised.
1472 if( ccgsl_pointer == 0 ){
1473 gsl_error( "matrix_complex is null", __FILE__, __LINE__, exception::GSL_EFAILED );
1474 return vector_complex();
1475 }
1476 vector_complex v(0);
1477 gsl_vector_complex_view w = gsl_matrix_complex_row( ccgsl_pointer, i );
1479 return v;
1480 }
1486 inline int
1488 return gsl_permute_matrix_complex( p.get(), get() ); }
1489 };
1490
1491 // Extra functions for vector_complex allocation from matrix_complex objects
1493 inline vector_complex vector_complex::alloc_row_from_matrix( matrix_complex& m, size_t const i ){
1494 return vector_complex ( gsl_vector_complex_alloc_row_from_matrix( m.get(), i ) ); }
1495 inline vector_complex vector_complex::alloc_col_from_matrix( matrix_complex& m, size_t const i ){
1496 return vector_complex ( gsl_vector_complex_alloc_col_from_matrix( m.get(), i ) ); }
1498}
1499#undef CCGSL_MTY
1500#endif
This class handles vector_complexs as shared handles.
This class can be used like a pointer for complex objects so that we can iterate over a vector (for e...
Definition: complex.hpp:691
This class handles complex numbers.
Definition: complex.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
A class template for the const iterators.
bool operator!=(iterator_t< reverse_t > const &i) const
Comparison with non-const iterator.
const_iterator_t()
The default constructor.
const_iterator_t< reverse_t > & operator=(const_iterator_t< reverse_t > const &i)
We can assign one output iterator from another.
bool operator!=(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
const_iterator_t< reverse_t > operator--(int)
The postfix – operator.
const_iterator_t< reverse_t > & operator++()
The prefix ++ operator.
const_iterator_t(matrix_complex const *v, size_t position)
This constructor allows vector to create non-default iterators.
const_iterator_t< reverse_t > & operator--()
The prefix – operator.
bool operator==(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
const_iterator_t(iterator_t< reverse_t > const &i)
A copy constructor.
bool operator==(iterator_t< reverse_t > const &i) const
Comparison with non-const iterator.
const_iterator_t< reverse_t > operator++(int)
The postfix ++ operator.
We create a suitable class for iterator types here.
vector_complex_ptr pointer
An iterator must have a pointer typea.
value_type reference
An iterator must have a reference type.
void increment()
Increment the iterator.
void decrement()
Decrement the iterator.
vector_complex value_type
An iterator must have a value type.
iterator_base()
The iterator is default constructible.
iterator_base(container *v, size_t position)
This constructor allows vector to create non-default iterators.
std::bidirectional_iterator_tag iterator_category
An iterator must have an iterator category.
reference operator*() const
Dereference the pointer.
bool operator==(iterator_base< container, content, reverse_t > const &i) const
The == operator.
size_t position
Mark position of iterator within matrix.
container * v
Store a pointer to a matrix we can iterate over: 0 if no matrix.
pointer operator->() const
Dereference the pointer.
bool operator!=(iterator_base< container, content, reverse_t > const &i) const
The != operator.
A class template for the two non-const iterators.
iterator_t< reverse_t > & operator--()
The prefix – operator.
iterator_t(matrix_complex *v, size_t position)
This constructor allows vector to create non-default iterators.
iterator_t< reverse_t > & operator=(iterator_t< reverse_t > const &i)
We can assign one output iterator from another.
iterator_t< reverse_t > & operator++()
The prefix ++ operator.
iterator_t< reverse_t > operator--(int)
The postfix – operator.
iterator_t< reverse_t > operator++(int)
The postfix ++ operator.
iterator_t()
The default constructor.
bool operator==(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
bool operator!=(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
This class handles matrix_complex objects as shared handles.
vector_complex column(size_t const j)
C++ version of gsl_matrix_complex_column().
static matrix_complex 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_complex_const_view_array_with_tda().
matrix_complex 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_complex_const_submatrix().
reverse_iterator rend()
Get iterator pointing beyond last vector_complex element.
vector_complex const const_diagonal() const
C++ version of gsl_matrix_complex_const_diagonal().
static matrix_complex const view_array(double const *base, size_t const n1, size_t const n2)
Another C++ version of gsl_matrix_complex_const_view_array().
void set_identity()
C++ version of gsl_matrix_complex_set_identity().
int mul_elements(matrix_complex const &b)
C++ version of gsl_matrix_complex_mul_elements().
complex_ptr const const_ptr(size_t const i, size_t const j) const
C++ version of gsl_matrix_complex_const_ptr().
vector_complex operator[](size_t const i)
This function allows us to use a matrix_complex like an array.
iterator_t< false > iterator
The iterator type.
matrix_complex & operator=(matrix_complex &&v)
Move operator.
int scale(complex const x)
C++ version of gsl_matrix_complex_scale().
vector_complex const subrow(size_t const i, size_t const offset, size_t const n) const
Another C++ version of gsl_matrix_complex_const_subrow().
matrix_complex(size_t const n1, size_t const n2)
The default constructor creates a new matrix_complex with n1 rows and n2 columns.
matrix_complex(gsl_matrix_complex *v)
Could construct from a gsl_matrix_complex.
int transpose_memcpy(matrix_complex const &src)
C++ version of gsl_matrix_complex_transpose_memcpy().
static matrix_complex view_array(double *base, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_view_array().
vector_complex row(size_t const i)
C++ version of gsl_matrix_complex_row().
gsl_matrix_complex * get()
Get the gsl_matrix_complex.
size_t size1() const
The number of rows of the matrix_complex.
complex get(size_t const i, size_t const j) const
C++ version of gsl_matrix_complex_get().
iterator begin()
Get iterator pointing to first vector_complex element.
void reset()
Stop sharing ownership of the shared pointer.
static matrix_complex view_vector(vector_complex &v, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_view_vector().
size_t * count
The shared reference count.
vector_complex const operator[](size_t const i) const
This function allows us to use a matrix_complex like an array.
size_t use_count() const
Find how many matrix_complex objects share this pointer.
matrix_complex submatrix(size_t const i, size_t const j, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_submatrix().
void set(size_t const i, size_t const j, complex x)
C++ version of gsl_matrix_complex_set().
vector_complex const subdiagonal(size_t const k) const
Another C++ version of gsl_matrix_complex_const_subdiagonal().
vector_complex const subcolumn(size_t const j, size_t const offset, size_t const n) const
Another C++ version of gsl_matrix_complex_const_subcolumn().
int fprintf(FILE *stream, char const *format) const
C++ version of gsl_matrix_complex_fprintf().
void swap(matrix_complex &m)
Swap two matrix_complex objects.
vector_complex superdiagonal(size_t const k)
C++ version of gsl_matrix_complex_superdiagonal().
vector_complex diagonal()
C++ version of gsl_matrix_complex_diagonal().
matrix_complex 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_complex_const_submatrix().
int swap_rows(size_t const i, size_t const j)
C++ version of gsl_matrix_complex_swap_rows().
iterator_t< true > reverse_iterator
The reverse_iterator type.
static matrix_complex view_vector_with_tda(vector_complex &v, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_complex_view_vector_with_tda().
~matrix_complex()
The destructor only deletes the pointers if count reaches zero.
const_iterator_t< true > const_reverse_iterator
The const_reverse_t type.
void set_zero()
C++ version of gsl_matrix_complex_set_zero().
const_iterator begin() const
Get iterator pointing to first vector_complex element.
vector_complex subdiagonal(size_t const k)
C++ version of gsl_matrix_complex_subdiagonal().
static matrix_complex const view_vector(vector_complex const &v, size_t const n1, size_t const n2)
Another C++ version of gsl_matrix_complex_const_view_vector().
const_iterator end() const
Get iterator pointing beyond last vector_complex element.
matrix_complex(block_complex &b, size_t const offset, size_t const n1, size_t const n2, size_t const d2)
C++ version of gsl_matrix_complex_alloc_from_block().
int add_constant(complex const x)
C++ version of gsl_matrix_complex_add_constant().
const_reverse_iterator rend() const
Get iterator pointing beyond last vector_complex element.
int conjtrans_memcpy(matrix_complex const &src)
C++ version of gsl_matrix_complex_conjtrans_memcpy().
int add(matrix_complex const &b)
C++ version of gsl_matrix_complex_add().
vector_complex const const_subdiagonal(size_t const k) const
C++ version of gsl_matrix_complex_const_subdiagonal().
vector_complex const diagonal() const
Another C++ version of gsl_matrix_complex_const_diagonal().
const_reverse_iterator rbegin() const
Get iterator pointing to first vector_complex element.
vector_complex const const_column(size_t const j) const
C++ version of gsl_matrix_complex_const_column().
static matrix_complex 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_complex_const_view_array_with_tda().
vector_complex const const_subcolumn(size_t const j, size_t const offset, size_t const n) const
C++ version of gsl_matrix_complex_const_subcolumn().
vector_complex const const_subrow(size_t const i, size_t const offset, size_t const n) const
C++ version of gsl_matrix_complex_const_subrow().
vector_complex subrow(size_t const i, size_t const offset, size_t const n)
C++ version of gsl_matrix_complex_subrow().
int swap_columns(size_t const i, size_t const j)
C++ version of gsl_matrix_complex_swap_columns().
int scale_columns(vector_complex const &x)
C++ version of gsl_matrix_complex_scale_columns().
int fread(FILE *stream)
C++ version of gsl_matrix_complex_fread().
int permute(permutation &p)
Permute the columns of this by permutation p.
matrix_complex(matrix_complex &&v)
Move constructor.
vector_complex const column(size_t const j) const
Another C++ version of gsl_matrix_complex_const_column().
matrix_complex(matrix_complex const &v)
The copy constructor.
int fscanf(FILE *stream)
C++ version of gsl_matrix_complex_fscanf().
size_t size_type
A container must have a size_type.
complex_ptr ptr(size_t const i, size_t const j)
C++ version of gsl_matrix_complex_ptr().
vector_complex const const_row(size_t const i) const
C++ version of gsl_matrix_complex_const_row().
int isneg() const
C++ version of gsl_matrix_complex_isneg().
static matrix_complex view_array_with_tda(double *base, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_complex_view_array_with_tda().
int scale_rows(vector_complex const &x)
C++ version of gsl_matrix_complex_scale_rows().
int add_diagonal(complex const x)
C++ version of gsl_matrix_complex_add_diagonal().
void tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex const &src)
Copy the upper or lower triangular part of matrix src to this.
vector_complex subcolumn(size_t const j, size_t const offset, size_t const n)
C++ version of gsl_matrix_complex_subcolumn().
vector_complex const row(size_t const i) const
Another C++ version of gsl_matrix_complex_const_row().
matrix_complex(matrix_complex &m, size_t const k1, size_t const k2, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_alloc_from_matrix().
int isnull() const
C++ version of gsl_matrix_complex_isnull().
int set_col(size_t const j, vector_complex const &v)
C++ version of gsl_matrix_complex_set_col().
static matrix_complex const const_view_vector_with_tda(vector_complex const &v, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_complex_const_view_vector_with_tda().
int div_elements(matrix_complex const &b)
C++ version of gsl_matrix_complex_div_elements().
static matrix_complex const const_view_vector(vector_complex const &v, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_const_view_vector().
void transpose_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex const &src)
Copy the upper or lower triangular part of matrix src to this.
const_iterator_t< false > const_iterator
The const_iterator type.
vector_complex const const_superdiagonal(size_t const k) const
C++ version of gsl_matrix_complex_const_superdiagonal().
static matrix_complex calloc(size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_calloc().
matrix_complex clone() const
The clone function.
matrix_complex & operator=(matrix_complex const &v)
The assignment operator.
iterator end()
Get iterator pointing beyond last vector_complex element.
int transpose()
C++ version of gsl_matrix_complex_transpose().
size_t size2() const
The number of columns of the matrix_complex.
vector_complex const superdiagonal(size_t const k) const
Another C++ version of gsl_matrix_complex_const_superdiagonal().
int swap_rowcol(size_t const i, size_t const j)
C++ version of gsl_matrix_complex_swap_rowcol().
static matrix_complex const view_vector_with_tda(vector_complex const &v, size_t const n1, size_t const n2, size_t const tda)
Another C++ version of gsl_matrix_complex_const_view_vector_with_tda().
int ispos() const
C++ version of gsl_matrix_complex_ispos().
matrix_complex()
The default constructor is only really useful for assigning to.
int set_row(size_t const i, vector_complex const &v)
C++ version of gsl_matrix_complex_set_row().
gsl_matrix_complex const * get() const
Get the gsl_matrix_complex.
reverse_iterator rbegin()
Get iterator pointing to first vector_complex element.
gsl_matrix_complex * ccgsl_pointer
The shared pointer.
int fwrite(FILE *stream) const
C++ version of gsl_matrix_complex_fwrite().
int memcpy(matrix_complex const &src)
C++ version of gsl_matrix_complex_memcpy().
int sub(matrix_complex const &b)
C++ version of gsl_matrix_complex_sub().
void set_all(complex x)
C++ version of gsl_matrix_complex_set_all().
bool unique() const
Find if this is the only object sharing the gsl_matrix_complex.
int get_row(vector_complex &v, size_t const i) const
C++ version of gsl_matrix_complex_get_row().
int get_col(vector_complex &v, size_t const j) const
C++ version of gsl_matrix_complex_get_col().
int isnonneg() const
C++ version of gsl_matrix_complex_isnonneg().
matrix_complex(std::initializer_list< std::initializer_list< std::complex< double > > > initializer_list)
Could construct from a std::initializer_list in C++11.
static matrix_complex const const_view_array(double const *base, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_const_view_array().
This class handles GSL permutation objects.
Definition: permutation.hpp:33
gsl_permutation * get()
Get the gsl_permutation.
This class handles vector_complex objects as shared handles.
static vector_complex alloc_row_from_matrix(matrix_complex &m, size_t const i)
C++ version of gsl_vector_complex_alloc_row_from_matrix().
static vector_complex alloc_col_from_matrix(matrix_complex &m, size_t const j)
C++ version of gsl_vector_complex_alloc_col_from_matrix().
gsl_vector_complex * get()
Get the gsl_vector_complex.
void wrap_gsl_vector_complex_without_ownership(gsl_vector_complex *v)
This function is intended mainly for internal use.
vector_complex()
The default constructor is only really useful for assigning to.
#define CCGSL_MTY
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.
vector_complex & operator*()
Dereference operator.
vector_complex_ptr(vector_complex const &v)
Typically we have to construct from a vector_complex.
vector_complex * operator->()
Dereference operator.