DOLFIN
DOLFIN C++ interface
|
#include <PETScMatrix.h>
Public Member Functions | |
PETScMatrix () | |
Create empty matrix (on MPI_COMM_WORLD) | |
PETScMatrix (MPI_Comm comm) | |
Create empty matrix. | |
PETScMatrix (Mat A) | |
PETScMatrix (const PETScMatrix &A) | |
Copy constructor. | |
virtual | ~PETScMatrix () |
Destructor. | |
void | init (const TensorLayout &tensor_layout) |
Initialize zero tensor using tensor 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::int64_t, std::int64_t > | local_range (std::size_t dim) const |
Return local ownership range. | |
std::size_t | nnz () const |
Return number of non-zero entries in matrix (collective) | |
virtual void | zero () |
Set all entries to zero and keep any sparse structure. | |
virtual void | apply (std::string mode) |
MPI_Comm | mpi_comm () const |
Return MPI communicator. | |
virtual std::string | str (bool verbose) const |
Return informal string representation (pretty-print) | |
virtual std::shared_ptr< GenericMatrix > | copy () const |
Return copy of matrix. | |
virtual void | init_vector (GenericVector &z, std::size_t dim) const |
virtual 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. | |
virtual 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 using global indices. | |
virtual void | set_local (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 using local indices. | |
virtual 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 using global indices. | |
virtual void | add_local (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 using local indices. | |
virtual 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. | |
virtual void | getrow (std::size_t row, std::vector< std::size_t > &columns, std::vector< double > &values) const |
Get non-zero values of given row. | |
virtual void | setrow (std::size_t row, const std::vector< std::size_t > &columns, const std::vector< double > &values) |
Set values for given row. | |
virtual void | zero (std::size_t m, const dolfin::la_index *rows) |
Set given rows (global row indices) to zero. | |
virtual void | zero_local (std::size_t m, const dolfin::la_index *rows) |
Set given rows (local row indices) to zero. | |
virtual void | ident (std::size_t m, const dolfin::la_index *rows) |
Set given rows (global row indices) to identity matrix. | |
virtual void | ident_local (std::size_t m, const dolfin::la_index *rows) |
Set given rows (local row indices) to identity matrix. | |
virtual void | mult (const GenericVector &x, GenericVector &y) const |
Compute matrix-vector product y = Ax. | |
virtual void | transpmult (const GenericVector &x, GenericVector &y) const |
virtual void | get_diagonal (GenericVector &x) const |
Get diagonal of a matrix. | |
virtual void | set_diagonal (const GenericVector &x) |
Set diagonal of a matrix. | |
virtual const PETScMatrix & | operator*= (double a) |
Multiply matrix by given number. | |
virtual const PETScMatrix & | operator/= (double a) |
Divide matrix by given number. | |
virtual const GenericMatrix & | operator= (const GenericMatrix &A) |
Assignment operator. | |
virtual bool | is_symmetric (double tol) const |
Test if matrix is symmetric. | |
virtual GenericLinearAlgebraFactory & | factory () const |
Return linear algebra backend factory. | |
void | set_options_prefix (std::string options_prefix) |
std::string | get_options_prefix () const |
void | set_from_options () |
Call PETSc function MatSetFromOptions on the PETSc Mat object. | |
const PETScMatrix & | operator= (const PETScMatrix &A) |
Assignment operator. | |
void | set_nullspace (const VectorSpaceBasis &nullspace) |
void | set_near_nullspace (const VectorSpaceBasis &nullspace) |
void | binary_dump (std::string file_name) const |
Dump matrix to PETSc binary format. | |
Public Member Functions inherited from dolfin::GenericMatrix | |
virtual | ~GenericMatrix () |
Destructor. | |
virtual std::size_t | rank () const |
Return tensor rank (number of dimensions) | |
virtual void | get (double *block, const dolfin::la_index *num_rows, const dolfin::la_index *const *rows) const |
Get block of values. | |
virtual void | set (const double *block, const dolfin::la_index *num_rows, const dolfin::la_index *const *rows) |
Set block of values using global indices. | |
virtual void | set_local (const double *block, const dolfin::la_index *num_rows, const dolfin::la_index *const *rows) |
Set block of values using local indices. | |
virtual void | add (const double *block, const dolfin::la_index *num_rows, const dolfin::la_index *const *rows) |
Add block of values using global indices. | |
virtual void | add_local (const double *block, const dolfin::la_index *num_rows, const dolfin::la_index *const *rows) |
Add block of values using local indices. | |
virtual void | add (const double *block, const std::vector< ArrayView< const dolfin::la_index >> &rows) |
Add block of values using global indices. | |
virtual void | add_local (const double *block, const std::vector< ArrayView< const dolfin::la_index >> &rows) |
Add block of values using local indices. | |
const GenericMatrix & | operator+= (const GenericMatrix &A) |
Add given matrix. | |
const GenericMatrix & | operator-= (const GenericMatrix &A) |
Subtract given matrix. | |
virtual double | operator() (dolfin::la_index i, dolfin::la_index j) const |
Get value of given entry. | |
virtual double | getitem (std::pair< dolfin::la_index, dolfin::la_index > ij) const |
Get value of given entry. | |
virtual void | setitem (std::pair< dolfin::la_index, dolfin::la_index > ij, double value) |
virtual void | ident_zeros (double tol=DOLFIN_EPS) |
Insert one on the diagonal for all zero rows. | |
Public Member Functions inherited from dolfin::GenericTensor | |
virtual | ~GenericTensor () |
Destructor. | |
Public Member Functions inherited from dolfin::LinearAlgebraObject | |
virtual const LinearAlgebraObject * | instance () const |
Return concrete instance / unwrap (const version) | |
virtual LinearAlgebraObject * | instance () |
Return concrete instance / unwrap (non-const version) | |
virtual std::shared_ptr< const LinearAlgebraObject > | shared_instance () const |
Return concrete shared ptr instance / unwrap (const version) | |
virtual std::shared_ptr< LinearAlgebraObject > | shared_instance () |
Return concrete shared ptr instance / unwrap (non-const version) | |
Public Member Functions inherited from dolfin::Variable | |
Variable () | |
Create unnamed variable. | |
Variable (const std::string name, const std::string label) | |
Create variable with given name and label. | |
Variable (const Variable &variable) | |
Copy constructor. | |
virtual | ~Variable () |
Destructor. | |
const Variable & | operator= (const Variable &variable) |
Assignment operator. | |
void | rename (const std::string name, const std::string label) |
Rename variable. | |
std::string | name () const |
Return name. | |
std::string | label () const |
Return label (description) | |
std::size_t | id () const |
Public Member Functions inherited from dolfin::PETScBaseMatrix | |
PETScBaseMatrix () | |
Constructor. | |
PETScBaseMatrix (Mat A) | |
Constructor. | |
PETScBaseMatrix (const PETScBaseMatrix &A) | |
Copy constructor. | |
~PETScBaseMatrix () | |
Destructor. | |
std::size_t | size (std::size_t dim) const |
Return number of rows (dim = 0) or columns (dim = 1) | |
std::pair< std::int64_t, std::int64_t > | size () const |
std::pair< std::int64_t, std::int64_t > | local_range (std::size_t dim) const |
Return local range along dimension dim. | |
void | init_vector (GenericVector &z, std::size_t dim) const |
Mat | mat () const |
Return PETSc Mat pointer. | |
MPI_Comm | mpi_comm () const |
Return the MPI communicator. | |
Public Member Functions inherited from dolfin::PETScObject | |
PETScObject () | |
Constructor. Ensures that PETSc has been initialised. | |
virtual | ~PETScObject () |
Destructor. | |
Additional Inherited Members | |
Static Public Member Functions inherited from dolfin::PETScObject | |
static void | petsc_error (int error_code, std::string filename, std::string petsc_function) |
Print error message for PETSc calls that return an error. | |
Public Attributes inherited from dolfin::Variable | |
Parameters | parameters |
Parameters. | |
Protected Member Functions inherited from dolfin::GenericLinearOperator | |
virtual void | init_layout (const GenericVector &x, const GenericVector &y, GenericLinearOperator *wrapper) |
Protected Attributes inherited from dolfin::PETScBaseMatrix | |
Mat | _matA |
This class provides a simple matrix class based on PETSc. It is a wrapper for a PETSc matrix pointer (Mat) implementing the GenericMatrix interface.
The interface is intentionally simple. For advanced usage, access the PETSc Mat pointer using the function mat() and use the standard PETSc interface.
|
explicit |
Create a wrapper around a PETSc Mat pointer. The Mat object should have been created, e.g. via PETSc MatCreate.
|
virtual |
Finalize assembly of tensor. The following values are recognized for the mode parameter:
add - corresponds to PETSc MatAssemblyBegin+End(MAT_FINAL_ASSEMBLY) insert - corresponds to PETSc MatAssemblyBegin+End(MAT_FINAL_ASSEMBLY) flush - corresponds to PETSc MatAssemblyBegin+End(MAT_FLUSH_ASSEMBLY)
Implements dolfin::GenericMatrix.
std::string PETScMatrix::get_options_prefix | ( | ) | const |
Returns the prefix used by PETSc when searching the options database
|
inlinevirtual |
Initialize vector z to be compatible with the matrix-vector product y = Ax. In the parallel case, both size and layout are important.
z | (GenericVector&) Vector to initialise |
dim | (std::size_t) The dimension (axis): dim = 0 –> z = y, dim = 1 –> z = x |
Implements dolfin::GenericMatrix.
void PETScMatrix::set_near_nullspace | ( | const VectorSpaceBasis & | nullspace | ) |
Attach near nullspace to matrix (used by preconditioners, such as smoothed aggregation algerbraic multigrid)
void PETScMatrix::set_nullspace | ( | const VectorSpaceBasis & | nullspace | ) |
Attach nullspace to matrix (typically used by Krylov solvers when solving singular systems)
void PETScMatrix::set_options_prefix | ( | std::string | options_prefix | ) |
Sets the prefix used by PETSc when searching the options database
|
virtual |
Matrix-vector product, y = A^T x. The y vector must either be zero-sized or have correct size and parallel layout.
Implements dolfin::GenericMatrix.