Reference documentation for deal.II version 9.2.0
\(\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\}}\)
dynamic_sparsity_pattern.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
17 #include <deal.II/base/utilities.h>
18 
21 
22 #include <algorithm>
23 #include <cmath>
24 #include <functional>
25 #include <numeric>
26 #include <set>
27 
29 
30 
31 
32 template <typename ForwardIterator>
33 void
35  ForwardIterator end,
36  const bool indices_are_sorted)
37 {
38  const int n_elements = end - begin;
39  if (n_elements <= 0)
40  return;
41 
42  const size_type stop_size = entries.size() + n_elements;
43 
44  if (indices_are_sorted == true && n_elements > 3)
45  {
46  // in debug mode, check whether the
47  // indices really are sorted.
48 #ifdef DEBUG
49  {
50  ForwardIterator test = begin, test1 = begin;
51  ++test1;
52  for (; test1 != end; ++test, ++test1)
53  Assert(*test1 > *test, ExcInternalError());
54  }
55 #endif
56 
57  if (entries.size() == 0 || entries.back() < *begin)
58  {
59  entries.insert(entries.end(), begin, end);
60  return;
61  }
62 
63  // find a possible insertion point for
64  // the first entry. check whether the
65  // first entry is a duplicate before
66  // actually doing something.
67  ForwardIterator my_it = begin;
68  size_type col = *my_it;
69  std::vector<size_type>::iterator it =
70  Utilities::lower_bound(entries.begin(), entries.end(), col);
71  while (*it == col)
72  {
73  ++my_it;
74  if (my_it == end)
75  break;
76  col = *my_it;
77  // check the very next entry in the
78  // current array
79  ++it;
80  if (it == entries.end())
81  break;
82  if (*it > col)
83  break;
84  if (*it == col)
85  continue;
86  // ok, it wasn't the very next one, do a
87  // binary search to find the insert point
88  it = Utilities::lower_bound(it, entries.end(), col);
89  if (it == entries.end())
90  break;
91  }
92  // all input entries were duplicates.
93  if (my_it == end)
94  return;
95 
96  // resize vector by just inserting the
97  // list
98  const size_type pos1 = it - entries.begin();
99  Assert(pos1 <= entries.size(), ExcInternalError());
100  entries.insert(it, my_it, end);
101  it = entries.begin() + pos1;
102  Assert(entries.size() >= static_cast<size_type>(it - entries.begin()),
103  ExcInternalError());
104 
105  // now merge the two lists.
106  std::vector<size_type>::iterator it2 = it + (end - my_it);
107 
108  // as long as there are indices both in
109  // the end of the entries list and in the
110  // input list
111  while (my_it != end && it2 != entries.end())
112  {
113  if (*my_it < *it2)
114  *it++ = *my_it++;
115  else if (*my_it == *it2)
116  {
117  *it++ = *it2++;
118  ++my_it;
119  }
120  else
121  *it++ = *it2++;
122  }
123  // in case there are indices left in the
124  // input list
125  while (my_it != end)
126  *it++ = *my_it++;
127 
128  // in case there are indices left in the
129  // end of entries
130  while (it2 != entries.end())
131  *it++ = *it2++;
132 
133  // resize and return
134  const size_type new_size = it - entries.begin();
135  Assert(new_size <= stop_size, ExcInternalError());
136  entries.resize(new_size);
137  return;
138  }
139 
140  // unsorted case or case with too few
141  // elements
142  ForwardIterator my_it = begin;
143 
144  // If necessary, increase the size of the
145  // array.
146  if (stop_size > entries.capacity())
147  entries.reserve(stop_size);
148 
149  size_type col = *my_it;
150  std::vector<size_type>::iterator it, it2;
151  // insert the first element as for one
152  // entry only first check the last
153  // element (or if line is still empty)
154  if ((entries.size() == 0) || (entries.back() < col))
155  {
156  entries.push_back(col);
157  it = entries.end() - 1;
158  }
159  else
160  {
161  // do a binary search to find the place
162  // where to insert:
163  it2 = Utilities::lower_bound(entries.begin(), entries.end(), col);
164 
165  // If this entry is a duplicate, continue
166  // immediately Insert at the right place
167  // in the vector. Vector grows
168  // automatically to fit elements. Always
169  // doubles its size.
170  if (*it2 != col)
171  it = entries.insert(it2, col);
172  else
173  it = it2;
174  }
175 
176  ++my_it;
177  // Now try to be smart and insert with
178  // bias in the direction we are
179  // walking. This has the advantage that
180  // for sorted lists, we always search in
181  // the right direction, what should
182  // decrease the work needed in here.
183  for (; my_it != end; ++my_it)
184  {
185  col = *my_it;
186  // need a special insertion command when
187  // we're at the end of the list
188  if (col > entries.back())
189  {
190  entries.push_back(col);
191  it = entries.end() - 1;
192  }
193  // search to the right (preferred search
194  // direction)
195  else if (col > *it)
196  {
197  it2 = Utilities::lower_bound(it++, entries.end(), col);
198  if (*it2 != col)
199  it = entries.insert(it2, col);
200  }
201  // search to the left
202  else if (col < *it)
203  {
204  it2 = Utilities::lower_bound(entries.begin(), it, col);
205  if (*it2 != col)
206  it = entries.insert(it2, col);
207  }
208  // if we're neither larger nor smaller,
209  // then this was a duplicate and we can
210  // just continue.
211  }
212 }
213 
214 
217 {
218  return entries.capacity() * sizeof(size_type) + sizeof(Line);
219 }
220 
221 
223  : have_entries(false)
224  , rows(0)
225  , cols(0)
226  , rowset(0)
227 {}
228 
229 
230 
232  : Subscriptor()
233  , have_entries(false)
234  , rows(0)
235  , cols(0)
236  , rowset(0)
237 {
238  (void)s;
239  Assert(s.rows == 0 && s.cols == 0,
240  ExcMessage(
241  "This constructor can only be called if the provided argument "
242  "is the sparsity pattern for an empty matrix. This constructor can "
243  "not be used to copy-construct a non-empty sparsity pattern."));
244 }
245 
246 
247 
249  const size_type n,
250  const IndexSet &rowset_)
251  : have_entries(false)
252  , rows(0)
253  , cols(0)
254  , rowset(0)
255 {
256  reinit(m, n, rowset_);
257 }
258 
259 
261  : have_entries(false)
262  , rows(0)
263  , cols(0)
264  , rowset(0)
265 {
266  reinit(rowset_.size(), rowset_.size(), rowset_);
267 }
268 
269 
271  : have_entries(false)
272  , rows(0)
273  , cols(0)
274  , rowset(0)
275 {
276  reinit(n, n);
277 }
278 
279 
280 
283 {
284  (void)s;
285  Assert(s.rows == 0 && s.cols == 0,
286  ExcMessage(
287  "This operator can only be called if the provided argument "
288  "is the sparsity pattern for an empty matrix. This operator can "
289  "not be used to copy a non-empty sparsity pattern."));
290 
291  Assert(rows == 0 && cols == 0,
292  ExcMessage("This operator can only be called if the current object is "
293  "empty."));
294 
295  return *this;
296 }
297 
298 
299 
300 void
302  const size_type n,
303  const IndexSet &rowset_)
304 {
305  have_entries = false;
306  rows = m;
307  cols = n;
308  rowset = rowset_;
309 
310  Assert(rowset.size() == 0 || rowset.size() == m,
311  ExcMessage(
312  "The IndexSet argument to this function needs to either "
313  "be empty (indicating the complete set of rows), or have size "
314  "equal to the desired number of rows as specified by the "
315  "first argument to this function. (Of course, the number "
316  "of indices in this IndexSet may be less than the number "
317  "of rows, but the *size* of the IndexSet must be equal.)"));
318 
319  std::vector<Line> new_lines(rowset.size() == 0 ? rows : rowset.n_elements());
320  lines.swap(new_lines);
321 }
322 
323 
324 
325 void
327 {}
328 
329 
330 
331 bool
333 {
334  return ((rows == 0) && (cols == 0));
335 }
336 
337 
338 
341 {
342  if (!have_entries)
343  return 0;
344 
345  size_type m = 0;
346  for (const auto &line : lines)
347  {
348  m = std::max(m, static_cast<size_type>(line.entries.size()));
349  }
350 
351  return m;
352 }
353 
354 
355 
356 bool
358 {
361  Assert(
362  rowset.size() == 0 || rowset.is_element(i),
363  ExcMessage(
364  "The row IndexSet does not contain the index i. This sparsity pattern "
365  "object cannot know whether the entry (i, j) exists or not."));
366 
367  // Avoid a segmentation fault in below code if the row index happens to
368  // not be present in the IndexSet rowset:
369  if (!(rowset.size() == 0 || rowset.is_element(i)))
370  return false;
371 
372  if (!have_entries)
373  return false;
374 
375  const size_type rowindex =
376  rowset.size() == 0 ? i : rowset.index_within_set(i);
377 
378  return std::binary_search(lines[rowindex].entries.begin(),
379  lines[rowindex].entries.end(),
380  j);
381 }
382 
383 
384 
385 void
387 {
389 
390  // loop over all elements presently in the sparsity pattern and add the
391  // transpose element. note:
392  //
393  // 1. that the sparsity pattern changes which we work on, but not the present
394  // row
395  //
396  // 2. that the @p{add} function can be called on elements that already exist
397  // without any harm
398  for (size_type row = 0; row < lines.size(); ++row)
399  {
400  const size_type rowindex =
401  rowset.size() == 0 ? row : rowset.nth_index_in_set(row);
402 
403  for (const size_type row_entry : lines[row].entries)
404  // add the transpose entry if this is not the diagonal
405  if (rowindex != row_entry)
406  add(row_entry, rowindex);
407  }
408 }
409 
410 
411 
412 void
414 {
415  AssertIndexRange(row, n_rows());
416  if (!have_entries)
417  return;
418 
419  if (rowset.size() > 0 && !rowset.is_element(row))
420  return;
421 
422  const size_type rowindex =
423  rowset.size() == 0 ? row : rowset.index_within_set(row);
424 
425  AssertIndexRange(rowindex, lines.size());
426  lines[rowindex].entries = std::vector<size_type>();
427 }
428 
429 
430 
433 {
435  view.reinit(rows.n_elements(), this->n_cols());
436  AssertDimension(rows.size(), this->n_rows());
437 
438  const auto end = rows.end();
440  for (auto it = rows.begin(); it != end; ++it, ++view_row)
441  {
442  const size_type rowindex =
443  rowset.size() == 0 ? *it : rowset.index_within_set(*it);
444 
445  view.lines[view_row].entries = lines[rowindex].entries;
446  view.have_entries |= (lines[rowindex].entries.size() > 0);
447  }
448  return view;
449 }
450 
451 
452 
453 template <typename SparsityPatternTypeLeft, typename SparsityPatternTypeRight>
454 void
456  const SparsityPatternTypeLeft & sp_A,
457  const SparsityPatternTypeRight &sp_B)
458 {
459  Assert(sp_A.n_rows() == sp_B.n_rows(),
460  ExcDimensionMismatch(sp_A.n_rows(), sp_B.n_rows()));
461 
462  this->reinit(sp_A.n_cols(), sp_B.n_cols());
463  // we will go through all the
464  // rows in the matrix A, and for each column in a row we add the whole
465  // row of matrix B with that row number. This means that we will insert
466  // a lot of entries to each row, which is best handled by the
467  // DynamicSparsityPattern class.
468 
469  std::vector<size_type> new_cols;
470  new_cols.reserve(sp_B.max_entries_per_row());
471 
472  // C_{kl} = A_{ik} B_{il}
473  for (size_type i = 0; i < sp_A.n_rows(); ++i)
474  {
475  // get all column numbers from sp_B in a temporary vector:
476  new_cols.resize(sp_B.row_length(i));
477  {
478  const auto last_il = sp_B.end(i);
479  auto * col_ptr = new_cols.data();
480  for (auto il = sp_B.begin(i); il != last_il; ++il)
481  *col_ptr++ = il->column();
482  }
483  std::sort(new_cols.begin(), new_cols.end());
484 
485  // now for each k, add new_cols to the target sparsity
486  const auto last_ik = sp_A.end(i);
487  for (auto ik = sp_A.begin(i); ik != last_ik; ++ik)
488  this->add_entries(ik->column(), new_cols.begin(), new_cols.end(), true);
489  }
490 }
491 
492 
493 
494 template <typename SparsityPatternTypeLeft, typename SparsityPatternTypeRight>
495 void
497  const SparsityPatternTypeLeft & left,
498  const SparsityPatternTypeRight &right)
499 {
500  Assert(left.n_cols() == right.n_rows(),
501  ExcDimensionMismatch(left.n_cols(), right.n_rows()));
502 
503  this->reinit(left.n_rows(), right.n_cols());
504 
505  typename SparsityPatternTypeLeft::iterator it_left = left.begin(),
506  end_left = left.end();
507  for (; it_left != end_left; ++it_left)
508  {
509  const unsigned int j = it_left->column();
510 
511  // We are sitting on entry (i,j) of the left sparsity pattern. We then
512  // need to add all entries (i,k) to the final sparsity pattern where (j,k)
513  // exists in the right sparsity pattern -- i.e., we need to iterate over
514  // row j.
515  typename SparsityPatternTypeRight::iterator it_right = right.begin(j),
516  end_right = right.end(j);
517  for (; it_right != end_right; ++it_right)
518  this->add(it_left->row(), it_right->column());
519  }
520 }
521 
522 
523 
524 void
525 DynamicSparsityPattern::print(std::ostream &out) const
526 {
527  for (size_type row = 0; row < lines.size(); ++row)
528  {
529  out << '[' << (rowset.size() == 0 ? row : rowset.nth_index_in_set(row));
530 
531  for (const auto entry : lines[row].entries)
532  out << ',' << entry;
533 
534  out << ']' << std::endl;
535  }
536 
537  AssertThrow(out, ExcIO());
538 }
539 
540 
541 
542 void
544 {
545  for (size_type row = 0; row < lines.size(); ++row)
546  {
547  const size_type rowindex =
548  rowset.size() == 0 ? row : rowset.nth_index_in_set(row);
549 
550  for (const auto entry : lines[row].entries)
551  // while matrix entries are usually
552  // written (i,j), with i vertical and
553  // j horizontal, gnuplot output is
554  // x-y, that is we have to exchange
555  // the order of output
556  out << entry << " " << -static_cast<signed int>(rowindex) << std::endl;
557  }
558 
559 
560  AssertThrow(out, ExcIO());
561 }
562 
563 
564 
567 {
568  size_type b = 0;
569  for (size_type row = 0; row < lines.size(); ++row)
570  {
571  const size_type rowindex =
572  rowset.size() == 0 ? row : rowset.nth_index_in_set(row);
573 
574  for (const auto entry : lines[row].entries)
575  if (static_cast<size_type>(
576  std::abs(static_cast<int>(rowindex - entry))) > b)
577  b = std::abs(static_cast<signed int>(rowindex - entry));
578  }
579 
580  return b;
581 }
582 
583 
584 
587 {
588  if (!have_entries)
589  return 0;
590 
591  size_type n = 0;
592  for (const auto &line : lines)
593  {
594  n += line.entries.size();
595  }
596 
597  return n;
598 }
599 
600 
601 
602 IndexSet
604 {
605  std::set<types::global_dof_index> cols;
606  for (const auto &line : lines)
607  cols.insert(line.entries.begin(), line.entries.end());
608 
609  IndexSet res(this->n_cols());
610  res.add_indices(cols.begin(), cols.end());
611  return res;
612 }
613 
614 
615 
616 IndexSet
618 {
619  const IndexSet all_rows = complete_index_set(this->n_rows());
620  const IndexSet &locally_stored_rows = rowset.size() == 0 ? all_rows : rowset;
621 
622  std::vector<types::global_dof_index> rows;
623  auto line = lines.begin();
624  AssertDimension(locally_stored_rows.n_elements(), lines.size());
625  for (const auto &row : locally_stored_rows)
626  {
627  if (line->entries.size() > 0)
628  rows.push_back(row);
629 
630  ++line;
631  }
632 
633  IndexSet res(this->n_rows());
634  res.add_indices(rows.begin(), rows.end());
635  return res;
636 }
637 
638 
639 
642 {
643  size_type mem = sizeof(DynamicSparsityPattern) +
645  sizeof(rowset);
646 
647  for (const auto &line : lines)
649 
650  return mem;
651 }
652 
653 
654 
658  const DynamicSparsityPattern::size_type col) const
659 {
660  AssertIndexRange(row, n_rows());
661  AssertIndexRange(col, n_cols());
662  Assert(rowset.size() == 0 || rowset.is_element(row), ExcInternalError());
663 
664  const DynamicSparsityPattern::size_type local_row =
665  rowset.size() ? rowset.index_within_set(row) : row;
666 
667  // now we need to do a binary search. Note that col indices are assumed to
668  // be sorted.
669  const auto &cols = lines[local_row].entries;
670  auto it = Utilities::lower_bound(cols.begin(), cols.end(), col);
671 
672  if ((it != cols.end()) && (*it == col))
673  return (it - cols.begin());
674  else
676 }
677 
678 
679 
680 // explicit instantiations
681 template void
683 template void
685  const size_type *,
686  const bool);
687 #ifndef DEAL_II_VECTOR_ITERATOR_IS_POINTER
688 template void
689 DynamicSparsityPattern::Line::add_entries(std::vector<size_type>::iterator,
690  std::vector<size_type>::iterator,
691  const bool);
692 template void
694  std::vector<size_type>::const_iterator,
695  std::vector<size_type>::const_iterator,
696  const bool);
697 #endif
698 
699 template void
701  const DynamicSparsityPattern &);
702 template void
704  const SparsityPattern &);
705 template void
707  const DynamicSparsityPattern &);
708 template void
710  const SparsityPattern &);
711 
712 template void
714  const SparsityPattern &);
715 template void
717  const SparsityPattern &);
718 template void
720  const DynamicSparsityPattern &);
721 template void
723  const DynamicSparsityPattern &);
724 
IndexSet
Definition: index_set.h:74
dynamic_sparsity_pattern.h
DynamicSparsityPattern::n_nonzero_elements
size_type n_nonzero_elements() const
Definition: dynamic_sparsity_pattern.cc:586
DynamicSparsityPattern::add_entries
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_unique_and_sorted=false)
Definition: dynamic_sparsity_pattern.h:1051
DynamicSparsityPattern::n_cols
size_type n_cols() const
Definition: dynamic_sparsity_pattern.h:1024
StandardExceptions::ExcIO
static ::ExceptionBase & ExcIO()
DynamicSparsityPattern::nonempty_rows
IndexSet nonempty_rows() const
Definition: dynamic_sparsity_pattern.cc:617
DynamicSparsityPattern::Line::add_entries
void add_entries(ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted)
Definition: dynamic_sparsity_pattern.cc:34
DynamicSparsityPattern::begin
iterator begin() const
Definition: dynamic_sparsity_pattern.h:1105
memory_consumption.h
IndexSet::size
size_type size() const
Definition: index_set.h:1635
utilities.h
DynamicSparsityPattern::get_view
DynamicSparsityPattern get_view(const IndexSet &rows) const
Definition: dynamic_sparsity_pattern.cc:432
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
DynamicSparsityPattern::max_entries_per_row
size_type max_entries_per_row() const
Definition: dynamic_sparsity_pattern.cc:340
DynamicSparsityPattern::nonempty_cols
IndexSet nonempty_cols() const
Definition: dynamic_sparsity_pattern.cc:603
DynamicSparsityPattern::print_gnuplot
void print_gnuplot(std::ostream &out) const
Definition: dynamic_sparsity_pattern.cc:543
DynamicSparsityPattern::Line::entries
std::vector< size_type > entries
Definition: dynamic_sparsity_pattern.h:708
LACExceptions::ExcNotQuadratic
static ::ExceptionBase & ExcNotQuadratic()
Subscriptor
Definition: subscriptor.h:62
IndexSet::n_elements
size_type n_elements() const
Definition: index_set.h:1833
DynamicSparsityPattern::empty
bool empty() const
Definition: dynamic_sparsity_pattern.cc:332
DynamicSparsityPattern::have_entries
bool have_entries
Definition: dynamic_sparsity_pattern.h:676
DynamicSparsityPattern::print
void print(std::ostream &out) const
Definition: dynamic_sparsity_pattern.cc:525
DynamicSparsityPattern::Line::memory_consumption
size_type memory_consumption() const
Definition: dynamic_sparsity_pattern.cc:216
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
DynamicSparsityPattern::column_index
size_type column_index(const size_type row, const size_type col) const
Definition: dynamic_sparsity_pattern.cc:656
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
DynamicSparsityPattern::exists
bool exists(const size_type i, const size_type j) const
Definition: dynamic_sparsity_pattern.cc:357
Utilities::lower_bound
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:1102
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
DynamicSparsityPattern::n_rows
size_type n_rows() const
Definition: dynamic_sparsity_pattern.h:1016
MemoryConsumption::memory_consumption
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Definition: memory_consumption.h:268
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
DynamicSparsityPattern::size_type
types::global_dof_index size_type
Definition: dynamic_sparsity_pattern.h:329
DynamicSparsityPattern::memory_consumption
size_type memory_consumption() const
Definition: dynamic_sparsity_pattern.cc:641
DynamicSparsityPattern
Definition: dynamic_sparsity_pattern.h:323
IndexSet::is_element
bool is_element(const size_type index) const
Definition: index_set.h:1766
SparsityPattern
Definition: sparsity_pattern.h:865
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
DynamicSparsityPattern::operator=
DynamicSparsityPattern & operator=(const DynamicSparsityPattern &)
Definition: dynamic_sparsity_pattern.cc:282
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
DynamicSparsityPattern::clear_row
void clear_row(const size_type row)
Definition: dynamic_sparsity_pattern.cc:413
sparsity_pattern.h
DynamicSparsityPattern::symmetrize
void symmetrize()
Definition: dynamic_sparsity_pattern.cc:386
complete_index_set
IndexSet complete_index_set(const IndexSet::size_type N)
Definition: index_set.h:1014
DynamicSparsityPattern::cols
size_type cols
Definition: dynamic_sparsity_pattern.h:686
unsigned int
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
LinearAlgebra::CUDAWrappers::kernel::size_type
types::global_dof_index size_type
Definition: cuda_kernels.h:45
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
IndexSet::nth_index_in_set
size_type nth_index_in_set(const size_type local_index) const
Definition: index_set.h:1881
DynamicSparsityPattern::end
iterator end() const
Definition: dynamic_sparsity_pattern.h:1115
DynamicSparsityPattern::Line
Definition: dynamic_sparsity_pattern.h:701
numbers::invalid_size_type
const types::global_dof_index invalid_size_type
Definition: types.h:200
DynamicSparsityPattern::compute_mmult_pattern
void compute_mmult_pattern(const SparsityPatternTypeLeft &left, const SparsityPatternTypeRight &right)
Definition: dynamic_sparsity_pattern.cc:496
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
IndexSet::index_within_set
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1922
DynamicSparsityPattern::rows
size_type rows
Definition: dynamic_sparsity_pattern.h:681
DynamicSparsityPattern::reinit
void reinit(const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
Definition: dynamic_sparsity_pattern.cc:301
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
DynamicSparsityPattern::rowset
IndexSet rowset
Definition: dynamic_sparsity_pattern.h:692
DynamicSparsityPattern::compress
void compress()
Definition: dynamic_sparsity_pattern.cc:326
IndexSet::add_indices
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1704
DynamicSparsityPattern::add
void add(const size_type i, const size_type j)
Definition: dynamic_sparsity_pattern.h:1032
DynamicSparsityPattern::DynamicSparsityPattern
DynamicSparsityPattern()
Definition: dynamic_sparsity_pattern.cc:222
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
DynamicSparsityPattern::bandwidth
size_type bandwidth() const
Definition: dynamic_sparsity_pattern.cc:566
DynamicSparsityPattern::lines
std::vector< Line > lines
Definition: dynamic_sparsity_pattern.h:736
DynamicSparsityPattern::compute_Tmmult_pattern
void compute_Tmmult_pattern(const SparsityPatternTypeLeft &left, const SparsityPatternTypeRight &right)
Definition: dynamic_sparsity_pattern.cc:455
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)