17 #include <deal.II/base/std_cxx14/memory.h> 18 #include <deal.II/base/vector_slice.h> 19 #include <deal.II/base/utilities.h> 20 #include <deal.II/lac/sparsity_pattern.h> 21 #include <deal.II/lac/sparsity_tools.h> 22 #include <deal.II/lac/full_matrix.h> 23 #include <deal.II/lac/dynamic_sparsity_pattern.h> 32 DEAL_II_NAMESPACE_OPEN
48 store_diagonal_first_in_row(false)
63 store_diagonal_first_in_row(false)
67 ExcMessage(
"This constructor can only be called if the provided argument " 68 "is the sparsity pattern for an empty matrix. This constructor can " 69 "not be used to copy-construct a non-empty sparsity pattern."));
78 const unsigned int max_per_row)
85 store_diagonal_first_in_row(m == n)
94 const std::vector<unsigned int> &row_lengths)
100 store_diagonal_first_in_row(m == n)
102 reinit (m, n, row_lengths);
108 const unsigned int max_per_row)
115 reinit (m, m, max_per_row);
121 const std::vector<unsigned int> &row_lengths)
128 reinit (m, m, row_lengths);
134 const unsigned int max_per_row,
171 original_last_before_side_diagonals
172 = (row > extra_off_diagonals ?
175 row-extra_off_diagonals) :
179 original_first_after_side_diagonals
180 = (row <
rows-extra_off_diagonals-1 ?
181 std::upper_bound (original_row_start,
183 row+extra_off_diagonals) :
191 next_free_slot = std::copy (original_row_start,
192 original_last_before_side_diagonals,
196 for (
size_type i=1; i<=std::min(row,extra_off_diagonals);
197 ++i, ++next_free_slot)
198 *next_free_slot = row-i;
199 for (
size_type i=1; i<=std::min(extra_off_diagonals,
rows-row-1);
200 ++i, ++next_free_slot)
201 *next_free_slot = row+i;
204 next_free_slot = std::copy (original_first_after_side_diagonals,
223 ExcMessage(
"This operator can only be called if the provided argument " 224 "is the sparsity pattern for an empty matrix. This operator can " 225 "not be used to copy a non-empty sparsity pattern."));
228 ExcMessage(
"This operator can only be called if the current object is " 239 const unsigned int max_per_row)
242 const std::vector<unsigned int> row_lengths (m, max_per_row);
243 reinit (m, n, row_lengths);
251 const VectorSlice<
const std::vector<unsigned int> > &row_lengths)
259 if ((m==0) || (n==0))
284 std::size_t vec_len = 0;
287 std::max(row_lengths[i], 1U) :
303 std::min (static_cast<size_type>(*std::max_element(row_lengths.begin(),
334 std::max(std::min(static_cast<size_type>(row_lengths[i-1]),n),
335 static_cast<size_type> (1U)) :
336 std::min(static_cast<size_type>(row_lengths[i-1]),n));
376 const std::size_t nonzero_elements
379 std::bind(std::not_equal_to<size_type>(), std::placeholders::_1,
invalid_entry));
381 std::unique_ptr<size_type[]> new_colnums (
new size_type[nonzero_elements]);
405 ? tmp_entries.begin()+1
406 : tmp_entries.begin(),
411 new_colnums[next_free_entry++] = tmp_entries[j];
415 next_row_start = next_free_entry;
433 (std::find (&new_colnums[
rowstart[line]+1],
434 &new_colnums[next_row_start],
436 &new_colnums[next_row_start]),
440 (std::adjacent_find(&new_colnums[
rowstart[line]+1],
441 &new_colnums[next_row_start]) ==
442 &new_colnums[next_row_start]),
447 Assert (next_free_entry == nonzero_elements,
455 colnums = std::move(new_colnums);
472 const bool do_diag_optimize = (sp.
n_rows() == sp.
n_cols());
473 std::vector<unsigned int> row_lengths (sp.
n_rows());
477 if (do_diag_optimize && !sp.
exists(i,i))
490 end_row = sp.
end (row);
492 for (; col_num != end_row; ++col_num)
495 if ((col!=row) || !do_diag_optimize)
516 const bool do_diag_optimize = (dsp.
n_rows() == dsp.
n_cols());
517 std::vector<unsigned int> row_lengths (dsp.
n_rows());
521 if (do_diag_optimize && !dsp.
exists(i,i))
531 for (
unsigned int index=0; index <
row_length; ++index)
534 if ((col!=row) || !do_diag_optimize)
547 template <
typename number>
553 const bool matrix_is_square = (matrix.m() == matrix.n());
555 std::vector<unsigned int> entries_per_row (matrix.m(), 0);
556 for (
size_type row=0; row<matrix.m(); ++row)
558 for (
size_type col=0; col<matrix.n(); ++col)
559 if (matrix(row,col) != 0)
560 ++entries_per_row[row];
563 (matrix(row,row) == 0))
564 ++entries_per_row[row];
567 reinit (matrix.m(), matrix.n(), entries_per_row);
574 std::vector<size_type> column_indices;
575 column_indices.reserve (entries_per_row.size() > 0
577 *std::max_element (entries_per_row.begin(),
578 entries_per_row.end())
581 for (
size_type row=0; row<matrix.m(); ++row)
583 column_indices.resize(entries_per_row[row]);
586 for (
size_type col=0; col<matrix.n(); ++col)
587 if (matrix(row,col) != 0)
589 column_indices[current_index] = col;
600 column_indices[current_index] = row;
608 add_entries(row, column_indices.begin(), column_indices.end(),
true);
619 const std::vector<unsigned int> &row_lengths)
621 reinit(m, n, make_slice(row_lengths));
697 = Utilities::lower_bound<const size_type *> (sorted_region_start,
736 template <
typename ForwardIterator>
739 ForwardIterator begin,
741 const bool indices_are_sorted)
743 if (indices_are_sorted ==
true)
747 ForwardIterator it =
begin;
748 bool has_larger_entries =
false;
756 has_larger_entries =
true;
759 if (has_larger_entries ==
false)
760 for ( ; it !=
end; ++it)
771 for (ForwardIterator p =
begin; p !=
end; ++p)
778 for (ForwardIterator it =
begin; it !=
end; ++it)
794 if (
colnums[k] == j)
return true;
818 std::pair<SparsityPattern::size_type, SparsityPattern::size_type>
838 return std::make_pair (row,col);
890 out <<
']' << std::endl;
911 out <<
colnums[j] <<
" " << -
static_cast<signed int>(i) << std::endl;
919 unsigned int m = this->
n_rows();
920 unsigned int n = this->
n_cols();
921 out <<
"<svg xmlns=\"http://www.w3.org/2000/svg\" version=\"1.1\" viewBox=\"0 0 " << n+2
922 <<
" " << m+2 <<
" \">\n" 923 "<style type=\"text/css\" >\n" 930 " <rect width=\"" << n+2 <<
"\" height=\"" << m+2 <<
"\" fill=\"rgb(128, 128, 128)\"/>\n" 931 " <rect x=\"1\" y=\"1\" width=\"" << n+0.1 <<
"\" height=\"" << m+0.1
932 <<
"\" fill=\"rgb(255, 255, 255)\"/>\n\n";
937 for (; it!=
end; ++it)
939 out <<
" <rect class=\"pixel\" x=\"" << it->
column()+1
940 <<
"\" y=\"" << it->
row()+1
941 <<
"\" width=\".9\" height=\".9\"/>\n";
943 out <<
"</svg>" << std::endl;
959 if (static_cast<size_type>(std::abs(static_cast<int>(i-
colnums[j]))) > b)
960 b = std::abs(static_cast<signed int>(i-
colnums[j]));
983 out.write (reinterpret_cast<const char *>(
rowstart.get()),
985 - reinterpret_cast<const char *>(
rowstart.get()));
987 out.write (reinterpret_cast<const char *>(
colnums.get()),
989 - reinterpret_cast<const char *>(
colnums.get()));
1025 in.read (reinterpret_cast<char *>(
rowstart.get()),
1027 - reinterpret_cast<char *>(
rowstart.get()));
1032 in.read (reinterpret_cast<char *>(
colnums.get()),
1034 - reinterpret_cast<char *>(
colnums.get()));
1055 template void SparsityPattern::add_entries<const SparsityPattern::size_type *> (
const size_type,
1059 #ifndef DEAL_II_VECTOR_ITERATOR_IS_POINTER 1060 template void SparsityPattern::add_entries<std::vector<SparsityPattern::size_type>::const_iterator>
1062 std::vector<size_type>::const_iterator,
1063 std::vector<size_type>::const_iterator,
1066 template void SparsityPattern::add_entries<std::vector<SparsityPattern::size_type>::iterator>
1068 std::vector<size_type>::iterator,
1069 std::vector<size_type>::iterator,
1072 DEAL_II_NAMESPACE_CLOSE
Iterator lower_bound(Iterator first, Iterator last, const T &val)
void block_write(std::ostream &out) const
const types::global_dof_index invalid_size_type
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false)
static ::ExceptionBase & ExcMatrixIsCompressed()
#define AssertDimension(dim1, dim2)
std::size_t n_nonzero_elements() const
size_type operator()(const size_type i, const size_type j) const
static ::ExceptionBase & ExcIO()
size_type column_number(const size_type row, const size_type index) const
static const size_type invalid_entry
std::size_t memory_consumption() const
std::pair< size_type, size_type > matrix_position(const std::size_t global_index) const
void print_gnuplot(std::ostream &out) const
bool exists(const size_type i, const size_type j) const
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
bool is_compressed() const
size_type bandwidth() const
std::unique_ptr< size_type[]> colnums
static ::ExceptionBase & ExcNotCompressed()
static ::ExceptionBase & ExcMessage(std::string arg1)
void add(const size_type i, const size_type j)
#define Assert(cond, exc)
void print_svg(std::ostream &out) const
void print(std::ostream &out) const
SparsityPattern & operator=(const SparsityPattern &)
size_type max_entries_per_row() const
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
void reinit(const size_type m, const size_type n, const unsigned int max_per_row)
unsigned int max_row_length
static ::ExceptionBase & ExcNotQuadratic()
types::global_dof_index size_type
static ::ExceptionBase & ExcNotEnoughSpace(int arg1, int arg2)
size_type row_length(const size_type row) const
bool store_diagonal_first_in_row
std::unique_ptr< std::size_t[]> rowstart
static ::ExceptionBase & ExcEmptyObject()
void block_read(std::istream &in)
bool exists(const size_type i, const size_type j) const
unsigned int row_length(const size_type row) const
size_type row_position(const size_type i, const size_type j) const
static ::ExceptionBase & ExcInternalError()