Class SUNMatrixWrapper

Class Documentation

class amici::SUNMatrixWrapper

A RAII wrapper for SUNMatrix structs.

This can create dense, sparse, or banded matrices using the respective constructor.

Public Functions

SUNMatrixWrapper() = default
SUNMatrixWrapper(sunindextype M, sunindextype N, sunindextype NNZ, int sparsetype)

Create sparse matrix. See SUNSparseMatrix in sunmatrix_sparse.h.

Parameters
  • M – Number of rows

  • N – Number of columns

  • NNZ – Number of nonzeros

  • sparsetype – Sparse type

SUNMatrixWrapper(sunindextype M, sunindextype N)

Create dense matrix. See SUNDenseMatrix in sunmatrix_dense.h.

Parameters
  • M – Number of rows

  • N – Number of columns

SUNMatrixWrapper(sunindextype M, sunindextype ubw, sunindextype lbw)

Create banded matrix. See SUNBandMatrix in sunmatrix_band.h.

Parameters
  • M – Number of rows and columns

  • ubw – Upper bandwidth

  • lbw – Lower bandwidth

SUNMatrixWrapper(const SUNMatrixWrapper &A, realtype droptol, int sparsetype)

Create sparse matrix from dense or banded matrix. See SUNSparseFromDenseMatrix and SUNSparseFromBandMatrix in sunmatrix_sparse.h.

Parameters
  • A – Wrapper for dense matrix

  • droptol – tolerance for dropping entries

  • sparsetype – Sparse type

explicit SUNMatrixWrapper(SUNMatrix mat)

Wrap existing SUNMatrix.

Parameters

mat

~SUNMatrixWrapper()
SUNMatrixWrapper(const SUNMatrixWrapper &other)

Copy constructor.

Parameters

other

SUNMatrixWrapper(SUNMatrixWrapper &&other)

Move constructor.

Parameters

other

SUNMatrixWrapper &operator=(const SUNMatrixWrapper &other)

Copy assignment.

Parameters

other

Returns

SUNMatrixWrapper &operator=(SUNMatrixWrapper &&other)

Move assignment.

Parameters

other

Returns

void reallocate(sunindextype nnz)

Reallocate space for sparse matrix according to specified nnz.

Parameters

nnz – new number of nonzero entries

void realloc()

Reallocate space for sparse matrix to used space according to last entry in indexptrs.

SUNMatrix get() const

Get the wrapped SUNMatrix.

Note

Even though the returned matrix_ pointer is const qualified, matrix_->content will not be const. This is a shortcoming in the underlying C library, which we cannot address and it is not intended that any of those values are modified externally. If matrix_->content is manipulated, cpp:meth:SUNMatrixWrapper:refresh needs to be called.

Returns

raw SunMatrix object

inline sunindextype rows() const

Get the number of rows.

Returns

number of rows

inline sunindextype columns() const

Get the number of columns.

Returns

number of columns

sunindextype num_nonzeros() const

Get the number of specified non-zero elements (sparse matrices only)

Note

value will be 0 before indexptrs are set.

Returns

number of nonzero entries

sunindextype num_indexptrs() const

Get the number of indexptrs that can be specified (sparse matrices only)

Returns

number of indexptrs

sunindextype capacity() const

Get the number of allocated data elements.

Returns

number of allocated entries

realtype *data()

Get raw data of a sparse matrix.

Returns

pointer to first data entry

const realtype *data() const

Get const raw data of a sparse matrix.

Returns

pointer to first data entry

inline realtype get_data(sunindextype idx) const

Get data of a sparse matrix.

Parameters

idx – data index

Returns

idx-th data entry

inline realtype get_data(sunindextype irow, sunindextype icol) const

Get data entry for a dense matrix.

Parameters
  • irow – row

  • icol – col

Returns

A(irow,icol)

inline void set_data(sunindextype idx, realtype data)

Set data entry for a sparse matrix.

Parameters
  • idx – data index

  • data – data for idx-th entry

inline void set_data(sunindextype irow, sunindextype icol, realtype data)

Set data entry for a dense matrix.

Parameters
  • irow – row

  • icol – col

  • data – data for idx-th entry

inline sunindextype get_indexval(sunindextype idx) const

Get the index value of a sparse matrix.

Parameters

idx – data index

Returns

row (CSC) or column (CSR) for idx-th data entry

inline void set_indexval(sunindextype idx, sunindextype val)

Set the index value of a sparse matrix.

Parameters
  • idx – data index

  • val – row (CSC) or column (CSR) for idx-th data entry

inline void set_indexvals(const gsl::span<const sunindextype> vals)

Set the index values of a sparse matrix.

Parameters

vals – rows (CSC) or columns (CSR) for data entries

inline sunindextype get_indexptr(sunindextype ptr_idx) const

Get the index pointer of a sparse matrix.

Parameters

ptr_idx – pointer index

Returns

index where the ptr_idx-th column (CSC) or row (CSR) starts

inline void set_indexptr(sunindextype ptr_idx, sunindextype ptr)

Set the index pointer of a sparse matrix.

Parameters
  • ptr_idx – pointer index

  • ptr – data-index where the ptr_idx-th column (CSC) or row (CSR) starts

inline void set_indexptrs(const gsl::span<const sunindextype> ptrs)

Set the index pointers of a sparse matrix.

Parameters

ptrs – starting data-indices where the columns (CSC) or rows (CSR) start

int sparsetype() const

Get the type of sparse matrix.

Returns

matrix type

void scale(realtype a)

multiply with a scalar (in-place)

Parameters

a – scalar value to multiply matrix

void multiply(N_Vector c, const_N_Vector b, realtype alpha = 1.0) const

N_Vector interface for multiply.

Parameters
  • c – output vector, may already contain values

  • b – multiplication vector

  • alpha – scalar coefficient for matrix

inline void multiply(AmiVector &c, AmiVector const &b, realtype alpha = 1.0) const

AmiVector interface for multiply.

Parameters
  • c – output vector, may already contain values

  • b – multiplication vector

  • alpha – scalar coefficient for matrix

void multiply(gsl::span<realtype> c, gsl::span<const realtype> b, const realtype alpha = 1.0) const

Perform matrix vector multiplication c += alpha * A*b.

Parameters
  • c – output vector, may already contain values

  • b – multiplication vector

  • alpha – scalar coefficient

void multiply(N_Vector c, const_N_Vector b, gsl::span<const int> cols, bool transpose) const

Perform reordered matrix vector multiplication c += A[:,cols]*b.

Parameters
  • c – output vector, may already contain values

  • b – multiplication vector

  • cols – int vector for column reordering

  • transpose – bool transpose A before multiplication

void multiply(gsl::span<realtype> c, gsl::span<const realtype> b, gsl::span<const int> cols, bool transpose) const

Perform reordered matrix vector multiplication c += A[:,cols]*b.

Parameters
  • c – output vector, may already contain values

  • b – multiplication vector

  • cols – int vector for column reordering

  • transpose – bool transpose A before multiplication

void sparse_multiply(SUNMatrixWrapper &C, const SUNMatrixWrapper &B) const

Perform matrix matrix multiplication C = A * B for sparse A, B, C.

Note

will overwrite existing data, indexptrs, indexvals for C, but will use preallocated space for these vars

Parameters
  • C – output matrix,

  • B – multiplication matrix

void sparse_add(const SUNMatrixWrapper &A, realtype alpha, const SUNMatrixWrapper &B, realtype beta)

Perform sparse matrix matrix addition C = alpha * A + beta * B.

Note

will overwrite existing data, indexptrs, indexvals for C, but will use preallocated space for these vars

Parameters
  • A – addition matrix

  • alpha – scalar A

  • B – addition matrix

  • beta – scalar B

void sparse_sum(const std::vector<SUNMatrixWrapper> &mats)

Perform matrix-matrix addition A = sum(mats(0)…mats(len(mats)))

Note

will overwrite existing data, indexptrs, indexvals for A, but will use preallocated space for these vars

Parameters

mats – vector of sparse matrices

sunindextype scatter(const sunindextype k, const realtype beta, sunindextype *w, gsl::span<realtype> x, const sunindextype mark, SUNMatrixWrapper *C, sunindextype nnz) const

Compute x = x + beta * A(:,k), where x is a dense vector and A(:,k) is sparse, and update the sparsity pattern for C(:,j) if applicable.

This function currently has two purposes:

  • perform parts of sparse matrix-matrix multiplication C(:,j)=A(:,k)*B(k,j) enabled by passing beta=B(k,j), x=C(:,j), C=C, w=sparsity of C(:,j) from B(k,0…j-1), nnz=nnz(C(:,0…j-1)

  • add the k-th column of the sparse matrix A multiplied by beta to the dense vector x. enabled by passing beta=*, x=x, C=nullptr, w=nullptr, nnz=*

Parameters
  • k – column index

  • beta – scaling factor

  • w – index workspace, (w[i]<mark) indicates non-zeroness of C(i,j) (dimension: m), if this is a nullptr, sparsity pattern of C will not be updated (if applicable).

  • x – dense output vector (dimension: m)

  • mark – marker for w to indicate nonzero pattern

  • C – sparse output matrix, if this is a nullptr, sparsity pattern of C will not be updated

  • nnz – number of nonzeros that were already written to C

Returns

updated number of nonzeros in C

void transpose(SUNMatrixWrapper &C, const realtype alpha, sunindextype blocksize) const

Compute transpose A’ of sparse matrix A and writes it to the matrix C = alpha * A’.

Parameters
  • C – output matrix (sparse or dense)

  • alpha – scalar multiplier

  • blocksize – blocksize for transposition. For full matrix transpose set to ncols/nrows

void to_dense(SUNMatrixWrapper &D) const

Writes a sparse matrix A to a dense matrix D.

Parameters

D – dense output matrix

void to_diag(N_Vector v) const

Writes the diagonal of sparse matrix A to a dense vector v.

Parameters

v – dense outut vector

void zero()

Set to 0.0, for sparse matrices also resets indexptr/indexvals.

inline SUNMatrix_ID matrix_id() const

Get matrix id.

Returns

SUNMatrix_ID

void refresh()

Update internal cache, needs to be called after external manipulation of matrix_->content.