ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
matrix_complex_float.hpp
Go to the documentation of this file.
1/*
2 * $Id: matrix_complex_float.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_FLOAT_HPP
21#define CCGSL_MATRIX_COMPLEX_FLOAT_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_float.h>
29
30#include"exception.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_float( size_t const n1, size_t const n2 ){
79 if( n1 > 0 and n2 > 0 ) ccgsl_pointer = gsl_matrix_complex_float_alloc( n1, n2 );
80 else { ccgsl_pointer = new gsl_matrix_complex_float; 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_float_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_float( gsl_matrix_complex_float* 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_float( std::initializer_list<std::initializer_list<std::complex<float> > > 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_float_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_float_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_float 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_float_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_float is now shared.
159 return *this;
160 }
161 // clone()
168 matrix_complex_float copy( get()->size1, get()->size2 );
169 // Now copy
170 gsl_matrix_complex_float_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_float_free( ccgsl_pointer );
183 else delete ccgsl_pointer; }
184 delete count;
185 }
186 }
190 void reset(){ matrix_complex_float().swap( *this ); }
191#ifdef __GXX_EXPERIMENTAL_CXX0X__
197 std::swap( count, v.count );
198 v.ccgsl_pointer = nullptr;
199 }
206 matrix_complex_float( std::move( v ) ).swap( *this );
207 return *this;
208 }
209#endif
210 private:
224 vector_complex_float& operator*(){ return *this; }
230 };
235 template<typename container,typename content,bool reverse_t> class iterator_base {
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_float_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
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_float,vector_complex_float,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:
460 // We need a constructor for matrix_complex_float
468 };
472 template<bool reverse_t> class const_iterator_t
473 : public iterator_base<matrix_complex_float const,vector_complex_float 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_float
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_float* 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_float const& src){
694 gsl_matrix_complex_float_tricpy( Uplo, Diag, get(), src.get() );
695 }
702 void
703 transpose_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex_float const& src){
704 gsl_matrix_complex_float_transpose_tricpy( Uplo, Diag, get(), src.get() );
705 }
706 // view operations
715 matrix_complex_float submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ){
716 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
717 *m = gsl_matrix_complex_float_submatrix( get(), i, j, n1, n2 ).matrix;
718 return matrix_complex_float( m );
719 }
725 vector_complex_float row( size_t const i ){
726 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
727 *w = gsl_matrix_complex_float_row( get(), i ).vector;
728 return vector_complex_float( w );
729 }
735 vector_complex_float column( size_t const j ){
736 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
737 *w = gsl_matrix_complex_float_column( get(), j ).vector;
738 return vector_complex_float( w );
739 }
745 diagonal(){ gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
746 *w = gsl_matrix_complex_float_diagonal( get() ).vector;
747 return vector_complex_float( w );
748 }
755 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
756 *w = gsl_matrix_complex_float_subdiagonal( get(), k ).vector;
757 return vector_complex_float( w );
758 }
765 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
766 *w = gsl_matrix_complex_float_superdiagonal( get(), k ).vector;
767 return vector_complex_float( w );
768 }
776 vector_complex_float subrow( size_t const i, size_t const offset, size_t const n ){
777 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
778 *w = gsl_matrix_complex_float_subrow( get(), i, offset, n ).vector;
779 return vector_complex_float( w );
780 }
788 vector_complex_float subcolumn( size_t const j, size_t const offset, size_t const n ){
789 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
790 *w = gsl_matrix_complex_float_subcolumn( get(), j, offset, n ).vector;
791 return vector_complex_float( w );
792 }
800 static matrix_complex_float view_array( float* base, size_t const n1, size_t const n2 ){
801 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
802 *m = gsl_matrix_complex_float_view_array( base, n1, n2 ).matrix;
803 return matrix_complex_float( m );
804 }
813 static matrix_complex_float view_array_with_tda( float* base, size_t const n1, size_t const n2, size_t const tda ){
814 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
815 *m = gsl_matrix_complex_float_view_array_with_tda( base, n1, n2, tda ).matrix;
816 return matrix_complex_float( m );
817 }
825 static matrix_complex_float view_vector( vector_complex_float& v, size_t const n1, size_t const n2 ){
826 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
827 *m = gsl_matrix_complex_float_view_vector( v.get(), n1, n2 ).matrix;
828 return matrix_complex_float( m );
829 }
838 static matrix_complex_float view_vector_with_tda( vector_complex_float& v, size_t const n1, size_t const n2, size_t const tda ){
839 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
840 *m = gsl_matrix_complex_float_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
841 return matrix_complex_float( m );
842 }
843 // const versions ...
852 matrix_complex_float const const_submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ) const {
853 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
854 *m = gsl_matrix_complex_float_const_submatrix( get(), i, j, n1, n2 ).matrix;
855 return matrix_complex_float( m );
856 }
862 vector_complex_float const const_row( size_t const i ) const {
863 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
864 *w = gsl_matrix_complex_float_const_row( get(), i ).vector;
865 return vector_complex_float( w );
866 }
872 vector_complex_float const const_column( size_t const j ) const {
873 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
874 *w = gsl_matrix_complex_float_const_column( get(), j ).vector;
875 return vector_complex_float( w );
876 }
882 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
883 *w = gsl_matrix_complex_float_const_diagonal( get() ).vector;
884 return vector_complex_float( w );
885 }
891 vector_complex_float const const_subdiagonal( size_t const k ) const {
892 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
893 *w = gsl_matrix_complex_float_const_subdiagonal( get(), k ).vector;
894 return vector_complex_float( w );
895 }
901 vector_complex_float const const_superdiagonal( size_t const k ) const {
902 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
903 *w = gsl_matrix_complex_float_const_superdiagonal( get(), k ).vector;
904 return vector_complex_float( w );
905 }
913 vector_complex_float const const_subrow( size_t const i, size_t const offset, size_t const n ) const {
914 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
915 *w = gsl_matrix_complex_float_const_subrow( get(), i, offset, n ).vector;
916 return vector_complex_float( w );
917 }
925 vector_complex_float const const_subcolumn( size_t const j, size_t const offset, size_t const n ) const {
926 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
927 *w = gsl_matrix_complex_float_const_subcolumn( get(), j, offset, n ).vector;
928 return vector_complex_float( w );
929 }
937 static matrix_complex_float const const_view_array( float const* base, size_t const n1, size_t const n2 ){
938 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
939 *m = gsl_matrix_complex_float_const_view_array( base, n1, n2 ).matrix;
940 return matrix_complex_float( m );
941 }
950 static matrix_complex_float const
951 const_view_array_with_tda( float const* base, size_t const n1, size_t const n2, size_t const tda ){
952 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
953 *m = gsl_matrix_complex_float_const_view_array_with_tda( base, n1, n2, tda ).matrix;
954 return matrix_complex_float( m );
955 }
963 static matrix_complex_float const const_view_vector( vector_complex_float const& v, size_t const n1, size_t const n2 ){
964 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
965 *m = gsl_matrix_complex_float_const_view_vector( v.get(), n1, n2 ).matrix;
966 return matrix_complex_float( m );
967 }
976 static matrix_complex_float const
977 const_view_vector_with_tda( vector_complex_float const& v, size_t const n1, size_t const n2, size_t const tda ){
978 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
979 *m = gsl_matrix_complex_float_const_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
980 return matrix_complex_float( m );
981 }
990 matrix_complex_float const submatrix( size_t const i, size_t const j, size_t const n1, size_t const n2 ) const {
991 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
992 *m = gsl_matrix_complex_float_const_submatrix( get(), i, j, n1, n2 ).matrix;
993 return matrix_complex_float( m );
994 }
1000 vector_complex_float const row( size_t const i ) const {
1001 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
1002 *w = gsl_matrix_complex_float_const_row( get(), i ).vector;
1003 return vector_complex_float( w );
1004 }
1010 vector_complex_float const column( size_t const j ) const {
1011 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
1012 *w = gsl_matrix_complex_float_const_column( get(), j ).vector;
1013 return vector_complex_float( w );
1014 }
1020 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
1021 *w = gsl_matrix_complex_float_const_diagonal( get() ).vector;
1022 return vector_complex_float( w );
1023 }
1029 vector_complex_float const subdiagonal( size_t const k ) const {
1030 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
1031 *w = gsl_matrix_complex_float_const_subdiagonal( get(), k ).vector;
1032 return vector_complex_float( w );
1033 }
1039 vector_complex_float const superdiagonal( size_t const k ) const {
1040 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
1041 *w = gsl_matrix_complex_float_const_superdiagonal( get(), k ).vector;
1042 return vector_complex_float( w );
1043 }
1051 vector_complex_float const subrow( size_t const i, size_t const offset, size_t const n ) const {
1052 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
1053 *w = gsl_matrix_complex_float_const_subrow( get(), i, offset, n ).vector;
1054 return vector_complex_float( w );
1055 }
1063 vector_complex_float const subcolumn( size_t const j, size_t const offset, size_t const n ) const {
1064 gsl_vector_complex_float* w = static_cast<gsl_vector_complex_float*>( malloc( sizeof( gsl_vector_complex_float ) ) );
1065 *w = gsl_matrix_complex_float_const_subcolumn( get(), j, offset, n ).vector;
1066 return vector_complex_float( w );
1067 }
1075 static matrix_complex_float const view_array( float const* base, size_t const n1, size_t const n2 ){
1076 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
1077 *m = gsl_matrix_complex_float_const_view_array( base, n1, n2 ).matrix;
1078 return matrix_complex_float( m );
1079 }
1088 static matrix_complex_float const
1089 view_array_with_tda( float const* base, size_t const n1, size_t const n2, size_t const tda ){
1090 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
1091 *m = gsl_matrix_complex_float_const_view_array_with_tda( base, n1, n2, tda ).matrix;
1092 return matrix_complex_float( m );
1093 }
1101 static matrix_complex_float const view_vector( vector_complex_float const& v, size_t const n1, size_t const n2 ){
1102 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
1103 *m = gsl_matrix_complex_float_const_view_vector( v.get(), n1, n2 ).matrix;
1104 return matrix_complex_float( m );
1105 }
1114 static matrix_complex_float const
1115 view_vector_with_tda( vector_complex_float const& v, size_t const n1, size_t const n2, size_t const tda ){
1116 gsl_matrix_complex_float* m = static_cast<gsl_matrix_complex_float*>( malloc( sizeof( gsl_matrix_complex_float ) ) );
1117 *m = gsl_matrix_complex_float_const_view_vector_with_tda( v.get(), n1, n2, tda ).matrix;
1118 return matrix_complex_float( m );
1119 }
1120 private:
1124 gsl_matrix_complex_float* ccgsl_pointer;
1128 size_t* count;
1129 public:
1130 // shared reference functions
1135 gsl_matrix_complex_float* get(){ return ccgsl_pointer; }
1140 gsl_matrix_complex_float const* get() const { return ccgsl_pointer; }
1146 bool unique() const { return count != 0 and *count == 1; }
1151 size_t use_count() const { return count == 0 ? 0 : *count; }
1157#ifdef __GXX_EXPERIMENTAL_CXX0X__
1158 explicit
1159#endif
1160 operator bool() const { return ccgsl_pointer != 0; }
1161 // GSL functions
1169 static matrix_complex_float calloc( size_t const n1, size_t const n2 ){ return matrix_complex_float( gsl_matrix_complex_float_calloc( n1, n2 ) ); }
1173 void set_zero(){ gsl_matrix_complex_float_set_zero( get() ); }
1178 void set_all( complex_float x ){ gsl_matrix_complex_float_set_all( get(), x ); }
1184 int memcpy( matrix_complex_float const& src ){ return gsl_matrix_complex_float_memcpy( get(), src.get() ); }
1190 int add( matrix_complex_float const& b ){ return gsl_matrix_complex_float_add( get(), b.get() ); }
1196 int sub( matrix_complex_float const& b ){ return gsl_matrix_complex_float_sub( get(), b.get() ); }
1202 int scale( complex_float const x ){ return gsl_matrix_complex_float_scale( get(), x ); }
1208 int add_constant( complex_float const x ){ return gsl_matrix_complex_float_add_constant( get(), x ); }
1213 int isnull() const { return gsl_matrix_complex_float_isnull( get() ); }
1218 int ispos() const { return gsl_matrix_complex_float_ispos( get() ); }
1223 int isneg() const { return gsl_matrix_complex_float_isneg( get() ); }
1228 int isnonneg() const { return gsl_matrix_complex_float_isnonneg( get() ); }
1235 complex_float get( size_t const i, size_t const j ) const { return gsl_matrix_complex_float_get( get(), i, j ); }
1242 void set( size_t const i, size_t const j, complex_float x ){ gsl_matrix_complex_float_set( get(), i, j, x ); }
1249 complex_float_ptr ptr( size_t const i, size_t const j ){
1250 if( i >= ccgsl_pointer->size1 )
1251 gsl_error( "Index out of range", __FILE__, __LINE__, exception::GSL_EINVAL );
1252 if( j >= ccgsl_pointer->size2 )
1253 gsl_error( "Index out of range", __FILE__, __LINE__, exception::GSL_EINVAL );
1254 return complex_float_ptr( get()->data + CCGSL_MTY * (i * get()->tda + j ) ); }
1261 complex_float_ptr const const_ptr( size_t const i, size_t const j ) const {
1262 if( i >= ccgsl_pointer->size1 )
1263 gsl_error( "Index out of range", __FILE__, __LINE__, exception::GSL_EINVAL );
1264 if( j >= ccgsl_pointer->size2 )
1265 gsl_error( "Index out of range", __FILE__, __LINE__, exception::GSL_EINVAL );
1266 return complex_float_ptr( get()->data + CCGSL_MTY * (i * get()->tda + j ) ); }
1272 int fread( FILE* stream ){ return gsl_matrix_complex_float_fread( stream, get() ); }
1278 int fwrite( FILE* stream ) const { return gsl_matrix_complex_float_fwrite( stream, get() ); }
1284 int fscanf( FILE* stream ){ return gsl_matrix_complex_float_fscanf( stream, get() ); }
1291 int fprintf( FILE* stream, char const* format ) const {
1292 return gsl_matrix_complex_float_fprintf( stream, get(), format ); }
1301 matrix_complex_float( block_complex_float& b, size_t const offset, size_t const n1, size_t const n2, size_t const d2 ){
1302 ccgsl_pointer = gsl_matrix_complex_float_alloc_from_block( b.get(), offset, n1, n2, d2 );
1303 // just plausibly we could allocate vector_complex_float but not count
1304 try { count = new size_t; } catch( std::bad_alloc& e ){
1305 // try to tidy up before rethrowing
1306 gsl_matrix_complex_float_free( ccgsl_pointer );
1307 throw e;
1308 }
1309 *count = 1; // initially there is just one reference to ccgsl_pointer
1310 }
1319 matrix_complex_float( matrix_complex_float& m, size_t const k1, size_t const k2, size_t const n1, size_t const n2 ){
1320 ccgsl_pointer = gsl_matrix_complex_float_alloc_from_matrix( m.get(), k1, k2, n1, n2 );
1321 // just plausibly we could allocate matrix_complex_float but not count
1322 try { count = new size_t; } catch( std::bad_alloc& e ){
1323 // try to tidy up before rethrowing
1324 gsl_matrix_complex_float_free( ccgsl_pointer );
1325 throw e;
1326 }
1327 *count = 1; // initially there is just one reference to ccgsl_pointer
1328 }
1329 // More functions
1333 void set_identity(){ gsl_matrix_complex_float_set_identity( get() ); }
1340 int swap_rows( size_t const i, size_t const j ){ return gsl_matrix_complex_float_swap_rows( get(), i, j ); }
1347 int swap_columns( size_t const i, size_t const j ){
1348 return gsl_matrix_complex_float_swap_columns( get(), i, j ); }
1356 int swap_rowcol( size_t const i, size_t const j ){ return gsl_matrix_complex_float_swap_rowcol( get(), i, j ); }
1361 int transpose(){ return gsl_matrix_complex_float_transpose( get() ); }
1368 return gsl_matrix_complex_float_transpose_memcpy( get(), src.get() ); }
1374 int
1376 return gsl_matrix_complex_float_mul_elements( get(), b.get() ); }
1384 return gsl_matrix_complex_float_div_elements( get(), b.get() ); }
1391 return gsl_matrix_complex_float_conjtrans_memcpy( get(), src.get() ); }
1398 return gsl_matrix_complex_float_scale_rows( get(), x.get() ); }
1405 return gsl_matrix_complex_float_scale_columns( get(), x.get() ); }
1412 return gsl_matrix_complex_float_add_diagonal( get(), x ); }
1419 int get_row( vector_complex_float& v, size_t const i ) const {
1420 return gsl_matrix_complex_float_get_row( v.get(), get(), i ); }
1427 int get_col( vector_complex_float& v, size_t const j ) const {
1428 return gsl_matrix_complex_float_get_col( v.get(), get(), j ); }
1435 int set_row( size_t const i, vector_complex_float const& v ){
1436 return gsl_matrix_complex_float_set_row( get(), i, v.get() ); }
1443 int set_col( size_t const j, vector_complex_float const& v ){
1444 return gsl_matrix_complex_float_set_col( get(), j, v.get() ); }
1445 // Extra functions for []
1454 // First check that iterator is initialised.
1455 if( ccgsl_pointer == 0 ){
1456 gsl_error( "matrix_complex_float is null", __FILE__, __LINE__, exception::GSL_EFAILED );
1457 return vector_complex_float();
1458 }
1460 gsl_vector_complex_float_view w = gsl_matrix_complex_float_row( ccgsl_pointer, i );
1462 return v;
1463 }
1471 vector_complex_float const operator[]( size_t const i ) const {
1472 // First check that iterator is initialised.
1473 if( ccgsl_pointer == 0 ){
1474 gsl_error( "matrix_complex_float is null", __FILE__, __LINE__, exception::GSL_EFAILED );
1475 return vector_complex_float();
1476 }
1478 gsl_vector_complex_float_view w = gsl_matrix_complex_float_row( ccgsl_pointer, i );
1480 return v;
1481 }
1487 inline int
1489 return gsl_permute_matrix_complex_float( p.get(), get() ); }
1490 };
1491
1492 // Extra functions for vector_complex_float allocation from matrix_complex_float objects
1494 inline vector_complex_float vector_complex_float::alloc_row_from_matrix( matrix_complex_float& m, size_t const i ){
1495 return vector_complex_float ( gsl_vector_complex_float_alloc_row_from_matrix( m.get(), i ) ); }
1496 inline vector_complex_float vector_complex_float::alloc_col_from_matrix( matrix_complex_float& m, size_t const i ){
1497 return vector_complex_float ( gsl_vector_complex_float_alloc_col_from_matrix( m.get(), i ) ); }
1499}
1500#undef CCGSL_MTY
1501#endif
This class handles vector_complex_floats as shared handles.
This class can be used like a pointer for complex_float objects so that we can iterate over a vector ...
This class handles complex_float numbers.
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.
const_iterator_t< reverse_t > operator--(int)
The postfix – 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.
bool operator!=(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
const_iterator_t(matrix_complex_float const *v, size_t position)
This constructor allows vector to create non-default iterators.
const_iterator_t< reverse_t > & operator=(const_iterator_t< reverse_t > const &i)
We can assign one output iterator from another.
const_iterator_t< reverse_t > & operator--()
The prefix – operator.
const_iterator_t< reverse_t > operator++(int)
The postfix ++ operator.
const_iterator_t< reverse_t > & operator++()
The prefix ++ operator.
bool operator==(iterator_t< reverse_t > const &i) const
Comparison with non-const iterator.
We create a suitable class for iterator types here.
bool operator!=(iterator_base< container, content, reverse_t > const &i) const
The != operator.
vector_complex_float_ptr pointer
An iterator must have a pointer typea.
pointer operator->() const
Dereference the pointer.
container * v
Store a pointer to a matrix we can iterate over: 0 if no matrix.
iterator_base(container *v, size_t position)
This constructor allows vector to create non-default iterators.
vector_complex_float value_type
An iterator must have a value type.
bool operator==(iterator_base< container, content, reverse_t > const &i) const
The == operator.
reference operator*() const
Dereference the pointer.
value_type reference
An iterator must have a reference type.
iterator_base()
The iterator is default constructible.
size_t position
Mark position of iterator within matrix.
std::bidirectional_iterator_tag iterator_category
An iterator must have an iterator category.
A class template for the two non-const iterators.
bool operator!=(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
iterator_t(matrix_complex_float *v, size_t position)
This constructor allows vector to create non-default iterators.
bool operator==(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
iterator_t< reverse_t > operator++(int)
The postfix ++ operator.
iterator_t< reverse_t > operator--(int)
The postfix – operator.
iterator_t< reverse_t > & operator--()
The prefix – operator.
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.
This class handles matrix_complex_float objects as shared handles.
int ispos() const
C++ version of gsl_matrix_complex_float_ispos().
vector_complex_float diagonal()
C++ version of gsl_matrix_complex_float_diagonal().
static matrix_complex_float calloc(size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_float_calloc().
iterator begin()
Get iterator pointing to first vector_complex_float element.
matrix_complex_float & operator=(matrix_complex_float &&v)
Move operator.
complex_float get(size_t const i, size_t const j) const
C++ version of gsl_matrix_complex_float_get().
int scale_rows(vector_complex_float const &x)
C++ version of gsl_matrix_complex_float_scale_rows().
void set_identity()
C++ version of gsl_matrix_complex_float_set_identity().
int isneg() const
C++ version of gsl_matrix_complex_float_isneg().
int set_col(size_t const j, vector_complex_float const &v)
C++ version of gsl_matrix_complex_float_set_col().
static matrix_complex_float const const_view_vector_with_tda(vector_complex_float const &v, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_complex_float_const_view_vector_with_tda().
vector_complex_float operator[](size_t const i)
This function allows us to use a matrix_complex_float like an array.
const_iterator end() const
Get iterator pointing beyond last vector_complex_float element.
void set(size_t const i, size_t const j, complex_float x)
C++ version of gsl_matrix_complex_float_set().
size_t size1() const
The number of rows of the matrix_complex_float.
vector_complex_float const subcolumn(size_t const j, size_t const offset, size_t const n) const
Another C++ version of gsl_matrix_complex_float_const_subcolumn().
size_t size_type
A container must have a size_type.
int fprintf(FILE *stream, char const *format) const
C++ version of gsl_matrix_complex_float_fprintf().
complex_float_ptr ptr(size_t const i, size_t const j)
C++ version of gsl_matrix_complex_float_ptr().
void set_zero()
C++ version of gsl_matrix_complex_float_set_zero().
static matrix_complex_float const view_vector_with_tda(vector_complex_float const &v, size_t const n1, size_t const n2, size_t const tda)
Another C++ version of gsl_matrix_complex_float_const_view_vector_with_tda().
vector_complex_float const operator[](size_t const i) const
This function allows us to use a matrix_complex_float like an array.
vector_complex_float subdiagonal(size_t const k)
C++ version of gsl_matrix_complex_float_subdiagonal().
complex_float_ptr const const_ptr(size_t const i, size_t const j) const
C++ version of gsl_matrix_complex_float_const_ptr().
matrix_complex_float & operator=(matrix_complex_float const &v)
The assignment operator.
matrix_complex_float()
The default constructor is only really useful for assigning to.
void set_all(complex_float x)
C++ version of gsl_matrix_complex_float_set_all().
vector_complex_float subrow(size_t const i, size_t const offset, size_t const n)
C++ version of gsl_matrix_complex_float_subrow().
bool unique() const
Find if this is the only object sharing the gsl_matrix_complex_float.
vector_complex_float const const_subdiagonal(size_t const k) const
C++ version of gsl_matrix_complex_float_const_subdiagonal().
vector_complex_float subcolumn(size_t const j, size_t const offset, size_t const n)
C++ version of gsl_matrix_complex_float_subcolumn().
reverse_iterator rbegin()
Get iterator pointing to first vector_complex_float element.
int mul_elements(matrix_complex_float const &b)
C++ version of gsl_matrix_complex_float_mul_elements().
void tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex_float const &src)
Copy the upper or lower triangular part of matrix src to this.
int swap_columns(size_t const i, size_t const j)
C++ version of gsl_matrix_complex_float_swap_columns().
vector_complex_float const diagonal() const
Another C++ version of gsl_matrix_complex_float_const_diagonal().
vector_complex_float const row(size_t const i) const
Another C++ version of gsl_matrix_complex_float_const_row().
iterator_t< true > reverse_iterator
The reverse_iterator type.
static matrix_complex_float const const_view_array_with_tda(float const *base, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_complex_float_const_view_array_with_tda().
void swap(matrix_complex_float &m)
Swap two matrix_complex_float objects.
vector_complex_float const const_subrow(size_t const i, size_t const offset, size_t const n) const
C++ version of gsl_matrix_complex_float_const_subrow().
matrix_complex_float(matrix_complex_float &&v)
Move constructor.
~matrix_complex_float()
The destructor only deletes the pointers if count reaches zero.
matrix_complex_float(size_t const n1, size_t const n2)
The default constructor creates a new matrix_complex_float with n1 rows and n2 columns.
vector_complex_float const const_subcolumn(size_t const j, size_t const offset, size_t const n) const
C++ version of gsl_matrix_complex_float_const_subcolumn().
size_t size2() const
The number of columns of the matrix_complex_float.
static matrix_complex_float view_vector_with_tda(vector_complex_float &v, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_complex_float_view_vector_with_tda().
vector_complex_float row(size_t const i)
C++ version of gsl_matrix_complex_float_row().
size_t use_count() const
Find how many matrix_complex_float objects share this pointer.
vector_complex_float const const_row(size_t const i) const
C++ version of gsl_matrix_complex_float_const_row().
int scale(complex_float const x)
C++ version of gsl_matrix_complex_float_scale().
vector_complex_float column(size_t const j)
C++ version of gsl_matrix_complex_float_column().
matrix_complex_float(gsl_matrix_complex_float *v)
Could construct from a gsl_matrix_complex_float.
int get_col(vector_complex_float &v, size_t const j) const
C++ version of gsl_matrix_complex_float_get_col().
iterator_t< false > iterator
The iterator type.
int add_constant(complex_float const x)
C++ version of gsl_matrix_complex_float_add_constant().
static matrix_complex_float view_array(float *base, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_float_view_array().
matrix_complex_float(matrix_complex_float const &v)
The copy constructor.
matrix_complex_float submatrix(size_t const i, size_t const j, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_float_submatrix().
matrix_complex_float clone() const
The clone function.
static matrix_complex_float const view_vector(vector_complex_float const &v, size_t const n1, size_t const n2)
Another C++ version of gsl_matrix_complex_float_const_view_vector().
int conjtrans_memcpy(matrix_complex_float const &src)
C++ version of gsl_matrix_complex_float_conjtrans_memcpy().
int fread(FILE *stream)
C++ version of gsl_matrix_complex_float_fread().
reverse_iterator rend()
Get iterator pointing beyond last vector_complex_float element.
matrix_complex_float(matrix_complex_float &m, size_t const k1, size_t const k2, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_float_alloc_from_matrix().
vector_complex_float const subrow(size_t const i, size_t const offset, size_t const n) const
Another C++ version of gsl_matrix_complex_float_const_subrow().
int transpose()
C++ version of gsl_matrix_complex_float_transpose().
int swap_rows(size_t const i, size_t const j)
C++ version of gsl_matrix_complex_float_swap_rows().
matrix_complex_float(std::initializer_list< std::initializer_list< std::complex< float > > > initializer_list)
Could construct from a std::initializer_list in C++11.
void reset()
Stop sharing ownership of the shared pointer.
const_reverse_iterator rend() const
Get iterator pointing beyond last vector_complex_float element.
vector_complex_float const const_diagonal() const
C++ version of gsl_matrix_complex_float_const_diagonal().
int scale_columns(vector_complex_float const &x)
C++ version of gsl_matrix_complex_float_scale_columns().
vector_complex_float const const_column(size_t const j) const
C++ version of gsl_matrix_complex_float_const_column().
const_reverse_iterator rbegin() const
Get iterator pointing to first vector_complex_float element.
vector_complex_float const const_superdiagonal(size_t const k) const
C++ version of gsl_matrix_complex_float_const_superdiagonal().
static matrix_complex_float const const_view_array(float const *base, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_float_const_view_array().
static matrix_complex_float view_array_with_tda(float *base, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_complex_float_view_array_with_tda().
matrix_complex_float 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_float_const_submatrix().
vector_complex_float const superdiagonal(size_t const k) const
Another C++ version of gsl_matrix_complex_float_const_superdiagonal().
gsl_matrix_complex_float * get()
Get the gsl_matrix_complex_float.
const_iterator begin() const
Get iterator pointing to first vector_complex_float element.
int permute(permutation &p)
Permute the columns of this by permutation p.
int swap_rowcol(size_t const i, size_t const j)
C++ version of gsl_matrix_complex_float_swap_rowcol().
gsl_matrix_complex_float const * get() const
Get the gsl_matrix_complex_float.
size_t * count
The shared reference count.
int isnull() const
C++ version of gsl_matrix_complex_float_isnull().
matrix_complex_float(block_complex_float &b, size_t const offset, size_t const n1, size_t const n2, size_t const d2)
C++ version of gsl_matrix_complex_float_alloc_from_block().
int fscanf(FILE *stream)
C++ version of gsl_matrix_complex_float_fscanf().
static matrix_complex_float const const_view_vector(vector_complex_float const &v, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_float_const_view_vector().
int transpose_memcpy(matrix_complex_float const &src)
C++ version of gsl_matrix_complex_float_transpose_memcpy().
int add_diagonal(complex_float const x)
C++ version of gsl_matrix_complex_float_add_diagonal().
int get_row(vector_complex_float &v, size_t const i) const
C++ version of gsl_matrix_complex_float_get_row().
int set_row(size_t const i, vector_complex_float const &v)
C++ version of gsl_matrix_complex_float_set_row().
const_iterator_t< false > const_iterator
The const_iterator type.
matrix_complex_float 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_float_const_submatrix().
const_iterator_t< true > const_reverse_iterator
The const_reverse_t type.
int memcpy(matrix_complex_float const &src)
C++ version of gsl_matrix_complex_float_memcpy().
gsl_matrix_complex_float * ccgsl_pointer
The shared pointer.
void transpose_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex_float const &src)
Copy the upper or lower triangular part of matrix src to this.
iterator end()
Get iterator pointing beyond last vector_complex_float element.
int div_elements(matrix_complex_float const &b)
C++ version of gsl_matrix_complex_float_div_elements().
vector_complex_float superdiagonal(size_t const k)
C++ version of gsl_matrix_complex_float_superdiagonal().
int fwrite(FILE *stream) const
C++ version of gsl_matrix_complex_float_fwrite().
vector_complex_float const column(size_t const j) const
Another C++ version of gsl_matrix_complex_float_const_column().
int sub(matrix_complex_float const &b)
C++ version of gsl_matrix_complex_float_sub().
int add(matrix_complex_float const &b)
C++ version of gsl_matrix_complex_float_add().
static matrix_complex_float view_vector(vector_complex_float &v, size_t const n1, size_t const n2)
C++ version of gsl_matrix_complex_float_view_vector().
int isnonneg() const
C++ version of gsl_matrix_complex_float_isnonneg().
static matrix_complex_float const view_array_with_tda(float const *base, size_t const n1, size_t const n2, size_t const tda)
Another C++ version of gsl_matrix_complex_float_const_view_array_with_tda().
vector_complex_float const subdiagonal(size_t const k) const
Another C++ version of gsl_matrix_complex_float_const_subdiagonal().
static matrix_complex_float const view_array(float const *base, size_t const n1, size_t const n2)
Another C++ version of gsl_matrix_complex_float_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_float objects as shared handles.
vector_complex_float()
The default constructor is only really useful for assigning to.
gsl_vector_complex_float * get()
Get the gsl_vector_complex_float.
void wrap_gsl_vector_complex_float_without_ownership(gsl_vector_complex_float *v)
This function is intended mainly for internal use.
static vector_complex_float alloc_row_from_matrix(matrix_complex_float &m, size_t const i)
C++ version of gsl_vector_complex_float_alloc_row_from_matrix().
static vector_complex_float alloc_col_from_matrix(matrix_complex_float &m, size_t const j)
C++ version of gsl_vector_complex_float_alloc_col_from_matrix().
#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_float & operator*()
Dereference operator.
vector_complex_float * operator->()
Dereference operator.
vector_complex_float_ptr(vector_complex_float const &v)
Typically we have to construct from a vector_complex_float.