Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/trilinos_tpetra_sparse_matrix.h>
Classes | |
struct | Traits |
Public Types | |
using | size_type = ::types::global_dof_index |
using | iterator |
using | const_iterator |
using | value_type = Number |
Public Member Functions | |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Constructors and initialization. | |
SparseMatrix () | |
SparseMatrix (const SparsityPattern< MemorySpace > &sparsity_pattern) | |
SparseMatrix (const size_type m, const size_type n, const unsigned int n_max_entries_per_row) | |
SparseMatrix (const size_type m, const size_type n, const std::vector< unsigned int > &n_entries_per_row) | |
SparseMatrix (SparseMatrix< Number, MemorySpace > &&other) noexcept | |
SparseMatrix (const SparseMatrix< Number, MemorySpace > &)=delete | |
SparseMatrix< Number, MemorySpace > & | operator= (const SparseMatrix< Number, MemorySpace > &)=delete |
SparseMatrix< Number, MemorySpace > & | operator= (SparseMatrix< Number, MemorySpace > &&other) noexcept |
virtual | ~SparseMatrix () override=default |
template<typename SparsityPatternType > | |
void | reinit (const SparsityPatternType &sparsity_pattern) |
void | reinit (const SparsityPattern< MemorySpace > &sparsity_pattern) |
Constructors and initialization using an IndexSet description | |
SparseMatrix (const IndexSet ¶llel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const unsigned int n_max_entries_per_row=0) | |
SparseMatrix (const IndexSet ¶llel_partitioning, const MPI_Comm communicator, const std::vector< unsigned int > &n_entries_per_row) | |
SparseMatrix (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const size_type n_max_entries_per_row=0) | |
SparseMatrix (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const MPI_Comm communicator, const std::vector< unsigned int > &n_entries_per_row) | |
template<typename SparsityPatternType > | |
std::enable_if_t< !std::is_same_v< SparsityPatternType, ::SparseMatrix< double > > > | reinit (const IndexSet ¶llel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false) |
template<typename SparsityPatternType > | |
std::enable_if_t< !std::is_same_v< SparsityPatternType, ::SparseMatrix< double > > > | reinit (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false) |
void | reinit (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const ::SparseMatrix< Number > &dealii_sparse_matrix, const MPI_Comm communicator=MPI_COMM_WORLD, const double drop_tolerance=1e-13, const bool copy_values=true, const ::SparsityPattern *use_this_sparsity=nullptr) |
Information on the matrix | |
size_type | m () const |
size_type | n () const |
unsigned int | local_size () const |
std::pair< size_type, size_type > | local_range () const |
bool | in_local_range (const size_type index) const |
size_t | n_nonzero_elements () const |
unsigned int | row_length (const size_type row) const |
bool | is_compressed () const |
MPI_Comm | get_mpi_communicator () const |
Modifying entries | |
SparseMatrix & | operator= (const double d) |
SparseMatrix & | operator*= (const Number factor) |
SparseMatrix & | operator/= (const Number factor) |
void | copy_from (const SparseMatrix< Number, MemorySpace > &source) |
void | add (const size_type i, const size_type j, const Number value) |
void | add (const size_type row, const size_type n_cols, const size_type *col_indices, const Number *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false) |
void | add (const Number factor, const SparseMatrix< Number, MemorySpace > &matrix) |
void | set (const size_type i, const size_type j, const Number value) |
void | set (const std::vector< size_type > &indices, const FullMatrix< Number > &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< Number > &full_matrix, const bool elide_zero_values=false) |
void | set (const size_type row, const std::vector< size_type > &col_indices, const std::vector< Number > &values, const bool elide_zero_values=false) |
void | set (const size_type row, const size_type n_cols, const size_type *col_indices, const Number *values, const bool elide_zero_values=false) |
void | clear_row (const size_type row, const Number new_diag_value=0) |
void | clear_rows (const ArrayView< const size_type > &rows, const Number new_diag_value=0) |
void | clear () |
Entry Access | |
Number | operator() (const size_type i, const size_type j) const |
Number | el (const size_type i, const size_type j) const |
Number | diag_element (const size_type i) const |
Multiplications | |
template<typename InputVectorType > | |
void | vmult (InputVectorType &dst, const InputVectorType &src) const |
template<typename InputVectorType > | |
void | Tvmult (InputVectorType &dst, const InputVectorType &src) const |
template<typename InputVectorType > | |
void | vmult_add (InputVectorType &dst, const InputVectorType &src) const |
template<typename InputVectorType > | |
void | Tvmult_add (InputVectorType &dst, const InputVectorType &src) const |
Number | matrix_norm_square (const Vector< Number, MemorySpace > &v) const |
Number | matrix_scalar_product (const Vector< Number, MemorySpace > &u, const Vector< Number, MemorySpace > &v) const |
Number | residual (Vector< Number, MemorySpace > &dst, const Vector< Number, MemorySpace > &x, const Vector< Number, MemorySpace > &b) const |
Matrix norms | |
Number | frobenius_norm () const |
Mixed Stuff | |
void | print (std::ostream &out, const bool print_detailed_trilinos_information=false) const |
void | compress (VectorOperation::values operation) |
void | resume_fill () |
const TpetraTypes::MatrixType< Number, MemorySpace > & | trilinos_matrix () const |
TpetraTypes::MatrixType< Number, MemorySpace > & | trilinos_matrix () |
Teuchos::RCP< const TpetraTypes::MatrixType< Number, MemorySpace > > | trilinos_rcp () const |
Teuchos::RCP< TpetraTypes::MatrixType< Number, MemorySpace > > | trilinos_rcp () |
Partitioners | |
IndexSet | locally_owned_domain_indices () const |
IndexSet | locally_owned_range_indices () const |
Iterators | |
const_iterator | begin () const |
iterator | begin () |
const_iterator | end () const |
iterator | end () |
const_iterator | begin (const size_type r) const |
iterator | begin (const size_type r) |
const_iterator | end (const size_type r) const |
iterator | end (const size_type r) |
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 & | ExcAccessToNonlocalRow (std::size_t arg1) |
static ::ExceptionBase & | ExcMatrixNotCompressed () |
static ::ExceptionBase & | ExcSourceEqualsDestination () |
static ::ExceptionBase & | ExcColMapMismatch () |
static ::ExceptionBase & | ExcDomainMapMismatch () |
static ::ExceptionBase & | ExcInvalidIndex (size_type arg1, size_type arg2) |
static ::ExceptionBase & | ExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4) |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Private Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
Number | element (const size_type i, const size_type j, const bool no_error) const |
void | prepare_add () |
void | prepare_set () |
void | check_no_subscribers () const noexcept |
Private Attributes | |
Teuchos::RCP< TpetraTypes::MapType< MemorySpace > > | column_space_map |
Teuchos::RCP< TpetraTypes::MatrixType< Number, MemorySpace > > | matrix |
bool | compressed |
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 |
Friends | |
class | BlockMatrixBase< SparseMatrix< Number, MemorySpace > > |
This class implements a wrapper to use the Trilinos distributed sparse matrix class Tpetra::CrsMatrix. This is precisely the kind of matrix we deal with all the time - we most likely get it from some assembly process, where also entries not locally owned might need to be written and hence need to be forwarded to the owner process. This class is designed to be used in a distributed memory architecture with an MPI compiler on the bottom, but it works equally well for serial processes. The only requirement for this class to work is that Trilinos has been installed with the same compiler as is used for generating deal.II.
Moreover, this class takes an optional template argument for Kokkos::Nodes, allowing the usage of different Kokkos::Nodes. Kokkos allows the writing of portable applications targeting, for example, CUDA, OpenMP, Serial, or Threads, as backends for the execution and memory spaces. The backend is chosen by choosing the corresponding Kokkos Node.
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. This class is templated and can be used with different scalar types. However, Trilinos need to be installed with complex support for usage with complex scalar types.
Definition at line 110 of file trilinos_tpetra_sparse_matrix.h.
using LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::size_type = ::types::global_dof_index |
Declare the type for container size.
Definition at line 116 of file trilinos_tpetra_sparse_matrix.h.
using LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::iterator |
Declare an alias for the iterator class.
Definition at line 147 of file trilinos_tpetra_sparse_matrix.h.
using LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::const_iterator |
Declare an alias for the const iterator class.
Definition at line 153 of file trilinos_tpetra_sparse_matrix.h.
using LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::value_type = Number |
Declare an alias for the type used to store matrix elements, in analogy to all the other container classes.
Definition at line 160 of file trilinos_tpetra_sparse_matrix.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.
LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::SparseMatrix | ( | ) |
Default constructor. Generates an empty (zero-size) matrix.
LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::SparseMatrix | ( | const SparsityPattern< MemorySpace > & | sparsity_pattern | ) |
Generate a matrix from a TpetraWrappers::SparsityPattern object.
LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::SparseMatrix | ( | const size_type | m, |
const size_type | n, | ||
const unsigned int | n_max_entries_per_row ) |
LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::SparseMatrix | ( | const size_type | m, |
const size_type | n, | ||
const std::vector< unsigned int > & | n_entries_per_row ) |
|
noexcept |
Move constructor. Create a new sparse matrix by stealing the internal data of the other
object.
|
delete |
Copy constructor is deleted.
|
overridevirtualdefault |
Destructor. Made virtual so that one can use pointers to objects of this class.
LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::SparseMatrix | ( | const IndexSet & | parallel_partitioning, |
const MPI_Comm | communicator = MPI_COMM_WORLD, | ||
const unsigned int | n_max_entries_per_row = 0 ) |
Constructor using an IndexSet and an MPI communicator to describe the parallel partitioning. The parameter n_max_entries_per_row
sets the number of nonzero entries in each row that will be allocated. Note that this number does not need to be exact, and it is even allowed that the actual matrix structure has more nonzero entries than specified in the constructor. However it is still advantageous to provide good estimates here since this will considerably increase the performance of the matrix setup. However, there is no effect in the performance of matrix-vector products, since Trilinos reorganizes the matrix memory prior to use (in the compress() step).
LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::SparseMatrix | ( | const IndexSet & | parallel_partitioning, |
const MPI_Comm | communicator, | ||
const std::vector< unsigned int > & | n_entries_per_row ) |
Same as before, but now set the number of non-zero entries in each matrix row separately. Since we know the number of elements in the matrix exactly in this case, we can already allocate the right amount of memory, which makes the creation process including the insertion of nonzero elements by the respective SparseMatrix::reinit call considerably faster.
LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::SparseMatrix | ( | const IndexSet & | row_parallel_partitioning, |
const IndexSet & | col_parallel_partitioning, | ||
const MPI_Comm | communicator = MPI_COMM_WORLD, | ||
const size_type | n_max_entries_per_row = 0 ) |
This constructor is similar to the one above, but it now takes two different IndexSet partitions for row and columns. This interface is meant to be used for generating rectangular matrices, where the first index set describes the parallel partitioning of the degrees of freedom associated with the matrix rows and the second one the partitioning of the matrix columns. The second index set specifies the partitioning of the vectors this matrix is to be multiplied with, not the distribution of the elements that actually appear in the matrix.
The parameter n_max_entries_per_row
defines how much memory will be allocated for each row. This number does not need to be accurate, as the structure is reorganized in the compress() call.
LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::SparseMatrix | ( | const IndexSet & | row_parallel_partitioning, |
const IndexSet & | col_parallel_partitioning, | ||
const MPI_Comm | communicator, | ||
const std::vector< unsigned int > & | n_entries_per_row ) |
Same as before, but now set the number of non-zero entries in each matrix row separately. Since we know the number of elements in the matrix exactly in this case, we can already allocate the right amount of memory, which makes the creation process including the insertion of nonzero elements by the respective SparseMatrix::reinit call considerably faster.
|
delete |
operator= is deleted.
|
noexcept |
Move assignment operator.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::reinit | ( | const SparsityPatternType & | sparsity_pattern | ) |
This function initializes the Trilinos matrix with a deal.II sparsity pattern, i.e. it makes the underlying Trilinos Tpetra::CrsMatrix know the position of nonzero entries according to the sparsity pattern. This function is meant for use in serial programs, where there is no need to specify how the matrix is going to be distributed among different processors. This function works in parallel, too, but it is recommended to manually specify the parallel partitioning of the matrix using a Tpetra::Map. When run in parallel, it is currently necessary that each processor holds the sparsity_pattern structure because each processor sets its rows.
This is a collective operation that needs to be called on all processors in order to avoid a dead lock.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::reinit | ( | const SparsityPattern< MemorySpace > & | sparsity_pattern | ) |
This function reinitializes the Trilinos sparse matrix from a (possibly distributed) Trilinos sparsity pattern. It also works in parallel. In that case, the partitioning of the Trilinos sparsity pattern is used.
This is a collective operation that needs to be called on all processors in order to avoid a dead lock.
std::enable_if_t< !std::is_same_v< SparsityPatternType, ::SparseMatrix< double > > > LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::reinit | ( | const IndexSet & | parallel_partitioning, |
const SparsityPatternType & | sparsity_pattern, | ||
const MPI_Comm | communicator = MPI_COMM_WORLD, | ||
const bool | exchange_data = false ) |
This function is initializes the Trilinos Tpetra matrix according to the specified sparsity_pattern
, and also reassigns the matrix rows to different processes according to the user-supplied index set parallel_partitioning
and parallel communicator. In programs following the style of the tutorial programs, this function (and the respective call for a rectangular matrix) are the natural way to initialize the matrix size, its distribution among the MPI processes (if run in parallel) as well as the location of non-zero elements. Trilinos stores the sparsity pattern internally, so it won't be needed any more after this call, in contrast to the deal.II own object. The optional argument exchange_data
can be used for reinitialization with a sparsity pattern that is not fully constructed. If the flag is not set, each processor just sets the elements in the sparsity pattern that belong to its rows.
This is a collective operation that needs to be called on all processors in order to avoid a dead lock.
std::enable_if_t< !std::is_same_v< SparsityPatternType, ::SparseMatrix< double > > > LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::reinit | ( | const IndexSet & | row_parallel_partitioning, |
const IndexSet & | col_parallel_partitioning, | ||
const SparsityPatternType & | sparsity_pattern, | ||
const MPI_Comm | communicator = MPI_COMM_WORLD, | ||
const bool | exchange_data = false ) |
This function is similar to the other initialization function above, but now also reassigns the matrix rows and columns according to two user-supplied index sets. To be used for rectangular matrices. The optional argument exchange_data
can be used for reinitialization with a sparsity pattern that is not fully constructed. This feature is only implemented for input sparsity patterns of type DynamicSparsityPattern.
This is a collective operation that needs to be called on all processors in order to avoid a dead lock.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::reinit | ( | const IndexSet & | row_parallel_partitioning, |
const IndexSet & | col_parallel_partitioning, | ||
const ::SparseMatrix< Number > & | dealii_sparse_matrix, | ||
const MPI_Comm | communicator = MPI_COMM_WORLD, | ||
const double | drop_tolerance = 1e-13, | ||
const bool | copy_values = true, | ||
const ::SparsityPattern * | use_this_sparsity = nullptr ) |
This function initializes the Trilinos matrix using the deal.II sparse matrix and the entries stored therein. It uses a threshold to copy only elements with modulus larger than the threshold (so zeros in the deal.II matrix can be filtered away). In contrast to the other reinit function with deal.II sparse matrix argument, this function takes a parallel partitioning specified by the user instead of internally generating it.
The optional parameter copy_values
decides whether only the sparsity structure of the input matrix should be used or the matrix entries should be copied, too.
This is a collective operation that needs to be called on all processors in order to avoid a dead lock.
|
inline |
Return the number of rows in this matrix.
Definition at line 1742 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return the number of columns in this matrix.
Definition at line 1751 of file trilinos_tpetra_sparse_matrix.h.
unsigned int LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::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().
std::pair< size_type, size_type > LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::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()
.
|
inline |
Return whether index
is in the local range or not, see also local_range().
Definition at line 1874 of file trilinos_tpetra_sparse_matrix.h.
size_t LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::n_nonzero_elements | ( | ) | const |
Return the total number of nonzero elements of this matrix (summed over all MPI processes).
unsigned int SparseMatrix< Number, MemorySpace >::row_length | ( | const size_type | row | ) | const |
Number of entries in a specific row.
Definition at line 1887 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return the state of the matrix, i.e., whether compress() needs to be called after an operation requiring data exchange. A call to compress() is also needed when the method set() has been called (even when working in serial).
Definition at line 1764 of file trilinos_tpetra_sparse_matrix.h.
MPI_Comm LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::get_mpi_communicator | ( | ) | const |
Return the underlying MPI communicator.
SparseMatrix & LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::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 keeps the sparsity pattern previously used.
SparseMatrix & LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::operator*= | ( | const Number | factor | ) |
Multiply the entire matrix by a fixed factor.
SparseMatrix & LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::operator/= | ( | const Number | factor | ) |
Divide the entire matrix by a fixed factor.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::copy_from | ( | const SparseMatrix< Number, MemorySpace > & | source | ) |
Copy the given (Trilinos) matrix (sparsity pattern and entries).
|
inline |
Add value
to the element (i,j). Just as the respective call in deal.II SparseMatrix<Number, MemorySpace> class. Moreover, if value
is not a finite number an exception is thrown.
Definition at line 1715 of file trilinos_tpetra_sparse_matrix.h.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::add | ( | const size_type | row, |
const size_type | n_cols, | ||
const size_type * | col_indices, | ||
const Number * | 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. Just as the respective call in deal.II SparseMatrix<Number, MemorySpace> class. 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 LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::add | ( | const Number | factor, |
const SparseMatrix< Number, MemorySpace > & | matrix ) |
Add matrix
scaled by factor
to this matrix, i.e. the matrix factor*matrix
is added to this
. If the sparsity pattern of the calling matrix does not contain all the elements in the sparsity pattern of the input matrix, this function will throw an exception.
|
inline |
Set the element (i,j) to value
.
This function is able to insert new elements into the matrix as long as compress() has not been called, so the sparsity pattern will be extended. When compress() is called for the first time (or in case the matrix is initialized from a sparsity pattern), no new elements can be added and an insertion of elements at positions which have not been initialized will throw an exception.
For the case that the matrix is constructed without a sparsity pattern and new matrix entries are added on demand, please note the following behavior imposed by the underlying Tpetra::CrsMatrix data structure: If the same matrix entry is inserted more than once, the matrix entries will be added upon calling compress() (since Tpetra does not track values to the same entry before the final compress() is called), even if VectorOperation::insert is specified as argument to compress(). In the case you cannot make sure that matrix entries are only set once, initialize the matrix with a sparsity pattern to fix the matrix structure before inserting elements.
Definition at line 1704 of file trilinos_tpetra_sparse_matrix.h.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::set | ( | const std::vector< size_type > & | indices, |
const FullMatrix< Number > & | 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.
This function is able to insert new elements into the matrix as long as compress() has not been called, so the sparsity pattern will be extended. After compress() has been called for the first time or the matrix has been initialized from a sparsity pattern, extending the sparsity pattern is no longer possible and an insertion of elements at positions which have not been initialized will throw an exception.
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.
For the case that the matrix is constructed without a sparsity pattern and new matrix entries are added on demand, please note the following behavior imposed by the underlying Tpetra::CrsMatrix data structure: If the same matrix entry is inserted more than once, the matrix entries will be added upon calling compress() (since Epetra does not track values to the same entry before the final compress() is called), even if VectorOperation::insert is specified as argument to compress(). In the case you cannot make sure that matrix entries are only set once, initialize the matrix with a sparsity pattern to fix the matrix structure before inserting elements.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::set | ( | const std::vector< size_type > & | row_indices, |
const std::vector< size_type > & | col_indices, | ||
const FullMatrix< Number > & | 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 LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::set | ( | const size_type | row, |
const std::vector< size_type > & | col_indices, | ||
const std::vector< Number > & | 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.
This function is able to insert new elements into the matrix as long as compress() has not been called, so the sparsity pattern will be extended. After compress() has been called for the first time or the matrix has been initialized from a sparsity pattern, extending the sparsity pattern is no longer possible and an insertion of elements at positions which have not been initialized will throw an exception.
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.
For the case that the matrix is constructed without a sparsity pattern and new matrix entries are added on demand, please note the following behavior imposed by the underlying Tpetra::CrsMatrix data structure: If the same matrix entry is inserted more than once, the matrix entries will be added upon calling compress() (since Epetra does not track values to the same entry before the final compress() is called), even if VectorOperation::insert is specified as argument to compress(). In the case you cannot make sure that matrix entries are only set once, initialize the matrix with a sparsity pattern to fix the matrix structure before inserting elements.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::set | ( | const size_type | row, |
const size_type | n_cols, | ||
const size_type * | col_indices, | ||
const Number * | 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.
This function is able to insert new elements into the matrix as long as compress() has not been called, so the sparsity pattern will be extended. After compress() has been called for the first time or the matrix has been initialized from a sparsity pattern, extending the sparsity pattern is no longer possible and an insertion of elements at positions which have not been initialized will throw an exception.
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.
For the case that the matrix is constructed without a sparsity pattern and new matrix entries are added on demand, please note the following behavior imposed by the underlying Tpetra::CrsMatrix data structure: If the same matrix entry is inserted more than once, the matrix entries will be added upon calling compress() (since Epetra does not track values to the same entry before the final compress() is called), even if VectorOperation::insert is specified as argument to compress(). In the case you cannot make sure that matrix entries are only set once, initialize the matrix with a sparsity pattern to fix the matrix structure before inserting elements.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::clear_row | ( | const size_type | row, |
const Number | 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 the entries to zero.
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.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::clear_rows | ( | const ArrayView< const size_type > & | rows, |
const Number | 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.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::clear | ( | ) |
Release all memory and return to a state just like after having called the default constructor.
This is a collective operation that needs to be called on all processors in order to avoid a dead lock.
Number LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::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. As in the deal.II sparse matrix class, we throw an exception if the respective entry doesn't exist in the sparsity pattern of this class, which is requested from Trilinos. Moreover, an exception will be thrown when the requested element is not saved on the calling process.
Number LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::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. On the other hand, if you want to be sure the entry exists, you should use operator() instead.
The lack of error checking in this function can also yield surprising results if you have a parallel matrix. In that case, just because you get a zero result from this function does not mean that either the entry does not exist in the sparsity pattern or that it does but has a value of zero. Rather, it could also be that it simply isn't stored on the current processor; in that case, it may be stored on a different processor, and possibly so with a nonzero value.
Number LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::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 and it also throws an error if (i,i) is not element of the local matrix.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::vmult | ( | InputVectorType & | dst, |
const InputVectorType & | src ) const |
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::Tvmult | ( | InputVectorType & | dst, |
const InputVectorType & | src ) const |
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::vmult_add | ( | InputVectorType & | dst, |
const InputVectorType & | 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.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::Tvmult_add | ( | InputVectorType & | dst, |
const InputVectorType & | 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.
Number LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::matrix_norm_square | ( | const Vector< Number, MemorySpace > & | 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 SparseMatrix
class used in deal.II (i.e. the original one, not the Trilinos wrapper class) since Trilinos doesn't support this operation and needs a temporary vector.
The vector has to be initialized with the same IndexSet the matrix was initialized with.
In case of a localized Vector, this function will only work when running on one processor, since the matrix object is inherently distributed. Otherwise, an exception will be thrown.
Number LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::matrix_scalar_product | ( | const Vector< Number, MemorySpace > & | u, |
const Vector< Number, MemorySpace > & | 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 SparseMatrix
class used in deal.II (i.e. the original one, not the Trilinos wrapper class) since Trilinos doesn't support this operation and needs a temporary vector.
The vector u
has to be initialized with the same IndexSet that was used for the row indices of the matrix and the vector v
has to be initialized with the same IndexSet that was used for the column indices of the matrix.
In case of a localized Vector, this function will only work when running on one processor, since the matrix object is inherently distributed. Otherwise, an exception will be thrown.
This function is only implemented for square matrices.
|
inline |
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.
The vectors dst
and b
have to be initialized with the same IndexSet that was used for the row indices of the matrix and the vector x
has to be initialized with the same IndexSet that was used for the column indices of the matrix.
In case of a localized Vector, this function will only work when running on one processor, since the matrix object is inherently distributed. Otherwise, an exception will be thrown.
Definition at line 1726 of file trilinos_tpetra_sparse_matrix.h.
Number LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::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.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::print | ( | std::ostream & | out, |
const bool | print_detailed_trilinos_information = false ) const |
Print the matrix to the given stream, using the format (line,col) value, i.e. one nonzero entry of the matrix per line. The optional flag outputs the sparsity pattern in Trilinos style, where the data is sorted according to the processor number when printed to the stream, as well as a summary of the matrix like the global size.
void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::compress | ( | VectorOperation::values | operation | ) |
This command does two things:
In both cases, this function compresses the data structures and allows the resulting matrix to be used in all other operations like matrix- vector products. This is a collective operation, i.e., it needs to be run on all processors when used in parallel.
See Compressing distributed objects for more information.
operation
can be safely omitted, as that parameter is not used at all and is only present to ensure compatibility with other SparseMatrix classes. void LinearAlgebra::TpetraWrappers::SparseMatrix< Number, MemorySpace >::resume_fill | ( | ) |
This function must be called to allow for changes to the structure of the matrix again after compress() was called. Once you are done modifying the matrix structure, you must call compress() again.
|
inline |
Return a const reference to the underlying Trilinos Tpetra::CrsMatrix class.
Definition at line 1919 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return a (modifiable) reference to the underlying Trilinos Tpetra::CrsMatrix class.
Definition at line 1928 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return a const Teuchos::RCP to the underlying Trilinos Tpetra::CrsMatrix class.
Definition at line 1937 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return a (modifiable) Teuchos::RCP to the underlying Trilinos Tpetra::CrsMatrix class.
Definition at line 1946 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return the partitioning of the domain space of this matrix, i.e., the partitioning of the vectors this matrix has to be multiplied with.
Definition at line 1955 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return the partitioning of the range space of this matrix, i.e., the partitioning of the vectors that are result from matrix-vector products.
Definition at line 1964 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return an iterator pointing to the first element of the matrix.
The elements accessed by iterators within each row are ordered in the way in which Trilinos stores them, though the implementation guarantees that all elements of one row are accessed before the elements of the next row. If your algorithm relies on visiting elements within one row, you will need to consult with the Trilinos documentation on the order in which it stores data. It is, however, generally not a good and long-term stable idea to rely on the order in which receive elements if you iterate over them.
When you iterate over the elements of a parallel matrix, you will only be able to access the locally owned rows. (You can access the other rows as well, but they will look empty.) In that case, you probably want to call the begin() function that takes the row as an argument to limit the range of elements to loop over.
Definition at line 1773 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Like the function above, but for non-const matrices.
Definition at line 1823 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return an iterator pointing the element past the last one of this matrix.
Definition at line 1782 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Like the function above, but for non-const matrices.
Definition at line 1832 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return an iterator pointing to the first element 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)
. The returned iterator may not be dereferenceable in that case if neither row r
nor any of the following rows contain any nonzero entries.
The elements accessed by iterators within each row are ordered in the way in which Trilinos stores them, though the implementation guarantees that all elements of one row are accessed before the elements of the next row. If your algorithm relies on visiting elements within one row, you will need to consult with the Trilinos documentation on the order in which it stores data. It is, however, generally not a good and long-term stable idea to rely on the order in which receive elements if you iterate over them.
Definition at line 1791 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Like the function above, but for non-const matrices.
Definition at line 1841 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Return an iterator pointing the element past the last one of row r
, or past the end of the entire sparsity pattern if none of the rows after r
contain any entries at all.
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.
Definition at line 1803 of file trilinos_tpetra_sparse_matrix.h.
|
inline |
Like the function above, but for non-const matrices.
Definition at line 1854 of file trilinos_tpetra_sparse_matrix.h.
|
private |
Helper function for operator() and el().
|
inlineprivate |
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.
This function is called from BlockMatrixBase.
Definition at line 1901 of file trilinos_tpetra_sparse_matrix.h.
|
inlineprivate |
Same as prepare_add() but prepare the matrix for setting elements if the representation of elements in this class requires such an operation.
This function is called from BlockMatrixBase.
Definition at line 1910 of file trilinos_tpetra_sparse_matrix.h.
|
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.
|
friend |
Definition at line 1274 of file trilinos_tpetra_sparse_matrix.h.
|
private |
Pointer to the user-supplied Tpetra Trilinos mapping of the matrix columns that assigns parts of the matrix to the individual processes.
Definition at line 1236 of file trilinos_tpetra_sparse_matrix.h.
|
private |
A sparse matrix object in Trilinos to be used for finite element based problems which allows for assembling into non-local elements. The actual type, a sparse matrix, is set in the constructor.
Definition at line 1243 of file trilinos_tpetra_sparse_matrix.h.
|
private |
A boolean variable to hold information on whether the matrix is fill complete or if the matrix is in compute mode.
Definition at line 1249 of file trilinos_tpetra_sparse_matrix.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.