linalg::matrix Class Reference

A wrapper class for GNU Scientific Library matrices. More...

#include <linalg.hpp>

Collaboration diagram for linalg::matrix:

Collaboration graph
[legend]

List of all members.

Public Member Functions

 matrix ()
 Allocate 1x1 matrix zero matrix.
 matrix (const size_t m, const size_t n, const double fillvalue=0)
 Allocate m by n matrix, filled with some value.
 matrix (gsl_matrix *M)
 Assign a matrix from a GSL matrix pointer.
 matrix (const matrix &M)
 Copy constructor.
matrix operator= (const matrix &M)
 Assignment.
 ~matrix ()
 Clears the memory allocated by the GSL.
size_t precision () const
 Number of decimal digits to output.
size_t rows () const
 Number of rows.
size_t cols () const
 Number of columns.
double & operator() (const size_t i, const size_t j)
 Fortran-style parenthetical indexing (hence Octave-style too) Indices start from 1.
const double & operator() (const size_t i, const size_t j) const
 Constant version of previous function.
vector_view operator() (const size_t i, const slice &b)
 Indexing for vectorisation. The GSL provides limited vectorisation routines which can be useful for operating on blocks of matrices all at once.
const vector_view operator() (const size_t i, const slice &b) const
 Constant version of previous function.
vector_view operator() (const slice &a, const size_t j)
 Indexing for vectorisation. The GSL provides limited vectorisation routines which can be useful for operating on blocks of matrices all at once.
const vector_view operator() (const slice &a, const size_t j) const
 Constant version of previous function.
matrix operator* (const double a) const
 Scale the matrix.
matrix operator+ (const matrix &N) const
 Matrix addition.
matrix operator* (const matrix &N) const
 Matrix multiplication.
matrix operator- (const matrix &N) const
 Matrix subtraction.
vector operator* (const vector &v) const
 Matrix vector product (gaxpy operation).
matrix inv () const
 Inverse.
matrix T () const
 Tranpose. Returns a copy of the original matrix, transposed. The original matrix is not modified.
double tr () const
 Trace. The trace of a square matrix.
double det () const
 Determinant. The determinant computed via an LU factorisation.
vector inv (const vector &w) const
 Solves Mv = w for v with LU factorisation. Given another vector, will return the solution vector v to the equation Mv = w.
double cond () const
 L2 condition number, using SVD.

Static Public Member Functions

static void set_precision (size_t p)
 Set the number of decimal digits to output.

Private Member Functions

LUmatrixLU () const
 LU decomposition, pivots in U.
void SVDfactor () const

Private Attributes

gsl_matrix * A
LUmatrixLUptr
bool LUfactored
double condition_number
gsl_vector * SVD
bool SVDfactored

Static Private Attributes

static size_t precsn = 4
 Precision.
static const double eps = std::numeric_limits<double>::epsilon()
 Machine epsilon.

Friends

class vector_view

Classes

class  LUmatrix
 A private matrix member class for matrices factored in LU form. More...


Detailed Description

A wrapper class for GNU Scientific Library matrices.

Constructor & Destructor Documentation

linalg::matrix::matrix (  ) 

Allocate 1x1 matrix zero matrix.

00022                 {
00023     A = gsl_matrix_calloc(1,1); //Allocate 1x1 matrix zero matrix.    
00024     LUfactored = false;
00025     SVDfactored = false;
00026   }

linalg::matrix::matrix ( const size_t  m,
const size_t  n,
const double  fillvalue = 0 
)

Allocate m by n matrix, filled with some value.

Parameters:
m - Number of rows.
n - Number of columns.
fillvalue - Value to fill matrix with, default is 0.
00028                                                     {
00029     A = gsl_matrix_alloc(m,n);
00030     gsl_matrix_set_all(A,fillvalue);
00031     LUfactored = false;
00032     SVDfactored = false;
00033   }

linalg::matrix::matrix ( gsl_matrix *  M  ) 

Assign a matrix from a GSL matrix pointer.

00035                              {
00036     A = M;
00037     LUfactored = false;
00038     SVDfactored = false;
00039   }

linalg::matrix::matrix ( const matrix M  ) 

Copy constructor.

00041                                {
00042     A = gsl_matrix_alloc(M.A->size1, M.A->size2);
00043     gsl_matrix_memcpy(A,M.A);
00044 
00045     LUfactored = M.LUfactored;
00046     if(LUfactored)
00047       LUptr = new LUmatrix(*M.LUptr);
00048     
00049     condition_number = M.condition_number;
00050 
00051     SVDfactored = M.SVDfactored;
00052     if(SVDfactored){
00053       SVD = gsl_vector_alloc(M.SVD->size);
00054       gsl_vector_memcpy(SVD,M.SVD);
00055     }       
00056   }

linalg::matrix::~matrix (  ) 

Clears the memory allocated by the GSL.

00085                  {
00086     if(A != 0) //Has subclass already deleted this matrix?
00087       gsl_matrix_free(A);
00088     if(LUfactored)
00089       delete LUptr;
00090     if(SVDfactored)
00091       gsl_vector_free(SVD);
00092   }


Member Function Documentation

matrix linalg::matrix::operator= ( const matrix M  ) 

Assignment.

00058                                          {
00059     if(this != &M){
00060       gsl_matrix_free(A);
00061       A = gsl_matrix_alloc(M.A->size1, M.A->size2);
00062       gsl_matrix_memcpy(A,M.A);
00063 
00064       if(LUfactored)
00065         delete LUptr;
00066 
00067       LUfactored = M.LUfactored;
00068       if(LUfactored){
00069         LUptr = new LUmatrix(*M.LUptr);
00070       }
00071       
00072       if(SVDfactored)
00073         gsl_vector_free(SVD);
00074 
00075       SVDfactored = M.SVDfactored;
00076       if(SVDfactored){
00077         condition_number = M.condition_number;
00078         SVD = gsl_vector_alloc(M.SVD->size);
00079         gsl_vector_memcpy(SVD,M.SVD);
00080       }       
00081     }
00082     return *this;
00083   }

size_t linalg::matrix::precision (  )  const

Number of decimal digits to output.

00098                                 {
00099     return precsn;
00100   }

void linalg::matrix::set_precision ( size_t  p  )  [static]

Set the number of decimal digits to output.

00102                                     {
00103     precsn = p;
00104   }

size_t linalg::matrix::rows (  )  const

Number of rows.

00106                            {
00107     return A -> size1;
00108   }

size_t linalg::matrix::cols (  )  const

Number of columns.

00110                            {
00111     return A -> size2;
00112   }

double & linalg::matrix::operator() ( const size_t  i,
const size_t  j 
) [inline]

Fortran-style parenthetical indexing (hence Octave-style too) Indices start from 1.

Parameters:
i - Row number.
j - Column number.
Returns:
A reference to the matrix element.
00414                                                                  {
00415     try{
00416       if(LUfactored)
00417         delete LUptr;
00418       if(SVDfactored)
00419         gsl_vector_free(SVD);
00420   
00421       LUfactored = false;
00422       SVDfactored = false;
00423       return(*gsl_matrix_ptr(A,i-1,j-1));
00424     }
00425     catch(indexOutOfRange& exc){
00426       exc.i = i;
00427       exc.j = j;
00428       exc.m = A -> size1;
00429       exc.n = A -> size2;
00430       throw exc;
00431     }
00432   }

const double & linalg::matrix::operator() ( const size_t  i,
const size_t  j 
) const [inline]

Constant version of previous function.

00434                                                               {
00435     try{
00436       return(*gsl_matrix_const_ptr(A,i-1,j-1));
00437     }
00438     catch(indexOutOfRange& exc){
00439       exc.i = i;
00440       exc.j = j;
00441       exc.m = A -> size1;
00442       exc.n = A -> size2;
00443       throw exc;
00444     }
00445   }

vector_view linalg::matrix::operator() ( const size_t  i,
const slice b 
)

Indexing for vectorisation. The GSL provides limited vectorisation routines which can be useful for operating on blocks of matrices all at once.

Parameters:
i - Row to be sliced.
b - Slice range.
Returns:
A vector_view of the sliced row.
00118                                                                {
00119     vector_view x_sub(A,i,b);
00120 
00121     if(LUfactored)
00122       delete LUptr;
00123     if(SVDfactored)
00124       gsl_vector_free(SVD);
00125 
00126     LUfactored = false;
00127     SVDfactored = false;
00128     return x_sub;
00129   }

const vector_view linalg::matrix::operator() ( const size_t  i,
const slice b 
) const

Constant version of previous function.

00130                                                                            {
00131     vector_view x_sub(A,i,b);
00132     return x_sub;
00133   }

vector_view linalg::matrix::operator() ( const slice a,
const size_t  j 
)

Indexing for vectorisation. The GSL provides limited vectorisation routines which can be useful for operating on blocks of matrices all at once.

Parameters:
a - Slice range
j - column to be sliced.
Returns:
A vector_view of the sliced column.
00135                                                               {
00136     vector_view x_sub(A,a,j);
00137 
00138     if(LUfactored)
00139       delete LUptr;
00140     if(SVDfactored)
00141       gsl_vector_free(SVD);
00142     
00143     LUfactored = false;
00144     SVDfactored = false;
00145     return x_sub;
00146   }

const vector_view linalg::matrix::operator() ( const slice a,
const size_t  j 
) const

Constant version of previous function.

00147                                                                           {
00148     vector_view x_sub(A,a,j);
00149     return x_sub;
00150   }

matrix linalg::matrix::operator* ( const double  a  )  const

Scale the matrix.

00153                                                {
00154     matrix Z = *this;
00155     gsl_matrix_scale(Z.A,a);
00156     return Z;
00157   }

matrix linalg::matrix::operator+ ( const matrix N  )  const

Matrix addition.

00159                                                {
00160     matrix Z = *this;
00161     try{
00162       gsl_matrix_add(Z.A,N.A);
00163     }
00164     catch(inconformantSizes& exc){
00165       exc.n_A = A->size1;
00166       exc.m_A = A->size2;
00167       exc.n_B = N.A -> size1;
00168       exc.m_B = N.A -> size2;
00169       throw exc;
00170     }
00171     return Z;
00172   }

matrix linalg::matrix::operator* ( const matrix N  )  const

Matrix multiplication.

00174                                               {
00175     matrix Z(A->size1,N.A->size2);
00176     try{
00177       gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1,  A, N.A, 0, Z.A);
00178     }
00179     catch(inconformantSizes& exc){
00180       exc.n_A = A->size1;
00181       exc.m_A = A->size2;
00182       exc.n_B = N.A -> size1;
00183       exc.m_B = N.A -> size2;
00184       throw exc;
00185     }
00186     return Z;
00187   }

matrix linalg::matrix::operator- ( const matrix N  )  const

Matrix subtraction.

00189                                                {
00190     matrix Z = *this;
00191     try{
00192       gsl_matrix_sub(Z.A,N.A);
00193     }
00194     catch(inconformantSizes& exc){
00195       exc.n_A = A->size1;
00196       exc.m_A = A->size2;
00197       exc.n_B = N.A -> size1;
00198       exc.m_B = N.A -> size2;
00199       throw exc;
00200     }
00201     return Z;
00202   }

vector linalg::matrix::operator* ( const vector v  )  const

Matrix vector product (gaxpy operation).

00204                                                {
00205     vector w(A -> size1);
00206     try{
00207       gsl_blas_dgemv(CblasNoTrans, 1, A, v.x, 0, w.x);
00208     }
00209     catch(inconformantSizes& exc){
00210       exc.n_A = A->size1;
00211       exc.m_A = A->size2;
00212       exc.n_B = v.x -> size;
00213       exc.m_B = 1;
00214       throw exc;
00215     }
00216     return w;
00217   }

matrix linalg::matrix::inv (  )  const

Inverse.

This function returns the inverse matrix, computed via an LU factorisation. If the matrix is ill-conditioned or the inverse does not exist, it will happily return a matrix full of NaN.

Returns:
The inverse matrix.
00240                           {
00241     if(A -> size1 != A -> size2){
00242       matrixNotSquare exc;
00243       exc.reason = "Cannot invert non-square matrices.";
00244       exc.file = __FILE__;
00245       exc.line = __LINE__;
00246       throw exc;
00247     }
00248     
00249     if(!LUfactored){
00250       LUptr =  LU();    
00251       LUfactored = true;
00252     }
00253     
00254     matrix Z(LUptr->matrix_ptr()->size1, LUptr->matrix_ptr() -> size2);
00255 
00256     gsl_linalg_LU_invert(LUptr->matrix_ptr(), LUptr->perm_ptr(), Z.A);
00257     return Z;
00258   }

Here is the call graph for this function:

matrix linalg::matrix::T (  )  const

Tranpose. Returns a copy of the original matrix, transposed. The original matrix is not modified.

Returns:
Transposed copy of matrix.
00260                          {
00261     matrix Z(A->size2, A->size1);
00262     for(size_t i = 1; i <= A->size1; i++)
00263       for(size_t j = 1; j <= A->size2; j++)
00264         Z(j,i) = gsl_matrix_get(A,i-1,j-1);
00265     return Z;
00266   }

double linalg::matrix::tr (  )  const

Trace. The trace of a square matrix.

Returns:
The trace.
00268                          {
00269     if(A -> size1 != A -> size2){
00270       matrixNotSquare exc;
00271       exc.reason = "Cannot find trace of non-square matrix.";
00272       exc.file = __FILE__;
00273       exc.line = __LINE__;
00274       throw exc;
00275     }
00276 
00277     double result=0;
00278     for(size_t i = 1; i <= A->size1; i++)
00279       result += gsl_matrix_get(A,i,i);
00280     return result;
00281   }

double linalg::matrix::det (  )  const

Determinant. The determinant computed via an LU factorisation.

Returns:
The determinant.
00283                           {
00284     if(A -> size1 != A -> size2){
00285       matrixNotSquare exc;
00286       exc.reason = "Cannot find determinant of non-square matrices.";
00287       exc.file = __FILE__;
00288       exc.line = __LINE__;
00289       throw exc;
00290     }
00291     
00292     if( !LUfactored ){
00293       LUptr = LU();
00294       LUfactored = true;
00295     }
00296 
00297     double out = gsl_linalg_LU_det(LUptr->matrix_ptr(),LUptr->sgn());
00298     return out;
00299   }

Here is the call graph for this function:

vector linalg::matrix::inv ( const vector w  )  const

Solves Mv = w for v with LU factorisation. Given another vector, will return the solution vector v to the equation Mv = w.

Parameters:
w - Another vector, must have the same size as number of rows in the matrix.
Returns:
The solution v to the equation Mv = w.
00301                                          {
00302     if(A -> size1 != A -> size2){
00303       matrixNotSquare exc;
00304       exc.reason = "Cannot invert non-square matrices.";
00305       exc.file = __FILE__;
00306       exc.line = __LINE__;
00307       throw exc;
00308     }
00309 
00310     vector v(w.size());
00311     if( !LUfactored){
00312       LUptr = LU();
00313       LUfactored = true;
00314     }
00315 
00316     gsl_linalg_LU_solve(LUptr->matrix_ptr(), LUptr->perm_ptr(), w.x, v.x);
00317     return v;
00318   }

Here is the call graph for this function:

double linalg::matrix::cond (  )  const

L2 condition number, using SVD.

This function computes the L2 condition number of a matrix, the largest singular value divided by the smallest. It directly uses the SVD decomposition (as provided by the GSL), and it may take a long time for large matrices.

Returns:
The L2 condition number.
00320                            {
00321     if(SVDfactored)
00322       return condition_number;
00323 
00324     SVDfactor();
00325     condition_number =  gsl_vector_get(SVD,0) 
00326                         /gsl_vector_get(SVD,std::max(A->size1, A->size2)-1);
00327     return condition_number;
00328   }

Here is the call graph for this function:

matrix::LUmatrix * linalg::matrix::LU (  )  const [private]

LU decomposition, pivots in U.

This function will compute the LU factorisation of a matrix, where the pivots are in the upper triangular matrix U and L has ones along the diagonal (not stored). For efficiency, L and U are both internally stored within a single matrix.

Returns:
A pointer to the private class LUmatrix.
00222                                   {
00223     if(A -> size1 != A -> size2){
00224       matrixNotSquare exc;
00225       exc.reason = "Cannot LU-factorise non-square matrices.";
00226       exc.file = __FILE__;
00227       exc.line = __LINE__;
00228       throw exc;
00229     }
00230 
00231     if(!LUfactored){
00232       LUptr = new LUmatrix(A);
00233       gsl_linalg_LU_decomp(LUptr->matrix_ptr(), LUptr->perm_ptr(), 
00234                            LUptr->sgn_ptr());
00235       LUfactored = true;
00236     }
00237     return LUptr;    
00238   }

Here is the call graph for this function:

void linalg::matrix::SVDfactor (  )  const [private]

00330                               {
00331     if(SVDfactored)
00332       return;
00333 
00334     gsl_matrix * M;
00335     if(A -> size1 < A -> size2){ //Transpose so that m>=n
00336       M = gsl_matrix_alloc(A -> size2, A -> size1);
00337       gsl_matrix_transpose_memcpy(M,A);
00338     }
00339     else{
00340       M = gsl_matrix_alloc(A -> size1, A -> size2);
00341       gsl_matrix_memcpy(M,A);
00342     }
00343     
00344     SVD = gsl_vector_alloc(M -> size2);
00345     gsl_matrix* V = gsl_matrix_alloc(M -> size2, M -> size2);
00346     gsl_vector* work = gsl_vector_alloc(M -> size2);
00347     gsl_linalg_SV_decomp(M,V,SVD,work);
00348     SVDfactored = true;
00349     gsl_vector_free(work);
00350     gsl_matrix_free(V);
00351     gsl_matrix_free(M);
00352   }


Friends And Related Function Documentation

friend class vector_view [friend]


Member Data Documentation

gsl_matrix* linalg::matrix::A [private]

size_t linalg::matrix::precsn = 4 [static, private]

Precision.

const double linalg::matrix::eps = std::numeric_limits<double>::epsilon() [static, private]

Machine epsilon.

LUmatrix* linalg::matrix::LUptr [mutable, private]

bool linalg::matrix::LUfactored [mutable, private]

double linalg::matrix::condition_number [mutable, private]

gsl_vector* linalg::matrix::SVD [mutable, private]

bool linalg::matrix::SVDfactored [mutable, private]


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

Generated on Fri Jun 6 17:28:29 2008 by  doxygen 1.5.6