Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
TrilinosWrappers::SparsityPattern Class Reference

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

Inheritance diagram for TrilinosWrappers::SparsityPattern:
Inheritance graph
[legend]

Public Types

using size_type = ::types::global_dof_index
 
using const_iterator = SparsityPatternIterators::Iterator
 

Public Member Functions

size_type n_rows () const
 
size_type n_cols () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Basic constructors and initialization
 SparsityPattern ()
 
 SparsityPattern (const size_type m, const size_type n, const size_type n_entries_per_row=0)
 
 SparsityPattern (const size_type m, const size_type n, const std::vector< size_type > &n_entries_per_row)
 
 SparsityPattern (SparsityPattern &&other) noexcept
 
 SparsityPattern (const SparsityPattern &input_sparsity_pattern)
 
virtual ~SparsityPattern () override=default
 
void reinit (const size_type m, const size_type n, const size_type n_entries_per_row=0)
 
void reinit (const size_type m, const size_type n, const std::vector< size_type > &n_entries_per_row)
 
void copy_from (const SparsityPattern &input_sparsity_pattern)
 
template<typename SparsityPatternType >
void copy_from (const SparsityPatternType &nontrilinos_sparsity_pattern)
 
SparsityPatternoperator= (const SparsityPattern &input_sparsity_pattern)
 
void clear ()
 
void compress ()
 
Constructors and initialization using an IndexSet description
 SparsityPattern (const IndexSet &parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const size_type n_entries_per_row=0)
 
 SparsityPattern (const IndexSet &parallel_partitioning, const MPI_Comm communicator, const std::vector< size_type > &n_entries_per_row)
 
 SparsityPattern (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const size_type n_entries_per_row=0)
 
 SparsityPattern (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const MPI_Comm communicator, const std::vector< size_type > &n_entries_per_row)
 
 SparsityPattern (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const IndexSet &writable_rows, const MPI_Comm communicator=MPI_COMM_WORLD, const size_type n_entries_per_row=0)
 
void reinit (const IndexSet &parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const size_type n_entries_per_row=0)
 
void reinit (const IndexSet &parallel_partitioning, const MPI_Comm communicator, const std::vector< size_type > &n_entries_per_row)
 
void reinit (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const size_type n_entries_per_row=0)
 
void reinit (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const IndexSet &writeable_rows, const MPI_Comm communicator=MPI_COMM_WORLD, const size_type n_entries_per_row=0)
 
void reinit (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const MPI_Comm communicator, const std::vector< size_type > &n_entries_per_row)
 
template<typename SparsityPatternType >
void reinit (const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const SparsityPatternType &nontrilinos_sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false)
 
template<typename SparsityPatternType >
void reinit (const IndexSet &parallel_partitioning, const SparsityPatternType &nontrilinos_sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false)
 
Information on the sparsity pattern
bool is_compressed () const
 
unsigned int max_entries_per_row () const
 
unsigned int local_size () const
 
std::pair< size_type, size_typelocal_range () const
 
bool in_local_range (const size_type index) const
 
std::uint64_t n_nonzero_elements () const
 
size_type row_length (const size_type row) const
 
size_type bandwidth () const
 
bool empty () const
 
bool exists (const size_type i, const size_type j) const
 
bool row_is_stored_locally (const size_type i) const
 
std::size_t memory_consumption () const
 
Adding entries
void add (const size_type i, const size_type j)
 
template<typename ForwardIterator >
void add_entries (const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false)
 
virtual void add_row_entries (const size_type &row, const ArrayView< const size_type > &columns, const bool indices_are_sorted=false) override
 
virtual void add_entries (const ArrayView< const std::pair< size_type, size_type > > &entries)
 
Access of underlying Trilinos data
const Epetra_FECrsGraph & trilinos_sparsity_pattern () const
 
const Epetra_Map & domain_partitioner () const
 
const Epetra_Map & range_partitioner () const
 
MPI_Comm get_mpi_communicator () const
 
Partitioners
IndexSet locally_owned_domain_indices () const
 
IndexSet locally_owned_range_indices () const
 
Iterators
const_iterator begin () const
 
const_iterator end () const
 
const_iterator begin (const size_type r) const
 
const_iterator end (const size_type r) const
 
Input/Output
void write_ascii ()
 
void print (std::ostream &out, const bool write_extended_trilinos_info=false) const
 
void print_gnuplot (std::ostream &out) const
 
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 ::ExceptionBaseExcTrilinosError (int arg1)
 
static ::ExceptionBaseExcInvalidIndex (size_type arg1, size_type arg2)
 
static ::ExceptionBaseExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4)
 
static ::ExceptionBaseExcAccessToNonPresentElement (size_type arg1, size_type arg2)
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Protected Member Functions

virtual void resize (const size_type rows, const size_type cols)
 

Protected Attributes

size_type rows
 
size_type cols
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void check_no_subscribers () const noexcept
 

Private Attributes

std::unique_ptr< Epetra_Map > column_space_map
 
std::unique_ptr< Epetra_FECrsGraph > graph
 
std::unique_ptr< Epetra_CrsGraph > nonlocal_graph
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Friends

class TrilinosWrappers::SparseMatrix
 
class SparsityPatternIterators::Accessor
 
class SparsityPatternIterators::Iterator
 

Detailed Description

This class implements a wrapper class to use the Trilinos distributed sparsity pattern class Epetra_FECrsGraph. This class is designed to be used for construction of parallel Trilinos matrices. The functionality of this class is modeled after the existing sparsity pattern classes, with the difference that this class can work fully in parallel according to a partitioning of the sparsity pattern rows.

This class has many similarities to the DynamicSparsityPattern, since it can dynamically add elements to the pattern without any memory being previously reserved for it. However, it also has a method SparsityPattern::compress(), that finalizes the pattern and enables its use with Trilinos sparse matrices.

Definition at line 273 of file trilinos_sparsity_pattern.h.

Member Typedef Documentation

◆ size_type

Declare type for container size.

Definition at line 279 of file trilinos_sparsity_pattern.h.

◆ const_iterator

Declare an alias for the iterator class.

Definition at line 284 of file trilinos_sparsity_pattern.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ SparsityPattern() [1/10]

SparsityPattern::SparsityPattern ( )

Default constructor. Generates an empty (zero-size) sparsity pattern.

Definition at line 85 of file trilinos_sparsity_pattern.cc.

◆ SparsityPattern() [2/10]

SparsityPattern::SparsityPattern ( const size_type  m,
const size_type  n,
const size_type  n_entries_per_row = 0 
)

Generate a sparsity pattern that is completely stored locally, having \(m\) rows and \(n\) columns. The resulting matrix will be completely stored locally, too.

It is possible to specify the number of columns entries per row using the optional n_entries_per_row argument. However, this value does not need to be accurate or even given at all, since one does usually not have this kind of information before building the sparsity pattern (the usual case when the function DoFTools::make_sparsity_pattern() is called). The entries are allocated dynamically in a similar manner as for the deal.II DynamicSparsityPattern classes. However, a good estimate will reduce the setup time of the sparsity pattern.

Definition at line 100 of file trilinos_sparsity_pattern.cc.

◆ SparsityPattern() [3/10]

SparsityPattern::SparsityPattern ( const size_type  m,
const size_type  n,
const std::vector< size_type > &  n_entries_per_row 
)

Generate a sparsity pattern that is completely stored locally, having \(m\) rows and \(n\) columns. The resulting matrix will be completely stored locally, too.

The vector n_entries_per_row specifies the number of entries in each row (an information usually not available, though).

Definition at line 109 of file trilinos_sparsity_pattern.cc.

◆ SparsityPattern() [4/10]

SparsityPattern::SparsityPattern ( SparsityPattern &&  other)
noexcept

Move constructor. Create a new sparse matrix by stealing the internal data.

Definition at line 119 of file trilinos_sparsity_pattern.cc.

◆ SparsityPattern() [5/10]

SparsityPattern::SparsityPattern ( const SparsityPattern input_sparsity_pattern)

Copy constructor. Sets the calling sparsity pattern to be the same as the input sparsity pattern.

Definition at line 130 of file trilinos_sparsity_pattern.cc.

◆ ~SparsityPattern()

virtual TrilinosWrappers::SparsityPattern::~SparsityPattern ( )
overridevirtualdefault

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

◆ SparsityPattern() [6/10]

SparsityPattern::SparsityPattern ( const IndexSet parallel_partitioning,
const MPI_Comm  communicator = MPI_COMM_WORLD,
const size_type  n_entries_per_row = 0 
)

Constructor for a square sparsity pattern using an IndexSet and an MPI communicator for the description of the parallel partitioning. Moreover, the number of nonzero entries in the rows of the sparsity pattern can be specified. Note that this number does not need to be exact, and it is even allowed that the actual sparsity structure has more nonzero entries than specified in the constructor. However it is still advantageous to provide good estimates here since a good value will avoid repeated allocation of memory, which considerably increases the performance when creating the sparsity pattern.

Definition at line 146 of file trilinos_sparsity_pattern.cc.

◆ SparsityPattern() [7/10]

SparsityPattern::SparsityPattern ( const IndexSet parallel_partitioning,
const MPI_Comm  communicator,
const std::vector< size_type > &  n_entries_per_row 
)

Same as before, but now use the exact number of nonzeros in each m row. Since we know the number of elements in the sparsity pattern exactly in this case, we can already allocate the right amount of memory, which makes the creation process by the respective SparsityPattern::reinit call considerably faster. However, this is a rather unusual situation, since knowing the number of entries in each row is usually connected to knowing the indices of nonzero entries, which the sparsity pattern is designed to describe.

Definition at line 158 of file trilinos_sparsity_pattern.cc.

◆ SparsityPattern() [8/10]

SparsityPattern::SparsityPattern ( const IndexSet row_parallel_partitioning,
const IndexSet col_parallel_partitioning,
const MPI_Comm  communicator = MPI_COMM_WORLD,
const size_type  n_entries_per_row = 0 
)

This constructor is similar to the one above, but it now takes two different index sets to describe the parallel partitioning of rows and columns. This interface is meant to be used for generating rectangular sparsity pattern. Note that there is no real parallelism along the columns – the processor that owns a certain row always owns all the column elements, no matter how far they might be spread out. The second Epetra_Map is only used to specify the number of columns and for internal arrangements when doing matrix-vector products with vectors based on that column map.

The number of columns entries per row is specified as the maximum number of entries argument.

Definition at line 171 of file trilinos_sparsity_pattern.cc.

◆ SparsityPattern() [9/10]

SparsityPattern::SparsityPattern ( const IndexSet row_parallel_partitioning,
const IndexSet col_parallel_partitioning,
const MPI_Comm  communicator,
const std::vector< size_type > &  n_entries_per_row 
)

This constructor is similar to the one above, but it now takes two different index sets for rows and columns. This interface is meant to be used for generating rectangular matrices, where one map specifies the parallel distribution of rows and the second one specifies the distribution of degrees of freedom associated with matrix columns. This second map is however not used for the distribution of the columns themselves – rather, all column elements of a row are stored on the same processor. The vector n_entries_per_row specifies the number of entries in each row of the newly generated matrix.

Definition at line 184 of file trilinos_sparsity_pattern.cc.

◆ SparsityPattern() [10/10]

SparsityPattern::SparsityPattern ( const IndexSet row_parallel_partitioning,
const IndexSet col_parallel_partitioning,
const IndexSet writable_rows,
const MPI_Comm  communicator = MPI_COMM_WORLD,
const size_type  n_entries_per_row = 0 
)

This constructor constructs general sparsity patterns, possible non-square ones. Constructing a sparsity pattern this way allows the user to explicitly specify the rows into which we are going to add elements. This set is required to be a superset of the first index set row_parallel_partitioning that includes also rows that are owned by another processor (ghost rows). Note that elements can only be added to rows specified by writable_rows.

This method is beneficial when the rows to which a processor is going to write can be determined before actually inserting elements into the matrix. For the typical parallel::distributed::Triangulation class used in deal.II, we know that a processor only will add row elements for what we call the locally relevant dofs (see DoFTools::extract_locally_relevant_dofs). The other constructors methods use general Trilinos facilities that allow to add elements to arbitrary rows (as done by all the other reinit functions). However, this flexibility come at a cost, the most prominent being that adding elements into the same matrix from multiple threads in shared memory is not safe whenever MPI is used. For these settings, the current method is the one to choose: It will store the off-processor data as an additional sparsity pattern (that is then passed to the Trilinos matrix via the reinit method) which can be organized in such a way that thread-safety can be ensured (as long as the user makes sure to never write into the same matrix row simultaneously, of course).

Definition at line 198 of file trilinos_sparsity_pattern.cc.

Member Function Documentation

◆ reinit() [1/9]

void SparsityPattern::reinit ( const size_type  m,
const size_type  n,
const size_type  n_entries_per_row = 0 
)

Initialize a sparsity pattern that is completely stored locally, having \(m\) rows and \(n\) columns. The resulting matrix will be completely stored locally.

The number of columns entries per row is specified as the maximum number of entries argument. This does not need to be an accurate number since the entries are allocated dynamically in a similar manner as for the deal.II DynamicSparsityPattern classes, but a good estimate will reduce the setup time of the sparsity pattern.

Definition at line 214 of file trilinos_sparsity_pattern.cc.

◆ reinit() [2/9]

void SparsityPattern::reinit ( const size_type  m,
const size_type  n,
const std::vector< size_type > &  n_entries_per_row 
)

Initialize a sparsity pattern that is completely stored locally, having \(m\) rows and \(n\) columns. The resulting matrix will be completely stored locally.

The vector n_entries_per_row specifies the number of entries in each row.

Definition at line 227 of file trilinos_sparsity_pattern.cc.

◆ copy_from() [1/2]

void SparsityPattern::copy_from ( const SparsityPattern input_sparsity_pattern)

Copy function. Sets the calling sparsity pattern to be the same as the input sparsity pattern.

Definition at line 646 of file trilinos_sparsity_pattern.cc.

◆ copy_from() [2/2]

template<typename SparsityPatternType >
void SparsityPattern::copy_from ( const SparsityPatternType &  nontrilinos_sparsity_pattern)

Copy function from one of the deal.II sparsity patterns. If used in parallel, this function uses an ad-hoc partitioning of the rows and columns.

Definition at line 662 of file trilinos_sparsity_pattern.cc.

◆ operator=()

SparsityPattern & SparsityPattern::operator= ( const SparsityPattern input_sparsity_pattern)

Copy operator. This operation is only allowed for empty objects, to avoid potentially very costly operations automatically synthesized by the compiler. Use copy_from() instead if you know that you really want to copy a sparsity pattern with non-trivial content.

Definition at line 637 of file trilinos_sparsity_pattern.cc.

◆ clear()

void SparsityPattern::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.

Definition at line 679 of file trilinos_sparsity_pattern.cc.

◆ compress()

void SparsityPattern::compress ( )

In analogy to our own SparsityPattern class, this function compresses the sparsity pattern and allows the resulting pattern to be used for actually generating a (Trilinos-based) matrix. This function also exchanges non-local data that might have accumulated during the addition of new elements. This function must therefore be called once the structure is fixed. This is a collective operation, i.e., it needs to be run on all processors when used in parallel.

Definition at line 701 of file trilinos_sparsity_pattern.cc.

◆ reinit() [3/9]

void SparsityPattern::reinit ( const IndexSet parallel_partitioning,
const MPI_Comm  communicator = MPI_COMM_WORLD,
const size_type  n_entries_per_row = 0 
)

Reinitialization function for generating a square sparsity pattern using an IndexSet and an MPI communicator for the description of the parallel partitioning and the number of nonzero entries in the rows of the sparsity pattern. Note that this number does not need to be exact, and it is even allowed that the actual sparsity 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 when creating the sparsity pattern.

This function does not create any entries by itself, but provides the correct data structures that can be used by the respective add() function.

Definition at line 464 of file trilinos_sparsity_pattern.cc.

◆ reinit() [4/9]

void SparsityPattern::reinit ( const IndexSet parallel_partitioning,
const MPI_Comm  communicator,
const std::vector< size_type > &  n_entries_per_row 
)

Same as before, but now use the exact number of nonzeros in each m row. Since we know the number of elements in the sparsity pattern exactly in this case, we can already allocate the right amount of memory, which makes process of adding entries to the sparsity pattern considerably faster. However, this is a rather unusual situation, since knowing the number of entries in each row is usually connected to knowing the indices of nonzero entries, which the sparsity pattern is designed to describe.

Definition at line 479 of file trilinos_sparsity_pattern.cc.

◆ reinit() [5/9]

void SparsityPattern::reinit ( const IndexSet row_parallel_partitioning,
const IndexSet col_parallel_partitioning,
const MPI_Comm  communicator = MPI_COMM_WORLD,
const size_type  n_entries_per_row = 0 
)

This reinit function is similar to the one above, but it now takes two different index sets for rows and columns. This interface is meant to be used for generating rectangular sparsity pattern, where one index set describes the parallel partitioning of the dofs associated with the sparsity pattern rows and the other one of the sparsity pattern columns. Note that there is no real parallelism along the columns – the processor that owns a certain row always owns all the column elements, no matter how far they might be spread out. The second IndexSet is only used to specify the number of columns and for internal arrangements when doing matrix-vector products with vectors based on an EpetraMap based on that IndexSet.

The number of columns entries per row is specified by the argument n_entries_per_row.

Definition at line 494 of file trilinos_sparsity_pattern.cc.

◆ reinit() [6/9]

void SparsityPattern::reinit ( const IndexSet row_parallel_partitioning,
const IndexSet col_parallel_partitioning,
const IndexSet writeable_rows,
const MPI_Comm  communicator = MPI_COMM_WORLD,
const size_type  n_entries_per_row = 0 
)

This reinit function is used to specify general matrices, possibly non-square ones. In addition to the arguments of the other reinit method above, it allows the user to explicitly specify the rows into which we are going to add elements. This set is a superset of the first index set row_parallel_partitioning that includes also rows that are owned by another processor (ghost rows).

This method is beneficial when the rows to which a processor is going to write can be determined before actually inserting elements into the matrix. For the typical parallel::distributed::Triangulation class used in deal.II, we know that a processor only will add row elements for what we call the locally relevant dofs (see DoFTools::extract_locally_relevant_dofs). Trilinos matrices allow to add elements to arbitrary rows (as done by all the other reinit functions) and this is what all the other reinit methods do, too. However, this flexibility come at a cost, the most prominent being that adding elements into the same matrix from multiple threads in shared memory is not safe whenever MPI is used. For these settings, the current method is the one to choose: It will store the off-processor data as an additional sparsity pattern (that is then passed to the Trilinos matrix via the reinit method) which can be organized in such a way that thread-safety can be ensured (as long as the user makes sure to never write into the same matrix row simultaneously, of course).

Definition at line 538 of file trilinos_sparsity_pattern.cc.

◆ reinit() [7/9]

void SparsityPattern::reinit ( const IndexSet row_parallel_partitioning,
const IndexSet col_parallel_partitioning,
const MPI_Comm  communicator,
const std::vector< size_type > &  n_entries_per_row 
)

Same as before, but now using a vector n_entries_per_row for specifying the number of entries in each row of the sparsity pattern.

Definition at line 516 of file trilinos_sparsity_pattern.cc.

◆ reinit() [8/9]

template<typename SparsityPatternType >
void SparsityPattern::reinit ( const IndexSet row_parallel_partitioning,
const IndexSet col_parallel_partitioning,
const SparsityPatternType &  nontrilinos_sparsity_pattern,
const MPI_Comm  communicator = MPI_COMM_WORLD,
const bool  exchange_data = false 
)

Reinit function. Takes one of the deal.II sparsity patterns and the parallel partitioning of the rows and columns specified by two index sets and a parallel communicator for initializing the current Trilinos sparsity pattern. 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.

Definition at line 585 of file trilinos_sparsity_pattern.cc.

◆ reinit() [9/9]

template<typename SparsityPatternType >
void SparsityPattern::reinit ( const IndexSet parallel_partitioning,
const SparsityPatternType &  nontrilinos_sparsity_pattern,
const MPI_Comm  communicator = MPI_COMM_WORLD,
const bool  exchange_data = false 
)

Reinit function. Takes one of the deal.II sparsity patterns and a parallel partitioning of the rows and columns for initializing the current Trilinos sparsity pattern. 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.

Definition at line 611 of file trilinos_sparsity_pattern.cc.

◆ is_compressed()

bool TrilinosWrappers::SparsityPattern::is_compressed ( ) const

Return the state of the sparsity pattern, i.e., whether compress() needs to be called after an operation requiring data exchange.

◆ max_entries_per_row()

unsigned int SparsityPattern::max_entries_per_row ( ) const

Return the maximum number of entries per row on the current processor.

Definition at line 932 of file trilinos_sparsity_pattern.cc.

◆ local_size()

unsigned int SparsityPattern::local_size ( ) const

Return the local dimension of the sparsity pattern, i.e. the number of rows stored on the present MPI process. In the sequential case, this number is the same as n_rows(), but for parallel matrices it may be smaller.

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

Definition at line 905 of file trilinos_sparsity_pattern.cc.

◆ local_range()

std::pair< SparsityPattern::size_type, SparsityPattern::size_type > SparsityPattern::local_range ( ) const

Return a pair of indices indicating which rows of this sparsity pattern 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,n_rows()), otherwise it will be a pair (i,i+n), where n=local_size().

Definition at line 913 of file trilinos_sparsity_pattern.cc.

◆ in_local_range()

bool TrilinosWrappers::SparsityPattern::in_local_range ( const size_type  index) const

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

◆ n_nonzero_elements()

std::uint64_t SparsityPattern::n_nonzero_elements ( ) const

Return the number of nonzero elements of this sparsity pattern.

Definition at line 924 of file trilinos_sparsity_pattern.cc.

◆ row_length()

SparsityPattern::size_type SparsityPattern::row_length ( const size_type  row) const

Return the number of entries in the given row.

In a parallel context, the row in question may of course not be stored on the current processor, and in that case it is not possible to query the number of entries in it. In that case, the returned value is static_cast<size_type>(-1).

Definition at line 940 of file trilinos_sparsity_pattern.cc.

◆ bandwidth()

SparsityPattern::size_type SparsityPattern::bandwidth ( ) const

Compute the bandwidth of the matrix represented by this structure. The bandwidth is the maximum of \(|i-j|\) for which the index pair \((i,j)\) represents a nonzero entry of the matrix. Consequently, the maximum bandwidth a \(n\times m\) matrix can have is \(\max\{n-1,m-1\}\).

Definition at line 878 of file trilinos_sparsity_pattern.cc.

◆ empty()

bool TrilinosWrappers::SparsityPattern::empty ( ) const

Return whether the object is empty. It is empty if no memory is allocated, which is the same as when both dimensions are zero.

◆ exists()

bool SparsityPattern::exists ( const size_type  i,
const size_type  j 
) const

Return whether the index (i,j) exists in the sparsity pattern (i.e., it may be non-zero) or not.

Definition at line 797 of file trilinos_sparsity_pattern.cc.

◆ row_is_stored_locally()

bool SparsityPattern::row_is_stored_locally ( const size_type  i) const

Return whether a given row is stored in the current object on this process.

Definition at line 788 of file trilinos_sparsity_pattern.cc.

◆ memory_consumption()

std::size_t SparsityPattern::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object. Currently not implemented for this class.

Definition at line 1072 of file trilinos_sparsity_pattern.cc.

◆ add()

void TrilinosWrappers::SparsityPattern::add ( const size_type  i,
const size_type  j 
)

Add the element (i,j) to the sparsity pattern.

◆ add_entries() [1/2]

template<typename ForwardIterator >
void TrilinosWrappers::SparsityPattern::add_entries ( const size_type  row,
ForwardIterator  begin,
ForwardIterator  end,
const bool  indices_are_sorted = false 
)

Add several elements in one row to the sparsity pattern.

◆ add_row_entries()

void SparsityPattern::add_row_entries ( const size_type row,
const ArrayView< const size_type > &  columns,
const bool  indices_are_sorted = false 
)
overridevirtual

Optimized function for adding new entries to a given row.

Implements SparsityPatternBase.

Definition at line 960 of file trilinos_sparsity_pattern.cc.

◆ trilinos_sparsity_pattern()

const Epetra_FECrsGraph & TrilinosWrappers::SparsityPattern::trilinos_sparsity_pattern ( ) const

Return a const reference to the underlying Trilinos Epetra_CrsGraph data that stores the sparsity pattern.

◆ domain_partitioner()

const Epetra_Map & SparsityPattern::domain_partitioner ( ) const

Return a const reference to the underlying Trilinos Epetra_Map that sets the parallel partitioning of the domain space of this sparsity pattern, i.e., the partitioning of the vectors matrices based on this sparsity pattern are multiplied with.

Definition at line 970 of file trilinos_sparsity_pattern.cc.

◆ range_partitioner()

const Epetra_Map & SparsityPattern::range_partitioner ( ) const

Return a const reference to the underlying Trilinos Epetra_Map that sets the partitioning of the range space of this sparsity pattern, i.e., the partitioning of the vectors that are result from matrix-vector products.

Definition at line 981 of file trilinos_sparsity_pattern.cc.

◆ get_mpi_communicator()

MPI_Comm SparsityPattern::get_mpi_communicator ( ) const

Return the underlying MPI communicator.

Definition at line 992 of file trilinos_sparsity_pattern.cc.

◆ locally_owned_domain_indices()

IndexSet TrilinosWrappers::SparsityPattern::locally_owned_domain_indices ( ) const

Return the partitioning of the domain space of this pattern, i.e., the partitioning of the vectors a matrix based on this sparsity pattern has to be multiplied with.

◆ locally_owned_range_indices()

IndexSet TrilinosWrappers::SparsityPattern::locally_owned_range_indices ( ) const

Return the partitioning of the range space of this pattern, i.e., the partitioning of the vectors that are the result from matrix-vector products from a matrix based on this pattern.

◆ begin() [1/2]

const_iterator TrilinosWrappers::SparsityPattern::begin ( ) const

Iterator starting at the first entry.

◆ end() [1/2]

const_iterator TrilinosWrappers::SparsityPattern::end ( ) const

Final iterator.

◆ begin() [2/2]

const_iterator TrilinosWrappers::SparsityPattern::begin ( const size_type  r) const

Iterator starting at the first entry of row r.

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

◆ end() [2/2]

const_iterator TrilinosWrappers::SparsityPattern::end ( const size_type  r) const

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

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

◆ write_ascii()

void SparsityPattern::write_ascii ( )

Abstract Trilinos object that helps view in ASCII other Trilinos objects. Currently this function is not implemented. TODO: Not implemented.

Definition at line 1003 of file trilinos_sparsity_pattern.cc.

◆ print()

void SparsityPattern::print ( std::ostream &  out,
const bool  write_extended_trilinos_info = false 
) const

Print (the locally owned part of) the sparsity pattern to the given stream, using the format (line,col). The optional flag outputs the sparsity pattern in Trilinos style, where even the according processor number is printed to the stream, as well as a summary before actually writing the entries.

Definition at line 1014 of file trilinos_sparsity_pattern.cc.

◆ print_gnuplot()

void SparsityPattern::print_gnuplot ( std::ostream &  out) const

Print the sparsity of the matrix in a format that gnuplot understands and which can be used to plot the sparsity pattern in a graphical way. The format consists of pairs i j of nonzero elements, each representing one entry of this matrix, one per line of the output file. Indices are counted from zero on, as usual. Since sparsity patterns are printed in the same way as matrices are displayed, we print the negative of the column index, which means that the (0,0) element is in the top left rather than in the bottom left corner.

Print the sparsity pattern in gnuplot by setting the data style to dots or points and use the plot command.

Definition at line 1041 of file trilinos_sparsity_pattern.cc.

◆ add_entries() [2/2]

void SparsityPatternBase::add_entries ( const ArrayView< const std::pair< size_type, size_type > > &  entries)
virtual

General function for adding new entries from an unsorted list of pairs.

This function is useful when multiple entries need to be added which do not correspond to a particular row, e.g., when assembling a flux sparsity pattern.

Reimplemented from SparsityPatternBase.

Definition at line 109 of file sparsity_pattern_base.cc.

◆ n_rows()

size_type SparsityPatternBase::n_rows ( ) const
inherited

Return number of rows of this matrix, which equals the dimension of the image space.

◆ n_cols()

size_type SparsityPatternBase::n_cols ( ) const
inherited

Return number of columns of this matrix, which equals the dimension of the range space.

◆ resize()

virtual void SparsityPatternBase::resize ( const size_type  rows,
const size_type  cols 
)
protectedvirtualinherited

Internal function for updating the stored size of the sparsity pattern.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
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.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
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.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
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.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
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.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ TrilinosWrappers::SparseMatrix

friend class TrilinosWrappers::SparseMatrix
friend

Definition at line 993 of file trilinos_sparsity_pattern.h.

◆ SparsityPatternIterators::Accessor

Definition at line 994 of file trilinos_sparsity_pattern.h.

◆ SparsityPatternIterators::Iterator

Definition at line 995 of file trilinos_sparsity_pattern.h.

Member Data Documentation

◆ column_space_map

std::unique_ptr<Epetra_Map> TrilinosWrappers::SparsityPattern::column_space_map
private

Pointer to the user-supplied Epetra Trilinos mapping of the matrix columns that assigns parts of the matrix to the individual processes.

Definition at line 976 of file trilinos_sparsity_pattern.h.

◆ graph

std::unique_ptr<Epetra_FECrsGraph> TrilinosWrappers::SparsityPattern::graph
private

A sparsity pattern object in Trilinos to be used for finite element based problems which allows for adding non-local elements to the pattern.

Definition at line 983 of file trilinos_sparsity_pattern.h.

◆ nonlocal_graph

std::unique_ptr<Epetra_CrsGraph> TrilinosWrappers::SparsityPattern::nonlocal_graph
private

A sparsity pattern object for the non-local part of the sparsity pattern that is going to be sent to the owning processor. Only used when the particular constructor or reinit method with writable_rows argument is set

Definition at line 991 of file trilinos_sparsity_pattern.h.

◆ rows

size_type SparsityPatternBase::rows
protectedinherited

Number of rows that this sparsity pattern shall represent.

Definition at line 121 of file sparsity_pattern_base.h.

◆ cols

size_type SparsityPatternBase::cols
protectedinherited

Number of columns that this sparsity pattern shall represent.

Definition at line 126 of file sparsity_pattern_base.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
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.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
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.

◆ object_info

const std::type_info* Subscriptor::object_info
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.

◆ mutex

std::mutex Subscriptor::mutex
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.


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