Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/petsc_matrix_free.h>
Public Types | |
using | const_iterator = MatrixIterators::const_iterator |
using | size_type = types::global_dof_index |
using | value_type = PetscScalar |
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 () |
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 |
void | reinit (Mat A) |
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) |
MatrixBase & | add (const PetscScalar factor, const MatrixBase &other) |
void | clear_row (const size_type row, const PetscScalar new_diag_value=0) |
void | clear_rows (const ArrayView< const size_type > &rows, const PetscScalar new_diag_value=0) |
void | clear_rows_columns (const std::vector< size_type > &row_and_column_indices, 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 | local_domain_size () const |
std::pair< size_type, size_type > | local_domain () const |
MPI_Comm | get_mpi_communicator () const |
std::uint64_t | 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) |
PetscScalar | residual (VectorBase &dst, const VectorBase &x, const VectorBase &b) const |
const_iterator | begin () const |
const_iterator | begin (const size_type r) const |
const_iterator | end () 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 |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcSourceEqualsDestination () |
static ::ExceptionBase & | ExcWrongMode (int arg1, int arg2) |
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 Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
void | do_reinit (const MPI_Comm comm, const unsigned int m, const unsigned int n, const unsigned int local_rows, const unsigned int local_columns) |
void | check_no_subscribers () const noexcept |
Static Private Member Functions | |
static int | matrix_free_mult (Mat A, Vec src, Vec dst) |
Private Attributes | |
std::vector< PetscInt > | column_indices |
std::vector< PetscScalar > | column_values |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
Static Private Attributes | |
static std::mutex | mutex |
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.
Definition at line 58 of file petsc_matrix_free.h.
|
inherited |
Declare an alias for the iterator class.
Definition at line 289 of file petsc_matrix_base.h.
|
inherited |
Declare type for container size.
Definition at line 294 of file petsc_matrix_base.h.
|
inherited |
Declare an alias in analogy to all the other container classes.
Definition at line 299 of file petsc_matrix_base.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 229 of file subscriptor.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 234 of file subscriptor.h.
MatrixFree< dim, Number, VectorizedArrayType >::MatrixFree | ( | ) |
Default constructor. Create an empty matrix object.
Definition at line 27 of file petsc_matrix_free.cc.
MatrixFree< dim, Number, VectorizedArrayType >::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 35 of file petsc_matrix_free.cc.
MatrixFree< dim, Number, VectorizedArrayType >::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 46 of file petsc_matrix_free.cc.
MatrixFree< dim, Number, VectorizedArrayType >::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 68 of file petsc_matrix_free.cc.
MatrixFree< dim, Number, VectorizedArrayType >::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 78 of file petsc_matrix_free.cc.
void MatrixFree< dim, Number, VectorizedArrayType >::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 100 of file petsc_matrix_free.cc.
void MatrixFree< dim, Number, VectorizedArrayType >::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 116 of file petsc_matrix_free.cc.
void MatrixFree< dim, Number, VectorizedArrayType >::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 141 of file petsc_matrix_free.cc.
void MatrixFree< dim, Number, VectorizedArrayType >::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 152 of file petsc_matrix_free.cc.
void MatrixFree< dim, Number, VectorizedArrayType >::clear | ( | ) |
Release all memory and return to a state just like after having called the default constructor.
Definition at line 169 of file petsc_matrix_free.cc.
|
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.
|
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.
|
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.
|
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.
|
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 181 of file petsc_matrix_free.cc.
|
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 194 of file petsc_matrix_free.cc.
|
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 213 of file petsc_matrix_free.cc.
|
inherited |
This method associates the PETSc Mat to the instance of the class. This is particularly useful when performing PETSc to Deal.II operations since it allows to reuse the Deal.II MatrixBase and the PETSc Mat without incurring in memory copies.
Definition at line 92 of file petsc_matrix_base.cc.
|
inherited |
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
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.
|
inherited |
Add the matrix other
scaled by the factor factor
to the current matrix.
Definition at line 516 of file petsc_matrix_base.cc.
|
inherited |
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 147 of file petsc_matrix_base.cc.
|
inherited |
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 155 of file petsc_matrix_base.cc.
|
inherited |
Same as clear_rows(), except that the function also zeros the columns.
Definition at line 184 of file petsc_matrix_base.cc.
|
inherited |
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 244 of file petsc_matrix_base.cc.
|
inherited |
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.
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 216 of file petsc_matrix_base.cc.
|
inherited |
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 232 of file petsc_matrix_base.cc.
|
inherited |
Return the number of rows in this matrix.
Definition at line 286 of file petsc_matrix_base.cc.
|
inherited |
Return the number of columns in this matrix.
Definition at line 299 of file petsc_matrix_base.cc.
|
inherited |
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 312 of file petsc_matrix_base.cc.
|
inherited |
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 325 of file petsc_matrix_base.cc.
Return whether index
is in the local range or not, see also local_range().
|
inherited |
Return the local number of columns stored on the present MPI process.
To figure out which elements exactly are stored locally, use local_domain().
Definition at line 339 of file petsc_matrix_base.cc.
|
inherited |
Return a pair of indices indicating which columns of this matrix are stored locally. The first number is the index of the first column stored, the second the index of the one past the last one that is stored locally.
Definition at line 352 of file petsc_matrix_base.cc.
|
inherited |
Return the underlying MPI communicator.
|
inherited |
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 368 of file petsc_matrix_base.cc.
|
inherited |
Number of entries in a specific row.
Definition at line 382 of file petsc_matrix_base.cc.
|
inherited |
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 418 of file petsc_matrix_base.cc.
|
inherited |
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 431 of file petsc_matrix_base.cc.
|
inherited |
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 444 of file petsc_matrix_base.cc.
|
inherited |
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 456 of file petsc_matrix_base.cc.
|
inherited |
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 467 of file petsc_matrix_base.cc.
|
inherited |
Return the trace of the matrix, i.e. the sum of all diagonal entries in the matrix.
Definition at line 480 of file petsc_matrix_base.cc.
|
inherited |
Multiply the entire matrix by a fixed factor.
Definition at line 493 of file petsc_matrix_base.cc.
|
inherited |
Divide the entire matrix by a fixed factor.
Definition at line 504 of file petsc_matrix_base.cc.
|
inherited |
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 660 of file petsc_matrix_base.cc.
|
inherited |
Iterator starting at the first entry. This can only be called on a processor owning the entire matrix. In all other cases refer to the version of begin() taking a row number as an argument.
|
inherited |
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 dereferenceable in that case.
|
inherited |
Final iterator. This can only be called on a processor owning the entire matrix. In all other cases refer to the version of end() taking a row number as an argument.
|
inherited |
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 dereferenceable. This is in particular the case if it is the end iterator for the last row of a matrix.
|
inherited |
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 676 of file petsc_matrix_base.cc.
|
inherited |
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 682 of file petsc_matrix_base.cc.
|
inherited |
Make an in-place transpose of a matrix.
Definition at line 688 of file petsc_matrix_base.cc.
|
inherited |
Test whether a matrix is symmetric. Default tolerance is \(1000\times32\)-bit machine precision.
Definition at line 700 of file petsc_matrix_base.cc.
|
inherited |
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 710 of file petsc_matrix_base.cc.
|
inherited |
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 722 of file petsc_matrix_base.cc.
|
inherited |
Print the elements of a matrix to the given output stream.
[in,out] | out | The output stream to which to write. |
[in] | alternative_output | This argument is ignored. It exists for compatibility with similar functions in other matrix classes. |
Definition at line 740 of file petsc_matrix_base.cc.
|
inherited |
Return the number bytes consumed by this matrix on this CPU.
Definition at line 787 of file petsc_matrix_base.cc.
|
protectedinherited |
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.
|
protectedinherited |
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.
|
protectedinherited |
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.
|
protectedinherited |
Same as prepare_add() but prepare the matrix for setting elements if the representation of elements in this class requires such an operation.
|
protectedinherited |
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 644 of file petsc_matrix_base.cc.
|
protectedinherited |
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 652 of file petsc_matrix_base.cc.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 135 of file subscriptor.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 155 of file subscriptor.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 300 of file subscriptor.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 317 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 203 of file subscriptor.cc.
|
inlineinherited |
Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.
This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.
Definition at line 309 of file subscriptor.h.
|
privatenoexceptinherited |
Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.
Definition at line 52 of file subscriptor.cc.
|
protectedinherited |
A generic matrix object in PETSc. The actual type, a sparse matrix, is set in the constructor.
Definition at line 1011 of file petsc_matrix_base.h.
|
protectedinherited |
Store whether the last action was a write or add operation.
Definition at line 1016 of file petsc_matrix_base.h.
|
mutableprivateinherited |
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 1106 of file petsc_matrix_base.h.
|
mutableprivateinherited |
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 1116 of file petsc_matrix_base.h.
|
mutableprivateinherited |
Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).
The creator (and owner) of an object is counted in the map below if HE manages to supply identification.
We use the mutable
keyword in order to allow subscription to constant objects also.
This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic
class template.
Definition at line 218 of file subscriptor.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 224 of file subscriptor.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.
Definition at line 240 of file subscriptor.h.
|
mutableprivateinherited |
Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.
Definition at line 248 of file subscriptor.h.
|
staticprivateinherited |
A mutex used to ensure data consistency when accessing the mutable
members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers()
.
Definition at line 271 of file subscriptor.h.