uBLASMatrix.h

Note

The documentation on this page was automatically extracted from the DOLFIN C++ code and may need to be edited or expanded.

class uBLASMatrix

Parent class(es)

This class provides a simple matrix class based on uBLAS. It is a simple wrapper for a uBLAS matrix implementing the GenericMatrix interface.

The interface is intentionally simple. For advanced usage, access the underlying uBLAS matrix and use the standard uBLAS interface which is documented at http://www.boost.org/libs/numeric/ublas/doc/index.htm.

Developer note: specialised member functions must be inlined to avoid link errors.

uBLASMatrix()

Create empty matrix

uBLASMatrix(std::size_t M, std::size_t N)

Create M x N matrix

uBLASMatrix(const uBLASMatrix &A)

Copy constructor

explicit uBLASMatrix(const ublas::matrix_expression<E> &A)

Create matrix from given uBLAS matrix expression

void init(const TensorLayout &tensor_layout)

Initialize zero tensor using tenor layout

bool empty() const

Return true if empty

std::size_t size(std::size_t dim) const

Return size of given dimension

std::pair<std::size_t, std::size_t> local_range(std::size_t dim) const

Return local ownership range

void zero()

Set all entries to zero and keep any sparse structure

void apply(std::string mode)

Finalize assembly of tensor

MPI_Comm mpi_comm() const

Return MPI communicator

std::string str(bool verbose) const

Return informal string representation (pretty-print)

std::shared_ptr<GenericMatrix> copy() const

Return copy of matrix

void resize(std::size_t M, std::size_t N)

Resize matrix to M x N

void init_vector(GenericVector &z, std::size_t dim) const

Intialixe vector z to be compatible with the matrix-vector product y = Ax. In the parallel case, both size and layout are important.

Arguments
dim (std::size_t)
The dimension (axis): dim = 0 –> z = y, dim = 1 –> z = x
void get(double *block, std::size_t m, const dolfin::la_index *rows, std::size_t n, const dolfin::la_index *cols) const

Get block of values

void set(const double *block, std::size_t m, const dolfin::la_index *rows, std::size_t n, const dolfin::la_index *cols)

Set block of values

void add(const double *block, std::size_t m, const dolfin::la_index *rows, std::size_t n, const dolfin::la_index *cols)

Add block of values

void axpy(double a, const GenericMatrix &A, bool same_nonzero_pattern)

Add multiple of given matrix (AXPY operation)

double norm(std::string norm_type) const

Return norm of matrix

void getrow(std::size_t row, std::vector<std::size_t> &columns, std::vector<double> &values) const

Get non-zero values of given row

void setrow(std::size_t row_idx, const std::vector<std::size_t> &columns, const std::vector<double> &values)

Set values for given row

void zero(std::size_t m, const dolfin::la_index *rows)

Set given rows to zero

void ident(std::size_t m, const dolfin::la_index *rows)

Set given rows to identity matrix

void mult(const GenericVector &x, GenericVector &y) const

Matrix-vector product, y = Ax

void transpmult(const GenericVector &x, GenericVector &y) const

Matrix-vector product, y = A^T x

void set_diagonal(const GenericVector &x)

Set diagonal of a matrix

const uBLASMatrix<Mat> &operator*=(double a)

Multiply matrix by given number

const uBLASMatrix<Mat> &operator/=(double a)

Divide matrix by given number

const GenericMatrix &operator=(const GenericMatrix &A)

Assignment operator

boost::tuples::tuple<const std::size_t *, const std::size_t *, const double *, int> data() const

Return pointers to underlying compresssed storage data See GenericMatrix for documentation.

GenericLinearAlgebraFactory &factory() const

Return linear algebra backend factory

const Mat &mat() const

Return reference to uBLAS matrix (const version)

Mat &mat()

Return reference to uBLAS matrix (non-const version)

void solve(uBLASVector &x, const uBLASVector &b) const

Solve Ax = b out-of-place using uBLAS (A is not destroyed)

void solve_in_place(uBLASVector &x, const uBLASVector &b)

Solve Ax = b in-place using uBLAS(A is destroyed)

void invert()

Compute inverse of matrix

void lump(uBLASVector &m) const

Lump matrix into vector m

void compress()

Compress matrix (eliminate all non-zeros from a sparse matrix)

double operator()(dolfin::la_index i, dolfin::la_index j) const

Access value of given entry

const uBLASMatrix<Mat> &operator=(const uBLASMatrix<Mat> &A)

Assignment operator