GammaLib  2.1.0.dev
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GMatrixSparse Class Reference

Sparse matrix class interface definition. More...

#include <GMatrixSparse.hpp>

Inheritance diagram for GMatrixSparse:
GMatrixBase GBase

Public Member Functions

 GMatrixSparse (void)
 Void matrix constructor. More...
 
 GMatrixSparse (const int &rows, const int &columns, const int &elements=0)
 Matrix constructor. More...
 
 GMatrixSparse (const GMatrix &matrix)
 GMatrix to GMatrixSparse storage class convertor. More...
 
 GMatrixSparse (const GMatrixSparse &matrix)
 Copy constructor. More...
 
 GMatrixSparse (const GMatrixSymmetric &matrix)
 GMatrixSymmetric to GMatrixSparse storage class convertor. More...
 
virtual ~GMatrixSparse (void)
 Destructor. More...
 
virtual double & operator() (const int &row, const int &column)
 Return reference to matrix element. More...
 
virtual const double & operator() (const int &row, const int &column) const
 Return reference to matrix element (const version) More...
 
virtual GVector operator* (const GVector &vector) const
 Vector multiplication. More...
 
virtual bool operator== (const GMatrixSparse &matrix) const
 Equalty operator. More...
 
virtual bool operator!= (const GMatrixSparse &matrix) const
 Non-equality operator. More...
 
virtual GMatrixSparseoperator= (const GMatrixSparse &matrix)
 Matrix assignment operator. More...
 
virtual GMatrixSparseoperator= (const double &value)
 Value assignment operator. More...
 
virtual GMatrixSparse operator+ (const GMatrixSparse &matrix) const
 Binary matrix addition. More...
 
virtual GMatrixSparse operator+ (const double &scalar) const
 Binary matrix scalar addition. More...
 
virtual GMatrixSparse operator- (const GMatrixSparse &matrix) const
 Binary matrix subtraction. More...
 
virtual GMatrixSparse operator- (const double &scalar) const
 Binary matrix scalar subtraction. More...
 
virtual GMatrixSparse operator* (const GMatrixSparse &matrix) const
 Binary matrix multiplication. More...
 
virtual GMatrixSparse operator- (void) const
 Negate matrix elements. More...
 
virtual GMatrixSparseoperator+= (const GMatrixSparse &matrix)
 Unary matrix addition operator. More...
 
virtual GMatrixSparseoperator+= (const double &scalar)
 Unary matrix scalar addition operator. More...
 
virtual GMatrixSparseoperator-= (const GMatrixSparse &matrix)
 Unary matrix subtraction operator. More...
 
virtual GMatrixSparseoperator-= (const double &scalar)
 Unary matrix scalar subtraction operator. More...
 
virtual GMatrixSparseoperator*= (const GMatrixSparse &matrix)
 Unary matrix multiplication operator. More...
 
virtual GMatrixSparseoperator*= (const double &scalar)
 Scale matrix elements. More...
 
virtual GMatrixSparseoperator/= (const double &scalar)
 Divide matrix elements. More...
 
virtual void clear (void)
 Clear matrix. More...
 
virtual GMatrixSparseclone (void) const
 Clone matrix. More...
 
virtual std::string classname (void) const
 Return class name. More...
 
virtual double & at (const int &row, const int &column)
 Return reference to matrix element. More...
 
virtual const double & at (const int &row, const int &column) const
 Return reference to matrix element (const version) More...
 
virtual GVector row (const int &row) const
 Extract row as vector from matrix. More...
 
virtual void row (const int &row, const GVector &vector)
 Set row in matrix. More...
 
virtual GVector column (const int &column) const
 Extract column as vector from matrix. More...
 
virtual void column (const int &column, const GVector &vector)
 Insert vector column into matrix. More...
 
virtual void add_to_row (const int &row, const GVector &vector)
 Add row to matrix elements. More...
 
virtual void add_to_column (const int &column, const GVector &vector)
 Add vector column into matrix. More...
 
virtual double fill (void) const
 Returns fill of matrix. More...
 
virtual double min (void) const
 Return minimum matrix element. More...
 
virtual double max (void) const
 Return maximum matrix element. More...
 
virtual double sum (void) const
 Sum matrix elements. More...
 
virtual std::string print (const GChatter &chatter=NORMAL) const
 Print matrix. More...
 
void column (const int &column, const double *values, const int *rows, int number)
 Insert compressed array into matrix column. More...
 
void add_to_column (const int &column, const double *values, const int *rows, int number)
 Add compressed array into matrix column. More...
 
GMatrixSparse transpose (void) const
 Return transposed matrix. More...
 
GMatrixSparse invert (void) const
 Return inverted matrix. More...
 
GVector solve (const GVector &vector) const
 Solves linear matrix equation. More...
 
GMatrixSparse abs (void) const
 Return absolute of matrix. More...
 
GMatrixSparse cholesky_decompose (const bool &compress=true) const
 Return Cholesky decomposition. More...
 
GVector cholesky_solver (const GVector &vector, const bool &compress=true) const
 Cholesky solver. More...
 
GMatrixSparse cholesky_invert (const bool &compress=true) const
 Invert matrix using a Cholesky decomposition. More...
 
void set_mem_block (const int &block)
 Set memory block size. More...
 
void stack_init (const int &size=0, const int &entries=0)
 Initialises matrix filling stack. More...
 
int stack_push_column (const GVector &vector, const int &col)
 Push a vector column on the matrix stack. More...
 
int stack_push_column (const double *values, const int *rows, const int &number, const int &col)
 Push a compressed array on the matrix stack. More...
 
void stack_flush (void)
 Flush matrix stack. More...
 
void stack_destroy (void)
 Destroy matrix stack. More...
 
- Public Member Functions inherited from GMatrixBase
 GMatrixBase (void)
 Void constructor. More...
 
 GMatrixBase (const GMatrixBase &matrix)
 Copy constructor. More...
 
virtual ~GMatrixBase (void)
 Destructor. More...
 
virtual GMatrixBaseoperator= (const GMatrixBase &matrix)
 Assignment operator. More...
 
virtual bool operator== (const GMatrixBase &matrix) const
 Equalty operator. More...
 
virtual bool operator!= (const GMatrixBase &matrix) const
 Non-equality operator. More...
 
bool is_empty (void) const
 Check if matrix is empty. More...
 
const int & size (void) const
 Return number of matrix elements. More...
 
const int & columns (void) const
 Return number of matrix columns. More...
 
const int & rows (void) const
 Return number of matrix rows. More...
 
- Public Member Functions inherited from GBase
virtual ~GBase (void)
 Destructor. More...
 

Private Member Functions

void init_members (void)
 Initialise class mambers. More...
 
void copy_members (const GMatrixSparse &m)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void alloc_members (const int &rows, const int &cols, const int &elements=0)
 Allocate matrix. More...
 
void init_stack_members (void)
 Initialise fill stack. More...
 
void free_stack_members (void)
 Delete fill-stack class members. More...
 
int get_index (const int &row, const int &column) const
 Determines element index for (row,column) More...
 
void fill_pending (void)
 Fills pending matrix element. More...
 
void alloc_elements (int start, const int &num)
 Allocate memory for new matrix elements. More...
 
void free_elements (const int &start, const int &num)
 Free memory for obsolete matrix elements. More...
 
void remove_zero_row_col (void)
 Remove rows and columns with zeros. More...
 
void insert_zero_row_col (const int &rows, const int &cols)
 Insert zero rows and columns. More...
 
void mix_column_prepare (const int *src1_row, int src1_num, const int *src2_row, int src2_num, int *num_1, int *num_2, int *num_mix)
 Prepare mix of sparse columns. More...
 
void mix_column (const double *src1_data, const int *src1_row, int src1_num, const double *src2_data, const int *src2_row, int src2_num, double *dst_data, int *dst_row, int *dst_num)
 Mix of sparse columns. More...
 
void insert_row (const int &row, const GVector &vector, const bool &add)
 Insert row in matrix. More...
 

Private Attributes

int * m_rowinx
 Row-indices of all elements. More...
 
double m_zero
 The zero element (needed for data access) More...
 
double m_fill_val
 Element to be filled. More...
 
int m_fill_row
 Row into which element needs to be filled. More...
 
int m_fill_col
 Column into which element needs to be filled. More...
 
int m_mem_block
 Memory block to be allocated at once. More...
 
GSparseSymbolicm_symbolic
 Holds GSparseSymbolic object after decomposition. More...
 
GSparseNumericm_numeric
 Holds GSparseNumeric object after decomposition. More...
 
int m_stack_max_entries
 Maximum number of entries in the stack. More...
 
int m_stack_size
 Maximum number of elements in the stack. More...
 
int m_stack_entries
 Number of entries in the stack. More...
 
int * m_stack_colinx
 Column index for each entry [m_stack_entries]. More...
 
int * m_stack_start
 Start in stack for each entry [m_stack_entries+1]. More...
 
double * m_stack_data
 Stack data [m_stack_size]. More...
 
int * m_stack_rowinx
 Stack row indices [m_stack_size]. More...
 
int * m_stack_work
 Stack flush integer working array [m_cols]. More...
 
int * m_stack_rows
 Stack push integer working array [m_cols]. More...
 
double * m_stack_values
 Stack push double buffer [m_cols]. More...
 

Friends

class GSparseSymbolic
 
class GSparseNumeric
 
GMatrixSparse cs_symperm (const GMatrixSparse &matrix, const int *pinv)
 cs_symperm More...
 
GMatrixSparse cs_transpose (const GMatrixSparse &matrix, int values)
 
double cs_cumsum (int *p, int *c, int n)
 cs_cumsum More...
 

Additional Inherited Members

- Protected Member Functions inherited from GMatrixBase
void init_members (void)
 Initialise class members. More...
 
void copy_members (const GMatrixBase &matrix)
 Copy class members. More...
 
void free_members (void)
 Delete class members. More...
 
void select_non_zero (void)
 Setup compressed matrix factorisation. More...
 
void scale_elements (const double &scalar)
 Scale all matrix elements with a scalar. More...
 
void set_all_elements (const double &value)
 Set all elements to a specific value. More...
 
double get_min_element (void) const
 Return minimum matrix element. More...
 
double get_max_element (void) const
 Returns maximum matrix element. More...
 
double get_element_sum (void) const
 Returns sum over all matrix elements. More...
 
std::string print_elements (const GChatter &chatter=NORMAL, const int &num=15) const
 Print all matrix elements. More...
 
std::string print_row_compression (const GChatter &chatter=NORMAL) const
 Print row compression scheme if it exists. More...
 
std::string print_col_compression (const GChatter &chatter=NORMAL) const
 Print column compression scheme if it exists. More...
 
- Protected Attributes inherited from GMatrixBase
int m_rows
 Number of rows. More...
 
int m_cols
 Number of columns. More...
 
int m_elements
 Number of elements stored in matrix. More...
 
int m_alloc
 Size of allocated matrix memory. More...
 
int m_num_rowsel
 Number of selected rows (for compressed decomposition) More...
 
int m_num_colsel
 Number of selected columns (for compressed decomposition) More...
 
int * m_colstart
 Column start indices (m_cols+1) More...
 
int * m_rowsel
 Row selection (for compressed decomposition) More...
 
int * m_colsel
 Column selection (for compressed decomposition) More...
 
double * m_data
 Matrix data. More...
 

Detailed Description

Sparse matrix class interface definition.

This class implements a sparse matrix. The class only stores non-zero elements, which can considerably reduce the memory requirements for large systems that are sparsely filled. Also the computations will be faster because only non-zero elements will be used in the operations.

For a description of the common matrix methods, please refer to the GMatrixBase class.

The class has been inspired from the code CSparse that can be downloaded from http://www.cise.ufl.edu/research/sparse/CSparse/ which has been written by Timothy A. Davis to accompany his book "Direct Methods for Sparse Linear Systems". In the CSparse code, a sparse matrix is implemented by the structure cs_sparse. The members of this structure map as follows to the members of the GMatrixSparse class:

CSparse GMatrixSparse Dimension Description
nzmax m_elements n.a. Maximum number of entries
m m_rows n.a. Number of rows
n m_cols n.a. Number of columns
*p *m_colstart m_cols+1 Column pointers (index of first element of each column in m_data array)
*i *m_rowinx m_elements Row indices for all non-zero elements
*x *m_data m_elements Numerical values of all non-zero elements

Matrix allocation is done using the constructors

GMatrixSparse sparsematrix(rows, columns);
GMatrixSparse sparsematrix(rows, columns, elements);
GMatrixSparse sparsematrix(matrix);
GMatrixSparse sparsematrix(sparsematrix);
GMatrixSparse sparsematrix(symmetricmatrix);

where rows and columns specify the number of rows and columns of the matrix and elements is the size of the memory that should be pre-allocated for matrix storage (optional). Storage conversion constructors exist that allow allocating a sparse matrix by copying from a general matrix of type GMatrix and a symmetric matrix of type GMatrixSymmetric.

Because the location of a specific matrix element is not necessarily known in advance, the class implements a so called pending element for matrix element setting. A pending element is an element that is scheduled for insertion in the matrix, but which has not yet been inserted. The matrix element access operators () and at() return the reference to this pending element if a memory location for a given row and column does not yet exist. The protected method fill_pending() is then called before all operations that are done on the matrix to insert the element before the operation is executed.

Another mechanism that has been implemented to speed up the fill of the sparse matrix is a "fill stack" that is used for insertion or addition of matrix columns by the

matrix.column(index, vector);
matrix.column(index, values, rows, number);
matrix.add_to_column(index, vector);
matrix.add_to_column(index, values, rows, number);

methods. The fill stack is a buffer that implements a queue for columns that are to be inserted into the matrix. This is more efficient than insertion of elements one by one, as every insertion of an element may require to shift all elements back by one memory location.

The fill stack is used as follows:

matrix.stack_init(size, entries);
...
matrix.column(index, vector); // (or one of the other column methods)
...
matrix.stack_flush();
...
matrix.stack_destroy();

The stack_init(size, entries) method initialises the fill stack, where size is the size of the allocated memory buffer and entries is the maximum number of entries (i.e. columns) that will be held by the buffer. If size is set to 0 (the default value), a default size value of 512 is used. If entries is set to 0 (the default value), the number of matrix columns is taken as default entries value. Note that a too large number of elements will produce some overhead due to fill stack management, hence entries should not exceed a value of the order of 10-100.

The stack_flush() method flushes the stack, which is mandatory before any usage of the matrix. Note that the fill stack IS NOT INSERTED AUTOMATICALLY before any matrix operation, hence manual stack flushing is needed to make all filled matrix elements available for usage. The stack_destroy() method will flush the stack and free all stack elements. This method should be called once no filling is required anymore. If stack_destroy() is called immediately after filling, no call to stack_flush() is needed as the stack_destroy() method flushes the stack before destroying it. The matrix stack is also destroyed by the sparse matrix destructor, hence manual stack destruction is not mandatory.

Except for *m_rowinx which is implemented on the level of GMatrixSparse, all other members are implemented by the base class GMatrixBase.

Definition at line 187 of file GMatrixSparse.hpp.

Constructor & Destructor Documentation

GMatrixSparse::GMatrixSparse ( void  )

Void matrix constructor.

Constructs empty sparse matrix. The number of rows and columns of the matrix will be zero.

Definition at line 101 of file GMatrixSparse.cpp.

References init_members().

Referenced by clone().

GMatrixSparse::GMatrixSparse ( const int &  rows,
const int &  columns,
const int &  elements = 0 
)
explicit

Matrix constructor.

Parameters
[in]rowsNumber of rows [>=0].
[in]columnsNumber of columns [>=0].
[in]elementsNumber of allocated elements.
Exceptions
GException::invalid_argumentNumber of rows or columns is negative.

Constructs sparse matrix of dimension rows times columns. The optional elements argument specifies the size of the physical memory that should be allocated. By default, no memory will be allocated and memory allocation will be performed on-the-fly. If the amount of required memory is larger than the size specified by elements, additional momeory will be allocated automatically on-the-fly.

Definition at line 129 of file GMatrixSparse.cpp.

References alloc_members(), G_CONSTRUCTOR, init_members(), and gammalib::str().

GMatrixSparse::GMatrixSparse ( const GMatrix matrix)

GMatrix to GMatrixSparse storage class convertor.

Parameters
[in]matrixGeneric matrix (GMatrix).

Constructs a sparse matrix by using the number of rows and columns of a generic matrix and by assigning the elements of the generic matrix to the sparse matrix.

Definition at line 170 of file GMatrixSparse.cpp.

References alloc_members(), GMatrix::column(), column(), GMatrixBase::columns(), init_members(), GMatrixBase::m_cols, and GMatrixBase::rows().

GMatrixSparse::GMatrixSparse ( const GMatrixSparse matrix)

Copy constructor.

Parameters
[in]matrixMatrix.

Constructs matrix by copying an existing matrix.

Definition at line 224 of file GMatrixSparse.cpp.

References copy_members(), and init_members().

GMatrixSparse::GMatrixSparse ( const GMatrixSymmetric matrix)

GMatrixSymmetric to GMatrixSparse storage class convertor.

Parameters
[in]matrixSymmetric matrix (GMatrixSymmetric).

Constructs a sparse matrix by using the number of rows and columns of a symmetric matrix and by assigning the elements of the symmetric matrix to the sparse matrix.

Definition at line 198 of file GMatrixSparse.cpp.

References alloc_members(), GMatrixSymmetric::column(), column(), GMatrixBase::columns(), init_members(), GMatrixBase::m_cols, and GMatrixBase::rows().

GMatrixSparse::~GMatrixSparse ( void  )
virtual

Destructor.

Definition at line 241 of file GMatrixSparse.cpp.

References free_members().

Member Function Documentation

GMatrixSparse GMatrixSparse::abs ( void  ) const

Return absolute of matrix.

Returns
Absolute of matrix

Returns matrix where all elements of the matrix have been replaced by their absolute values.

Definition at line 1618 of file GMatrixSparse.cpp.

References abs(), fill_pending(), GMatrixBase::m_data, and GMatrixBase::m_elements.

void GMatrixSparse::add_to_column ( const int &  column,
const GVector vector 
)
virtual

Add vector column into matrix.

Parameters
[in]columnMatrix column [0,...,columns()[.
[in]vectorVector.
Exceptions
GException::out_of_rangeInvalid matrix column specified.
GException::invalid_argumentMatrix dimension mismatches the vector size.

Adds the contect of a vector to a matrix column.

This is the main driver routine to add data to a matrix. It handles both normal and stack-based filled. Note that there is another instance of this method that takes a compressed array.

Implements GMatrixBase.

Definition at line 1270 of file GMatrixSparse.cpp.

References column(), G_ADD_TO_COLUMN, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, m_fill_col, m_fill_row, m_fill_val, m_rowinx, GMatrixBase::m_rows, m_stack_data, GVector::non_zeros(), GVector::size(), stack_push_column(), and gammalib::str().

Referenced by GCTAOnOffObservation::likelihood_cstat(), GObservation::likelihood_gaussian_binned(), GObservation::likelihood_poisson_binned(), GObservation::likelihood_poisson_unbinned(), and GCTAOnOffObservation::likelihood_wstat().

void GMatrixSparse::add_to_column ( const int &  column,
const double *  values,
const int *  rows,
int  number 
)

Add compressed array into matrix column.

Parameters
[in]columnMatrix column [0,...,columns()[.
[in]valuesCompressed array.
[in]rowsRow indices of array.
[in]numberNumber of elements in array.
Exceptions
GException::out_of_rangeInvalid matrix column specified.
GException::invalid_argumentMatrix dimension mismatches the vector size.

Adds the content of a compressed array into a matrix column.

This is the main driver routine to add data to a matrix. It handles both normal and stack-based filled. Note that there is another instance of this method that takes a vector.

Definition at line 1393 of file GMatrixSparse.cpp.

References column(), fill_pending(), G_ADD_TO_COLUMN2, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, m_rowinx, GMatrixBase::m_rows, m_stack_data, mix_column(), gammalib::number(), stack_push_column(), and gammalib::str().

void GMatrixSparse::add_to_row ( const int &  row,
const GVector vector 
)
virtual

Add row to matrix elements.

Parameters
[in]rowMatrix row [0,...,row()[.
[in]vectorVector.

Implements GMatrixBase.

Definition at line 1243 of file GMatrixSparse.cpp.

References insert_row().

void GMatrixSparse::alloc_elements ( int  start,
const int &  num 
)
private

Allocate memory for new matrix elements.

Parameters
[in]startIndex of first allocated element.
[in]numNumber of elements to be allocated.

Inserts a memory allocation for 'num' elements at the index 'start'. The new elements are filled with 0.0 and the corresponding row indices are set to 0.

NOTE: This method does not take care of updating 'm_colstart'. This has to be done by the client.

Definition at line 3097 of file GMatrixSparse.cpp.

References GMatrixBase::m_alloc, GMatrixBase::m_data, GMatrixBase::m_elements, m_mem_block, and m_rowinx.

Referenced by alloc_members(), cholesky_decompose(), column(), GSparseSymbolic::cs_amd(), and fill_pending().

void GMatrixSparse::alloc_members ( const int &  rows,
const int &  columns,
const int &  elements = 0 
)
private

Allocate matrix.

Parameters
[in]rowsNumber of rows.
[in]columnsNumber of columns.
[in]elementsNumber of matrix elements to be physically allocated.

This is the main constructor code that allocates and initialises memory for matrix elements. The method only allocates elements if both rows and columns are positive. Otherwise the method does nothing and will set the m_rows and m_cols attributes to zero.

Definition at line 2833 of file GMatrixSparse.cpp.

References alloc_elements(), GMatrixBase::columns(), GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_rows, and GMatrixBase::rows().

Referenced by GMatrixSparse().

double & GMatrixSparse::at ( const int &  row,
const int &  column 
)
virtual

Return reference to matrix element.

Parameters
[in]rowMatrix row [0,...,rows()[.
[in]columnMatrix column [0,...,columns()[.
Returns
Reference to matrix element.
Exceptions
GException::out_of_rangeMatrix row or column out of range.

Implements GMatrixBase.

Definition at line 785 of file GMatrixSparse.cpp.

References G_AT, GMatrixBase::m_cols, and GMatrixBase::m_rows.

const double & GMatrixSparse::at ( const int &  row,
const int &  column 
) const
virtual

Return reference to matrix element (const version)

Parameters
[in]rowMatrix row [0,...,rows()[.
[in]columnMatrix column [0,...,columns()[.
Returns
Reference to matrix element.
Exceptions
GException::out_of_rangeMatrix row or column out of range.

Implements GMatrixBase.

Definition at line 810 of file GMatrixSparse.cpp.

References G_AT, GMatrixBase::m_cols, and GMatrixBase::m_rows.

GMatrixSparse GMatrixSparse::cholesky_decompose ( const bool &  compress = true) const

Return Cholesky decomposition.

Parameters
[in]compressUse zero-row/column compression (defaults to true).
Returns
Cholesky decomposition of matrix

Returns the Cholesky decomposition of a sparse matrix. The decomposition is stored within a GMatrixSparse object.

Definition at line 1748 of file GMatrixSparse.cpp.

References alloc_elements(), GSparseNumeric::cholesky_numeric_analysis(), GSparseSymbolic::cholesky_symbolic_analysis(), fill_pending(), free_elements(), GSparseSymbolic, insert_zero_row_col(), GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, GSparseNumeric::m_L, m_numeric, m_rowinx, GMatrixBase::m_rows, m_symbolic, and remove_zero_row_col().

Referenced by cholesky_invert(), GOptimizerLM::errors(), GObservations::errors_hessian(), and solve().

GMatrixSparse GMatrixSparse::cholesky_invert ( const bool &  compress = true) const

Invert matrix using a Cholesky decomposition.

Parameters
[in]compressUse zero-row/column compression.
Returns
Inverted matrix.

Inverts the matrix using a Cholesky decomposition.

Definition at line 2072 of file GMatrixSparse.cpp.

References cholesky_decompose(), cholesky_solver(), column(), GMatrixBase::m_cols, and GMatrixBase::m_rows.

Referenced by invert().

GVector GMatrixSparse::cholesky_solver ( const GVector vector,
const bool &  compress = true 
) const

Cholesky solver.

Parameters
[in]vectorSolution vector.
[in]compressRequest matrix compression.
Exceptions
GException::invalid_argumentMatrix and vector do not match.
GException::invalid_valueMatrix has not been factorised.

Solves the linear equation A*x=b using a Cholesky decomposition of A. This function is to be applied on a GMatrixSparse matrix for which a Choleksy factorization has been produced using 'cholesky_decompose'.

Definition at line 1832 of file GMatrixSparse.cpp.

References G_CHOL_SOLVE, iperm(), GMatrixBase::m_cols, GMatrixBase::m_colsel, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_num_colsel, GMatrixBase::m_num_rowsel, GSparseSymbolic::m_pinv, m_rowinx, GMatrixBase::m_rows, GMatrixBase::m_rowsel, m_symbolic, perm(), row(), GVector::size(), and gammalib::str().

Referenced by cholesky_invert(), GOptimizerLM::errors(), GObservations::errors_hessian(), and solve().

std::string GMatrixSparse::classname ( void  ) const
inlinevirtual

Return class name.

Returns
String containing the class name ("GMatrixSparse").

Implements GMatrixBase.

Definition at line 327 of file GMatrixSparse.hpp.

void GMatrixSparse::clear ( void  )
virtual

Clear matrix.

Implements GMatrixBase.

Definition at line 750 of file GMatrixSparse.cpp.

References free_members(), and init_members().

Referenced by GCTAEdispRmf::init_members(), and GRmf::init_members().

GMatrixSparse * GMatrixSparse::clone ( void  ) const
virtual

Clone matrix.

Returns
Pointer to deep copy of matrix.

Implements GMatrixBase.

Definition at line 768 of file GMatrixSparse.cpp.

References GMatrixSparse().

GVector GMatrixSparse::column ( const int &  column) const
virtual

Extract column as vector from matrix.

Parameters
[in]columnMatrix column [0,...,columns()[.
Returns
Vector of matrix column elements (rows() elements).
Exceptions
GException::out_of_rangeInvalid matrix column specified.

This method extracts a matrix column into a vector.

Implements GMatrixBase.

Definition at line 915 of file GMatrixSparse.cpp.

References G_EXTRACT_COLUMN, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, m_fill_col, m_fill_row, m_fill_val, m_rowinx, and GMatrixBase::m_rows.

Referenced by add_to_column(), cholesky_invert(), column(), GResponse::convolve(), GModelData::eval(), GResponse::eval_probs(), get_index(), GMatrixSparse(), GObservation::likelihood_gaussian_binned(), GObservation::likelihood_poisson_binned(), GObservation::likelihood_poisson_unbinned(), GObservation::model(), operator()(), operator+=(), operator-=(), and operator=().

void GMatrixSparse::column ( const int &  column,
const GVector vector 
)
virtual

Insert vector column into matrix.

Parameters
[in]columnMatrix column [0,...,columns()[.
[in]vectorVector.
Exceptions
GException::out_of_rangeInvalid matrix column specified.
GException::invalid_argumentVector size does not match number of matrix rows.

Inserts the content of a vector into a matrix column. Any previous content in the matrix column will be overwritted.

This is the main driver routine to insert data into a matrix. Note that there is another instance of this function that takes a compressed array.

Implements GMatrixBase.

Definition at line 964 of file GMatrixSparse.cpp.

References alloc_elements(), column(), free_elements(), G_SET_COLUMN, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, m_fill_col, m_fill_row, m_fill_val, m_rowinx, GMatrixBase::m_rows, row(), GVector::size(), and gammalib::str().

void GMatrixSparse::column ( const int &  column,
const double *  values,
const int *  rows,
int  number 
)

Insert compressed array into matrix column.

Parameters
[in]columnMatrix column [0,...,columns()[.
[in]valuesCompressed array.
[in]rowsRow numbers of array.
[in]numberNumber of elements in array.
Exceptions
GException::out_of_rangeInvalid matrix column specified.
GException::invalid_argumentRow numbers are not smaller than number of matrix rows

Inserts the content of a copressed array into a matrix column. Any previous content in the matrix column will be overwritted.

This is the main driver routine to insert data into a matrix. Note that there is another instance of this function that takes a vector.

Definition at line 1107 of file GMatrixSparse.cpp.

References alloc_elements(), column(), free_elements(), G_SET_COLUMN2, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, m_fill_col, m_fill_row, m_fill_val, m_rowinx, GMatrixBase::m_rows, gammalib::number(), row(), and gammalib::str().

void GMatrixSparse::copy_members ( const GMatrixSparse matrix)
private

Copy class members.

Parameters
[in]matrixMatrix to be copied.

Definition at line 2764 of file GMatrixSparse.cpp.

References GSparseNumeric, GSparseSymbolic, GMatrixBase::m_alloc, GMatrixBase::m_elements, m_fill_col, m_fill_row, m_fill_val, m_mem_block, m_numeric, m_rowinx, m_symbolic, and m_zero.

Referenced by GMatrixSparse(), and operator=().

double GMatrixSparse::fill ( void  ) const
virtual

Returns fill of matrix.

The fill of a matrix is defined as the number non-zero elements devided by the number of total elements. By definition, the fill is comprised in the interval [0,..,1]. The fill of an undefined matrix is defined to be 0.

Implements GMatrixBase.

Definition at line 1645 of file GMatrixSparse.cpp.

References GMatrixBase::m_cols, GMatrixBase::m_elements, m_fill_val, GMatrixBase::m_rows, and GMatrixBase::size().

Referenced by print().

void GMatrixSparse::fill_pending ( void  )
private

Fills pending matrix element.

If 'm_fill_val' is non-zero a pending matrix element exists that should be filled into (row,col)=(m_fill_row,m_fill_col). This routine performs the filling of the matrix with this element and resets 'm_fill_val' to zero. This routine allows for element-by-element filling of a sparse matrix. This is, of course, very time consuming and should in general be avoided. However, it allows to design a sparse matrix class that hides the matrix sparsity completely to the user.

Definition at line 3011 of file GMatrixSparse.cpp.

References alloc_elements(), GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, m_fill_col, m_fill_row, m_fill_val, and m_rowinx.

Referenced by abs(), add_to_column(), cholesky_decompose(), insert_row(), insert_zero_row_col(), invert(), operator()(), operator*=(), operator-(), operator=(), remove_zero_row_col(), stack_flush(), and transpose().

void GMatrixSparse::free_elements ( const int &  start,
const int &  num 
)
private

Free memory for obsolete matrix elements.

Parameters
[in]startIndex of first freed element.
[in]numNumber of elements to be freed.

Free memory for 'num' elements starting from index 'start'.

NOTE: This method does not take care of updating 'm_colstart'. This has to be done by the client.

Definition at line 3202 of file GMatrixSparse.cpp.

References GMatrixBase::m_alloc, GMatrixBase::m_data, GMatrixBase::m_elements, m_mem_block, and m_rowinx.

Referenced by cholesky_decompose(), column(), GSparseSymbolic::cs_fkeep(), and cs_symperm().

void GMatrixSparse::free_members ( void  )
private

Delete class members.

Definition at line 2801 of file GMatrixSparse.cpp.

References free_stack_members(), m_numeric, m_rowinx, and m_symbolic.

Referenced by clear(), operator=(), and ~GMatrixSparse().

void GMatrixSparse::free_stack_members ( void  )
private

Delete fill-stack class members.

Deletes all memory that has been allocated using the init_stack_members() method.

Definition at line 2901 of file GMatrixSparse.cpp.

References m_stack_colinx, m_stack_data, m_stack_rowinx, m_stack_rows, m_stack_start, m_stack_values, and m_stack_work.

Referenced by free_members(), stack_destroy(), and stack_init().

int GMatrixSparse::get_index ( const int &  row,
const int &  column 
) const
private

Determines element index for (row,column)

Parameters
[in]rowMatrix row [0,...,m_rows[.
[in]columnMatrix column [0,...,m_cols[.
Returns
Element index.
Exceptions
GException::out_of_rangeMatrix row or column out of range

Returns the index in the compressed array for (row,col). The following special results exist:

  -1: Requested index does not exist in the matrix.
  m_elements: Requested index is the pending element.

Definition at line 2942 of file GMatrixSparse.cpp.

References column(), G_GET_INDEX, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_elements, m_fill_col, m_fill_row, m_fill_val, m_rowinx, and GMatrixBase::m_rows.

Referenced by operator()().

void GMatrixSparse::init_members ( void  )
private
void GMatrixSparse::init_stack_members ( void  )
private

Initialise fill stack.

The fill stack is used to fill the sparse matrix without any prior know- ledge about the number of location of the non-zero matrix elements. Values are first filled into the stack and flushed into the sparse matrix once the stack is full.

Definition at line 2876 of file GMatrixSparse.cpp.

References m_stack_colinx, m_stack_data, m_stack_entries, m_stack_max_entries, m_stack_rowinx, m_stack_rows, m_stack_size, m_stack_start, m_stack_values, and m_stack_work.

Referenced by init_members(), stack_destroy(), and stack_init().

void GMatrixSparse::insert_row ( const int &  row,
const GVector vector,
const bool &  add 
)
private

Insert row in matrix.

Parameters
[in]rowMatrix row [0,...,row()[.
[in]vectorVector.
[in]addAdd vector to existing elements?
Exceptions
GException::out_of_rangeInvalid matrix row specified.
GException::invalid_argumentVector size incompatible with number of matrix columns.
Todo:
Remove elements that are empty after addition

Definition at line 3675 of file GMatrixSparse.cpp.

References fill_pending(), G_INSERT_ROW, G_SET_ROW, GMatrixBase::m_alloc, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, m_mem_block, m_rowinx, GMatrixBase::m_rows, row(), GVector::size(), and gammalib::str().

Referenced by add_to_row(), and row().

void GMatrixSparse::insert_zero_row_col ( const int &  rows,
const int &  cols 
)
private

Insert zero rows and columns.

Parameters
[in]rowsNumber of rows.
[in]colsNumber of columns.

Insert zero rows and columns into matrix. Since for a sparse matrix this does not require any allocation of additional memory, the data are not moved by this function, but the pointers are re-arranged.

Definition at line 3404 of file GMatrixSparse.cpp.

References fill_pending(), GMatrixBase::m_cols, GMatrixBase::m_colsel, GMatrixBase::m_colstart, GMatrixBase::m_elements, GMatrixBase::m_num_colsel, m_rowinx, GMatrixBase::m_rows, GMatrixBase::m_rowsel, and GMatrixBase::rows().

Referenced by cholesky_decompose().

GMatrixSparse GMatrixSparse::invert ( void  ) const

Return inverted matrix.

Returns
Inverted matrix.

Returns inverse of matrix. Inversion is done for the moment using Cholesky decomposition. This does not work on any kind of matrix.

Todo:
Specify in documentation for which kind of matrix the method works.

Definition at line 1554 of file GMatrixSparse.cpp.

References cholesky_invert(), and fill_pending().

Referenced by GObservations::likelihood::covariance().

double GMatrixSparse::max ( void  ) const
virtual

Return maximum matrix element.

Implements GMatrixBase.

Definition at line 1698 of file GMatrixSparse.cpp.

References GMatrixBase::m_cols, GMatrixBase::m_data, GMatrixBase::m_elements, m_fill_val, and GMatrixBase::m_rows.

double GMatrixSparse::min ( void  ) const
virtual

Return minimum matrix element.

Implements GMatrixBase.

Definition at line 1672 of file GMatrixSparse.cpp.

References GMatrixBase::m_cols, GMatrixBase::m_data, GMatrixBase::m_elements, m_fill_val, and GMatrixBase::m_rows.

void GMatrixSparse::mix_column ( const double *  src1_data,
const int *  src1_row,
int  src1_num,
const double *  src2_data,
const int *  src2_row,
int  src2_num,
double *  dst_data,
int *  dst_row,
int *  dst_num 
)
private

Mix of sparse columns.

Parameters
[in]src1_dataData array [0...src1_num-1] of first column.
[in]src1_rowRow index array [0...src1_num-1] of first column.
[in]src1_numNumber of elements in first column.
[in]src2_dataData array [0...src2_num-1] of second column.
[in]src2_rowRow index array [0...src2_num-1] of second column.
[in]src2_numNumber of elements in second column.
[in]dst_dataData array [0...dst_num-1] of mixed column.
[in]dst_rowRow index array [0...dst_num-1] of mixed column.
[out]dst_numNumber of elements in mixed column.

This method mixes two sparse matrix columns into a single column.

Definition at line 3581 of file GMatrixSparse.cpp.

Referenced by add_to_column(), stack_flush(), and stack_push_column().

void GMatrixSparse::mix_column_prepare ( const int *  src1_row,
int  src1_num,
const int *  src2_row,
int  src2_num,
int *  num_1,
int *  num_2,
int *  num_mix 
)
private

Prepare mix of sparse columns.

Parameters
[in]src1_rowRow index array [0...src1_num-1] of first column.
[in]src1_numNumber of elements in first column.
[in]src2_rowRow index array [0...src2_num-1] of second column.
[in]src2_numNumber of elements in second column.
[out]num_1Number of elements only found in column 1.
[out]num_2Number of elements only found in column 2.
[out]num_mixNumber of elements found in both columns.

This method prepares the mix of two sparse matrix columns into a single column. It counts the number of elements in two columns that have the same row and puts this number into num_mix. The number of elements that are only present in the first column will be put in num_1, the number of elements that are only present in the second column will be put in num_2.

Definition at line 3496 of file GMatrixSparse.cpp.

Referenced by stack_flush(), and stack_push_column().

bool GMatrixSparse::operator!= ( const GMatrixSparse matrix) const
virtual

Non-equality operator.

Parameters
[in]matrixMatrix.

This operator checks if two matrices are not identical. Two matrices are considered not identical if they differ in their dimensions or if at least one element differs.

Definition at line 548 of file GMatrixSparse.cpp.

References operator==().

double & GMatrixSparse::operator() ( const int &  row,
const int &  column 
)
virtual

Return reference to matrix element.

Parameters
[in]rowMatrix row [0,...,rows()-1].
[in]columnMatrix column [0,...,columns()-1].
Returns
Reference to matrix element.

Returns the reference to the matrix element at row and column. If the matrix element does not yet exist, a reference to the pending element with a value of 0.0 is returned.

If the matrix row or column are invalid an exception is thrown by the get_index() method.

Implements GMatrixBase.

Definition at line 352 of file GMatrixSparse.cpp.

References column(), fill_pending(), get_index(), GMatrixBase::m_data, m_fill_col, m_fill_row, m_fill_val, and row().

const double & GMatrixSparse::operator() ( const int &  row,
const int &  column 
) const
virtual

Return reference to matrix element (const version)

Parameters
[in]rowMatrix row [0,...,rows()-1].
[in]columnMatrix column [0,...,columns()-1].
Returns
Const reference to matrix element.

Returns a const reference to the matrix element at row and column. If the matrix element does not yet exist, a reference to the zero element is returned. If the matrix element corresponds to the pending element, a reference to the pending element is returned. Otherwise, a reference to the matrix elements is returned.

If the matrix row or column are invalid an exception is thrown by the get_index() method.

Implements GMatrixBase.

Definition at line 391 of file GMatrixSparse.cpp.

References get_index(), GMatrixBase::m_data, GMatrixBase::m_elements, m_fill_val, and m_zero.

GVector GMatrixSparse::operator* ( const GVector vector) const
virtual

Vector multiplication.

Parameters
[in]vectorVector.
Exceptions
GException::invalid_argumentVector size differs from number of columns in matrix.

This method performs a vector multiplication of a matrix. The vector multiplication will produce a vector. The matrix multiplication can only be performed when the number of matrix columns is equal to the length of the vector.

The method includes any pending fill.

Implements GMatrixBase.

Definition at line 429 of file GMatrixSparse.cpp.

References G_OP_MUL_VEC, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, m_fill_col, m_fill_row, m_fill_val, m_rowinx, GMatrixBase::m_rows, GVector::size(), and gammalib::str().

GMatrixSparse GMatrixSparse::operator* ( const GMatrixSparse matrix) const
inlinevirtual

Binary matrix multiplication.

Parameters
[in]matrixMatrix.
Returns
Result of matrix multiplication.

Returns the product of two matrices. The method makes use of the unary multiplication operator.

Definition at line 429 of file GMatrixSparse.hpp.

GMatrixSparse & GMatrixSparse::operator*= ( const GMatrixSparse matrix)
virtual

Unary matrix multiplication operator.

Parameters
[in]matrixMatrix.
Exceptions
GException::invalid_argumentNumber of rows in matrix is different from number of matrix columns.

This method performs a matrix multiplication. The operation can only succeed when the dimensions of both matrices are compatible.

Todo:
Implement native sparse code

Definition at line 702 of file GMatrixSparse.cpp.

References G_OP_MAT_MUL, GMatrixBase::m_cols, GMatrixBase::m_elements, GMatrixBase::m_rows, row(), gammalib::str(), and sum().

GMatrixSparse & GMatrixSparse::operator*= ( const double &  scalar)
inlinevirtual

Scale matrix elements.

Parameters
[in]scalarScale factor.
Returns
Matrix with elements multiplied by scalar.

Returns a matrix where all elements have been multiplied by the specified scalar value.

Definition at line 447 of file GMatrixSparse.hpp.

References fill_pending(), and GMatrixBase::scale_elements().

GMatrixSparse GMatrixSparse::operator+ ( const GMatrixSparse matrix) const
inlinevirtual

Binary matrix addition.

Parameters
[in]matrixMatrix.
Returns
Result of matrix addition.

Returns the sum of two matrices. The method makes use of the unary addition operator.

Definition at line 358 of file GMatrixSparse.hpp.

GMatrixSparse GMatrixSparse::operator+ ( const double &  scalar) const
inlinevirtual

Binary matrix scalar addition.

Parameters
[in]scalarScalar.
Returns
Matrix with scalar added.

Returns a matrix where a scalar has been added to each matrix element.

Definition at line 375 of file GMatrixSparse.hpp.

GMatrixSparse & GMatrixSparse::operator+= ( const GMatrixSparse matrix)
virtual

Unary matrix addition operator.

Parameters
[in]matrixMatrix.
Exceptions
GException::invalid_argumentIncompatible number of matrix rows or columns.

This method performs a matrix addition. The operation can only succeed when the dimensions of both matrices are identical.

Todo:
Implement native sparse code

Definition at line 571 of file GMatrixSparse.cpp.

References column(), G_OP_ADD, GMatrixBase::m_cols, GMatrixBase::m_rows, and gammalib::str().

GMatrixSparse & GMatrixSparse::operator+= ( const double &  scalar)
virtual

Unary matrix scalar addition operator.

Parameters
[in]scalarScalar.

Adds a scalar to each matrix element.

Definition at line 609 of file GMatrixSparse.cpp.

References column(), and GMatrixBase::m_cols.

GMatrixSparse GMatrixSparse::operator- ( const GMatrixSparse matrix) const
inlinevirtual

Binary matrix subtraction.

Parameters
[in]matrixMatrix.
Returns
Result of matrix subtraction.

Returns the difference between two matrices. The method makes use of the unary subtraction operator.

Definition at line 393 of file GMatrixSparse.hpp.

GMatrixSparse GMatrixSparse::operator- ( const double &  scalar) const
inlinevirtual

Binary matrix scalar subtraction.

Parameters
[in]scalarScalar.
Returns
Matrix with scalar subtracted.

Returns a matrix where a scalar has been subtracted from each matrix element.

Definition at line 411 of file GMatrixSparse.hpp.

GMatrixSparse GMatrixSparse::operator- ( void  ) const
virtual

Negate matrix elements.

Returns
Matrix with negated elements.

Definition at line 479 of file GMatrixSparse.cpp.

References fill_pending(), GMatrixBase::m_data, and GMatrixBase::m_elements.

GMatrixSparse & GMatrixSparse::operator-= ( const GMatrixSparse matrix)
virtual

Unary matrix subtraction operator.

Parameters
[in]matrixMatrix.
Exceptions
GException::invalid_argumentIncompatible matrix size.

This method performs a matrix addition. The operation can only succeed when the dimensions of both matrices are identical.

Todo:
Implement native sparse code

Definition at line 636 of file GMatrixSparse.cpp.

References column(), G_OP_SUB, GMatrixBase::m_cols, GMatrixBase::m_rows, and gammalib::str().

GMatrixSparse & GMatrixSparse::operator-= ( const double &  scalar)
virtual

Unary matrix scalar subtraction operator.

Parameters
[in]scalarScalar.

Subtracts a scalar from each matrix element.

Definition at line 674 of file GMatrixSparse.cpp.

References column(), and GMatrixBase::m_cols.

GMatrixSparse & GMatrixSparse::operator/= ( const double &  scalar)
inlinevirtual

Divide matrix elements.

Parameters
[in]scalarScalar.
Returns
Matrix with elements divided by scalar.

Returns a matrix where all elements have been divided by the specified scalar value.

Definition at line 465 of file GMatrixSparse.hpp.

GMatrixSparse & GMatrixSparse::operator= ( const GMatrixSparse matrix)
virtual

Matrix assignment operator.

Parameters
[in]matrixMatrix.
Returns
Matrix.

Assigns the content of another matrix to the actual matrix instance.

Definition at line 265 of file GMatrixSparse.cpp.

References copy_members(), free_members(), init_members(), and GMatrixBase::operator=().

GMatrixSparse & GMatrixSparse::operator= ( const double &  value)
virtual

Value assignment operator.

Parameters
[in]valueValue.
Returns
Matrix.

Assigns the specified value to all elements of the matrix.

Definition at line 299 of file GMatrixSparse.cpp.

References column(), fill_pending(), GMatrixBase::is_empty(), GMatrixBase::m_cols, GMatrixBase::m_colstart, and GMatrixBase::m_rows.

bool GMatrixSparse::operator== ( const GMatrixSparse matrix) const
virtual

Equalty operator.

Parameters
[in]matrixMatrix.

This operator checks if two matrices are identical. Two matrices are considered identical if they have the same dimensions and identicial elements.

Todo:
Implement native sparse code

Definition at line 509 of file GMatrixSparse.cpp.

References GMatrixBase::m_cols, GMatrixBase::m_rows, and row().

Referenced by operator!=().

void GMatrixSparse::remove_zero_row_col ( void  )
private

Remove rows and columns with zeros.

Exceptions
GException::runtime_errorAll matrix elements are zero.

Remove all rows and columns that contain only zeros from matrix. This function is needed for compressed matrix factorisation. The resulting matrix has reduced size (number of rows and columns) and dimensions (number of elements). Note that the physical memory is not reduced by the method.

Definition at line 3298 of file GMatrixSparse.cpp.

References fill_pending(), G_REMOVE_ZERO, GMatrixBase::m_cols, GMatrixBase::m_colsel, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_num_colsel, GMatrixBase::m_num_rowsel, m_rowinx, GMatrixBase::m_rows, GMatrixBase::m_rowsel, row(), and GMatrixBase::select_non_zero().

Referenced by cholesky_decompose().

GVector GMatrixSparse::row ( const int &  row) const
virtual

Extract row as vector from matrix.

Parameters
[in]rowMatrix row [0,...,rows()[.
Returns
Vector of matrix row elements (columns() elements).
Exceptions
GException::out_of_rangeInvalid matrix row specified.

This method extracts a matrix row into a vector.

Implements GMatrixBase.

Definition at line 836 of file GMatrixSparse.cpp.

References G_EXTRACT_ROW, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, m_fill_col, m_fill_row, m_fill_val, m_rowinx, and GMatrixBase::m_rows.

Referenced by cholesky_solver(), column(), insert_row(), operator()(), operator*=(), operator==(), and remove_zero_row_col().

void GMatrixSparse::row ( const int &  row,
const GVector vector 
)
virtual

Set row in matrix.

Parameters
[in]rowMatrix row [0,...,rows()[.
[in]vectorVector of matrix row elements (columns() elements).

Implements GMatrixBase.

Definition at line 894 of file GMatrixSparse.cpp.

References insert_row().

void GMatrixSparse::set_mem_block ( const int &  block)
inline

Set memory block size.

Parameters
[in]blockMemory block size.

Sets the size of the memory block that will be allocated at once.

Definition at line 341 of file GMatrixSparse.hpp.

References m_mem_block.

GVector GMatrixSparse::solve ( const GVector vector) const

Solves linear matrix equation.

Parameters
[in]vectorSolution vector.
Returns
Solutions of matrix equation.

Solves the linear equation

\[M \times {\tt solution} = {\tt vector} \]

where \(M\) is the matrix, \({\tt vector}\) is the result, and \({\tt solution}\) is the solution. Solving is done using Cholesky decomposition. This does not work on any kind of matrix.

If the matrix is empty and the vector has a zero length the method returns an empty vector.

Todo:
Specify in documentation for which kind of matrix the method works.

Definition at line 1589 of file GMatrixSparse.cpp.

References cholesky_decompose(), cholesky_solver(), GMatrixBase::is_empty(), and GVector::size().

Referenced by GOptimizerLM::iteration().

void GMatrixSparse::stack_destroy ( void  )

Destroy matrix stack.

Flush and destroy matrix stack

Definition at line 2714 of file GMatrixSparse.cpp.

References free_stack_members(), init_stack_members(), and stack_flush().

Referenced by GObservations::likelihood::eval().

void GMatrixSparse::stack_flush ( void  )

Flush matrix stack.

Adds the content of the stack to the actual matrix. First, the total number of matrix elements is determined. Then new memory is allocated to hold all elements. Finally, all elements are filled into the new memory that is then replacing the old matrix memory.

The method uses the internal working array m_stack_work.

NOTE: The present algorithm assumes that each column occurs only once in the stack!

Definition at line 2485 of file GMatrixSparse.cpp.

References fill_pending(), GMatrixBase::is_empty(), GMatrixBase::m_alloc, GMatrixBase::m_cols, GMatrixBase::m_colstart, GMatrixBase::m_data, GMatrixBase::m_elements, m_mem_block, m_rowinx, m_stack_colinx, m_stack_data, m_stack_entries, m_stack_rowinx, m_stack_start, m_stack_work, mix_column(), and mix_column_prepare().

Referenced by stack_destroy(), and stack_push_column().

void GMatrixSparse::stack_init ( const int &  size = 0,
const int &  entries = 0 
)

Initialises matrix filling stack.

Parameters
[in]sizeStack size.
[in]entriesMaximum number of entries.

The matrix filling stack is used to allow for a fast column-wise filling of a sparse matrix. Columns are successively appended to a stack which is regularily flushed when it is full. This reduces memory copies and movements and increases filling speed.

Definition at line 2189 of file GMatrixSparse.cpp.

References free_stack_members(), G_SPARSE_MATRIX_DEFAULT_STACK_SIZE, init_stack_members(), GMatrixBase::m_cols, m_stack_colinx, m_stack_data, m_stack_max_entries, m_stack_rowinx, m_stack_rows, m_stack_size, m_stack_start, m_stack_values, and m_stack_work.

Referenced by GObservations::likelihood::eval(), and GObservation::model().

int GMatrixSparse::stack_push_column ( const GVector vector,
const int &  col 
)

Push a vector column on the matrix stack.

Parameters
[in]vectorVector.
[in]colMatrix column [0,...,columns()[.
Exceptions
GException::out_of_rangeInvalid matrix column specified.
GException::invalid_argumentMatrix dimension mismatches the vector size.

Adds the contect of a vector to the matrix stack. This method is identical to the GMatrixSparse::stack_push_column method that uses a compressed array, yet it takes a full vector. Internal working arrays are used to convert the full column vector in a compressed array and to hand it over to the compressed array version.

Definition at line 2233 of file GMatrixSparse.cpp.

References G_STACK_PUSH1, GMatrixBase::m_cols, GMatrixBase::m_rows, m_stack_rows, m_stack_values, GVector::size(), and gammalib::str().

Referenced by add_to_column().

int GMatrixSparse::stack_push_column ( const double *  values,
const int *  rows,
const int &  number,
const int &  col 
)

Push a compressed array on the matrix stack.

Parameters
[in]valuesCompressed array.
[in]rowsMatrix rows of array elements.
[in]numberNumber of elements in array.
[in]colMatrix column [0,...,columns()[.
Exceptions
GException::out_of_rangeInvalid matrix column specified.
GException::invalid_argumentMatrix dimension mismatches the vector size.

The method puts new data on the stack while assuring that each column is present only once. If an already existing column is encountered the data from the existing and new column are mixed and put into a new entry one the stack; the old entry is signalled as obsolete.

On return the method indicates the number of elements in the input array that have not been put onto the stack (due to memory limits). This value can be either '0' (i.e. all elements have been put on the stack) or 'number' (i.e. none of the elements have been put on the stack). Columns are not partially put onto the stack.

Definition at line 2296 of file GMatrixSparse.cpp.

References G_STACK_PUSH2, GMatrixBase::m_cols, GMatrixBase::m_rows, m_stack_colinx, m_stack_data, m_stack_entries, m_stack_max_entries, m_stack_rowinx, m_stack_size, m_stack_start, mix_column(), mix_column_prepare(), gammalib::number(), stack_flush(), and gammalib::str().

double GMatrixSparse::sum ( void  ) const
virtual

Sum matrix elements.

Implements GMatrixBase.

Definition at line 1724 of file GMatrixSparse.cpp.

References GMatrixBase::m_data, GMatrixBase::m_elements, and m_fill_val.

Referenced by operator*=().

GMatrixSparse GMatrixSparse::transpose ( void  ) const

Return transposed matrix.

Returns
Transposed matrix.

Returns transposed matrix of the matrix.

Definition at line 1528 of file GMatrixSparse.cpp.

References cs_transpose, and fill_pending().

Referenced by GObservation::likelihood_gaussian_binned(), GObservation::likelihood_poisson_binned(), and GObservation::likelihood_poisson_unbinned().

Friends And Related Function Documentation

double cs_cumsum ( int *  p,
int *  c,
int  n 
)
friend

cs_cumsum

Parameters
[out]pInteger array (n+1 elements).
[in]cInteger array (n elements).
[in]nNumber of elements.

Evaluate p[0..n] = cumulative sum of c[0..n-1]

Definition at line 3988 of file GMatrixSparse.cpp.

GMatrixSparse cs_symperm ( const GMatrixSparse matrix,
const int *  pinv 
)
friend

cs_symperm

Parameters
[in]matrixMatrix.
[in]pinvTBD.

Returns matrix(p,p) where matrix and matrix(p,p) are symmetric the upper part stored.

Definition at line 3835 of file GMatrixSparse.cpp.

GMatrixSparse cs_transpose ( const GMatrixSparse matrix,
int  values 
)
friend

Definition at line 3921 of file GMatrixSparse.cpp.

Referenced by transpose().

friend class GSparseNumeric
friend

Definition at line 191 of file GMatrixSparse.hpp.

Referenced by copy_members().

friend class GSparseSymbolic
friend

Definition at line 190 of file GMatrixSparse.hpp.

Referenced by cholesky_decompose(), and copy_members().

Member Data Documentation

int GMatrixSparse::m_fill_col
private

Column into which element needs to be filled.

Definition at line 302 of file GMatrixSparse.hpp.

Referenced by add_to_column(), column(), copy_members(), fill_pending(), get_index(), init_members(), operator()(), operator*(), print(), and row().

int GMatrixSparse::m_fill_row
private

Row into which element needs to be filled.

Definition at line 301 of file GMatrixSparse.hpp.

Referenced by add_to_column(), column(), copy_members(), fill_pending(), get_index(), init_members(), operator()(), operator*(), print(), and row().

double GMatrixSparse::m_fill_val
private
int GMatrixSparse::m_mem_block
private

Memory block to be allocated at once.

Definition at line 303 of file GMatrixSparse.hpp.

Referenced by alloc_elements(), copy_members(), free_elements(), init_members(), insert_row(), print(), set_mem_block(), and stack_flush().

GSparseNumeric* GMatrixSparse::m_numeric
private

Holds GSparseNumeric object after decomposition.

Definition at line 305 of file GMatrixSparse.hpp.

Referenced by cholesky_decompose(), copy_members(), free_members(), and init_members().

int* GMatrixSparse::m_stack_colinx
private

Column index for each entry [m_stack_entries].

Definition at line 311 of file GMatrixSparse.hpp.

Referenced by free_stack_members(), init_stack_members(), stack_flush(), stack_init(), and stack_push_column().

double* GMatrixSparse::m_stack_data
private

Stack data [m_stack_size].

Definition at line 313 of file GMatrixSparse.hpp.

Referenced by add_to_column(), free_stack_members(), init_stack_members(), print(), stack_flush(), stack_init(), and stack_push_column().

int GMatrixSparse::m_stack_entries
private

Number of entries in the stack.

Definition at line 310 of file GMatrixSparse.hpp.

Referenced by init_stack_members(), stack_flush(), and stack_push_column().

int GMatrixSparse::m_stack_max_entries
private

Maximum number of entries in the stack.

Definition at line 308 of file GMatrixSparse.hpp.

Referenced by init_stack_members(), stack_init(), and stack_push_column().

int* GMatrixSparse::m_stack_rowinx
private

Stack row indices [m_stack_size].

Definition at line 314 of file GMatrixSparse.hpp.

Referenced by free_stack_members(), init_stack_members(), stack_flush(), stack_init(), and stack_push_column().

int* GMatrixSparse::m_stack_rows
private

Stack push integer working array [m_cols].

Definition at line 316 of file GMatrixSparse.hpp.

Referenced by free_stack_members(), init_stack_members(), stack_init(), and stack_push_column().

int GMatrixSparse::m_stack_size
private

Maximum number of elements in the stack.

Definition at line 309 of file GMatrixSparse.hpp.

Referenced by init_stack_members(), print(), stack_init(), and stack_push_column().

int* GMatrixSparse::m_stack_start
private

Start in stack for each entry [m_stack_entries+1].

Definition at line 312 of file GMatrixSparse.hpp.

Referenced by free_stack_members(), init_stack_members(), stack_flush(), stack_init(), and stack_push_column().

double* GMatrixSparse::m_stack_values
private

Stack push double buffer [m_cols].

Definition at line 317 of file GMatrixSparse.hpp.

Referenced by free_stack_members(), init_stack_members(), stack_init(), and stack_push_column().

int* GMatrixSparse::m_stack_work
private

Stack flush integer working array [m_cols].

Definition at line 315 of file GMatrixSparse.hpp.

Referenced by free_stack_members(), init_stack_members(), stack_flush(), and stack_init().

GSparseSymbolic* GMatrixSparse::m_symbolic
private

Holds GSparseSymbolic object after decomposition.

Definition at line 304 of file GMatrixSparse.hpp.

Referenced by cholesky_decompose(), cholesky_solver(), copy_members(), free_members(), and init_members().

double GMatrixSparse::m_zero
private

The zero element (needed for data access)

Definition at line 299 of file GMatrixSparse.hpp.

Referenced by copy_members(), init_members(), and operator()().


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