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

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

Inheritance diagram for PETScWrappers::MatrixBase:
[legend]

Public Types

using const_iterator = MatrixIterators::const_iterator
 
using size_type = types::global_dof_index
 
using value_type = PetscScalar
 

Public Member Functions

 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
 
virtual const MPI_Comm & get_mpi_communicator () const =0
 
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)
 

Static Public Member Functions

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

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

Mat matrix
 
VectorOperation::values last_action
 

Private Attributes

std::vector< PetscInt > column_indices
 
std::vector< PetscScalar > column_values
 

Friends

template<class >
class ::BlockMatrixBase
 

Detailed Description

Base class for all matrix classes that are implemented on top of the PETSc matrix types. Since in PETSc all matrix types (i.e. sequential and parallel, sparse, blocked, etc.) are built by filling the contents of an abstract object that is only referenced through a pointer of a type that is independent of the actual matrix type, we can implement almost all functionality of matrices in this base class. Derived classes will then only have to provide the functionality to create one or the other kind of matrix.

The interface of this class is modeled after the existing SparseMatrix class in deal.II. It has almost the same member functions, and is often exchangeable. However, since PETSc only supports a single scalar type (either double, float, or a complex data type), it is not templated, and only works with whatever your PETSc installation has defined the data type PetscScalar to.

Note that PETSc only guarantees that operations do what you expect if the functions MatAssemblyBegin and MatAssemblyEnd have been called after matrix assembly. Therefore, you need to call SparseMatrix::compress() before you actually use the matrix. This also calls MatCompress that compresses the storage format for sparse matrices by discarding unused elements. PETSc allows to continue with assembling the matrix after calls to these functions, but since there are no more free entries available after that any more, it is better to only call SparseMatrix::compress() once at the end of the assembly stage and before the matrix is actively used.

Author
Wolfgang Bangerth, 2004

Definition at line 283 of file petsc_matrix_base.h.

Member Typedef Documentation

Declare an alias for the iterator class.

Definition at line 289 of file petsc_matrix_base.h.

Declare type for container size.

Definition at line 294 of file petsc_matrix_base.h.

Declare an alias in analogy to all the other container classes.

Definition at line 299 of file petsc_matrix_base.h.

Constructor & Destructor Documentation

PETScWrappers::MatrixBase::MatrixBase ( )

Default constructor.

Definition at line 75 of file petsc_matrix_base.cc.

PETScWrappers::MatrixBase::MatrixBase ( const MatrixBase )
delete

Copy constructor. It is deleted as copying this base class without knowing the concrete kind of matrix stored may both miss important details and be expensive if the matrix is large.

PETScWrappers::MatrixBase::~MatrixBase ( )
overridevirtual

Destructor. Made virtual so that one can use pointers to this class.

Definition at line 82 of file petsc_matrix_base.cc.

Member Function Documentation

MatrixBase& PETScWrappers::MatrixBase::operator= ( const MatrixBase )
delete

Copy operator. It is deleted as copying this base class without knowing the concrete kind of matrix stored may both miss important details and be expensive if the matrix is large.

MatrixBase & PETScWrappers::MatrixBase::operator= ( const value_type  d)

This operator assigns a scalar to a matrix. Since this does usually not make much sense (should we set all matrix entries to this value? Only the nonzero entries of the sparsity pattern?), this operation is only allowed if the actual value to be assigned is zero. This operator only exists to allow for the obvious notation matrix=0, which sets all elements of the matrix to zero, but keeps the sparsity pattern previously used.

Definition at line 109 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::clear ( )

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

Definition at line 90 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::set ( const size_type  i,
const size_type  j,
const PetscScalar  value 
)

Set the element (i,j) to value.

If the present object (from a derived class of this one) happens to be a sparse matrix, then this function adds a new entry to the matrix if it didn't exist before, very much in contrast to the SparseMatrix class which throws an error if the entry does not exist. If value is not a finite number an exception is thrown.

void PETScWrappers::MatrixBase::set ( const std::vector< size_type > &  indices,
const FullMatrix< PetscScalar > &  full_matrix,
const bool  elide_zero_values = false 
)

Set all elements given in a FullMatrix<double> into the sparse matrix locations given by indices. In other words, this function writes the elements in full_matrix into the calling matrix, using the local-to-global indexing specified by indices for both the rows and the columns of the matrix. This function assumes a quadratic sparse matrix and a quadratic full_matrix, the usual situation in FE calculations.

If the present object (from a derived class of this one) happens to be a sparse matrix, then this function adds some new entries to the matrix if they didn't exist before, very much in contrast to the SparseMatrix class which throws an error if the entry does not exist.

The optional parameter elide_zero_values can be used to specify whether zero values should be inserted anyway or they should be filtered away. The default value is false, i.e., even zero values are inserted/replaced.

void PETScWrappers::MatrixBase::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 
)

Same function as before, but now including the possibility to use rectangular full_matrices and different local-to-global indexing on rows and columns, respectively.

void PETScWrappers::MatrixBase::set ( const size_type  row,
const std::vector< size_type > &  col_indices,
const std::vector< PetscScalar > &  values,
const bool  elide_zero_values = false 
)

Set several elements in the specified row of the matrix with column indices as given by col_indices to the respective value.

If the present object (from a derived class of this one) happens to be a sparse matrix, then this function adds some new entries to the matrix if they didn't exist before, very much in contrast to the SparseMatrix class which throws an error if the entry does not exist.

The optional parameter elide_zero_values can be used to specify whether zero values should be inserted anyway or they should be filtered away. The default value is false, i.e., even zero values are inserted/replaced.

void PETScWrappers::MatrixBase::set ( const size_type  row,
const size_type  n_cols,
const size_type col_indices,
const PetscScalar *  values,
const bool  elide_zero_values = false 
)

Set several elements to values given by values in a given row in columns given by col_indices into the sparse matrix.

If the present object (from a derived class of this one) happens to be a sparse matrix, then this function adds some new entries to the matrix if they didn't exist before, very much in contrast to the SparseMatrix class which throws an error if the entry does not exist.

The optional parameter elide_zero_values can be used to specify whether zero values should be inserted anyway or they should be filtered away. The default value is false, i.e., even zero values are inserted/replaced.

void PETScWrappers::MatrixBase::add ( const size_type  i,
const size_type  j,
const PetscScalar  value 
)

Add value to the element (i,j).

If the present object (from a derived class of this one) happens to be a sparse matrix, then this function adds a new entry to the matrix if it didn't exist before, very much in contrast to the SparseMatrix class which throws an error if the entry does not exist. If value is not a finite number an exception is thrown.

void PETScWrappers::MatrixBase::add ( const std::vector< size_type > &  indices,
const FullMatrix< PetscScalar > &  full_matrix,
const bool  elide_zero_values = true 
)

Add all elements given in a FullMatrix<double> into sparse matrix locations given by indices. In other words, this function adds the elements in full_matrix to the respective entries in calling matrix, using the local-to-global indexing specified by indices for both the rows and the columns of the matrix. This function assumes a quadratic sparse matrix and a quadratic full_matrix, the usual situation in FE calculations.

If the present object (from a derived class of this one) happens to be a sparse matrix, then this function adds some new entries to the matrix if they didn't exist before, very much in contrast to the SparseMatrix class which throws an error if the entry does not exist.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

void PETScWrappers::MatrixBase::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 
)

Same function as before, but now including the possibility to use rectangular full_matrices and different local-to-global indexing on rows and columns, respectively.

void PETScWrappers::MatrixBase::add ( const size_type  row,
const std::vector< size_type > &  col_indices,
const std::vector< PetscScalar > &  values,
const bool  elide_zero_values = true 
)

Set several elements in the specified row of the matrix with column indices as given by col_indices to the respective value.

If the present object (from a derived class of this one) happens to be a sparse matrix, then this function adds some new entries to the matrix if they didn't exist before, very much in contrast to the SparseMatrix class which throws an error if the entry does not exist.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

void PETScWrappers::MatrixBase::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 
)

Add an array of values given by values in the given global matrix row at columns specified by col_indices in the sparse matrix.

If the present object (from a derived class of this one) happens to be a sparse matrix, then this function adds some new entries to the matrix if they didn't exist before, very much in contrast to the SparseMatrix class which throws an error if the entry does not exist.

The optional parameter elide_zero_values can be used to specify whether zero values should be added anyway or these should be filtered away and only non-zero data is added. The default value is true, i.e., zero values won't be added into the matrix.

void PETScWrappers::MatrixBase::clear_row ( const size_type  row,
const PetscScalar  new_diag_value = 0 
)

Remove all elements from this row by setting them to zero. The function does not modify the number of allocated nonzero entries, it only sets some entries to zero. It may drop them from the sparsity pattern, though (but retains the allocated memory in case new entries are again added later).

This operation is used in eliminating constraints (e.g. due to hanging nodes) and makes sure that we can write this modification to the matrix without having to read entries (such as the locations of non-zero elements) from it – without this operation, removing constraints on parallel matrices is a rather complicated procedure.

The second parameter can be used to set the diagonal entry of this row to a value different from zero. The default is to set it to zero.

Definition at line 125 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::clear_rows ( const std::vector< size_type > &  rows,
const PetscScalar  new_diag_value = 0 
)

Same as clear_row(), except that it works on a number of rows at once.

The second parameter can be used to set the diagonal entries of all cleared rows to something different from zero. Note that all of these diagonal entries get the same value – if you want different values for the diagonal entries, you have to set them by hand.

Definition at line 134 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::compress ( const VectorOperation::values  operation)

PETSc matrices store their own sparsity patterns. So, in analogy to our own SparsityPattern class, this function compresses the sparsity pattern and allows the resulting matrix to be used in all other operations where before only assembly functions were allowed. This function must therefore be called once you have assembled the matrix.

See Compressing distributed objects for more information.

Definition at line 191 of file petsc_matrix_base.cc.

PetscScalar PETScWrappers::MatrixBase::operator() ( const size_type  i,
const size_type  j 
) const

Return the value of the entry (i,j). This may be an expensive operation and you should always take care where to call this function. In contrast to the respective function in the MatrixBase class, we don't throw an exception if the respective entry doesn't exist in the sparsity pattern of this class, since PETSc does not transmit this information.

This function is therefore exactly equivalent to the el() function.

PetscScalar PETScWrappers::MatrixBase::el ( const size_type  i,
const size_type  j 
) const

Return the value of the matrix entry (i,j). If this entry does not exist in the sparsity pattern, then zero is returned. While this may be convenient in some cases, note that it is simple to write algorithms that are slow compared to an optimal solution, since the sparsity of the matrix is not used.

Definition at line 163 of file petsc_matrix_base.cc.

PetscScalar PETScWrappers::MatrixBase::diag_element ( const size_type  i) const

Return the main diagonal element in the ith row. This function throws an error if the matrix is not quadratic.

Since we do not have direct access to the underlying data structure, this function is no faster than the elementwise access using the el() function. However, we provide this function for compatibility with the SparseMatrix class.

Definition at line 179 of file petsc_matrix_base.cc.

MatrixBase::size_type PETScWrappers::MatrixBase::m ( ) const

Return the number of rows in this matrix.

Definition at line 235 of file petsc_matrix_base.cc.

MatrixBase::size_type PETScWrappers::MatrixBase::n ( ) const

Return the number of columns in this matrix.

Definition at line 248 of file petsc_matrix_base.cc.

MatrixBase::size_type PETScWrappers::MatrixBase::local_size ( ) const

Return the local dimension of the matrix, i.e. the number of rows stored on the present MPI process. For sequential matrices, this number is the same as m(), but for parallel matrices it may be smaller.

To figure out which elements exactly are stored locally, use local_range().

Definition at line 261 of file petsc_matrix_base.cc.

std::pair< MatrixBase::size_type, MatrixBase::size_type > PETScWrappers::MatrixBase::local_range ( ) const

Return a pair of indices indicating which rows of this matrix are stored locally. The first number is the index of the first row stored, the second the index of the one past the last one that is stored locally. If this is a sequential matrix, then the result will be the pair (0,m()), otherwise it will be a pair (i,i+n), where n=local_size().

Definition at line 274 of file petsc_matrix_base.cc.

bool PETScWrappers::MatrixBase::in_local_range ( const size_type  index) const

Return whether index is in the local range or not, see also local_range().

virtual const MPI_Comm& PETScWrappers::MatrixBase::get_mpi_communicator ( ) const
pure virtual

Return a reference to the MPI communicator object in use with this matrix. This function has to be implemented in derived classes.

Implemented in PETScWrappers::MPI::SparseMatrix, PETScWrappers::SparseMatrix, PETScWrappers::MatrixFree, and PETScWrappers::FullMatrix.

MatrixBase::size_type PETScWrappers::MatrixBase::n_nonzero_elements ( ) const

Return the number of nonzero elements of this matrix. Actually, it returns the number of entries in the sparsity pattern; if any of the entries should happen to be zero, it is counted anyway.

Definition at line 288 of file petsc_matrix_base.cc.

MatrixBase::size_type PETScWrappers::MatrixBase::row_length ( const size_type  row) const

Number of entries in a specific row.

Definition at line 300 of file petsc_matrix_base.cc.

PetscReal PETScWrappers::MatrixBase::l1_norm ( ) const

Return the l1-norm of the matrix, that is \(|M|_1=max_{all columns j}\sum_{all rows i} |M_ij|\), (max. sum of columns). This is the natural matrix norm that is compatible to the l1-norm for vectors, i.e. \(|Mv|_1\leq |M|_1 |v|_1\). (cf. Haemmerlin-Hoffmann: Numerische Mathematik)

Definition at line 336 of file petsc_matrix_base.cc.

PetscReal PETScWrappers::MatrixBase::linfty_norm ( ) const

Return the linfty-norm of the matrix, that is \(|M|_infty=max_{all rows i}\sum_{all columns j} |M_ij|\), (max. sum of rows). This is the natural matrix norm that is compatible to the linfty-norm of vectors, i.e. \(|Mv|_infty \leq |M|_infty |v|_infty\). (cf. Haemmerlin-Hoffmann: Numerische Mathematik)

Definition at line 349 of file petsc_matrix_base.cc.

PetscReal PETScWrappers::MatrixBase::frobenius_norm ( ) const

Return the frobenius norm of the matrix, i.e. the square root of the sum of squares of all entries in the matrix.

Definition at line 362 of file petsc_matrix_base.cc.

PetscScalar PETScWrappers::MatrixBase::matrix_norm_square ( const VectorBase v) const

Return the square of the norm of the vector \(v\) with respect to the norm induced by this matrix, i.e. \(\left(v,Mv\right)\). This is useful, e.g. in the finite element context, where the \(L_2\) norm of a function equals the matrix norm with respect to the mass matrix of the vector representing the nodal values of the finite element function.

Obviously, the matrix needs to be quadratic for this operation.

The implementation of this function is not as efficient as the one in the MatrixBase class used in deal.II (i.e. the original one, not the PETSc wrapper class) since PETSc doesn't support this operation and needs a temporary vector.

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

Definition at line 374 of file petsc_matrix_base.cc.

PetscScalar PETScWrappers::MatrixBase::matrix_scalar_product ( const VectorBase u,
const VectorBase v 
) const

Compute the matrix scalar product \(\left(u,Mv\right)\).

The implementation of this function is not as efficient as the one in the MatrixBase class used in deal.II (i.e. the original one, not the PETSc wrapper class) since PETSc doesn't support this operation and needs a temporary 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.

Definition at line 383 of file petsc_matrix_base.cc.

PetscScalar PETScWrappers::MatrixBase::trace ( ) const

Return the trace of the matrix, i.e. the sum of all diagonal entries in the matrix.

Definition at line 393 of file petsc_matrix_base.cc.

MatrixBase & PETScWrappers::MatrixBase::operator*= ( const PetscScalar  factor)

Multiply the entire matrix by a fixed factor.

Definition at line 406 of file petsc_matrix_base.cc.

MatrixBase & PETScWrappers::MatrixBase::operator/= ( const PetscScalar  factor)

Divide the entire matrix by a fixed factor.

Definition at line 417 of file petsc_matrix_base.cc.

MatrixBase & PETScWrappers::MatrixBase::add ( const PetscScalar  factor,
const MatrixBase other 
)

Add the matrix other scaled by the factor factor to the current matrix.

Definition at line 428 of file petsc_matrix_base.cc.

MatrixBase & PETScWrappers::MatrixBase::add ( const MatrixBase other,
const PetscScalar  factor 
)

Add the matrix other scaled by the factor factor to the current matrix.

Deprecated:
Use the function with order of arguments reversed instead.

Definition at line 440 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::vmult ( VectorBase dst,
const VectorBase src 
) const

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.

Definition at line 447 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::Tvmult ( VectorBase dst,
const VectorBase src 
) const

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 (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.

Definition at line 458 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::vmult_add ( VectorBase dst,
const VectorBase src 
) const

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 (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.

Definition at line 469 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::Tvmult_add ( VectorBase dst,
const VectorBase src 
) const

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 (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.

Definition at line 480 of file petsc_matrix_base.cc.

PetscScalar PETScWrappers::MatrixBase::residual ( VectorBase dst,
const VectorBase x,
const VectorBase b 
) const

Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst. The l2 norm of the residual vector is returned.

Source x and destination dst must not be the same vector.

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

Definition at line 578 of file petsc_matrix_base.cc.

const_iterator PETScWrappers::MatrixBase::begin ( ) const

Iterator starting at the first entry.

const_iterator PETScWrappers::MatrixBase::end ( ) const

Final iterator.

const_iterator PETScWrappers::MatrixBase::begin ( const size_type  r) const

Iterator starting at the first entry of row r.

Note that if the given row is empty, i.e. does not contain any nonzero entries, then the iterator returned by this function equals end(r). Note also that the iterator may not be dereferencable in that case.

const_iterator PETScWrappers::MatrixBase::end ( const size_type  r) const

Final iterator of row r. It points to the first element past the end of line r, or past the end of the entire sparsity pattern.

Note that the end iterator is not necessarily dereferencable. This is in particular the case if it is the end iterator for the last row of a matrix.

PETScWrappers::MatrixBase::operator Mat ( ) const

Conversion operator to gain access to the underlying PETSc type. If you do this, you cut this class off some information it may need, so this conversion operator should only be used if you know what you do. In particular, it should only be used for read-only operations into the matrix.

Definition at line 594 of file petsc_matrix_base.cc.

Mat & PETScWrappers::MatrixBase::petsc_matrix ( )

Return a reference to the underlying PETSc type. It can be used to modify the underlying data, so use it only when you know what you are doing.

Definition at line 600 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::transpose ( )

Make an in-place transpose of a matrix.

Definition at line 606 of file petsc_matrix_base.cc.

PetscBool PETScWrappers::MatrixBase::is_symmetric ( const double  tolerance = 1.e-12)

Test whether a matrix is symmetric. Default tolerance is \(1000\times32\)-bit machine precision.

Definition at line 613 of file petsc_matrix_base.cc.

PetscBool PETScWrappers::MatrixBase::is_hermitian ( const double  tolerance = 1.e-12)

Test whether a matrix is Hermitian, i.e. it is the complex conjugate of its transpose. Default tolerance is \(1000\times32\)-bit machine precision.

Definition at line 623 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::write_ascii ( const PetscViewerFormat  format = PETSC_VIEWER_DEFAULT)

Print the PETSc matrix object values using PETSc internal matrix viewer function MatView. The default format prints the non- zero matrix elements. For other valid view formats, consult http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatView.html

Definition at line 635 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::print ( std::ostream &  out,
const bool  alternative_output = false 
) const

Print the elements of a matrix to the given output stream.

Parameters
[in,out]outThe output stream to which to write.
[in]alternative_outputThis argument is ignored. It exists for compatibility with similar functions in other matrix classes.

Definition at line 650 of file petsc_matrix_base.cc.

std::size_t PETScWrappers::MatrixBase::memory_consumption ( ) const

Return the number bytes consumed by this matrix on this CPU.

Definition at line 681 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::prepare_action ( const VectorOperation::values  new_action)
protected

Ensure that the add/set mode that is required for actions following this call is compatible with the current mode. Should be called from all internal functions accessing matrix elements.

void PETScWrappers::MatrixBase::assert_is_compressed ( )
protected

Internal function that checks that there are no pending insert/add operations. Throws an exception otherwise. Useful before calling any PETSc internal functions modifying the matrix.

void PETScWrappers::MatrixBase::prepare_add ( )
protected

For some matrix storage formats, in particular for the PETSc distributed blockmatrices, set and add operations on individual elements can not be freely mixed. Rather, one has to synchronize operations when one wants to switch from setting elements to adding to elements. BlockMatrixBase automatically synchronizes the access by calling this helper function for each block. This function ensures that the matrix is in a state that allows adding elements; if it previously already was in this state, the function does nothing.

void PETScWrappers::MatrixBase::prepare_set ( )
protected

Same as prepare_add() but prepare the matrix for setting elements if the representation of elements in this class requires such an operation.

void PETScWrappers::MatrixBase::mmult ( MatrixBase C,
const MatrixBase B,
const VectorBase V 
) const
protected

Base function to perform the matrix-matrix multiplication \(C = AB\), or, if a vector \(V\) whose size is compatible with B is given, \(C = A \text{diag}(V) B\), where \(\text{diag}(V)\) defines a diagonal matrix with the vector entries.

This function assumes that the calling matrix \(A\) and \(B\) have compatible sizes. The size of \(C\) will be set within this function.

The content as well as the sparsity pattern of the matrix \(C\) will be reset by this function, so make sure that the sparsity pattern is not used somewhere else in your program. This is an expensive operation, so think twice before you use this function.

Definition at line 562 of file petsc_matrix_base.cc.

void PETScWrappers::MatrixBase::Tmmult ( MatrixBase C,
const MatrixBase B,
const VectorBase V 
) const
protected

Base function to perform the matrix-matrix multiplication with the transpose of this, i.e., \(C = A^T B\), or, if an optional vector \(V\) whose size is compatible with \(B\) is given, \(C = A^T \text{diag}(V) B\), where \(\text{diag}(V)\) defines a diagonal matrix with the vector entries.

This function assumes that the calling matrix \(A\) and \(B\) have compatible sizes. The size of \(C\) will be set within this function.

The content as well as the sparsity pattern of the matrix \(C\) will be changed by this function, so make sure that the sparsity pattern is not used somewhere else in your program. This is an expensive operation, so think twice before you use this function.

Definition at line 570 of file petsc_matrix_base.cc.

Friends And Related Function Documentation

template<class >
friend class ::BlockMatrixBase
friend

To allow calling protected prepare_add() and prepare_set().

Definition at line 1075 of file petsc_matrix_base.h.

Member Data Documentation

Mat PETScWrappers::MatrixBase::matrix
protected

A generic matrix object in PETSc. The actual type, a sparse matrix, is set in the constructor.

Definition at line 963 of file petsc_matrix_base.h.

VectorOperation::values PETScWrappers::MatrixBase::last_action
protected

Store whether the last action was a write or add operation.

Definition at line 968 of file petsc_matrix_base.h.

std::vector<PetscInt> PETScWrappers::MatrixBase::column_indices
mutableprivate

An internal array of integer values that is used to store the column indices when adding/inserting local data into the (large) sparse matrix.

This variable does not store any "state" of the matrix object. Rather, it is only used as a temporary buffer by some of the member functions of this class. As with all mutable member variables, the use of this variable is not thread-safe unless guarded by a mutex. However, since PETSc matrix operations are not thread-safe anyway, there is no need to attempt to make things thread-safe, and so there is no mutex associated with this variable.

Definition at line 1058 of file petsc_matrix_base.h.

std::vector<PetscScalar> PETScWrappers::MatrixBase::column_values
mutableprivate

An internal array of double values that is used to store the column indices when adding/inserting local data into the (large) sparse matrix.

The same comment as for the column_indices variable above applies.

Definition at line 1068 of file petsc_matrix_base.h.


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