16#ifndef dealii_block_sparsity_pattern_h
17#define dealii_block_sparsity_pattern_h
36template <
typename number>
77template <
typename SparsityPatternType>
160 SparsityPatternType &
168 const SparsityPatternType &
243 template <
typename ForwardIterator>
246 ForwardIterator
begin,
248 const bool indices_are_sorted =
false);
298 print(std::ostream &out)
const;
331 <<
"The blocks [" << arg1 <<
',' << arg2 <<
"] and [" << arg3
332 <<
',' << arg4 <<
"] have differing row numbers.");
341 <<
"The blocks [" << arg1 <<
',' << arg2 <<
"] and [" << arg3
342 <<
',' << arg4 <<
"] have differing column numbers.");
391 template <
typename number>
444 const std::vector<std::vector<unsigned int>> &row_lengths);
550 const std::vector<size_type> &col_block_sizes);
579 reinit(
const std::vector<size_type> &row_block_sizes,
580 const std::vector<size_type> &col_block_sizes);
587 reinit(
const std::vector<IndexSet> &partitioning);
613#ifdef DEAL_II_WITH_TRILINOS
664 const std::vector<size_type> &col_block_sizes);
674 const MPI_Comm &communicator = MPI_COMM_WORLD);
688 const std::vector<IndexSet> &row_parallel_partitioning,
689 const std::vector<IndexSet> &column_parallel_partitioning,
690 const std::vector<IndexSet> &writeable_rows,
691 const MPI_Comm & communicator = MPI_COMM_WORLD);
703 reinit(
const std::vector<size_type> &row_block_sizes,
704 const std::vector<size_type> &col_block_sizes);
711 reinit(
const std::vector<IndexSet> ¶llel_partitioning,
712 const MPI_Comm & communicator = MPI_COMM_WORLD);
720 reinit(
const std::vector<IndexSet> &row_parallel_partitioning,
721 const std::vector<IndexSet> &column_parallel_partitioning,
722 const MPI_Comm & communicator = MPI_COMM_WORLD);
733 reinit(
const std::vector<IndexSet> &row_parallel_partitioning,
734 const std::vector<IndexSet> &column_parallel_partitioning,
735 const std::vector<IndexSet> &writeable_rows,
736 const MPI_Comm & communicator = MPI_COMM_WORLD);
754template <
typename SparsityPatternType>
755inline SparsityPatternType &
761 return *sub_objects[row][column];
766template <
typename SparsityPatternType>
767inline const SparsityPatternType &
774 return *sub_objects[row][column];
779template <
typename SparsityPatternType>
788template <
typename SparsityPatternType>
792 return column_indices;
797template <
typename SparsityPatternType>
805 const std::pair<size_type, size_type> row_index =
806 row_indices.global_to_local(i),
808 column_indices.global_to_local(j);
809 sub_objects[row_index.first][col_index.first]->add(row_index.second,
815template <
typename SparsityPatternType>
816template <
typename ForwardIterator>
820 ForwardIterator
begin,
822 const bool indices_are_sorted)
825 if (block_column_indices.size() < this->n_block_cols())
827 block_column_indices.resize(this->n_block_cols());
828 counter_within_block.resize(this->n_block_cols());
843 if (block_column_indices[0].size() < n_cols)
844 for (
size_type i = 0; i < this->n_block_cols(); ++i)
845 block_column_indices[i].resize(n_cols);
849 for (
size_type i = 0; i < this->n_block_cols(); ++i)
850 counter_within_block[i] = 0;
861 for (ForwardIterator it =
begin; it !=
end; ++it)
865 const std::pair<size_type, size_type> col_index =
866 this->column_indices.global_to_local(col);
868 const size_type local_index = counter_within_block[col_index.first]++;
870 block_column_indices[col_index.first][local_index] = col_index.second;
877 for (
size_type i = 0; i < this->n_block_cols(); ++i)
878 length += counter_within_block[i];
887 const std::pair<size_type, size_type> row_index =
888 this->row_indices.global_to_local(row);
889 for (
size_type block_col = 0; block_col < n_block_cols(); ++block_col)
891 if (counter_within_block[block_col] == 0)
893 sub_objects[row_index.first][block_col]->add_entries(
895 block_column_indices[block_col].begin(),
896 block_column_indices[block_col].begin() +
897 counter_within_block[block_col],
904template <
typename SparsityPatternType>
912 const std::pair<size_type, size_type> row_index =
913 row_indices.global_to_local(i),
915 column_indices.global_to_local(j);
916 return sub_objects[row_index.first][col_index.first]->exists(
917 row_index.second, col_index.second);
922template <
typename SparsityPatternType>
927 const std::pair<size_type, size_type> row_index =
928 row_indices.global_to_local(row);
933 c += sub_objects[row_index.first][
b]->row_length(row_index.second);
940template <
typename SparsityPatternType>
949template <
typename SparsityPatternType>
959 const unsigned int index)
const
962 const std::pair<size_type, size_type> row_index =
971 unsigned int rowlen =
972 sub_objects[row_index.first][
b]->row_length(row_index.second);
973 if (index < c + rowlen)
974 return block_columns +
975 sub_objects[row_index.first][
b]->column_number(row_index.second,
978 block_columns +=
sub_objects[row_index.first][
b]->n_cols();
std::pair< unsigned int, size_type > global_to_local(const size_type i) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
void reinit(const std::vector< size_type > &row_block_sizes, const std::vector< size_type > &col_block_sizes)
std::vector< size_type > counter_within_block
BlockDynamicSparsityPattern()=default
size_type column_number(const size_type row, const unsigned int index) const
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
std::vector< std::vector< size_type > > block_column_indices
static ::ExceptionBase & ExcIncompatibleColNumbers(int arg1, int arg2, int arg3, int arg4)
static ::ExceptionBase & ExcIncompatibleRowNumbers(int arg1, int arg2, int arg3, int arg4)
#define Assert(cond, exc)
void copy_from(const BlockDynamicSparsityPattern &dsp)
Table< 2, SmartPointer< SparsityPatternType, BlockSparsityPatternBase< SparsityPatternType > > > sub_objects
void reinit(const size_type n_block_rows, const size_type n_block_columns)
#define AssertIndexRange(index, range)
BlockIndices column_indices
static ::ExceptionBase & ExcInternalError()
BlockSparsityPattern()=default
bool is_compressed() const
std::size_t memory_consumption() const
BlockSparsityPatternBase()
static const size_type invalid_entry
const SparsityPatternType & block(const size_type row, const size_type column) const
size_type n_block_rows() const
void print_gnuplot(std::ostream &out) const
void print(std::ostream &out) const
SparsityPatternType & block(const size_type row, const size_type column)
types::global_dof_index size_type
~BlockSparsityPatternBase() override
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false)
const BlockIndices & get_column_indices() const
size_type n_nonzero_elements() const
static const size_type invalid_entry
size_type max_entries_per_row() const
size_type n_block_cols() const
const BlockIndices & get_row_indices() const
unsigned int row_length(const size_type row) const
void add(const size_type i, const size_type j)
BlockSparsityPatternBase & operator=(const BlockSparsityPatternBase &)
void print_svg(std::ostream &out) const
bool exists(const size_type i, const size_type j) const
void reinit(const size_type n_block_rows, const size_type n_block_columns)
BlockSparsityPattern()=default
void reinit(const std::vector< size_type > &row_block_sizes, const std::vector< size_type > &col_block_sizes)
types::global_dof_index size_type
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
unsigned int global_dof_index