32template <
typename ForwardIterator>
36 const bool indices_are_sorted)
44 if (indices_are_sorted ==
true && n_elements > 3)
52 for (; test1 !=
end; ++test, ++test1)
67 ForwardIterator my_it =
begin;
69 std::vector<size_type>::iterator it =
106 std::vector<size_type>::iterator it2 = it + (
end - my_it);
115 else if (*my_it == *it2)
142 ForwardIterator my_it =
begin;
146 if (stop_size >
entries.capacity())
150 std::vector<size_type>::iterator it, it2;
183 for (; my_it !=
end; ++my_it)
218 return entries.capacity() *
sizeof(
size_type) +
sizeof(
Line);
232 , have_entries(false)
238 "This constructor can only be called if the provided argument "
239 "is the sparsity pattern for an empty matrix. This constructor can "
240 "not be used to copy-construct a non-empty sparsity pattern."));
249 , have_entries(false)
263 , have_entries(false)
277 "This operator can only be called if the provided argument "
278 "is the sparsity pattern for an empty matrix. This operator can "
279 "not be used to copy a non-empty sparsity pattern."));
282 ExcMessage(
"This operator can only be called if the current object is "
301 "The IndexSet argument to this function needs to either "
302 "be empty (indicating the complete set of rows), or have size "
303 "equal to the desired number of rows as specified by the "
304 "first argument to this function. (Of course, the number "
305 "of indices in this IndexSet may be less than the number "
306 "of rows, but the *size* of the IndexSet must be equal.)"));
310 lines.swap(new_lines);
324 return ((
rows == 0) && (
cols == 0));
336 for (
const auto &line :
lines)
350 const bool indices_are_sorted)
365 "The row IndexSet does not contain the index i. This sparsity pattern "
366 "object cannot know whether the entry (i, j) exists or not."));
379 return std::binary_search(
lines[rowindex].entries.begin(),
380 lines[rowindex].entries.end(),
406 if (rowindex != row_entry)
407 add(row_entry, rowindex);
427 lines[rowindex].entries = std::vector<size_type>();
436 view.
reinit(
rows.n_elements(), this->n_cols());
441 for (
auto it =
rows.begin(); it !=
end; ++it, ++view_row)
446 view.
lines[view_row].entries =
lines[rowindex].entries;
454template <
typename SparsityPatternTypeLeft,
typename SparsityPatternTypeRight>
457 const SparsityPatternTypeLeft & sp_A,
458 const SparsityPatternTypeRight &sp_B)
460 Assert(sp_A.n_rows() == sp_B.n_rows(),
463 this->
reinit(sp_A.n_cols(), sp_B.n_cols());
470 std::vector<size_type> new_cols;
471 new_cols.reserve(sp_B.max_entries_per_row());
474 for (
size_type i = 0; i < sp_A.n_rows(); ++i)
477 new_cols.resize(sp_B.row_length(i));
479 const auto last_il = sp_B.end(i);
480 auto * col_ptr = new_cols.data();
481 for (
auto il = sp_B.begin(i); il != last_il; ++il)
482 *col_ptr++ = il->column();
484 std::sort(new_cols.begin(), new_cols.end());
487 const auto last_ik = sp_A.end(i);
488 for (
auto ik = sp_A.begin(i); ik != last_ik; ++ik)
489 this->
add_entries(ik->column(), new_cols.begin(), new_cols.end(),
true);
495template <
typename SparsityPatternTypeLeft,
typename SparsityPatternTypeRight>
498 const SparsityPatternTypeLeft & left,
499 const SparsityPatternTypeRight &right)
501 Assert(left.n_cols() == right.n_rows(),
504 this->
reinit(left.n_rows(), right.n_cols());
506 typename SparsityPatternTypeLeft::iterator it_left = left.begin(),
507 end_left = left.end();
508 for (; it_left != end_left; ++it_left)
510 const unsigned int j = it_left->column();
516 typename SparsityPatternTypeRight::iterator it_right = right.begin(j),
517 end_right = right.end(j);
518 for (; it_right != end_right; ++it_right)
519 this->
add(it_left->row(), it_right->column());
532 for (
const auto entry :
lines[row].entries)
535 out <<
']' << std::endl;
551 for (
const auto entry :
lines[row].entries)
557 out << entry <<
" " << -
static_cast<signed int>(rowindex) << std::endl;
575 for (
const auto entry :
lines[row].entries)
577 std::abs(
static_cast<int>(rowindex - entry))) > b)
578 b =
std::abs(
static_cast<signed int>(rowindex - entry));
593 for (
const auto &line :
lines)
595 n += line.entries.size();
606 std::set<types::global_dof_index>
cols;
607 for (
const auto &line :
lines)
608 cols.insert(line.entries.begin(), line.entries.end());
623 std::vector<types::global_dof_index>
rows;
624 auto line =
lines.begin();
626 for (
const auto row : locally_stored_rows)
628 if (line->entries.size() > 0)
648 for (
const auto &line :
lines)
670 const auto &
cols =
lines[local_row].entries;
673 if ((it !=
cols.end()) && (*it == col))
674 return (it -
cols.begin());
688#ifndef DEAL_II_VECTOR_ITERATOR_IS_POINTER
691 std::vector<size_type>::iterator,
695 std::vector<size_type>::const_iterator,
696 std::vector<size_type>::const_iterator,
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_unique_and_sorted=false)
IndexSet nonempty_cols() const
DynamicSparsityPattern get_view(const IndexSet &rows) const
types::global_dof_index size_type
size_type memory_consumption() const
void compute_mmult_pattern(const SparsityPatternTypeLeft &left, const SparsityPatternTypeRight &right)
virtual void add_row_entries(const size_type &row, const ArrayView< const size_type > &columns, const bool indices_are_sorted=false) override
size_type column_index(const size_type row, const size_type col) const
DynamicSparsityPattern & operator=(const DynamicSparsityPattern &)
std::vector< Line > lines
void print(std::ostream &out) const
void reinit(const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
void clear_row(const size_type row)
size_type max_entries_per_row() const
void print_gnuplot(std::ostream &out) const
size_type bandwidth() const
void compute_Tmmult_pattern(const SparsityPatternTypeLeft &left, const SparsityPatternTypeRight &right)
size_type n_nonzero_elements() const
bool exists(const size_type i, const size_type j) const
void add(const size_type i, const size_type j)
IndexSet nonempty_rows() const
size_type index_within_set(const size_type global_index) const
size_type n_elements() const
bool is_element(const size_type index) const
size_type nth_index_in_set(const size_type local_index) const
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
virtual void resize(const size_type rows, const size_type cols)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIO()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
IndexSet complete_index_set(const IndexSet::size_type N)
std::enable_if_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
Iterator lower_bound(Iterator first, Iterator last, const T &val)
const types::global_dof_index invalid_size_type
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
void add_entries(ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted)
std::vector< size_type > entries
size_type memory_consumption() const