reconciliation_lib 0.1
Library for data reconciliation algorithms
Loading...
Searching...
No Matches
dvrlib::matrix Class Reference

RAII wrapper around a heap-allocated gsl_matrix. More...

#include <gsl_wrapper.h>

Public Member Functions

 matrix (int n1, int n2, bool id=false, const double *diag=nullptr)
 Allocate an n1 × n2 matrix.
 
 matrix (int n1, int n2, const double *x)
 Allocate an n1 × n2 matrix and copy values from the row-major array x.
 
template<int n>
 matrix (int n1, int n2, const double(*x)[n])
 Allocate an n1 × n2 matrix from a 2-D C array.
 
 matrix (const matrix &src)
 Copy constructor.
 
 matrix (const gsl_matrix *src)
 Construct from an existing (non-owned) gsl_matrix pointer.
 
 ~matrix ()
 Destructor — frees the underlying GSL allocation.
 
gsl_matrixgsl_internal ()
 Return a pointer to the underlying gsl_matrix (mutable).
 
const gsl_matrixgsl_internal () const
 Return a pointer to the underlying gsl_matrix (read-only).
 
int size1 () const
 Return the number of rows.
 
int size2 () const
 Return the number of columns.
 
void set (int i, int j, double val)
 Set element (i, j) to val.
 
double get (int i, int j) const
 Return element (i, j).
 
vector_view operator[] (int i)
 Return a mutable view of row i.
 
const vector_view operator[] (int i) const
 Return a read-only view of row i.
 
matrixoperator= (const matrix &src)
 Copy-assign from another matrix.
 
matrixoperator= (const matrix_view &src)
 Copy-assign from a view.
 
matrix operator+ (const matrix &src) const
 
matrixoperator+= (const matrix &src)
 
matrix operator- (const matrix &src) const
 
matrixoperator-= (const matrix &src)
 
matrix operator- () const
 Unary negation — returns -(*this).
 
vector operator* (const vector &src) const
 Matrix–vector product.
 
matrix operator* (const matrix &src) const
 Matrix–matrix product.
 
matrix operator* (double d) const
 Return a copy scaled by d.
 
matrixoperator*= (double d)
 Scale all elements by d in place.
 
matrix transpose () const
 Return the transpose.
 
matrix inverse () const
 Return the inverse (via LU decomposition).
 
vector linsolve (const vector &b) const
 Solve the linear system A*x = b and return x.
 
void svd (matrix &U, matrix &V, vector &s) const
 Compute the full singular value decomposition A = U * diag(s) * V^T.
 
vector svd () const
 Compute and return only the singular values.
 
matrix_view submatrix (int k1, int k2, int n1, int n2)
 Return a view of the n1 × n2 sub-matrix starting at (k1, k2).
 
const matrix_view submatrix (int k1, int k2, int n1, int n2) const
 Return a view of the n1 × n2 sub-matrix starting at (k1, k2).
 

Private Attributes

gsl_matrixm
 

Friends

class matrix_view
 

Detailed Description

RAII wrapper around a heap-allocated gsl_matrix.

Owns the underlying GSL allocation and frees it on destruction. Row access via operator[] returns a vector_view into the matrix data.

Constructor & Destructor Documentation

◆ matrix() [1/5]

dvrlib::matrix::matrix ( int  n1,
int  n2,
bool  id = false,
const double diag = nullptr 
)

Allocate an n1 × n2 matrix.

Parameters
[in]n1Number of rows.
[in]n2Number of columns.
[in]idIf true, initialise as identity (requires n1 == n2 when diag is null).
[in]diagOptional array of n1 diagonal values; if non-null overrides id.

◆ matrix() [2/5]

dvrlib::matrix::matrix ( int  n1,
int  n2,
const double x 
)

Allocate an n1 × n2 matrix and copy values from the row-major array x.

◆ matrix() [3/5]

template<int n>
dvrlib::matrix::matrix ( int  n1,
int  n2,
const double(*)  x[n] 
)

Allocate an n1 × n2 matrix from a 2-D C array.

The template parameter n must match n2.

◆ matrix() [4/5]

dvrlib::matrix::matrix ( const matrix src)

Copy constructor.

◆ matrix() [5/5]

dvrlib::matrix::matrix ( const gsl_matrix src)

Construct from an existing (non-owned) gsl_matrix pointer.

◆ ~matrix()

dvrlib::matrix::~matrix ( )

Destructor — frees the underlying GSL allocation.

Member Function Documentation

◆ get()

double dvrlib::matrix::get ( int  i,
int  j 
) const

Return element (i, j).

◆ gsl_internal() [1/2]

gsl_matrix * dvrlib::matrix::gsl_internal ( )

Return a pointer to the underlying gsl_matrix (mutable).

◆ gsl_internal() [2/2]

const gsl_matrix * dvrlib::matrix::gsl_internal ( ) const

Return a pointer to the underlying gsl_matrix (read-only).

◆ inverse()

matrix dvrlib::matrix::inverse ( ) const

Return the inverse (via LU decomposition).

◆ linsolve()

vector dvrlib::matrix::linsolve ( const vector b) const

Solve the linear system A*x = b and return x.

Parameters
[in]bRight-hand side vector (length size1()).

◆ operator*() [1/3]

matrix dvrlib::matrix::operator* ( const matrix src) const

Matrix–matrix product.

◆ operator*() [2/3]

vector dvrlib::matrix::operator* ( const vector src) const

Matrix–vector product.

◆ operator*() [3/3]

matrix dvrlib::matrix::operator* ( double  d) const

Return a copy scaled by d.

◆ operator*=()

matrix & dvrlib::matrix::operator*= ( double  d)

Scale all elements by d in place.

◆ operator+()

matrix dvrlib::matrix::operator+ ( const matrix src) const

◆ operator+=()

matrix & dvrlib::matrix::operator+= ( const matrix src)

◆ operator-() [1/2]

matrix dvrlib::matrix::operator- ( ) const

Unary negation — returns -(*this).

◆ operator-() [2/2]

matrix dvrlib::matrix::operator- ( const matrix src) const

◆ operator-=()

matrix & dvrlib::matrix::operator-= ( const matrix src)

◆ operator=() [1/2]

matrix & dvrlib::matrix::operator= ( const matrix src)

Copy-assign from another matrix.

◆ operator=() [2/2]

matrix & dvrlib::matrix::operator= ( const matrix_view src)

Copy-assign from a view.

◆ operator[]() [1/2]

vector_view dvrlib::matrix::operator[] ( int  i)

Return a mutable view of row i.

◆ operator[]() [2/2]

const vector_view dvrlib::matrix::operator[] ( int  i) const

Return a read-only view of row i.

◆ set()

void dvrlib::matrix::set ( int  i,
int  j,
double  val 
)

Set element (i, j) to val.

◆ size1()

int dvrlib::matrix::size1 ( ) const

Return the number of rows.

◆ size2()

int dvrlib::matrix::size2 ( ) const

Return the number of columns.

◆ submatrix() [1/2]

matrix_view dvrlib::matrix::submatrix ( int  k1,
int  k2,
int  n1,
int  n2 
)

Return a view of the n1 × n2 sub-matrix starting at (k1, k2).

Parameters
[in]k1Starting row index (0-based).
[in]k2Starting column index (0-based).
[in]n1Number of rows in the sub-matrix.
[in]n2Number of columns in the sub-matrix.

◆ submatrix() [2/2]

const matrix_view dvrlib::matrix::submatrix ( int  k1,
int  k2,
int  n1,
int  n2 
) const

Return a view of the n1 × n2 sub-matrix starting at (k1, k2).

Parameters
[in]k1Starting row index (0-based).
[in]k2Starting column index (0-based).
[in]n1Number of rows in the sub-matrix.
[in]n2Number of columns in the sub-matrix.

◆ svd() [1/2]

vector dvrlib::matrix::svd ( ) const

Compute and return only the singular values.

◆ svd() [2/2]

void dvrlib::matrix::svd ( matrix U,
matrix V,
vector s 
) const

Compute the full singular value decomposition A = U * diag(s) * V^T.

Parameters
[out]ULeft singular vectors (same size as *this).
[out]VRight singular vectors (size2() × size2()).
[out]sSingular values (length size2()).

◆ transpose()

matrix dvrlib::matrix::transpose ( ) const

Return the transpose.

Friends And Related Symbol Documentation

◆ matrix_view

Member Data Documentation

◆ m

gsl_matrix* dvrlib::matrix::m
private

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