Reference documentation for deal.II version 9.1.0-pre
Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
PETScWrappers::MatrixFree Class Referenceabstract

#include <deal.II/lac/petsc_matrix_free.h>

Inheritance diagram for PETScWrappers::MatrixFree:
[legend]

Public Member Functions

 MatrixFree ()
 
 MatrixFree (const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns)
 
 MatrixFree (const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const std::vector< unsigned int > &local_rows_per_process, const std::vector< unsigned int > &local_columns_per_process, const unsigned int this_process)
 
 MatrixFree (const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns)
 
 MatrixFree (const unsigned int m, const unsigned int n, const std::vector< unsigned int > &local_rows_per_process, const std::vector< unsigned int > &local_columns_per_process, const unsigned int this_process)
 
void reinit (const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns)
 
void reinit (const MPI_Comm &communicator, const unsigned int m, const unsigned int n, const std::vector< unsigned int > &local_rows_per_process, const std::vector< unsigned int > &local_columns_per_process, const unsigned int this_process)
 
void reinit (const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns)
 
void reinit (const unsigned int m, const unsigned int n, const std::vector< unsigned int > &local_rows_per_process, const std::vector< unsigned int > &local_columns_per_process, const unsigned int this_process)
 
void clear ()
 
const MPI_Comm & get_mpi_communicator () const override
 
virtual void vmult (VectorBase &dst, const VectorBase &src) const =0
 
virtual void Tvmult (VectorBase &dst, const VectorBase &src) const =0
 
virtual void vmult_add (VectorBase &dst, const VectorBase &src) const =0
 
virtual void Tvmult_add (VectorBase &dst, const VectorBase &src) const =0
 
virtual void vmult (Vec &dst, const Vec &src) const
 
- Public Member Functions inherited from PETScWrappers::MatrixBase
 MatrixBase ()
 
 MatrixBase (const MatrixBase &)=delete
 
MatrixBaseoperator= (const MatrixBase &)=delete
 
virtual ~MatrixBase () override
 
MatrixBaseoperator= (const value_type d)
 
void clear ()
 
void set (const size_type i, const size_type j, const PetscScalar value)
 
void set (const std::vector< size_type > &indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=false)
 
void set (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=false)
 
void set (const size_type row, const std::vector< size_type > &col_indices, const std::vector< PetscScalar > &values, const bool elide_zero_values=false)
 
void set (const size_type row, const size_type n_cols, const size_type *col_indices, const PetscScalar *values, const bool elide_zero_values=false)
 
void add (const size_type i, const size_type j, const PetscScalar value)
 
void add (const std::vector< size_type > &indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=true)
 
void add (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=true)
 
void add (const size_type row, const std::vector< size_type > &col_indices, const std::vector< PetscScalar > &values, const bool elide_zero_values=true)
 
void add (const size_type row, const size_type n_cols, const size_type *col_indices, const PetscScalar *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
 
void clear_row (const size_type row, const PetscScalar new_diag_value=0)
 
void clear_rows (const std::vector< size_type > &rows, const PetscScalar new_diag_value=0)
 
void compress (const VectorOperation::values operation)
 
PetscScalar operator() (const size_type i, const size_type j) const
 
PetscScalar el (const size_type i, const size_type j) const
 
PetscScalar diag_element (const size_type i) const
 
size_type m () const
 
size_type n () const
 
size_type local_size () const
 
std::pair< size_type, size_typelocal_range () const
 
bool in_local_range (const size_type index) const
 
size_type n_nonzero_elements () const
 
size_type row_length (const size_type row) const
 
PetscReal l1_norm () const
 
PetscReal linfty_norm () const
 
PetscReal frobenius_norm () const
 
PetscScalar matrix_norm_square (const VectorBase &v) const
 
PetscScalar matrix_scalar_product (const VectorBase &u, const VectorBase &v) const
 
PetscScalar trace () const
 
MatrixBaseoperator*= (const PetscScalar factor)
 
MatrixBaseoperator/= (const PetscScalar factor)
 
MatrixBaseadd (const PetscScalar factor, const MatrixBase &other)
 
MatrixBaseadd (const MatrixBase &other, const PetscScalar factor)
 
void vmult (VectorBase &dst, const VectorBase &src) const
 
void Tvmult (VectorBase &dst, const VectorBase &src) const
 
void vmult_add (VectorBase &dst, const VectorBase &src) const
 
void Tvmult_add (VectorBase &dst, const VectorBase &src) const
 
PetscScalar residual (VectorBase &dst, const VectorBase &x, const VectorBase &b) const
 
const_iterator begin () const
 
const_iterator end () const
 
const_iterator begin (const size_type r) const
 
const_iterator end (const size_type r) const
 
 operator Mat () const
 
Mat & petsc_matrix ()
 
void transpose ()
 
PetscBool is_symmetric (const double tolerance=1.e-12)
 
PetscBool is_hermitian (const double tolerance=1.e-12)
 
void write_ascii (const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
 
void print (std::ostream &out, const bool alternative_output=false) const
 
std::size_t memory_consumption () const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Private Member Functions

void do_reinit (const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns)
 

Static Private Member Functions

static int matrix_free_mult (Mat A, Vec src, Vec dst)
 

Private Attributes

MPI_Comm communicator
 

Additional Inherited Members

- Public Types inherited from PETScWrappers::MatrixBase
using const_iterator = MatrixIterators::const_iterator
 
using size_type = types::global_dof_index
 
using value_type = PetscScalar
 
- Static Public Member Functions inherited from PETScWrappers::MatrixBase
static::ExceptionBase & ExcSourceEqualsDestination ()
 
static::ExceptionBase & ExcWrongMode (int arg1, int arg2)
 
- Static Public Member Functions inherited from Subscriptor
static::ExceptionBase & ExcInUse (int arg1, std::string arg2, std::string arg3)
 
static::ExceptionBase & ExcNoSubscriber (std::string arg1, std::string arg2)
 
- Protected Member Functions inherited from PETScWrappers::MatrixBase
void prepare_action (const VectorOperation::values new_action)
 
void assert_is_compressed ()
 
void prepare_add ()
 
void prepare_set ()
 
void mmult (MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
 
void Tmmult (MatrixBase &C, const MatrixBase &B, const VectorBase &V) const
 
- Protected Attributes inherited from PETScWrappers::MatrixBase
Mat matrix
 
VectorOperation::values last_action
 

Detailed Description

Implementation of a parallel matrix class based on PETSc MatShell matrix-type. This base class implements only the interface to the PETSc matrix object, while all the functionality is contained in the matrix-vector multiplication which must be reimplemented in derived classes.

This interface is an addition to the MatrixFree class to realize user-defined matrix-classes together with PETSc solvers and functionalities. See also the documentation of MatrixFree class and step-37 and step-48.

Similar to other matrix classes in namespaces PETScWrappers and PETScWrappers::MPI, the MatrixFree class provides the usual matrix-vector multiplication vmult(VectorBase &dst, const VectorBase &src) which is pure virtual and must be reimplemented in derived classes. Besides the usual interface, this class has a matrix-vector multiplication vmult(Vec &dst, const Vec &src) taking PETSc Vec objects, which will be called by matrix_free_mult(Mat A, Vec src, Vec dst) registered as matrix-vector multiplication of this PETSc matrix object. The default implementation of the vmult function in the base class wraps the given PETSc vectors with the PETScWrappers::VectorBase class and then calls the usual vmult function with the usual interface.

Author
Wolfgang Bangerth, Martin Steigemann, 2012

Definition at line 60 of file petsc_matrix_free.h.

Constructor & Destructor Documentation

MatrixFree< dim, Number >::MatrixFree ( )

Default constructor. Create an empty matrix object.

Definition at line 28 of file petsc_matrix_free.cc.

MatrixFree< dim, Number >::MatrixFree ( const MPI_Comm &  communicator,
const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns 
)

Create a matrix object of dimensions m times n with communication happening over the provided communicator.

For the meaning of the local_rows and local_columns parameters, see the PETScWrappers::MPI::SparseMatrix class documentation.

As other PETSc matrices, also the matrix-free object needs to have a size and to perform matrix vector multiplications efficiently in parallel also local_rows and local_columns. But in contrast to PETSc::SparseMatrix classes a PETSc matrix-free object does not need any estimation of non_zero entries and has no option is_symmetric.

Definition at line 37 of file petsc_matrix_free.cc.

MatrixFree< dim, Number >::MatrixFree ( const MPI_Comm &  communicator,
const unsigned int  m,
const unsigned int  n,
const std::vector< unsigned int > &  local_rows_per_process,
const std::vector< unsigned int > &  local_columns_per_process,
const unsigned int  this_process 
)

Create a matrix object of dimensions m times n with communication happening over the provided communicator.

As other PETSc matrices, also the matrix-free object needs to have a size and to perform matrix vector multiplications efficiently in parallel also local_rows and local_columns. But in contrast to PETSc::SparseMatrix classes a PETSc matrix-free object does not need any estimation of non_zero entries and has no option is_symmetric.

Definition at line 49 of file petsc_matrix_free.cc.

MatrixFree< dim, Number >::MatrixFree ( const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns 
)

Constructor for the serial case: Same function as MatrixFree(), see above, with communicator = MPI_COMM_WORLD.

Definition at line 71 of file petsc_matrix_free.cc.

MatrixFree< dim, Number >::MatrixFree ( const unsigned int  m,
const unsigned int  n,
const std::vector< unsigned int > &  local_rows_per_process,
const std::vector< unsigned int > &  local_columns_per_process,
const unsigned int  this_process 
)

Constructor for the serial case: Same function as MatrixFree(), see above, with communicator = MPI_COMM_WORLD.

Definition at line 82 of file petsc_matrix_free.cc.

Member Function Documentation

void MatrixFree< dim, Number >::reinit ( const MPI_Comm &  communicator,
const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns 
)

Throw away the present matrix and generate one that has the same properties as if it were created by the constructor of this class with the same argument list as the present function.

Definition at line 104 of file petsc_matrix_free.cc.

void MatrixFree< dim, Number >::reinit ( const MPI_Comm &  communicator,
const unsigned int  m,
const unsigned int  n,
const std::vector< unsigned int > &  local_rows_per_process,
const std::vector< unsigned int > &  local_columns_per_process,
const unsigned int  this_process 
)

Throw away the present matrix and generate one that has the same properties as if it were created by the constructor of this class with the same argument list as the present function.

Definition at line 122 of file petsc_matrix_free.cc.

void MatrixFree< dim, Number >::reinit ( const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns 
)

Call the reinit() function above with communicator = MPI_COMM_WORLD.

Definition at line 147 of file petsc_matrix_free.cc.

void MatrixFree< dim, Number >::reinit ( const unsigned int  m,
const unsigned int  n,
const std::vector< unsigned int > &  local_rows_per_process,
const std::vector< unsigned int > &  local_columns_per_process,
const unsigned int  this_process 
)

Call the reinit() function above with communicator = MPI_COMM_WORLD.

Definition at line 158 of file petsc_matrix_free.cc.

void MatrixFree< dim, Number >::clear ( )

Release all memory and return to a state just like after having called the default constructor.

Definition at line 175 of file petsc_matrix_free.cc.

const MPI_Comm & MatrixFree< dim, Number >::get_mpi_communicator ( ) const
inlineoverridevirtual

Return a reference to the MPI communicator object in use with this matrix.

Implements PETScWrappers::MatrixBase.

Definition at line 294 of file petsc_matrix_free.h.

virtual void PETScWrappers::MatrixFree::vmult ( VectorBase dst,
const VectorBase src 
) const
pure virtual

Matrix-vector multiplication: let dst = M*src with M being this matrix.

Source and destination must not be the same vector.

Note that if the current object represents a parallel distributed matrix (of type PETScWrappers::MPI::SparseMatrix), then both vectors have to be distributed vectors as well. Conversely, if the matrix is not distributed, then neither of the vectors may be.

virtual void PETScWrappers::MatrixFree::Tvmult ( VectorBase dst,
const VectorBase src 
) const
pure virtual

Matrix-vector multiplication: let dst = MT*src with M being this matrix. This function does the same as vmult() but takes the transposed matrix.

Source and destination must not be the same vector.

Note that if the current object represents a parallel distributed matrix then both vectors have to be distributed vectors as well. Conversely, if the matrix is not distributed, then neither of the vectors may be.

virtual void PETScWrappers::MatrixFree::vmult_add ( VectorBase dst,
const VectorBase src 
) const
pure virtual

Adding Matrix-vector multiplication. Add M*src on dst with M being this matrix.

Source and destination must not be the same vector.

Note that if the current object represents a parallel distributed matrix then both vectors have to be distributed vectors as well. Conversely, if the matrix is not distributed, then neither of the vectors may be.

virtual void PETScWrappers::MatrixFree::Tvmult_add ( VectorBase dst,
const VectorBase src 
) const
pure virtual

Adding Matrix-vector multiplication. Add MT*src to dst with M being this matrix. This function does the same as vmult_add() but takes the transposed matrix.

Source and destination must not be the same vector.

Note that if the current object represents a parallel distributed matrix then both vectors have to be distributed vectors as well. Conversely, if the matrix is not distributed, then neither of the vectors may be.

void MatrixFree< dim, Number >::vmult ( Vec &  dst,
const Vec &  src 
) const
virtual

The matrix-vector multiplication called by matrix_free_mult(). This function can be reimplemented in derived classes for efficiency. The default implementation copies the given vectors into PETScWrappers::*Vector and calls vmult(VectorBase &dst, const VectorBase &src) which is purely virtual and must be reimplemented in derived classes.

Definition at line 187 of file petsc_matrix_free.cc.

int MatrixFree< dim, Number >::matrix_free_mult ( Mat  A,
Vec  src,
Vec  dst 
)
staticprivate

Callback-function registered as the matrix-vector multiplication of this matrix-free object called by PETSc routines. This function must be static and takes a PETSc matrix A, and vectors src and dst, where dst = A*src

Source and destination must not be the same vector.

This function calls vmult(Vec &dst, const Vec &src) which should be reimplemented in derived classes.

Definition at line 200 of file petsc_matrix_free.cc.

void MatrixFree< dim, Number >::do_reinit ( const unsigned int  m,
const unsigned int  n,
const unsigned int  local_rows,
const unsigned int  local_columns 
)
private

Do the actual work for the respective reinit() function and the matching constructor, i.e. create a matrix object. Getting rid of the previous matrix is left to the caller.

Definition at line 219 of file petsc_matrix_free.cc.

Member Data Documentation

MPI_Comm PETScWrappers::MatrixFree::communicator
private

Copy of the communicator object to be used for this parallel matrix- free object.

Definition at line 261 of file petsc_matrix_free.h.


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