Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/petsc_sparse_matrix.h>
Classes | |
struct | Traits |
Public Member Functions | |
SparseMatrix () | |
SparseMatrix (const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false) | |
SparseMatrix (const size_type m, const size_type n, const std::vector< size_type > &row_lengths, const bool is_symmetric=false) | |
template<typename SparsityPatternType > | |
SparseMatrix (const SparsityPatternType &sparsity_pattern, const bool preset_nonzero_locations=true) | |
SparseMatrix & | operator= (const double d) |
void | reinit (const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false) |
void | reinit (const size_type m, const size_type n, const std::vector< size_type > &row_lengths, const bool is_symmetric=false) |
template<typename SparsityPatternType > | |
void | reinit (const SparsityPatternType &sparsity_pattern, const bool preset_nonzero_locations=true) |
virtual const MPI_Comm & | get_mpi_communicator () const |
size_t | m () const |
size_t | n () const |
void | mmult (SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const |
void | Tmmult (SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const |
Public Member Functions inherited from PETScWrappers::MatrixBase | |
MatrixBase () | |
MatrixBase (const MatrixBase &)=delete | |
MatrixBase & | operator= (const MatrixBase &)=delete |
virtual | ~MatrixBase () |
MatrixBase & | operator= (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_type > | local_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 |
MatrixBase & | operator*= (const PetscScalar factor) |
MatrixBase & | operator/= (const PetscScalar factor) |
MatrixBase & | add (const PetscScalar factor, const MatrixBase &other) |
MatrixBase & | add (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 () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Member Functions | |
SparseMatrix (const SparseMatrix &)=delete | |
SparseMatrix & | operator= (const SparseMatrix &)=delete |
void | do_reinit (const size_type m, const size_type n, const size_type n_nonzero_per_row, const bool is_symmetric=false) |
void | do_reinit (const size_type m, const size_type n, const std::vector< size_type > &row_lengths, const bool is_symmetric=false) |
template<typename SparsityPatternType > | |
void | do_reinit (const SparsityPatternType &sparsity_pattern, const bool preset_nonzero_locations) |
Friends | |
class | BlockMatrixBase< SparseMatrix > |
Additional Inherited Members | |
Public Types inherited from PETScWrappers::MatrixBase | |
typedef MatrixIterators::const_iterator | const_iterator |
typedef types::global_dof_index | size_type |
typedef PetscScalar | value_type |
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 |
Implementation of a sequential sparse matrix class based on PETSc. All the functionality is actually in the base class, except for the calls to generate a sequential sparse matrix. This is possible since PETSc only works on an abstract matrix type and internally distributes to functions that do the actual work depending on the actual matrix type (much like using virtual functions). Only the functions creating a matrix of specific type differ, and are implemented in this particular class.
Definition at line 49 of file petsc_sparse_matrix.h.
SparseMatrix< number >::SparseMatrix | ( | ) |
Default constructor. Create an empty matrix.
Definition at line 31 of file petsc_sparse_matrix.cc.
SparseMatrix< number >::SparseMatrix | ( | const size_type | m, |
const size_type | n, | ||
const size_type | n_nonzero_per_row, | ||
const bool | is_symmetric = false |
||
) |
Create a sparse matrix of dimensions m
times n
, with an initial guess of n_nonzero_per_row
nonzero elements per row. PETSc is able to cope with the situation that more than this number of elements is later allocated for a row, but this involves copying data, and is thus expensive.
The is_symmetric
flag determines whether we should tell PETSc that the matrix is going to be symmetric (as indicated by the call MatSetOption(mat, MAT_SYMMETRIC)
. Note that the PETSc documentation states that one cannot form an ILU decomposition of a matrix for which this flag has been set to true
, only an ICC. The default value of this flag is false
.
Definition at line 41 of file petsc_sparse_matrix.cc.
SparseMatrix< number >::SparseMatrix | ( | const size_type | m, |
const size_type | n, | ||
const std::vector< size_type > & | row_lengths, | ||
const bool | is_symmetric = false |
||
) |
Initialize a rectangular matrix with m
rows and n
columns. The maximal number of nonzero entries for each row separately is given by the row_lengths
array.
Just as for the other constructors: PETSc is able to cope with the situation that more than this number of elements is later allocated for a row, but this involves copying data, and is thus expensive.
The is_symmetric
flag determines whether we should tell PETSc that the matrix is going to be symmetric (as indicated by the call MatSetOption(mat, MAT_SYMMETRIC)
. Note that the PETSc documentation states that one cannot form an ILU decomposition of a matrix for which this flag has been set to true
, only an ICC. The default value of this flag is false
.
Definition at line 51 of file petsc_sparse_matrix.cc.
|
explicit |
Initialize a sparse matrix using the given sparsity pattern.
Note that PETSc can be very slow if you do not provide it with a good estimate of the lengths of rows. Using the present function is a very efficient way to do this, as it uses the exact number of nonzero entries for each row of the matrix by using the given sparsity pattern argument. If the preset_nonzero_locations
flag is true
, this function in addition not only sets the correct row sizes up front, but also pre-allocated the correct nonzero entries in the matrix.
PETsc allows to later add additional nonzero entries to a matrix, by simply writing to these elements. However, this will then lead to additional memory allocations which are very inefficient and will greatly slow down your program. It is therefore significantly more efficient to get memory allocation right from the start.
Definition at line 63 of file petsc_sparse_matrix.cc.
|
privatedelete |
Purposefully not implemented
SparseMatrix & SparseMatrix< number >::operator= | ( | const double | 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 keep the sparsity pattern previously used.
Definition at line 72 of file petsc_sparse_matrix.cc.
void SparseMatrix< number >::reinit | ( | const size_type | m, |
const size_type | n, | ||
const size_type | n_nonzero_per_row, | ||
const bool | is_symmetric = false |
||
) |
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 81 of file petsc_sparse_matrix.cc.
void SparseMatrix< number >::reinit | ( | const size_type | m, |
const size_type | n, | ||
const std::vector< size_type > & | row_lengths, | ||
const bool | is_symmetric = false |
||
) |
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 97 of file petsc_sparse_matrix.cc.
void SparseMatrix< SparsityPatternType >::reinit | ( | const SparsityPatternType & | sparsity_pattern, |
const bool | preset_nonzero_locations = true |
||
) |
Initialize a sparse matrix using the given sparsity pattern.
Note that PETSc can be very slow if you do not provide it with a good estimate of the lengths of rows. Using the present function is a very efficient way to do this, as it uses the exact number of nonzero entries for each row of the matrix by using the given sparsity pattern argument. If the preset_nonzero_locations
flag is true
, this function in addition not only sets the correct row sizes up front, but also pre-allocated the correct nonzero entries in the matrix.
PETsc allows to later add additional nonzero entries to a matrix, by simply writing to these elements. However, this will then lead to additional memory allocations which are very inefficient and will greatly slow down your program. It is therefore significantly more efficient to get memory allocation right from the start.
Despite the fact that it would seem to be an obvious win, setting the preset_nonzero_locations
flag to true
doesn't seem to accelerate program. Rather on the contrary, it seems to be able to slow down entire programs somewhat. This is surprising, since we can use efficient function calls into PETSc that allow to create multiple entries at once; nevertheless, given the fact that it is inefficient, the respective flag has a default value equal to false
.
Definition at line 115 of file petsc_sparse_matrix.cc.
|
virtual |
Return a reference to the MPI communicator object in use with this matrix. Since this is a sequential matrix, it returns the MPI_COMM_SELF communicator.
Implements PETScWrappers::MatrixBase.
Definition at line 129 of file petsc_sparse_matrix.cc.
size_t SparseMatrix< number >::m | ( | ) | const |
Return the number of rows of this matrix.
Definition at line 246 of file petsc_sparse_matrix.cc.
size_t SparseMatrix< number >::n | ( | ) | const |
Return the number of columns of this matrix.
Definition at line 256 of file petsc_sparse_matrix.cc.
void SparseMatrix< number >::mmult | ( | SparseMatrix & | C, |
const SparseMatrix & | B, | ||
const MPI::Vector & | V = MPI::Vector() |
||
) | const |
Perform the matrix-matrix multiplication \(C = AB\), or, \(C = A \text{diag}(V) B\) given a compatible vector \(V\).
This function calls MatrixBase::mmult() to do the actual work.
Definition at line 266 of file petsc_sparse_matrix.cc.
void SparseMatrix< number >::Tmmult | ( | SparseMatrix & | C, |
const SparseMatrix & | B, | ||
const MPI::Vector & | V = MPI::Vector() |
||
) | const |
Perform the matrix-matrix multiplication with the transpose of this
, i.e., \(C = A^T B\), or, \(C = A^T \text{diag}(V) B\) given a compatible vector \(V\).
This function calls MatrixBase::Tmmult() to do the actual work.
Definition at line 277 of file petsc_sparse_matrix.cc.
|
privatedelete |
Purposefully not implemented
|
private |
Do the actual work for the respective reinit() function and the matching constructor, i.e. create a matrix. Getting rid of the previous matrix is left to the caller.
Definition at line 140 of file petsc_sparse_matrix.cc.
|
private |
Same as previous function.
Definition at line 163 of file petsc_sparse_matrix.cc.
|
private |
Same as previous function.
Definition at line 197 of file petsc_sparse_matrix.cc.
|
friend |
To allow calling protected prepare_add() and prepare_set().
Definition at line 271 of file petsc_sparse_matrix.h.