ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
gsl::matrix Class Reference

This class handles matrix objects as shared handles. More...

#include <matrix.hpp>

Collaboration diagram for gsl::matrix:
Collaboration graph

Classes

class  const_iterator_t
 A class template for the const iterators. More...
 
class  iterator_base
 We create a suitable class for iterator types here. More...
 
class  iterator_t
 A class template for the two non-const iterators. More...
 
struct  vector_ptr
 This is a pointer-like type for iterator return values. More...
 

Public Types

typedef const_iterator_t< false > const_iterator
 The const_iterator type. More...
 
typedef iterator_t< false > iterator
 The iterator type. More...
 
typedef const_iterator_t< true > const_reverse_iterator
 The const_reverse_t type. More...
 
typedef iterator_t< true > reverse_iterator
 The reverse_iterator type. More...
 
typedef size_t size_type
 A container must have a size_type. More...
 

Public Member Functions

 matrix ()
 The default constructor is only really useful for assigning to. More...
 
 matrix (size_t const n1, size_t const n2)
 This constructor creates a new matrix with n1 rows and n2 columns. More...
 
 matrix (gsl_matrix *v)
 Could construct from a gsl_matrix. More...
 
 matrix (std::initializer_list< std::initializer_list< double > > initializer_list)
 Could construct from a std::initializer_list in C++11. More...
 
 matrix (matrix const &v)
 The copy constructor. More...
 
matrixoperator= (matrix const &v)
 The assignment operator. More...
 
matrix clone () const
 The clone function. More...
 
 ~matrix ()
 The destructor only deletes the pointers if count reaches zero. More...
 
void wrap_gsl_matrix_without_ownership (gsl_matrix *v)
 This function is intended mainly for internal use. More...
 
void reset ()
 Stop sharing ownership of the shared pointer. More...
 
 matrix (matrix &&v)
 Move constructor. More...
 
matrixoperator= (matrix &&v)
 Move operator. More...
 
iterator begin ()
 Get iterator pointing to first vector element. More...
 
const_iterator begin () const
 Get iterator pointing to first vector element. More...
 
iterator end ()
 Get iterator pointing beyond last vector element. More...
 
const_iterator end () const
 Get iterator pointing beyond last vector element. More...
 
reverse_iterator rbegin ()
 Get iterator pointing to first vector element. More...
 
const_reverse_iterator rbegin () const
 Get iterator pointing to first vector element. More...
 
reverse_iterator rend ()
 Get iterator pointing beyond last vector element. More...
 
const_reverse_iterator rend () const
 Get iterator pointing beyond last vector element. More...
 
size_t size1 () const
 The number of rows of the matrix. More...
 
size_t size2 () const
 The number of columns of the matrix. More...
 
double * data ()
 Give access to the data block. More...
 
double const * data () const
 Give access to the data block. More...
 
void swap (matrix &m)
 Swap two matrix objects. More...
 
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. More...
 
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. More...
 
matrix submatrix (size_t const i, size_t const j, size_t const n1, size_t const n2)
 C++ version of gsl_matrix_submatrix(). More...
 
vector row (size_t const i)
 C++ version of gsl_matrix_row(). More...
 
vector column (size_t const j)
 C++ version of gsl_matrix_column(). More...
 
vector diagonal ()
 C++ version of gsl_matrix_diagonal(). More...
 
vector subdiagonal (size_t const k)
 C++ version of gsl_matrix_subdiagonal(). More...
 
vector superdiagonal (size_t const k)
 C++ version of gsl_matrix_superdiagonal(). More...
 
vector subrow (size_t const i, size_t const offset, size_t const n)
 C++ version of gsl_matrix_subrow(). More...
 
vector subcolumn (size_t const j, size_t const offset, size_t const n)
 C++ version of gsl_matrix_subcolumn(). More...
 
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(). More...
 
vector const const_row (size_t const i) const
 C++ version of gsl_matrix_const_row(). More...
 
vector const const_column (size_t const j) const
 C++ version of gsl_matrix_const_column(). More...
 
vector const const_diagonal () const
 C++ version of gsl_matrix_const_diagonal(). More...
 
vector const const_subdiagonal (size_t const k) const
 C++ version of gsl_matrix_const_subdiagonal(). More...
 
vector const const_superdiagonal (size_t const k) const
 C++ version of gsl_matrix_const_superdiagonal(). More...
 
vector const const_subrow (size_t const i, size_t const offset, size_t const n) const
 C++ version of gsl_matrix_const_subrow(). More...
 
vector const const_subcolumn (size_t const j, size_t const offset, size_t const n) const
 C++ version of gsl_matrix_const_subcolumn(). More...
 
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(). More...
 
vector const row (size_t const i) const
 Another C++ version of gsl_matrix_const_row(). More...
 
vector const column (size_t const j) const
 Another C++ version of gsl_matrix_const_column(). More...
 
vector const diagonal () const
 Another C++ version of gsl_matrix_const_diagonal(). More...
 
vector const subdiagonal (size_t const k) const
 Another C++ version of gsl_matrix_const_subdiagonal(). More...
 
vector const superdiagonal (size_t const k) const
 Another C++ version of gsl_matrix_const_superdiagonal(). More...
 
vector const subrow (size_t const i, size_t const offset, size_t const n) const
 Another C++ version of gsl_matrix_const_subrow(). More...
 
vector const subcolumn (size_t const j, size_t const offset, size_t const n) const
 Another C++ version of gsl_matrix_const_subcolumn(). More...
 
gsl_matrix * get ()
 Get the gsl_matrix. More...
 
gsl_matrix const * get () const
 Get the gsl_matrix. More...
 
bool unique () const
 Find if this is the only object sharing the gsl_matrix. More...
 
size_t use_count () const
 Find how many matrix objects share this pointer. More...
 
 operator bool () const
 Allow conversion to bool. More...
 
void set_zero ()
 C++ version of gsl_matrix_set_zero(). More...
 
void set_all (double x)
 C++ version of gsl_matrix_set_all(). More...
 
int memcpy (matrix const &src)
 C++ version of gsl_matrix_memcpy(). More...
 
double max () const
 C++ version of gsl_matrix_max(). More...
 
double min () const
 C++ version of gsl_matrix_min(). More...
 
void minmax (double *min_out, double *max_out) const
 C++ version of gsl_matrix_minmax(). More...
 
void minmax (double &min_out, double &max_out) const
 C++ version of gsl_matrix_minmax(). More...
 
int add (matrix const &b)
 C++ version of gsl_matrix_add(). More...
 
int sub (matrix const &b)
 C++ version of gsl_matrix_sub(). More...
 
int scale (double const x)
 C++ version of gsl_matrix_scale(). More...
 
int add_constant (double const x)
 C++ version of gsl_matrix_add_constant(). More...
 
int isnull () const
 C++ version of gsl_matrix_isnull(). More...
 
int ispos () const
 C++ version of gsl_matrix_ispos(). More...
 
int isneg () const
 C++ version of gsl_matrix_isneg(). More...
 
int isnonneg () const
 C++ version of gsl_matrix_isnonneg(). More...
 
double get (size_t const i, size_t const j) const
 C++ version of gsl_matrix_get(). More...
 
void set (size_t const i, size_t const j, double x)
 C++ version of gsl_matrix_set(). More...
 
double * ptr (size_t const i, size_t const j)
 C++ version of gsl_matrix_ptr(). More...
 
double const * const_ptr (size_t const i, size_t const j) const
 C++ version of gsl_matrix_const_ptr(). More...
 
int fread (FILE *stream)
 C++ version of gsl_matrix_fread(). More...
 
int fwrite (FILE *stream) const
 C++ version of gsl_matrix_fwrite(). More...
 
int fscanf (FILE *stream)
 C++ version of gsl_matrix_fscanf(). More...
 
int fprintf (FILE *stream, char const *format) const
 C++ version of gsl_matrix_fprintf(). More...
 
 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(). More...
 
 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(). More...
 
void set_identity ()
 C++ version of gsl_matrix_set_identity(). More...
 
int swap_rows (size_t const i, size_t const j)
 C++ version of gsl_matrix_swap_rows(). More...
 
int swap_columns (size_t const i, size_t const j)
 C++ version of gsl_matrix_swap_columns(). More...
 
int swap_rowcol (size_t const i, size_t const j)
 C++ version of gsl_matrix_swap_rowcol(). More...
 
int transpose ()
 C++ version of gsl_matrix_transpose(). More...
 
int transpose_memcpy (matrix const &src)
 C++ version of gsl_matrix_transpose_memcpy(). More...
 
void max_index (size_t &imax, size_t &jmax) const
 C++ version of gsl_matrix_max_index(). More...
 
void min_index (size_t &imin, size_t &jmin) const
 C++ version of gsl_matrix_min_index(). More...
 
void minmax_index (size_t &imin, size_t &jmin, size_t &imax, size_t &jmax) const
 C++ version of gsl_matrix_minmax_index(). More...
 
int mul_elements (matrix const &b)
 C++ version of gsl_matrix_mul_elements(). More...
 
int div_elements (matrix const &b)
 C++ version of gsl_matrix_div_elements(). More...
 
double norm1 () const
 C++ version of gsl_matrix_norm1(). More...
 
int scale_rows (vector const &x)
 C++ version of gsl_matrix_scale_rows(). More...
 
int scale_columns (vector const &x)
 C++ version of gsl_matrix_scale_columns(). More...
 
int add_diagonal (double const x)
 C++ version of gsl_matrix_add_diagonal(). More...
 
int get_row (vector &v, size_t const i) const
 C++ version of gsl_matrix_get_row(). More...
 
int get_col (vector &v, size_t const j) const
 C++ version of gsl_matrix_get_col(). More...
 
int set_row (size_t const i, vector const &v)
 C++ version of gsl_matrix_set_row(). More...
 
int set_col (size_t const j, vector const &v)
 C++ version of gsl_matrix_set_col(). More...
 
vector operator[] (size_t const i)
 This function allows us to use a matrix like an array. More...
 
vector const operator[] (size_t const i) const
 This function allows us to use a matrix like an array. More...
 
int permute (permutation &p)
 Permute the columns of this by permutation p. More...
 

Static Public Member Functions

static matrix view_array (double *base, size_t const n1, size_t const n2)
 C++ version of gsl_matrix_view_array(). More...
 
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(). More...
 
static matrix view_vector (vector &v, size_t const n1, size_t const n2)
 C++ version of gsl_matrix_view_vector(). More...
 
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(). More...
 
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(). More...
 
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(). More...
 
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(). More...
 
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(). More...
 
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(). More...
 
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(). More...
 
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(). More...
 
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(). More...
 
static matrix calloc (size_t const n1, size_t const n2)
 C++ version of gsl_matrix_calloc(). More...
 

Private Attributes

bool owns_data
 Used to allow a vector that does not own its data. More...
 
gsl_matrix * ccgsl_pointer
 The shared pointer. More...
 
size_t * count
 The shared reference count. More...
 

Friends

class multiroot::function_fdf
 
class multifit::nlinear::function_fdf
 
class multilarge::nlinear::function_fdf
 

Detailed Description

This class handles matrix objects as shared handles.

It models a random access container so that STL functions work with matrix.

Note that matrix_views are implemented as matrix objects here.

Note that in C++11 it is possible to iterate over the rows of a matrix using

for( auto row : matrix ){ ... }
This class handles matrix objects as shared handles.
Definition: matrix.hpp:72
vector row(size_t const i)
C++ version of gsl_matrix_row().
Definition: matrix.hpp:793

Otherwise,

for( matrix::iterator i = matrix.begin(); i != matrix.end(); ++i ){
vector row = *i;
...
}
A class template for the two non-const iterators.
Definition: matrix.hpp:426
iterator end()
Get iterator pointing beyond last vector element.
Definition: matrix.hpp:663
iterator begin()
Get iterator pointing to first vector element.
Definition: matrix.hpp:648
This class handles vector objects as shared handles.
Definition: vector.hpp:74

will achieve the same effect. But if the element pointed to by i is used more than once, do not use

i-> ...

because each call of operator->() creates a new and different vector.

Definition at line 72 of file matrix.hpp.

Member Typedef Documentation

◆ const_iterator

The const_iterator type.

Definition at line 626 of file matrix.hpp.

◆ const_reverse_iterator

The const_reverse_t type.

Definition at line 634 of file matrix.hpp.

◆ iterator

The iterator type.

Definition at line 630 of file matrix.hpp.

◆ reverse_iterator

The reverse_iterator type.

Definition at line 638 of file matrix.hpp.

◆ size_type

typedef size_t gsl::matrix::size_type

A container must have a size_type.

Definition at line 642 of file matrix.hpp.

Constructor & Destructor Documentation

◆ matrix() [1/8]

gsl::matrix::matrix ( )
inline

◆ matrix() [2/8]

gsl::matrix::matrix ( size_t const  n1,
size_t const  n2 
)
inlineexplicit

This constructor creates a new matrix with n1 rows and n2 columns.

Parameters
n1The number of rows in the matrix
n2The number of columns in the matrix

Definition at line 92 of file matrix.hpp.

References ccgsl_pointer, and count.

◆ matrix() [3/8]

gsl::matrix::matrix ( gsl_matrix *  v)
inlineexplicit

Could construct from a gsl_matrix.

This is not usually a good idea. In this case you should not use gsl_matrix_free() to deallocate the memory.

Parameters
vThe matrix

Definition at line 110 of file matrix.hpp.

References ccgsl_pointer, and count.

◆ matrix() [4/8]

gsl::matrix::matrix ( std::initializer_list< std::initializer_list< double > >  initializer_list)
inline

Could construct from a std::initializer_list in C++11.

Parameters
initializer_listThe initializer_list.

Definition at line 121 of file matrix.hpp.

References ccgsl_pointer, count, gsl::exception::GSL_EBADLEN, row(), and set().

◆ matrix() [5/8]

gsl::matrix::matrix ( matrix const &  v)
inline

The copy constructor.

This shares the matrix. Use clone() if you want a full copy.

Parameters
vThe matrix to copy.

Definition at line 154 of file matrix.hpp.

References count.

◆ ~matrix()

gsl::matrix::~matrix ( )
inline

The destructor only deletes the pointers if count reaches zero.

Definition at line 194 of file matrix.hpp.

References ccgsl_pointer, count, and owns_data.

◆ matrix() [6/8]

gsl::matrix::matrix ( matrix &&  v)
inline

Move constructor.

Parameters
vThe matrix to move.

Definition at line 236 of file matrix.hpp.

References count.

◆ matrix() [7/8]

gsl::matrix::matrix ( block b,
size_t const  offset,
size_t const  n1,
size_t const  n2,
size_t const  d2 
)
inline

C++ version of gsl_matrix_alloc_from_block().

Parameters
bThe block
offsetThe offset within the block
n1The number of rows in the matrix
n2The number of columns in the matrix
d2undocumented

Definition at line 1388 of file matrix.hpp.

References gsl::sf::mathieu::b(), ccgsl_pointer, and count.

◆ matrix() [8/8]

gsl::matrix::matrix ( matrix m,
size_t const  k1,
size_t const  k2,
size_t const  n1,
size_t const  n2 
)
inline

C++ version of gsl_matrix_alloc_from_matrix().

Parameters
mThe matrix
k1the row of m to take as row zero
k2the column of m to take as column zero
n1The number of rows in the matrix
n2The number of columns in the matrix

Definition at line 1406 of file matrix.hpp.

References ccgsl_pointer, count, and get().

Member Function Documentation

◆ add()

int gsl::matrix::add ( matrix const &  b)
inline

C++ version of gsl_matrix_add().

Parameters
bmatrix to add to this
Returns
error code on failure

Definition at line 1286 of file matrix.hpp.

References gsl::sf::mathieu::b(), and get().

◆ add_constant()

int gsl::matrix::add_constant ( double const  x)
inline

C++ version of gsl_matrix_add_constant().

Parameters
xconstant to add to each element of this
Returns
error code on failure

Definition at line 1304 of file matrix.hpp.

References get().

◆ add_diagonal()

int gsl::matrix::add_diagonal ( double const  x)
inline

C++ version of gsl_matrix_add_diagonal().

Parameters
xA constant
Returns
error code on failure

Definition at line 1544 of file matrix.hpp.

References get().

◆ begin() [1/2]

iterator gsl::matrix::begin ( )
inline

Get iterator pointing to first vector element.

Returns
iterator pointing to first vector element

Definition at line 648 of file matrix.hpp.

◆ begin() [2/2]

const_iterator gsl::matrix::begin ( ) const
inline

Get iterator pointing to first vector element.

Returns
iterator pointing to first vector element

Definition at line 655 of file matrix.hpp.

◆ calloc()

static matrix gsl::matrix::calloc ( size_t const  n1,
size_t const  n2 
)
inlinestatic

C++ version of gsl_matrix_calloc().

This constructs a matrix object with entries initialised to zero.

Parameters
n1The number of rows in the matrix
n2The number of columns in the matrix
Returns
A matrix initialised to zero

Definition at line 1241 of file matrix.hpp.

References matrix().

◆ clone()

matrix gsl::matrix::clone ( ) const
inline

The clone function.

Use this if you want a copy of the block that does not share the underlying data.

Returns
a new copy of this.

Definition at line 183 of file matrix.hpp.

References get(), size1(), and size2().

◆ column() [1/2]

vector gsl::matrix::column ( size_t const  j)
inline

C++ version of gsl_matrix_column().

Parameters
jA column index
Returns
A column as a vector

Definition at line 803 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ column() [2/2]

vector const gsl::matrix::column ( size_t const  j) const
inline

Another C++ version of gsl_matrix_const_column().

Parameters
jA column index
Returns
A column as a vector

Definition at line 1026 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ const_column()

vector const gsl::matrix::const_column ( size_t const  j) const
inline

C++ version of gsl_matrix_const_column().

Parameters
jA column index
Returns
A column as a vector

Definition at line 940 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ const_diagonal()

vector const gsl::matrix::const_diagonal ( ) const
inline

C++ version of gsl_matrix_const_diagonal().

Returns
The principal diagonal as a vector

Definition at line 949 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ const_ptr()

double const * gsl::matrix::const_ptr ( size_t const  i,
size_t const  j 
) const
inline

C++ version of gsl_matrix_const_ptr().

Parameters
iindex of row
jindex of column
Returns
pointer to element

Definition at line 1352 of file matrix.hpp.

References get().

◆ const_row()

vector const gsl::matrix::const_row ( size_t const  i) const
inline

C++ version of gsl_matrix_const_row().

Parameters
iA row index
Returns
A row as a vector

Definition at line 930 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ const_subcolumn()

vector const gsl::matrix::const_subcolumn ( size_t const  j,
size_t const  offset,
size_t const  n 
) const
inline

C++ version of gsl_matrix_const_subcolumn().

Parameters
jA column index
offsetA row offset
nThe number of elements
Returns
A subcolumn as a vector

Definition at line 993 of file matrix.hpp.

References get(), gsl::rstat::n(), and gsl::vector::vector().

◆ const_subdiagonal()

vector const gsl::matrix::const_subdiagonal ( size_t const  k) const
inline

C++ version of gsl_matrix_const_subdiagonal().

Parameters
kAn index
Returns
Subdiagonal k as a vector

Definition at line 959 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ const_submatrix()

matrix const gsl::matrix::const_submatrix ( size_t const  i,
size_t const  j,
size_t const  n1,
size_t const  n2 
) const
inline

C++ version of gsl_matrix_const_submatrix().

Parameters
iIndex in this of first row of submatrix
jIndex in this of first column of submatrix
n1Number of rows of submatrix
n2Number of columns of submatrix
Returns
The submatrix

Definition at line 920 of file matrix.hpp.

References get(), and matrix().

◆ const_subrow()

vector const gsl::matrix::const_subrow ( size_t const  i,
size_t const  offset,
size_t const  n 
) const
inline

C++ version of gsl_matrix_const_subrow().

Parameters
iA row index
offsetA column offset
nThe number of elements
Returns
A subrow as a vector

Definition at line 981 of file matrix.hpp.

References get(), gsl::rstat::n(), and gsl::vector::vector().

◆ const_superdiagonal()

vector const gsl::matrix::const_superdiagonal ( size_t const  k) const
inline

C++ version of gsl_matrix_const_superdiagonal().

Parameters
kAn index
Returns
Subdiagonal k as a vector

Definition at line 969 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ const_view_array()

static matrix const gsl::matrix::const_view_array ( double const *  base,
size_t const  n1,
size_t const  n2 
)
inlinestatic

C++ version of gsl_matrix_const_view_array().

Parameters
baseAn array of type double
n1The number of rows
n2The number of columns
Returns
A matrix

Definition at line 1091 of file matrix.hpp.

References matrix().

◆ const_view_array_with_tda()

static matrix const gsl::matrix::const_view_array_with_tda ( double const *  base,
size_t const  n1,
size_t const  n2,
size_t const  tda 
)
inlinestatic

C++ version of gsl_matrix_const_view_array_with_tda().

Parameters
baseAn array of type double
n1The number of rows
n2The number of columns
tdaThe number of columns in memory
Returns
A matrix

Definition at line 1105 of file matrix.hpp.

References matrix().

◆ const_view_vector()

static matrix const gsl::matrix::const_view_vector ( vector const &  v,
size_t const  n1,
size_t const  n2 
)
inlinestatic

C++ version of gsl_matrix_const_view_vector().

Parameters
vA vector
n1The number of rows
n2The number of columns
Returns
A matrix

Definition at line 1117 of file matrix.hpp.

References gsl::vector::get(), and matrix().

◆ const_view_vector_with_tda()

static matrix const gsl::matrix::const_view_vector_with_tda ( vector const &  v,
size_t const  n1,
size_t const  n2,
size_t const  tda 
)
inlinestatic

C++ version of gsl_matrix_const_view_vector_with_tda().

Parameters
vA vector
n1The number of rows
n2The number of columns
tdaThe number of columns in memory
Returns
A matrix

Definition at line 1131 of file matrix.hpp.

References gsl::vector::get(), and matrix().

◆ data() [1/2]

double * gsl::matrix::data ( )
inline

Give access to the data block.

The data() and size() functions mimic the functions of std::array<T> and std::vector<T>. This function can throw an exception or produce a GSL error if the matrix size2 and tda are not equal. However, new matrix objects are always initialised to have tda and size2 equal.

Returns
The data block.

Definition at line 726 of file matrix.hpp.

References ccgsl_pointer.

◆ data() [2/2]

double const * gsl::matrix::data ( ) const
inline

Give access to the data block.

The data() and size() functions mimic the functions of std::array<T> and std::vector<T>. This function can throw an exception or produce a GSL error if the matrix size1 and tda are not equal. However, new matrix objects are always initialised to have tda and size2 equal.

Returns
The data block.

Definition at line 740 of file matrix.hpp.

References ccgsl_pointer.

◆ diagonal() [1/2]

vector gsl::matrix::diagonal ( )
inline

C++ version of gsl_matrix_diagonal().

Returns
The principal diagonal as a vector

Definition at line 813 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ diagonal() [2/2]

vector const gsl::matrix::diagonal ( ) const
inline

Another C++ version of gsl_matrix_const_diagonal().

Returns
The principal diagonal as a vector

Definition at line 1035 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ div_elements()

int gsl::matrix::div_elements ( matrix const &  b)
inline

C++ version of gsl_matrix_div_elements().

Divide each element of this by the corrsponding element of b

Parameters
bAnother matrix
Returns
error code on failure

Definition at line 1517 of file matrix.hpp.

References gsl::sf::mathieu::b(), and get().

◆ end() [1/2]

iterator gsl::matrix::end ( )
inline

Get iterator pointing beyond last vector element.

Returns
iterator pointing beyond last vector element

Definition at line 663 of file matrix.hpp.

References ccgsl_pointer, and size1().

◆ end() [2/2]

const_iterator gsl::matrix::end ( ) const
inline

Get iterator pointing beyond last vector element.

Returns
iterator pointing beyond last vector element

Definition at line 671 of file matrix.hpp.

References ccgsl_pointer, and size1().

◆ fprintf()

int gsl::matrix::fprintf ( FILE *  stream,
char const *  format 
) const
inline

C++ version of gsl_matrix_fprintf().

Parameters
streamA C file stream
formatd, e, f or g
Returns
error code on failure

Definition at line 1378 of file matrix.hpp.

References get().

◆ fread()

int gsl::matrix::fread ( FILE *  stream)
inline

C++ version of gsl_matrix_fread().

Parameters
streamA C file stream
Returns
error code on failure

Definition at line 1359 of file matrix.hpp.

References get().

◆ fscanf()

int gsl::matrix::fscanf ( FILE *  stream)
inline

C++ version of gsl_matrix_fscanf().

Parameters
streamA C file stream
Returns
error code on failure

Definition at line 1371 of file matrix.hpp.

References get().

◆ fwrite()

int gsl::matrix::fwrite ( FILE *  stream) const
inline

C++ version of gsl_matrix_fwrite().

Parameters
streamA C file stream
Returns
error code on failure

Definition at line 1365 of file matrix.hpp.

References get().

◆ get() [1/3]

gsl_matrix * gsl::matrix::get ( )
inline

Get the gsl_matrix.

Returns
the gsl_matrix

Definition at line 1207 of file matrix.hpp.

References ccgsl_pointer.

Referenced by gsl::multilarge::linear::accumulate(), add(), add_constant(), add_diagonal(), gsl::linalg::balance_accum(), gsl::linalg::balance_columns(), gsl::linalg::balance_matrix(), gsl::linalg::bidiag_decomp(), gsl::linalg::bidiag_unpack(), gsl::linalg::bidiag_unpack2(), gsl::linalg::bidiag_unpack_B(), gsl::linalg::cholesky_band_decomp(), gsl::linalg::cholesky_band_invert(), gsl::linalg::cholesky_band_rcond(), gsl::linalg::cholesky_band_solve(), gsl::linalg::cholesky_band_svx(), gsl::linalg::cholesky_band_unpack(), gsl::linalg::cholesky_decomp(), gsl::linalg::cholesky_decomp1(), gsl::linalg::cholesky_decomp2(), gsl::linalg::cholesky_decomp_unit(), gsl::linalg::cholesky_invert(), gsl::linalg::cholesky_rcond(), gsl::linalg::cholesky_scale(), gsl::linalg::cholesky_scale_apply(), gsl::linalg::cholesky_solve(), gsl::linalg::cholesky_solve2(), gsl::linalg::cholesky_solve_mat(), gsl::linalg::cholesky_svx(), gsl::linalg::cholesky_svx2(), gsl::linalg::cholesky_svx_mat(), clone(), gsl::linalg::COD_decomp(), gsl::linalg::COD_decomp_e(), gsl::linalg::COD_lssolve(), gsl::linalg::COD_lssolve2(), gsl::linalg::COD_matZ(), gsl::linalg::COD_unpack(), column(), const_column(), const_diagonal(), const_ptr(), const_row(), const_subcolumn(), const_subdiagonal(), const_submatrix(), const_subrow(), const_superdiagonal(), gsl::multifit::nlinear::covar(), gsl::spmatrix::d2sp(), gsl::bspline::deriv_eval(), gsl::bspline::deriv_eval_nonzero(), gsl::multilarge::nlinear::df(), gsl::multifit::nlinear::df(), gsl::blas::dgemm(), gsl::blas::dgemv(), gsl::blas::dger(), diagonal(), div_elements(), gsl::blas::dsymm(), gsl::blas::dsymv(), gsl::blas::dsyr(), gsl::blas::dsyr2(), gsl::blas::dsyr2k(), gsl::blas::dsyrk(), gsl::blas::dtrmm(), gsl::blas::dtrmv(), gsl::blas::dtrsm(), gsl::blas::dtrsv(), gsl::multilarge::nlinear::eval_df(), gsl::multifit::nlinear::eval_fvv(), gsl::multifit::nlinear::fdfvv(), gsl::multilarge::nlinear::fdfvv(), gsl::multiroot::fdjacobian(), fprintf(), fread(), fscanf(), fwrite(), gsl::eigen::gen(), gsl::eigen::gen_QZ(), gsl::multilarge::linear::genform2(), gsl::eigen::gensymm(), gsl::eigen::gensymmv(), gsl::eigen::gensymmv_sort(), gsl::eigen::genv(), gsl::eigen::genv_QZ(), get_col(), get_row(), gsl::linalg::hessenberg_decomp(), gsl::linalg::hessenberg_set_zero(), gsl::linalg::hessenberg_submatrix(), gsl::linalg::hessenberg_unpack(), gsl::linalg::hessenberg_unpack_accum(), gsl::linalg::hesstri_decomp(), gsl::linalg::HH_solve(), gsl::linalg::HH_svx(), gsl::linalg::householder_hm(), gsl::linalg::householder_hm1(), gsl::linalg::householder_left(), gsl::linalg::householder_mh(), gsl::linalg::householder_right(), isneg(), isnonneg(), isnull(), ispos(), gsl::multilarge::linear::L_decomp(), gsl::linalg::L_solve_T(), gsl::linalg::ldlt_band_decomp(), gsl::linalg::ldlt_band_rcond(), gsl::linalg::ldlt_band_solve(), gsl::linalg::ldlt_band_svx(), gsl::linalg::ldlt_band_unpack(), gsl::linalg::ldlt_decomp(), gsl::linalg::ldlt_rcond(), gsl::linalg::ldlt_solve(), gsl::linalg::ldlt_svx(), gsl::multifit::linear(), gsl::multifit::linear_applyW(), gsl::multifit::linear_bsvd(), gsl::multifit::linear_est(), gsl::multifit::linear_genform2(), gsl::multifit::linear_L_decomp(), gsl::multifit::linear_Lk(), gsl::multifit::linear_Lsobolev(), gsl::multifit::linear_residuals(), gsl::multifit::linear_solve(), gsl::multifit::linear_stdform1(), gsl::multifit::linear_stdform2(), gsl::multifit::linear_svd(), gsl::multifit::linear_wgenform2(), gsl::multifit::linear_wstdform1(), gsl::multifit::linear_wstdform2(), gsl::linalg::LQ_decomp(), gsl::linalg::LQ_LQsolve(), gsl::linalg::LQ_Lsolve_T(), gsl::linalg::LQ_lssolve(), gsl::linalg::LQ_lssolve_T(), gsl::linalg::LQ_Lsvx_T(), gsl::linalg::LQ_QTvec(), gsl::linalg::LQ_solve_T(), gsl::linalg::LQ_svx_T(), gsl::linalg::LQ_unpack(), gsl::linalg::LQ_update(), gsl::linalg::LQ_vecQ(), gsl::linalg::LQ_vecQT(), gsl::linalg::LU_band_decomp(), gsl::linalg::LU_band_solve(), gsl::linalg::LU_band_svx(), gsl::linalg::LU_band_unpack(), gsl::linalg::LU_decomp(), gsl::linalg::LU_det(), gsl::linalg::LU_invert(), gsl::linalg::LU_invx(), gsl::linalg::LU_lndet(), gsl::linalg::LU_refine(), gsl::linalg::LU_sgndet(), gsl::linalg::LU_solve(), gsl::linalg::LU_svx(), matrix(), max(), max_index(), gsl::linalg::mcholesky_decomp(), gsl::linalg::mcholesky_invert(), gsl::linalg::mcholesky_rcond(), gsl::linalg::mcholesky_solve(), gsl::linalg::mcholesky_svx(), memcpy(), min(), min_index(), minmax(), minmax_index(), mul_elements(), gsl::eigen::nonsymm(), gsl::eigen::nonsymm_Z(), gsl::eigen::nonsymmv(), gsl::eigen::nonsymmv_Z(), norm1(), gsl::linalg::pcholesky_decomp(), gsl::linalg::pcholesky_decomp2(), gsl::linalg::pcholesky_invert(), gsl::linalg::pcholesky_rcond(), gsl::linalg::pcholesky_solve(), gsl::linalg::pcholesky_solve2(), gsl::linalg::pcholesky_svx(), gsl::linalg::pcholesky_svx2(), permute(), gsl::linalg::PTLQ_decomp(), gsl::linalg::PTLQ_decomp2(), gsl::linalg::PTLQ_LQsolve_T(), gsl::linalg::PTLQ_Lsolve_T(), gsl::linalg::PTLQ_Lsvx_T(), gsl::linalg::PTLQ_solve_T(), gsl::linalg::PTLQ_svx_T(), gsl::linalg::PTLQ_update(), ptr(), gsl::linalg::QL_decomp(), gsl::linalg::QL_unpack(), gsl::linalg::QR_band_decomp_L2(), gsl::linalg::QR_band_unpack_L2(), gsl::linalg::QR_decomp(), gsl::linalg::QR_decomp_old(), gsl::linalg::QR_lssolve(), gsl::linalg::QR_matQ(), gsl::linalg::QR_QRsolve(), gsl::linalg::QR_QTmat(), gsl::linalg::QR_QTmat_r(), gsl::linalg::QR_QTvec(), gsl::linalg::QR_Qvec(), gsl::linalg::QR_rcond(), gsl::linalg::QR_Rsolve(), gsl::linalg::QR_Rsvx(), gsl::linalg::QR_solve(), gsl::linalg::QR_svx(), gsl::linalg::QR_UD_decomp(), gsl::linalg::QR_UD_lssolve(), gsl::linalg::QR_unpack(), gsl::linalg::QR_unpack_r(), gsl::linalg::QR_update(), gsl::linalg::QR_UR_decomp(), gsl::linalg::QR_UU_decomp(), gsl::linalg::QR_UU_lssolve(), gsl::linalg::QR_UU_QTvec(), gsl::linalg::QR_UZ_decomp(), gsl::linalg::QRPT_decomp(), gsl::linalg::QRPT_decomp2(), gsl::linalg::QRPT_lssolve(), gsl::linalg::QRPT_lssolve2(), gsl::linalg::QRPT_QRsolve(), gsl::linalg::QRPT_rank(), gsl::linalg::QRPT_rcond(), gsl::linalg::QRPT_Rsolve(), gsl::linalg::QRPT_Rsvx(), gsl::linalg::QRPT_solve(), gsl::linalg::QRPT_svx(), gsl::linalg::QRPT_update(), gsl::linalg::R_solve(), gsl::linalg::R_svx(), gsl::multifit::robust(), gsl::multifit::robust_est(), gsl::multifit::robust_residuals(), row(), scale(), scale_columns(), scale_rows(), set(), set_all(), set_col(), set_identity(), set_row(), set_zero(), gsl::spmatrix::sp2d(), gsl::multilarge::linear::stdform1(), gsl::multilarge::linear::stdform2(), sub(), subcolumn(), subdiagonal(), submatrix(), subrow(), superdiagonal(), gsl::linalg::SV_decomp(), gsl::linalg::SV_decomp_jacobi(), gsl::linalg::SV_decomp_mod(), gsl::linalg::SV_leverage(), gsl::linalg::SV_solve(), swap_columns(), swap_rowcol(), swap_rows(), gsl::eigen::symm(), gsl::linalg::symmtd_decomp(), gsl::linalg::symmtd_unpack(), gsl::linalg::symmtd_unpack_T(), gsl::eigen::symmv(), gsl::eigen::symmv_sort(), transpose(), transpose_memcpy(), transpose_tricpy(), gsl::linalg::tri_invert(), gsl::linalg::tri_lower_invert(), gsl::linalg::tri_lower_rcond(), gsl::linalg::tri_lower_unit_invert(), gsl::linalg::tri_LTL(), gsl::linalg::tri_rcond(), gsl::linalg::tri_UL(), gsl::linalg::tri_upper_invert(), gsl::linalg::tri_upper_rcond(), gsl::linalg::tri_upper_unit_invert(), tricpy(), gsl::multifit::wlinear(), gsl::multifit::wlinear_svd(), gsl::multifit::wlinear_tsvd(), gsl::multifit::wlinear_usvd(), gsl::multilarge::linear::wstdform1(), and gsl::multilarge::linear::wstdform2().

◆ get() [2/3]

gsl_matrix const * gsl::matrix::get ( ) const
inline

Get the gsl_matrix.

Returns
the gsl_matrix

Definition at line 1212 of file matrix.hpp.

References ccgsl_pointer.

◆ get() [3/3]

double gsl::matrix::get ( size_t const  i,
size_t const  j 
) const
inline

C++ version of gsl_matrix_get().

Parameters
iindex of row
jindex of column
Returns
value of element

Definition at line 1331 of file matrix.hpp.

References get().

Referenced by get().

◆ get_col()

int gsl::matrix::get_col ( vector v,
size_t const  j 
) const
inline

C++ version of gsl_matrix_get_col().

Parameters
vA vector
jThe index of the column
Returns
error code on failure

Definition at line 1560 of file matrix.hpp.

References get(), and gsl::vector::get().

◆ get_row()

int gsl::matrix::get_row ( vector v,
size_t const  i 
) const
inline

C++ version of gsl_matrix_get_row().

Parameters
vA vector
iThe index of the row
Returns
error code on failure

Definition at line 1552 of file matrix.hpp.

References get(), and gsl::vector::get().

◆ isneg()

int gsl::matrix::isneg ( ) const
inline

C++ version of gsl_matrix_isneg().

Returns
+1 or 0 according as elements are all negative or not

Definition at line 1319 of file matrix.hpp.

References get().

◆ isnonneg()

int gsl::matrix::isnonneg ( ) const
inline

C++ version of gsl_matrix_isnonneg().

Returns
+1 or 0 according as elements are all nonnegative or not

Definition at line 1324 of file matrix.hpp.

References get().

◆ isnull()

int gsl::matrix::isnull ( ) const
inline

C++ version of gsl_matrix_isnull().

Returns
+1 or 0 according as elements are all zero or not

Definition at line 1309 of file matrix.hpp.

References get().

◆ ispos()

int gsl::matrix::ispos ( ) const
inline

C++ version of gsl_matrix_ispos().

Returns
+1 or 0 according as elements are all positive or not

Definition at line 1314 of file matrix.hpp.

References get().

◆ max()

double gsl::matrix::max ( ) const
inline

C++ version of gsl_matrix_max().

Returns
maximum element of matrix

Definition at line 1261 of file matrix.hpp.

References get().

◆ max_index()

void gsl::matrix::max_index ( size_t &  imax,
size_t &  jmax 
) const
inline

C++ version of gsl_matrix_max_index().

Parameters
imaxrow index of the first maximum element in the matrix
jmaxcolumn index of the first maximum element in the matrix

Definition at line 1485 of file matrix.hpp.

References get().

◆ memcpy()

int gsl::matrix::memcpy ( matrix const &  src)
inline

C++ version of gsl_matrix_memcpy().

Parameters
srcsource matrix
Returns
error code on failure

Definition at line 1256 of file matrix.hpp.

References get().

◆ min()

double gsl::matrix::min ( ) const
inline

C++ version of gsl_matrix_min().

Returns
minimum element of matrix

Definition at line 1266 of file matrix.hpp.

References get().

◆ min_index()

void gsl::matrix::min_index ( size_t &  imin,
size_t &  jmin 
) const
inline

C++ version of gsl_matrix_min_index().

Parameters
iminrow index of the first minimum element in the matrix
jmincolumn index of the first minimum element in the matrix

Definition at line 1492 of file matrix.hpp.

References get().

◆ minmax() [1/2]

void gsl::matrix::minmax ( double &  min_out,
double &  max_out 
) const
inline

C++ version of gsl_matrix_minmax().

Parameters
min_outminimum element of matrix
max_outmaximum element of matrix

Definition at line 1279 of file matrix.hpp.

References get().

◆ minmax() [2/2]

void gsl::matrix::minmax ( double *  min_out,
double *  max_out 
) const
inline

C++ version of gsl_matrix_minmax().

Parameters
min_outminimum element of matrix
max_outmaximum element of matrix

Definition at line 1272 of file matrix.hpp.

References get().

◆ minmax_index()

void gsl::matrix::minmax_index ( size_t &  imin,
size_t &  jmin,
size_t &  imax,
size_t &  jmax 
) const
inline

C++ version of gsl_matrix_minmax_index().

Parameters
iminrow index of the first minimum element in the matrix
jmincolumn index of the first minimum element in the matrix
imaxrow index of the first maximum element in the matrix
jmaxcolumn index of the first maximum element in the matrix

Definition at line 1501 of file matrix.hpp.

References get().

◆ mul_elements()

int gsl::matrix::mul_elements ( matrix const &  b)
inline

C++ version of gsl_matrix_mul_elements().

Multiply matrices elementwise.

Parameters
bAnother matrix
Returns
error code on failure

Definition at line 1509 of file matrix.hpp.

References gsl::sf::mathieu::b(), and get().

◆ norm1()

double gsl::matrix::norm1 ( ) const
inline

C++ version of gsl_matrix_norm1().

Returns
1-norm of matrix

Definition at line 1523 of file matrix.hpp.

References get().

◆ operator bool()

gsl::matrix::operator bool ( ) const
inlineexplicit

Allow conversion to bool.

Returns
true or false according as this contains a pointer to a gsl_matrix

Definition at line 1232 of file matrix.hpp.

References ccgsl_pointer.

◆ operator=() [1/2]

matrix & gsl::matrix::operator= ( matrix &&  v)
inline

Move operator.

Parameters
vThe matrix to move.
Returns
A reference to this.

Definition at line 246 of file matrix.hpp.

References matrix().

◆ operator=() [2/2]

matrix & gsl::matrix::operator= ( matrix const &  v)
inline

The assignment operator.

This makes a shared copy.

Parameters
vThe matrix to copy

Definition at line 163 of file matrix.hpp.

References ccgsl_pointer, count, and owns_data.

◆ operator[]() [1/2]

vector gsl::matrix::operator[] ( size_t const  i)
inline

This function allows us to use a matrix like an array.

Use with caution. Although matrix[i][j] is possible, it is less efficient than matrix::set(). The effect is the same as row().

Parameters
iThe index of the row
Returns
A vector representing a row

Definition at line 1586 of file matrix.hpp.

References ccgsl_pointer, gsl::exception::GSL_EFAULT, gsl::exception::GSL_EINVAL, size1(), and gsl::vector::wrap_gsl_vector_without_ownership().

◆ operator[]() [2/2]

vector const gsl::matrix::operator[] ( size_t const  i) const
inline

This function allows us to use a matrix like an array.

Use with caution. Although matrix[i][j] is possible, it is much less efficient than matrix::set(). The effect is the same as row()

Parameters
iThe index of the row
Returns
A vector representing a row

Definition at line 1620 of file matrix.hpp.

References ccgsl_pointer, gsl::exception::GSL_EFAULT, and gsl::vector::wrap_gsl_vector_without_ownership().

◆ permute()

int gsl::matrix::permute ( permutation p)
inline

Permute the columns of this by permutation p.

Parameters
pThe permutation
Returns
Error code on failure

Definition at line 1637 of file matrix.hpp.

References get(), and gsl::permutation::get().

◆ ptr()

double * gsl::matrix::ptr ( size_t const  i,
size_t const  j 
)
inline

C++ version of gsl_matrix_ptr().

Parameters
iindex of row
jindex of column
Returns
pointer to element

Definition at line 1345 of file matrix.hpp.

References get().

Referenced by gsl::matrix::iterator_base< container, content, reverse_t >::operator->().

◆ rbegin() [1/2]

reverse_iterator gsl::matrix::rbegin ( )
inline

Get iterator pointing to first vector element.

Returns
iterator pointing to first vector element

Definition at line 680 of file matrix.hpp.

References ccgsl_pointer, and size1().

◆ rbegin() [2/2]

const_reverse_iterator gsl::matrix::rbegin ( ) const
inline

Get iterator pointing to first vector element.

Returns
iterator pointing to first vector element

Definition at line 688 of file matrix.hpp.

References ccgsl_pointer, and size1().

◆ rend() [1/2]

reverse_iterator gsl::matrix::rend ( )
inline

Get iterator pointing beyond last vector element.

Returns
iterator pointing beyond last vector element

Definition at line 697 of file matrix.hpp.

◆ rend() [2/2]

const_reverse_iterator gsl::matrix::rend ( ) const
inline

Get iterator pointing beyond last vector element.

Returns
iterator pointing beyond last vector element

Definition at line 704 of file matrix.hpp.

◆ reset()

void gsl::matrix::reset ( )
inline

Stop sharing ownership of the shared pointer.

Definition at line 230 of file matrix.hpp.

References matrix().

◆ row() [1/2]

vector gsl::matrix::row ( size_t const  i)
inline

C++ version of gsl_matrix_row().

Parameters
iA row index
Returns
A row as a vector

Definition at line 793 of file matrix.hpp.

References get(), and gsl::vector::vector().

Referenced by matrix().

◆ row() [2/2]

vector const gsl::matrix::row ( size_t const  i) const
inline

Another C++ version of gsl_matrix_const_row().

Parameters
iA row index
Returns
A row as a vector

Definition at line 1016 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ scale()

int gsl::matrix::scale ( double const  x)
inline

C++ version of gsl_matrix_scale().

Parameters
xconstant to multiply this by
Returns
error code on failure

Definition at line 1298 of file matrix.hpp.

References get().

◆ scale_columns()

int gsl::matrix::scale_columns ( vector const &  x)
inline

C++ version of gsl_matrix_scale_columns().

Parameters
xA scalar
Returns
Error code on failure

Definition at line 1537 of file matrix.hpp.

References get(), and gsl::vector::get().

◆ scale_rows()

int gsl::matrix::scale_rows ( vector const &  x)
inline

C++ version of gsl_matrix_scale_rows().

Parameters
xA scalar
Returns
Error code on failure

Definition at line 1530 of file matrix.hpp.

References get(), and gsl::vector::get().

◆ set()

void gsl::matrix::set ( size_t const  i,
size_t const  j,
double  x 
)
inline

C++ version of gsl_matrix_set().

Parameters
iindex of row
jindex of column
xnew value for element

Definition at line 1338 of file matrix.hpp.

References get().

Referenced by matrix().

◆ set_all()

void gsl::matrix::set_all ( double  x)
inline

C++ version of gsl_matrix_set_all().

Parameters
xThe value to which all elements are set

Definition at line 1250 of file matrix.hpp.

References get().

◆ set_col()

int gsl::matrix::set_col ( size_t const  j,
vector const &  v 
)
inline

C++ version of gsl_matrix_set_col().

Parameters
jThe index of the column
vA vector
Returns
error code on failure

Definition at line 1576 of file matrix.hpp.

References get(), and gsl::vector::get().

◆ set_identity()

void gsl::matrix::set_identity ( )
inline

C++ version of gsl_matrix_set_identity().

Definition at line 1420 of file matrix.hpp.

References get().

◆ set_row()

int gsl::matrix::set_row ( size_t const  i,
vector const &  v 
)
inline

C++ version of gsl_matrix_set_row().

Parameters
iThe index of the row
vA vector
Returns
error code on failure

Definition at line 1568 of file matrix.hpp.

References get(), and gsl::vector::get().

◆ set_zero()

void gsl::matrix::set_zero ( )
inline

C++ version of gsl_matrix_set_zero().

Definition at line 1245 of file matrix.hpp.

References get().

◆ size1()

size_t gsl::matrix::size1 ( ) const
inline

The number of rows of the matrix.

Returns
The number of rows of the matrix

Definition at line 713 of file matrix.hpp.

References ccgsl_pointer.

Referenced by clone(), gsl::matrix::iterator_base< container, content, reverse_t >::decrement(), end(), gsl::matrix::iterator_base< container, content, reverse_t >::increment(), operator[](), and rbegin().

◆ size2()

size_t gsl::matrix::size2 ( ) const
inline

The number of columns of the matrix.

Returns
The number of columns of the matrix

Definition at line 718 of file matrix.hpp.

References ccgsl_pointer.

Referenced by clone().

◆ sub()

int gsl::matrix::sub ( matrix const &  b)
inline

C++ version of gsl_matrix_sub().

Parameters
bmatrix to subtract from this
Returns
error code on failure

Definition at line 1292 of file matrix.hpp.

References gsl::sf::mathieu::b(), and get().

◆ subcolumn() [1/2]

vector gsl::matrix::subcolumn ( size_t const  j,
size_t const  offset,
size_t const  n 
)
inline

C++ version of gsl_matrix_subcolumn().

Parameters
jA column index
offsetA row offset
nThe number of elements
Returns
A subcolumn as a vector

Definition at line 856 of file matrix.hpp.

References get(), gsl::rstat::n(), and gsl::vector::vector().

◆ subcolumn() [2/2]

vector const gsl::matrix::subcolumn ( size_t const  j,
size_t const  offset,
size_t const  n 
) const
inline

Another C++ version of gsl_matrix_const_subcolumn().

Parameters
jA column index
offsetA row offset
nThe number of elements
Returns
A subcolumn as a vector

Definition at line 1079 of file matrix.hpp.

References get(), gsl::rstat::n(), and gsl::vector::vector().

◆ subdiagonal() [1/2]

vector gsl::matrix::subdiagonal ( size_t const  k)
inline

C++ version of gsl_matrix_subdiagonal().

Parameters
kAn index
Returns
Subdiagonal k as a vector

Definition at line 822 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ subdiagonal() [2/2]

vector const gsl::matrix::subdiagonal ( size_t const  k) const
inline

Another C++ version of gsl_matrix_const_subdiagonal().

Parameters
kAn index
Returns
Subdiagonal k as a vector

Definition at line 1045 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ submatrix() [1/2]

matrix gsl::matrix::submatrix ( size_t const  i,
size_t const  j,
size_t const  n1,
size_t const  n2 
)
inline

C++ version of gsl_matrix_submatrix().

Parameters
iIndex in this of first row of submatrix
jIndex in this of first column of submatrix
n1Number of rows of submatrix
n2Number of columns of submatrix
Returns
The submatrix

Definition at line 783 of file matrix.hpp.

References get(), and matrix().

◆ submatrix() [2/2]

matrix const gsl::matrix::submatrix ( size_t const  i,
size_t const  j,
size_t const  n1,
size_t const  n2 
) const
inline

Another C++ version of gsl_matrix_const_submatrix().

Parameters
iIndex in this of first row of submatrix
jIndex in this of first column of submatrix
n1Number of rows of submatrix
n2Number of columns of submatrix
Returns
The submatrix

Definition at line 1006 of file matrix.hpp.

References get(), and matrix().

◆ subrow() [1/2]

vector gsl::matrix::subrow ( size_t const  i,
size_t const  offset,
size_t const  n 
)
inline

C++ version of gsl_matrix_subrow().

Parameters
iA row index
offsetA column offset
nThe number of elements
Returns
A subrow as a vector

Definition at line 844 of file matrix.hpp.

References get(), gsl::rstat::n(), and gsl::vector::vector().

◆ subrow() [2/2]

vector const gsl::matrix::subrow ( size_t const  i,
size_t const  offset,
size_t const  n 
) const
inline

Another C++ version of gsl_matrix_const_subrow().

Parameters
iA row index
offsetA column offset
nThe number of elements
Returns
A subrow as a vector

Definition at line 1067 of file matrix.hpp.

References get(), gsl::rstat::n(), and gsl::vector::vector().

◆ superdiagonal() [1/2]

vector gsl::matrix::superdiagonal ( size_t const  k)
inline

C++ version of gsl_matrix_superdiagonal().

Parameters
kAn index
Returns
Subdiagonal k as a vector

Definition at line 832 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ superdiagonal() [2/2]

vector const gsl::matrix::superdiagonal ( size_t const  k) const
inline

Another C++ version of gsl_matrix_const_superdiagonal().

Parameters
kAn index
Returns
Subdiagonal k as a vector

Definition at line 1055 of file matrix.hpp.

References get(), and gsl::vector::vector().

◆ swap()

void gsl::matrix::swap ( matrix m)
inline

Swap two matrix objects.

This works even if the matrix objects have different sizes because it swaps pointers.

Parameters
mThe matrix to swap with this.

Definition at line 752 of file matrix.hpp.

References ccgsl_pointer, and count.

◆ swap_columns()

int gsl::matrix::swap_columns ( size_t const  i,
size_t const  j 
)
inline

C++ version of gsl_matrix_swap_columns().

Parameters
iIndex of first column
jIndex of second column
Returns
error code on failure

Definition at line 1434 of file matrix.hpp.

References get().

◆ swap_rowcol()

int gsl::matrix::swap_rowcol ( size_t const  i,
size_t const  j 
)
inline

C++ version of gsl_matrix_swap_rowcol().

Swap row and column in place. Matrix must be square.

Parameters
iindex of row
jindex of column
Returns
error code on failure

Definition at line 1442 of file matrix.hpp.

References get().

◆ swap_rows()

int gsl::matrix::swap_rows ( size_t const  i,
size_t const  j 
)
inline

C++ version of gsl_matrix_swap_rows().

Parameters
iIndex of first row
jIndex of second row
Returns
error code on failure

Definition at line 1427 of file matrix.hpp.

References get().

◆ transpose()

int gsl::matrix::transpose ( )
inline

C++ version of gsl_matrix_transpose().

Returns
error code on failure.

Definition at line 1447 of file matrix.hpp.

References get().

◆ transpose_memcpy()

int gsl::matrix::transpose_memcpy ( matrix const &  src)
inline

C++ version of gsl_matrix_transpose_memcpy().

Parameters
srcmatrix whose transpose it to be copied to this
Returns
error code on failure

Definition at line 1453 of file matrix.hpp.

References get().

◆ transpose_tricpy()

void gsl::matrix::transpose_tricpy ( CBLAS_UPLO_t  Uplo,
CBLAS_DIAG_t  Diag,
matrix const &  src 
)
inline

Copy the upper or lower triangular part of matrix src to this.

Parameters
UploUpper or lower triangle: CblasUpper or CBlasLower
DiagDiagonal type
srcThe matrix to copy from

Definition at line 771 of file matrix.hpp.

References get().

◆ tricpy()

void gsl::matrix::tricpy ( CBLAS_UPLO_t  Uplo,
CBLAS_DIAG_t  Diag,
matrix const &  src 
)
inline

Copy the upper or lower triangular part of matrix src to this.

Parameters
UploUpper or lower triangle: CblasUpper or CBlasLower
DiagDiagonal type
srcThe matrix to copy from

Definition at line 762 of file matrix.hpp.

References get().

◆ unique()

bool gsl::matrix::unique ( ) const
inline

Find if this is the only object sharing the gsl_matrix.

Returns
true or falses according as this is the only matrix object sharing the gsl_matrix

Definition at line 1218 of file matrix.hpp.

References count.

◆ use_count()

size_t gsl::matrix::use_count ( ) const
inline

Find how many matrix objects share this pointer.

Returns
the number of matrix objects that share this pointer

Definition at line 1223 of file matrix.hpp.

References count.

◆ view_array() [1/2]

static matrix gsl::matrix::view_array ( double *  base,
size_t const  n1,
size_t const  n2 
)
inlinestatic

C++ version of gsl_matrix_view_array().

Parameters
baseAn array of type double
n1The number of rows
n2The number of columns
Returns
A matrix

Definition at line 868 of file matrix.hpp.

References matrix().

◆ view_array() [2/2]

static matrix const gsl::matrix::view_array ( double const *  base,
size_t const  n1,
size_t const  n2 
)
inlinestatic

Another C++ version of gsl_matrix_const_view_array().

Parameters
baseAn array of type double
n1The number of rows
n2The number of columns
Returns
A matrix

Definition at line 1143 of file matrix.hpp.

References matrix().

◆ view_array_with_tda() [1/2]

static matrix gsl::matrix::view_array_with_tda ( double *  base,
size_t const  n1,
size_t const  n2,
size_t const  tda 
)
inlinestatic

C++ version of gsl_matrix_view_array_with_tda().

Parameters
baseAn array of type double
n1The number of rows
n2The number of columns
tdaThe number of columns in memory
Returns
A matrix

Definition at line 881 of file matrix.hpp.

References matrix().

◆ view_array_with_tda() [2/2]

static matrix const gsl::matrix::view_array_with_tda ( double const *  base,
size_t const  n1,
size_t const  n2,
size_t const  tda 
)
inlinestatic

Another C++ version of gsl_matrix_const_view_array_with_tda().

Parameters
baseAn array of type double
n1The number of rows
n2The number of columns
tdaThe number of columns in memory
Returns
A matrix

Definition at line 1157 of file matrix.hpp.

References matrix().

◆ view_vector() [1/2]

static matrix gsl::matrix::view_vector ( vector v,
size_t const  n1,
size_t const  n2 
)
inlinestatic

C++ version of gsl_matrix_view_vector().

Parameters
vA vector
n1The number of rows
n2The number of columns
Returns
A matrix

Definition at line 893 of file matrix.hpp.

References gsl::vector::get(), and matrix().

◆ view_vector() [2/2]

static matrix const gsl::matrix::view_vector ( vector const &  v,
size_t const  n1,
size_t const  n2 
)
inlinestatic

Another C++ version of gsl_matrix_const_view_vector().

Parameters
vA vector
n1The number of rows
n2The number of columns
Returns
A matrix

Definition at line 1169 of file matrix.hpp.

References gsl::vector::get(), and matrix().

◆ view_vector_with_tda() [1/2]

static matrix gsl::matrix::view_vector_with_tda ( vector v,
size_t const  n1,
size_t const  n2,
size_t const  tda 
)
inlinestatic

C++ version of gsl_matrix_view_vector_with_tda().

Parameters
vA vector
n1The number of rows
n2The number of columns
tdaThe number of columns in memory
Returns
A matrix

Definition at line 906 of file matrix.hpp.

References gsl::vector::get(), and matrix().

◆ view_vector_with_tda() [2/2]

static matrix const gsl::matrix::view_vector_with_tda ( vector const &  v,
size_t const  n1,
size_t const  n2,
size_t const  tda 
)
inlinestatic

Another C++ version of gsl_matrix_const_view_vector_with_tda().

Parameters
vA vector
n1The number of rows
n2The number of columns
tdaThe number of columns in memory
Returns
A matrix

Definition at line 1183 of file matrix.hpp.

References gsl::vector::get(), and matrix().

◆ wrap_gsl_matrix_without_ownership()

void gsl::matrix::wrap_gsl_matrix_without_ownership ( gsl_matrix *  v)
inline

This function is intended mainly for internal use.

It allows this to point to a gsl_matrix without the possibility deleting it when this is no longer in scope. It is the responsibility of the programmer to delete v. The function is used internally for converting a function that takes a gsl::matrix* argument to one that takes a gsl_matrix* argument.

Parameters
vThe gsl_matrix

Definition at line 214 of file matrix.hpp.

References ccgsl_pointer, count, and owns_data.

Friends And Related Function Documentation

◆ multifit::nlinear::function_fdf

friend class multifit::nlinear::function_fdf
friend

Definition at line 74 of file matrix.hpp.

◆ multilarge::nlinear::function_fdf

friend class multilarge::nlinear::function_fdf
friend

Definition at line 75 of file matrix.hpp.

◆ multiroot::function_fdf

friend class multiroot::function_fdf
friend

Definition at line 73 of file matrix.hpp.

Member Data Documentation

◆ ccgsl_pointer

gsl_matrix* gsl::matrix::ccgsl_pointer
private

◆ count

size_t* gsl::matrix::count
private

The shared reference count.

Definition at line 1200 of file matrix.hpp.

Referenced by matrix(), operator=(), swap(), unique(), use_count(), wrap_gsl_matrix_without_ownership(), and ~matrix().

◆ owns_data

bool gsl::matrix::owns_data
private

Used to allow a vector that does not own its data.

Definition at line 1192 of file matrix.hpp.

Referenced by matrix(), operator=(), wrap_gsl_matrix_without_ownership(), and ~matrix().


The documentation for this class was generated from the following file: