20#ifndef CCGSL_MATRIX_COMPLEX_HPP
21#define CCGSL_MATRIX_COMPLEX_HPP
23#include<gsl/gsl_matrix.h>
25#ifdef __GXX_EXPERIMENTAL_CXX0X__
28#include<gsl/gsl_permute_matrix_complex_double.h>
79 if( n1 > 0 and n2 > 0 )
ccgsl_pointer = gsl_matrix_complex_alloc( n1, n2 );
83 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
85 if( n1 > 0 and n2 > 0 ) gsl_matrix_complex_free(
ccgsl_pointer );
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__,
119 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
121 if( n1 > 0 and n2 > 0 ) gsl_matrix_complex_free(
ccgsl_pointer );
128 for(
auto row : initializer_list ){
130 for(
auto x :
row ){
set( r, c, x ); ++c; }
170 gsl_matrix_complex_memcpy( copy.
get(),
get() );
191#ifdef __GXX_EXPERIMENTAL_CXX0X__
197 std::swap(
count, v.count );
198 v.ccgsl_pointer =
nullptr;
235 template<
typename container,
typename content,
bool reverse_t>
class iterator_base {
262 static content something;
267 }
else if(
v->ccgsl_pointer == 0 ){
286 static content something_base;
292 }
else if(
v->ccgsl_pointer == 0 ){
303#ifdef __GXX_EXPERIMENTAL_CXX0X__
304 return std::move(
ptr );
335 }
else if(
v->ccgsl_pointer == 0 ){
352 }
else if(
v->ccgsl_pointer == 0 ){
473 :
public iterator_base<matrix_complex const,vector_complex const,reverse_t>{
694 gsl_matrix_complex_tricpy( Uplo, Diag,
get(), src.
get() );
703 gsl_matrix_complex_transpose_tricpy( Uplo, Diag,
get(), src.
get() );
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;
725 gsl_vector_complex* w =
static_cast<gsl_vector_complex*
>( malloc(
sizeof( gsl_vector_complex ) ) );
726 *w = gsl_matrix_complex_row(
get(), i ).vector;
735 gsl_vector_complex* w =
static_cast<gsl_vector_complex*
>( malloc(
sizeof( gsl_vector_complex ) ) );
736 *w = gsl_matrix_complex_column(
get(), j ).vector;
744 diagonal(){ gsl_vector_complex* w =
static_cast<gsl_vector_complex*
>( malloc(
sizeof( gsl_vector_complex ) ) );
745 *w = gsl_matrix_complex_diagonal(
get() ).vector;
754 gsl_vector_complex* w =
static_cast<gsl_vector_complex*
>( malloc(
sizeof( gsl_vector_complex ) ) );
755 *w = gsl_matrix_complex_subdiagonal(
get(), k ).vector;
764 gsl_vector_complex* w =
static_cast<gsl_vector_complex*
>( malloc(
sizeof( gsl_vector_complex ) ) );
765 *w = gsl_matrix_complex_superdiagonal(
get(), k ).vector;
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;
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;
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;
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;
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;
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;
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;
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;
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;
881 gsl_vector_complex* w =
static_cast<gsl_vector_complex*
>( malloc(
sizeof( gsl_vector_complex ) ) );
882 *w = gsl_matrix_complex_const_diagonal(
get() ).vector;
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;
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;
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;
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;
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;
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;
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;
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;
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;
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;
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;
1019 gsl_vector_complex* w =
static_cast<gsl_vector_complex*
>( malloc(
sizeof( gsl_vector_complex ) ) );
1020 *w = gsl_matrix_complex_const_diagonal(
get() ).vector;
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;
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;
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;
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;
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;
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;
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;
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;
1156#ifdef __GXX_EXPERIMENTAL_CXX0X__
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() ); }
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 ); }
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 ); }
1301 ccgsl_pointer = gsl_matrix_complex_alloc_from_block(
b.get(), offset, n1, n2, d2 );
1303 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
1319 ccgsl_pointer = gsl_matrix_complex_alloc_from_matrix( m.
get(), k1, k2, n1, n2 );
1321 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
1339 int swap_rows(
size_t const i,
size_t const j ){
return gsl_matrix_complex_swap_rows(
get(), i, 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 ); }
1367 return gsl_matrix_complex_transpose_memcpy(
get(), src.
get() ); }
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() ); }
1411 return gsl_matrix_complex_add_diagonal(
get(), x ); }
1419 return gsl_matrix_complex_get_row( v.
get(),
get(), i ); }
1427 return gsl_matrix_complex_get_col( v.
get(),
get(), j ); }
1435 return gsl_matrix_complex_set_row(
get(), i, v.
get() ); }
1443 return gsl_matrix_complex_set_col(
get(), j, v.
get() ); }
1459 gsl_vector_complex_view w = gsl_matrix_complex_row(
ccgsl_pointer, i );
1477 gsl_vector_complex_view w = gsl_matrix_complex_row(
ccgsl_pointer, i );
1488 return gsl_permute_matrix_complex( p.
get(),
get() ); }
1494 return vector_complex ( gsl_vector_complex_alloc_row_from_matrix( m.get(), i ) ); }
1496 return vector_complex ( gsl_vector_complex_alloc_col_from_matrix( m.get(), i ) ); }
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...
This class handles complex numbers.
This class is used to handle gsl exceptions so that gsl can use these rather than the GSL error handl...
@ GSL_EFAILED
generic failure
@ GSL_EINVAL
invalid argument supplied by user
@ GSL_EBADLEN
matrix, vector lengths are not conformant
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.
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.
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
double b(int order, double qq)
C++ version of gsl_sf_mathieu_b().
gsl_sf_result result
Typedef for gsl_sf_result.
The gsl package creates an interface to the GNU Scientific Library for C++.
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.