20#ifndef CCGSL_MATRIX_HPP
21#define CCGSL_MATRIX_HPP
23#include<gsl/gsl_matrix.h>
24#include<gsl/gsl_permute_matrix_double.h>
38 namespace multilarge {
93 if( n1 > 0 and n2 > 0 )
ccgsl_pointer = gsl_matrix_alloc( n1, n2 );
97 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
116#ifdef __GXX_EXPERIMENTAL_CXX0X__
121 matrix( std::initializer_list<std::initializer_list<double> > initializer_list ) :
owns_data(true){
122 size_t const n1 = initializer_list.size();
123 size_t const n2 = initializer_list.begin()->size();
124 for(
auto l : initializer_list ){
125 if( l.size() != n2 ){
126 gsl::exception e(
"matrix rows have unequal sizes", __FILE__, __LINE__,
133 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
142 for(
auto row : initializer_list ){
144 for(
auto x :
row ){
set( r, c, x ); ++c; }
186 gsl_matrix_memcpy( copy.
get(),
get() );
231#ifdef __GXX_EXPERIMENTAL_CXX0X__
238 std::swap(
count, v.count );
239 v.ccgsl_pointer =
nullptr;
247 matrix( std::move( v ) ).swap( *
this );
276 template<
typename container,
typename content,
bool reverse_t>
class iterator_base {
303 static content something;
308 }
else if(
v->ccgsl_pointer == 0 ){
327 static content something_base;
328 static vector_ptr something( something_base );
333 }
else if(
v->ccgsl_pointer == 0 ){
344#ifdef __GXX_EXPERIMENTAL_CXX0X__
345 return std::move(
ptr );
376 }
else if(
v->ccgsl_pointer == 0 ){
393 }
else if(
v->ccgsl_pointer == 0 ){
727 if(
ccgsl_pointer == 0 ) gsl_error(
"null vector", __FILE__, __LINE__, GSL_EFAULT );
728#ifndef GSL_RANGE_CHECK_OFF
730 gsl_error(
"matrix size2 and tda do not match", __FILE__, __LINE__, GSL_EBADLEN );
741 if(
ccgsl_pointer == 0 ) gsl_error(
"null vector", __FILE__, __LINE__, GSL_EFAULT );
742#ifndef GSL_RANGE_CHECK_OFF
744 gsl_error(
"matrix size2 and tda do not match", __FILE__, __LINE__, GSL_EBADLEN );
763 gsl_matrix_tricpy( Uplo, Diag,
get(), src.
get() );
772 gsl_matrix_transpose_tricpy( Uplo, Diag,
get(), src.
get() );
784 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
785 *m = gsl_matrix_submatrix(
get(), i, j, n1, n2 ).
matrix;
794 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
804 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
805 *w = gsl_matrix_column(
get(), j ).
vector;
813 diagonal(){ gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
814 *w = gsl_matrix_diagonal(
get() ).
vector;
823 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
824 *w = gsl_matrix_subdiagonal(
get(), k ).
vector;
833 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
834 *w = gsl_matrix_superdiagonal(
get(), k ).
vector;
845 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
846 *w = gsl_matrix_subrow(
get(), i, offset,
n ).
vector;
857 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
858 *w = gsl_matrix_subcolumn(
get(), j, offset,
n ).
vector;
869 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
870 *m = gsl_matrix_view_array( base, n1, n2 ).
matrix;
882 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
883 *m = gsl_matrix_view_array_with_tda( base, n1, n2, tda ).
matrix;
894 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
895 *m = gsl_matrix_view_vector( v.
get(), n1, n2 ).
matrix;
907 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
908 *m = gsl_matrix_view_vector_with_tda( v.
get(), n1, n2, tda ).
matrix;
921 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
922 *m = gsl_matrix_const_submatrix(
get(), i, j, n1, n2 ).
matrix;
931 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
932 *w = gsl_matrix_const_row(
get(), i ).
vector;
941 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
942 *w = gsl_matrix_const_column(
get(), j ).
vector;
950 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
951 *w = gsl_matrix_const_diagonal(
get() ).
vector;
960 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
961 *w = gsl_matrix_const_subdiagonal(
get(), k ).
vector;
970 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
971 *w = gsl_matrix_const_superdiagonal(
get(), k ).
vector;
982 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
983 *w = gsl_matrix_const_subrow(
get(), i, offset,
n ).
vector;
994 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
995 *w = gsl_matrix_const_subcolumn(
get(), j, offset,
n ).
vector;
1006 matrix const submatrix(
size_t const i,
size_t const j,
size_t const n1,
size_t const n2 )
const {
1007 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1008 *m = gsl_matrix_const_submatrix(
get(), i, j, n1, n2 ).
matrix;
1017 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
1018 *w = gsl_matrix_const_row(
get(), i ).
vector;
1027 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
1028 *w = gsl_matrix_const_column(
get(), j ).
vector;
1036 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
1037 *w = gsl_matrix_const_diagonal(
get() ).
vector;
1046 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
1047 *w = gsl_matrix_const_subdiagonal(
get(), k ).
vector;
1056 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
1057 *w = gsl_matrix_const_superdiagonal(
get(), k ).
vector;
1067 vector const subrow(
size_t const i,
size_t const offset,
size_t const n )
const {
1068 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
1069 *w = gsl_matrix_const_subrow(
get(), i, offset,
n ).
vector;
1080 gsl_vector* w =
static_cast<gsl_vector*
>( malloc(
sizeof( gsl_vector ) ) );
1081 *w = gsl_matrix_const_subcolumn(
get(), j, offset,
n ).
vector;
1092 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1093 *m = gsl_matrix_const_view_array( base, n1, n2 ).
matrix;
1106 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1107 *m = gsl_matrix_const_view_array_with_tda( base, n1, n2, tda ).
matrix;
1118 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1119 *m = gsl_matrix_const_view_vector( v.
get(), n1, n2 ).
matrix;
1132 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1133 *m = gsl_matrix_const_view_vector_with_tda( v.
get(), n1, n2, tda ).
matrix;
1144 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1145 *m = gsl_matrix_const_view_array( base, n1, n2 ).
matrix;
1158 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1159 *m = gsl_matrix_const_view_array_with_tda( base, n1, n2, tda ).
matrix;
1170 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1171 *m = gsl_matrix_const_view_vector( v.
get(), n1, n2 ).
matrix;
1184 gsl_matrix* m =
static_cast<gsl_matrix*
>( malloc(
sizeof( gsl_matrix ) ) );
1185 *m = gsl_matrix_const_view_vector_with_tda( v.
get(), n1, n2, tda ).
matrix;
1229#ifdef __GXX_EXPERIMENTAL_CXX0X__
1241 static matrix calloc(
size_t const n1,
size_t const n2 ){
return matrix( gsl_matrix_calloc( n1, n2 ) ); }
1261 double max()
const {
return gsl_matrix_max(
get() ); }
1266 double min()
const {
return gsl_matrix_min(
get() ); }
1272 void minmax(
double* min_out,
double* max_out )
const {
1273 gsl_matrix_minmax(
get(), min_out, max_out ); }
1279 void minmax(
double& min_out,
double& max_out )
const {
1280 gsl_matrix_minmax(
get(), &min_out, &max_out ); }
1298 int scale(
double const x ){
return gsl_matrix_scale(
get(), x ); }
1314 int ispos()
const {
return gsl_matrix_ispos(
get() ); }
1319 int isneg()
const {
return gsl_matrix_isneg(
get() ); }
1331 double get(
size_t const i,
size_t const j )
const {
return gsl_matrix_get(
get(), i, j ); }
1338 void set(
size_t const i,
size_t const j,
double x ){ gsl_matrix_set(
get(), i, j, x ); }
1345 double*
ptr(
size_t const i,
size_t const j ){
return gsl_matrix_ptr(
get(), i, j ); }
1352 double const*
const_ptr(
size_t const i,
size_t const j )
const {
1353 return gsl_matrix_const_ptr(
get(), i, j ); }
1359 int fread( FILE* stream ){
return gsl_matrix_fread( stream,
get() ); }
1365 int fwrite( FILE* stream )
const {
return gsl_matrix_fwrite( stream,
get() ); }
1371 int fscanf( FILE* stream ){
return gsl_matrix_fscanf( stream,
get() ); }
1378 int fprintf( FILE* stream,
char const* format )
const {
1379 return gsl_matrix_fprintf( stream,
get(), format ); }
1388 matrix(
block&
b,
size_t const offset,
size_t const n1,
size_t const n2,
size_t const d2 ){
1389 ccgsl_pointer = gsl_matrix_alloc_from_block(
b.get(), offset, n1, n2, d2 );
1391 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
1406 matrix(
matrix& m,
size_t const k1,
size_t const k2,
size_t const n1,
size_t const n2 ){
1409 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
1427 int swap_rows(
size_t const i,
size_t const j ){
return gsl_matrix_swap_rows(
get(), i, j ); }
1435 return gsl_matrix_swap_columns(
get(), i, j ); }
1442 int swap_rowcol(
size_t const i,
size_t const j ){
return gsl_matrix_swap_rowcol(
get(), i, j ); }
1454 return gsl_matrix_transpose_memcpy(
get(), src.
get() ); }
1461 void max_index(
size_t* imax,
size_t* jmax )
const {
1462 gsl_matrix_max_index(
get(), imax, jmax ); }
1468 void min_index(
size_t* imin,
size_t* jmin )
const {
1469 gsl_matrix_min_index(
get(), imin, jmin ); }
1477 void minmax_index(
size_t* imin,
size_t* jmin,
size_t* imax,
size_t* jmax )
const {
1478 gsl_matrix_minmax_index(
get(), imin, jmin, imax, jmax ); }
1486 gsl_matrix_max_index(
get(), &imax, &jmax ); }
1493 gsl_matrix_min_index(
get(), &imin, &jmin ); }
1501 void minmax_index(
size_t& imin,
size_t& jmin,
size_t& imax,
size_t& jmax )
const {
1502 gsl_matrix_minmax_index(
get(), &imin, &jmin, &imax, &jmax ); }
1510 return gsl_matrix_mul_elements(
get(),
b.get() ); }
1518 return gsl_matrix_div_elements(
get(),
b.get() ); }
1524 return gsl_matrix_norm1(
get() ); }
1531 return gsl_matrix_scale_rows(
get(), x.
get() ); }
1538 return gsl_matrix_scale_columns(
get(), x.
get() ); }
1545 return gsl_matrix_add_diagonal(
get(), x ); }
1553 return gsl_matrix_get_row( v.
get(),
get(), i ); }
1561 return gsl_matrix_get_col( v.
get(),
get(), j ); }
1569 return gsl_matrix_set_row(
get(), i, v.
get() ); }
1577 return gsl_matrix_set_col(
get(), j, v.
get() ); }
1592#ifndef GSL_RANGE_CHECK_OFF
1593#ifndef __GXX_EXPERIMENTAL_CXX0X__
1598 gsl_error(
"trying to read beyond last row of matrix",
1600#ifdef __GXX_EXPERIMENTAL_CXX0X__
1638 return gsl_permute_matrix( p.
get(),
get() ); }
1644 return vector ( gsl_vector_alloc_row_from_matrix( m.get(), i ) ); }
1646 return vector ( gsl_vector_alloc_col_from_matrix( m.get(), i ) ); }
This class handles vectors as shared handles.
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
@ GSL_EFAULT
invalid pointer
Class that extends gsl_function_fdf so that it can be constructed from arbitrary function objects.
A class template for the const iterators.
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()
The default constructor.
bool operator!=(iterator_t< reverse_t > const &i) const
Comparison with non-const iterator.
const_iterator_t< reverse_t > & operator--()
The prefix – operator.
const_iterator_t< reverse_t > operator--(int)
The postfix – operator.
const_iterator_t(iterator_t< reverse_t > const &i)
A copy constructor.
const_iterator_t< reverse_t > & operator++()
The prefix ++ operator.
const_iterator_t< reverse_t > & operator=(const_iterator_t< reverse_t > const &i)
We can assign one output iterator from another.
const_iterator_t(matrix const *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.
const_iterator_t< reverse_t > operator++(int)
The postfix ++ operator.
We create a suitable class for iterator types here.
void decrement()
Decrement the iterator.
reference operator*() const
Dereference the pointer.
iterator_base()
The iterator is default constructible.
bool operator==(iterator_base< container, content, reverse_t > const &i) const
The == operator.
vector_ptr pointer
An iterator must have a pointer typea.
bool operator!=(iterator_base< container, content, reverse_t > const &i) const
The != operator.
void increment()
Increment the iterator.
vector value_type
An iterator must have a value type.
container * v
Store a pointer to a matrix we can iterate over: 0 if no matrix.
std::bidirectional_iterator_tag iterator_category
An iterator must have an iterator category.
value_type reference
An iterator must have a reference type.
pointer operator->() const
Dereference the pointer.
iterator_base(container *v, size_t position)
This constructor allows vector to create non-default iterators.
size_t position
Mark position of iterator within matrix.
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()
The default constructor.
iterator_t< reverse_t > & operator++()
The prefix ++ operator.
iterator_t< reverse_t > operator++(int)
The postfix ++ operator.
iterator_t< reverse_t > & operator--()
The prefix – operator.
iterator_t< reverse_t > operator--(int)
The postfix – operator.
bool operator==(const_iterator_t< reverse_t > const &i) const
Comparison with const iterator.
iterator_t(matrix *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.
This class handles matrix objects as shared handles.
int transpose_memcpy(matrix const &src)
C++ version of gsl_matrix_transpose_memcpy().
double * ptr(size_t const i, size_t const j)
C++ version of gsl_matrix_ptr().
int sub(matrix const &b)
C++ version of gsl_matrix_sub().
vector superdiagonal(size_t const k)
C++ version of gsl_matrix_superdiagonal().
int add(matrix const &b)
C++ version of gsl_matrix_add().
vector subdiagonal(size_t const k)
C++ version of gsl_matrix_subdiagonal().
void reset()
Stop sharing ownership of the shared pointer.
static matrix const const_view_array(double const *base, size_t const n1, size_t const n2)
C++ version of gsl_matrix_const_view_array().
~matrix()
The destructor only deletes the pointers if count reaches zero.
matrix(matrix &&v)
Move constructor.
iterator_t< true > reverse_iterator
The reverse_iterator type.
int isneg() const
C++ version of gsl_matrix_isneg().
static matrix const const_view_array_with_tda(double const *base, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_const_view_array_with_tda().
int add_diagonal(double const x)
C++ version of gsl_matrix_add_diagonal().
const_iterator_t< false > const_iterator
The const_iterator type.
iterator end()
Get iterator pointing beyond last vector element.
const_iterator begin() const
Get iterator pointing to first vector element.
matrix & operator=(matrix &&v)
Move operator.
reverse_iterator rend()
Get iterator pointing beyond last vector element.
int transpose()
C++ version of gsl_matrix_transpose().
static matrix view_array_with_tda(double *base, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_view_array_with_tda().
const_iterator end() const
Get iterator pointing beyond last vector element.
static matrix view_array(double *base, size_t const n1, size_t const n2)
C++ version of gsl_matrix_view_array().
double get(size_t const i, size_t const j) const
C++ version of gsl_matrix_get().
double norm1() const
C++ version of gsl_matrix_norm1().
vector const diagonal() const
Another C++ version of gsl_matrix_const_diagonal().
int set_col(size_t const j, vector const &v)
C++ version of gsl_matrix_set_col().
vector subrow(size_t const i, size_t const offset, size_t const n)
C++ version of gsl_matrix_subrow().
double * data()
Give access to the data block.
matrix(gsl_matrix *v)
Could construct from a gsl_matrix.
void minmax_index(size_t &imin, size_t &jmin, size_t &imax, size_t &jmax) const
C++ version of gsl_matrix_minmax_index().
matrix(matrix const &v)
The copy constructor.
static matrix const view_array(double const *base, size_t const n1, size_t const n2)
Another C++ version of gsl_matrix_const_view_array().
int swap_rowcol(size_t const i, size_t const j)
C++ version of gsl_matrix_swap_rowcol().
void swap(matrix &m)
Swap two matrix objects.
vector const column(size_t const j) const
Another C++ version of gsl_matrix_const_column().
int memcpy(matrix const &src)
C++ version of gsl_matrix_memcpy().
static matrix const view_vector(vector const &v, size_t const n1, size_t const n2)
Another C++ version of gsl_matrix_const_view_vector().
int fscanf(FILE *stream)
C++ version of gsl_matrix_fscanf().
const_iterator_t< true > const_reverse_iterator
The const_reverse_t type.
void set_all(double x)
C++ version of gsl_matrix_set_all().
int get_col(vector &v, size_t const j) const
C++ version of gsl_matrix_get_col().
vector const const_subdiagonal(size_t const k) const
C++ version of gsl_matrix_const_subdiagonal().
static matrix view_vector(vector &v, size_t const n1, size_t const n2)
C++ version of gsl_matrix_view_vector().
size_t use_count() const
Find how many matrix objects share this pointer.
void wrap_gsl_matrix_without_ownership(gsl_matrix *v)
This function is intended mainly for internal use.
int fwrite(FILE *stream) const
C++ version of gsl_matrix_fwrite().
vector const row(size_t const i) const
Another C++ version of gsl_matrix_const_row().
double max() const
C++ version of gsl_matrix_max().
int scale(double const x)
C++ version of gsl_matrix_scale().
int permute(permutation &p)
Permute the columns of this by permutation p.
iterator begin()
Get iterator pointing to first vector element.
vector const superdiagonal(size_t const k) const
Another C++ version of gsl_matrix_const_superdiagonal().
int set_row(size_t const i, vector const &v)
C++ version of gsl_matrix_set_row().
static matrix const view_vector_with_tda(vector const &v, size_t const n1, size_t const n2, size_t const tda)
Another C++ version of gsl_matrix_const_view_vector_with_tda().
int swap_columns(size_t const i, size_t const j)
C++ version of gsl_matrix_swap_columns().
vector const operator[](size_t const i) const
This function allows us to use a matrix like an array.
vector operator[](size_t const i)
This function allows us to use a matrix like an array.
static matrix const const_view_vector(vector const &v, size_t const n1, size_t const n2)
C++ version of gsl_matrix_const_view_vector().
int ispos() const
C++ version of gsl_matrix_ispos().
matrix submatrix(size_t const i, size_t const j, size_t const n1, size_t const n2)
C++ version of gsl_matrix_submatrix().
vector subcolumn(size_t const j, size_t const offset, size_t const n)
C++ version of gsl_matrix_subcolumn().
vector const const_row(size_t const i) const
C++ version of gsl_matrix_const_row().
vector const subrow(size_t const i, size_t const offset, size_t const n) const
Another C++ version of gsl_matrix_const_subrow().
vector row(size_t const i)
C++ version of gsl_matrix_row().
void set_identity()
C++ version of gsl_matrix_set_identity().
gsl_matrix * ccgsl_pointer
The shared pointer.
static matrix view_vector_with_tda(vector &v, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_view_vector_with_tda().
static matrix calloc(size_t const n1, size_t const n2)
C++ version of gsl_matrix_calloc().
vector const subcolumn(size_t const j, size_t const offset, size_t const n) const
Another C++ version of gsl_matrix_const_subcolumn().
const_reverse_iterator rend() const
Get iterator pointing beyond last vector element.
matrix(size_t const n1, size_t const n2)
This constructor creates a new matrix with n1 rows and n2 columns.
int div_elements(matrix const &b)
C++ version of gsl_matrix_div_elements().
int scale_rows(vector const &x)
C++ version of gsl_matrix_scale_rows().
double const * const_ptr(size_t const i, size_t const j) const
C++ version of gsl_matrix_const_ptr().
static matrix const const_view_vector_with_tda(vector const &v, size_t const n1, size_t const n2, size_t const tda)
C++ version of gsl_matrix_const_view_vector_with_tda().
void tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix const &src)
Copy the upper or lower triangular part of matrix src to this.
void transpose_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix const &src)
Copy the upper or lower triangular part of matrix src to this.
matrix const submatrix(size_t const i, size_t const j, size_t const n1, size_t const n2) const
Another C++ version of gsl_matrix_const_submatrix().
iterator_t< false > iterator
The iterator type.
size_t * count
The shared reference count.
void minmax(double *min_out, double *max_out) const
C++ version of gsl_matrix_minmax().
matrix & operator=(matrix const &v)
The assignment operator.
matrix(block &b, size_t const offset, size_t const n1, size_t const n2, size_t const d2)
C++ version of gsl_matrix_alloc_from_block().
int fprintf(FILE *stream, char const *format) const
C++ version of gsl_matrix_fprintf().
double min() const
C++ version of gsl_matrix_min().
vector const subdiagonal(size_t const k) const
Another C++ version of gsl_matrix_const_subdiagonal().
vector const const_superdiagonal(size_t const k) const
C++ version of gsl_matrix_const_superdiagonal().
bool unique() const
Find if this is the only object sharing the gsl_matrix.
static matrix const view_array_with_tda(double const *base, size_t const n1, size_t const n2, size_t const tda)
Another C++ version of gsl_matrix_const_view_array_with_tda().
reverse_iterator rbegin()
Get iterator pointing to first vector element.
int fread(FILE *stream)
C++ version of gsl_matrix_fread().
bool owns_data
Used to allow a vector that does not own its data.
void min_index(size_t &imin, size_t &jmin) const
C++ version of gsl_matrix_min_index().
vector diagonal()
C++ version of gsl_matrix_diagonal().
vector const const_subrow(size_t const i, size_t const offset, size_t const n) const
C++ version of gsl_matrix_const_subrow().
int add_constant(double const x)
C++ version of gsl_matrix_add_constant().
vector const const_diagonal() const
C++ version of gsl_matrix_const_diagonal().
const_reverse_iterator rbegin() const
Get iterator pointing to first vector element.
int isnonneg() const
C++ version of gsl_matrix_isnonneg().
size_t size1() const
The number of rows of the matrix.
void set(size_t const i, size_t const j, double x)
C++ version of gsl_matrix_set().
int swap_rows(size_t const i, size_t const j)
C++ version of gsl_matrix_swap_rows().
int scale_columns(vector const &x)
C++ version of gsl_matrix_scale_columns().
gsl_matrix * get()
Get the gsl_matrix.
void minmax(double &min_out, double &max_out) const
C++ version of gsl_matrix_minmax().
size_t size_type
A container must have a size_type.
int get_row(vector &v, size_t const i) const
C++ version of gsl_matrix_get_row().
matrix(std::initializer_list< std::initializer_list< double > > initializer_list)
Could construct from a std::initializer_list in C++11.
int mul_elements(matrix const &b)
C++ version of gsl_matrix_mul_elements().
gsl_matrix const * get() const
Get the gsl_matrix.
size_t size2() const
The number of columns of the matrix.
vector const const_column(size_t const j) const
C++ version of gsl_matrix_const_column().
void max_index(size_t &imax, size_t &jmax) const
C++ version of gsl_matrix_max_index().
void set_zero()
C++ version of gsl_matrix_set_zero().
int isnull() const
C++ version of gsl_matrix_isnull().
matrix(matrix &m, size_t const k1, size_t const k2, size_t const n1, size_t const n2)
C++ version of gsl_matrix_alloc_from_matrix().
matrix()
The default constructor is only really useful for assigning to.
matrix clone() const
The clone function.
matrix const const_submatrix(size_t const i, size_t const j, size_t const n1, size_t const n2) const
C++ version of gsl_matrix_const_submatrix().
double const * data() const
Give access to the data block.
vector const const_subcolumn(size_t const j, size_t const offset, size_t const n) const
C++ version of gsl_matrix_const_subcolumn().
vector column(size_t const j)
C++ version of gsl_matrix_column().
Class that extends gsl_multifit_nlinear_fdf so that it can be constructed from arbitrary function obj...
Class that extends gsl_multilarge_nlinear_fdf so that it can be constructed from arbitrary function o...
Class that extends gsl_multiroot_function_fdf so that it can be constructed from arbitrary function o...
This class handles GSL permutation objects.
gsl_permutation * get()
Get the gsl_permutation.
This class handles vector objects as shared handles.
gsl_vector * get()
Get the gsl_vector.
static vector alloc_row_from_matrix(matrix &m, size_t const i)
C++ version of gsl_vector_alloc_row_from_matrix().
void wrap_gsl_vector_without_ownership(gsl_vector *v)
This function is intended mainly for internal use.
vector()
The default constructor is only really useful for assigning to.
static vector alloc_col_from_matrix(matrix &m, size_t const j)
C++ version of gsl_vector_alloc_col_from_matrix().
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
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_ptr(vector const &v)
Typically we have to construct from a vector.
vector * operator->()
Dereference operator.
vector & operator*()
Dereference operator.