Reference documentation for deal.II version 9.0.0
sparsity_pattern.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2018 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
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>
24 
25 #include <iostream>
26 #include <iomanip>
27 #include <algorithm>
28 #include <cmath>
29 #include <numeric>
30 #include <functional>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 #ifdef DEAL_II_MSVC
35 __declspec(selectany) // Weak extern binding due to multiple link error
36 #endif
38 
39 
40 
42  :
43  max_dim(0),
44  max_vec_len(0),
45  rowstart(nullptr),
46  colnums(nullptr),
47  compressed(false),
48  store_diagonal_first_in_row(false)
49 {
50  reinit (0,0,0);
51 }
52 
53 
54 
56  :
57  Subscriptor(),
58  max_dim(0),
59  max_vec_len(0),
60  rowstart(nullptr),
61  colnums(nullptr),
62  compressed(false),
63  store_diagonal_first_in_row(false)
64 {
65  (void)s;
66  Assert (s.empty(),
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."));
70 
71  reinit (0,0,0);
72 }
73 
74 
75 
77  const size_type n,
78  const unsigned int max_per_row)
79  :
80  max_dim(0),
81  max_vec_len(0),
82  rowstart(nullptr),
83  colnums(nullptr),
84  compressed(false),
85  store_diagonal_first_in_row(m == n)
86 {
87  reinit (m,n,max_per_row);
88 }
89 
90 
91 
93  const size_type n,
94  const std::vector<unsigned int> &row_lengths)
95  :
96  max_dim(0),
97  max_vec_len(0),
98  rowstart(nullptr),
99  colnums(nullptr),
100  store_diagonal_first_in_row(m == n)
101 {
102  reinit (m, n, row_lengths);
103 }
104 
105 
106 
108  const unsigned int max_per_row)
109  :
110  max_dim(0),
111  max_vec_len(0),
112  rowstart(nullptr),
113  colnums(nullptr)
114 {
115  reinit (m, m, max_per_row);
116 }
117 
118 
119 
121  const std::vector<unsigned int> &row_lengths)
122  :
123  max_dim(0),
124  max_vec_len(0),
125  rowstart(nullptr),
126  colnums(nullptr)
127 {
128  reinit (m, m, row_lengths);
129 }
130 
131 
132 
134  const unsigned int max_per_row,
135  const size_type extra_off_diagonals)
136  :
137  max_dim(0),
138  max_vec_len(0),
139  rowstart(nullptr),
140  colnums(nullptr)
141 {
142  Assert (original.rows==original.cols, ExcNotQuadratic());
143  Assert (original.is_compressed(), ExcNotCompressed());
144 
145  reinit (original.rows, original.cols, max_per_row);
146 
147  // now copy the entries from the other object
148  for (size_type row=0; row<original.rows; ++row)
149  {
150  // copy the elements of this row of the other object
151  //
152  // note that the first object actually is the main-diagonal element,
153  // which we need not copy
154  //
155  // we do the copying in two steps: first we note that the elements in
156  // @p{original} are sorted, so we may first copy all the elements up to
157  // the first side-diagonal one which is to be filled in. then we insert
158  // the side-diagonals, finally copy the rest from that element onwards
159  // which is not a side-diagonal any more.
160  const size_type *const
161  original_row_start = &original.colnums[original.rowstart[row]] + 1;
162  // the following requires that @p{original} be compressed since
163  // otherwise there might be invalid_entry's
164  const size_type *const
165  original_row_end = &original.colnums[original.rowstart[row+1]];
166 
167  // find pointers before and after extra off-diagonals. if at top or
168  // bottom of matrix, then set these pointers such that no copying is
169  // necessary (see the @p{copy} commands)
170  const size_type *const
171  original_last_before_side_diagonals
172  = (row > extra_off_diagonals ?
173  Utilities::lower_bound (original_row_start,
174  original_row_end,
175  row-extra_off_diagonals) :
176  original_row_start);
177 
178  const size_type *const
179  original_first_after_side_diagonals
180  = (row < rows-extra_off_diagonals-1 ?
181  std::upper_bound (original_row_start,
182  original_row_end,
183  row+extra_off_diagonals) :
184  original_row_end);
185 
186  // find first free slot. the first slot in each row is the diagonal
187  // element
188  size_type *next_free_slot = &colnums[rowstart[row]] + 1;
189 
190  // copy elements before side-diagonals
191  next_free_slot = std::copy (original_row_start,
192  original_last_before_side_diagonals,
193  next_free_slot);
194 
195  // insert left and right 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;
202 
203  // copy rest
204  next_free_slot = std::copy (original_first_after_side_diagonals,
205  original_row_end,
206  next_free_slot);
207 
208  // this error may happen if the sum of previous elements per row and
209  // those of the new diagonals exceeds the maximum number of elements per
210  // row given to this constructor
211  Assert (next_free_slot <= &colnums[rowstart[row+1]],
212  ExcNotEnoughSpace (0,rowstart[row+1]-rowstart[row]));
213  };
214 }
215 
216 
217 
220 {
221  (void)s;
222  Assert (s.empty(),
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."));
226 
227  Assert (this->empty(),
228  ExcMessage("This operator can only be called if the current object is "
229  "empty."));
230 
231  return *this;
232 }
233 
234 
235 
236 void
238  const size_type n,
239  const unsigned int max_per_row)
240 {
241  // simply map this function to the other @p{reinit} function
242  const std::vector<unsigned int> row_lengths (m, max_per_row);
243  reinit (m, n, row_lengths);
244 }
245 
246 
247 
248 void
250  const size_type n,
251  const VectorSlice<const std::vector<unsigned int> > &row_lengths)
252 {
253  AssertDimension (row_lengths.size(), m);
254 
255  rows = m;
256  cols = n;
257 
258  // delete empty matrices
259  if ((m==0) || (n==0))
260  {
261  rowstart.reset();
262  colnums.reset();
263 
264  max_vec_len = max_dim = rows = cols = 0;
265  // if dimension is zero: ignore max_per_row
266  max_row_length = 0;
267  compressed = false;
268 
269  return;
270  }
271 
272  // first, if the matrix is quadratic, we will have to make sure that each
273  // row has at least one entry for the diagonal element. make this more
274  // obvious by having a variable which we can query
275  store_diagonal_first_in_row = (m == n);
276 
277  // find out how many entries we need in the @p{colnums} array. if this
278  // number is larger than @p{max_vec_len}, then we will need to reallocate
279  // memory
280  //
281  // note that the number of elements per row is bounded by the number of
282  // columns
283  //
284  std::size_t vec_len = 0;
285  for (size_type i=0; i<m; ++i)
286  vec_len += std::min(static_cast<size_type>(store_diagonal_first_in_row ?
287  std::max(row_lengths[i], 1U) :
288  row_lengths[i]),
289  n);
290 
291  // sometimes, no entries are requested in the matrix (this most often
292  // happens when blocks in a block matrix are simply zero). in that case,
293  // allocate exactly one element, to have a valid pointer to some memory
294  if (vec_len == 0)
295  {
296  vec_len = 1;
297  max_vec_len = vec_len;
298  colnums = std_cxx14::make_unique<size_type[]> (max_vec_len);
299  }
300 
301  max_row_length = (row_lengths.size() == 0 ?
302  0 :
303  std::min (static_cast<size_type>(*std::max_element(row_lengths.begin(),
304  row_lengths.end())),
305  n));
306 
307  if (store_diagonal_first_in_row && (max_row_length==0) && (m!=0))
308  max_row_length = 1;
309 
310  // allocate memory for the rowstart values, if necessary. even though we
311  // re-set the pointers again immediately after deleting their old content,
312  // set them to zero in between because the allocation might fail, in which
313  // case we get an exception and the destructor of this object will be called
314  // -- where we look at the non-nullness of the (now invalid) pointer again
315  // and try to delete the memory a second time.
316  if (rows > max_dim)
317  {
318  max_dim = rows;
319  rowstart = std_cxx14::make_unique<std::size_t[]> (max_dim+1);
320  }
321 
322  // allocate memory for the column numbers if necessary
323  if (vec_len > max_vec_len)
324  {
325  max_vec_len = vec_len;
326  colnums = std_cxx14::make_unique<size_type[]> (max_vec_len);
327  }
328 
329  // set the rowstart array
330  rowstart[0] = 0;
331  for (size_type i=1; i<=rows; ++i)
332  rowstart[i] = rowstart[i-1] +
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));
337  Assert ((rowstart[rows]==vec_len)
338  ||
339  ((vec_len == 1) && (rowstart[rows] == 0)),
340  ExcInternalError());
341 
342  // preset the column numbers by a value indicating it is not in use
343  std::fill_n (colnums.get(), vec_len, invalid_entry);
344 
345  // if diagonal elements are special: let the first entry in each row be the
346  // diagonal value
348  for (size_type i=0; i<rows; i++)
349  colnums[rowstart[i]] = i;
350 
351  compressed = false;
352 }
353 
354 
355 
356 void
358 {
359  // nothing to do if the object corresponds to an empty matrix
360  if ((rowstart==nullptr) && (colnums==nullptr))
361  {
362  compressed = true;
363  return;
364  }
365 
366  // do nothing if already compressed
367  if (compressed)
368  return;
369 
370  size_type next_free_entry = 0,
371  next_row_start = 0,
372  row_length = 0;
373 
374  // first find out how many non-zero elements there are, in order to allocate
375  // the right amount of memory
376  const std::size_t nonzero_elements
377  = std::count_if (&colnums[rowstart[0]],
378  &colnums[rowstart[rows]],
379  std::bind(std::not_equal_to<size_type>(), std::placeholders::_1, invalid_entry));
380  // now allocate the respective memory
381  std::unique_ptr<size_type[]> new_colnums (new size_type[nonzero_elements]);
382 
383 
384  // reserve temporary storage to store the entries of one row
385  std::vector<size_type> tmp_entries (max_row_length);
386 
387  // Traverse all rows
388  for (size_type line=0; line<rows; ++line)
389  {
390  // copy used entries, break if first unused entry is reached
391  row_length = 0;
392  for (size_type j=rowstart[line]; j<rowstart[line+1]; ++j,++row_length)
393  if (colnums[j] != invalid_entry)
394  tmp_entries[row_length] = colnums[j];
395  else
396  break;
397  // now @p{rowstart} is the number of entries in this line
398 
399  // Sort only beginning at the second entry, if optimized storage of
400  // diagonal entries is on.
401 
402  // if this line is empty or has only one entry, don't sort
403  if (row_length > 1)
404  std::sort ((store_diagonal_first_in_row)
405  ? tmp_entries.begin()+1
406  : tmp_entries.begin(),
407  tmp_entries.begin()+row_length);
408 
409  // insert column numbers into the new field
410  for (size_type j=0; j<row_length; ++j)
411  new_colnums[next_free_entry++] = tmp_entries[j];
412 
413  // note new start of this and the next row
414  rowstart[line] = next_row_start;
415  next_row_start = next_free_entry;
416 
417  // some internal checks: either the matrix is not quadratic, or if it
418  // is, then the first element of this row must be the diagonal element
419  // (i.e. with column index==line number)
420  // this test only makes sense if we have written to the index
421  // rowstart_line in new_colnums which is the case if row_length is not 0,
422  // so check this first
424  (row_length != 0 && new_colnums[rowstart[line]] == line),
425  ExcInternalError());
426  // assert that the first entry does not show up in the remaining ones
427  // and that the remaining ones are unique among themselves (this handles
428  // both cases, quadratic and rectangular matrices)
429  //
430  // the only exception here is if the row contains no entries at all
431  Assert ((rowstart[line] == next_row_start)
432  ||
433  (std::find (&new_colnums[rowstart[line]+1],
434  &new_colnums[next_row_start],
435  new_colnums[rowstart[line]]) ==
436  &new_colnums[next_row_start]),
437  ExcInternalError());
438  Assert ((rowstart[line] == next_row_start)
439  ||
440  (std::adjacent_find(&new_colnums[rowstart[line]+1],
441  &new_colnums[next_row_start]) ==
442  &new_colnums[next_row_start]),
443  ExcInternalError());
444  };
445 
446  // assert that we have used all allocated space, no more and no less
447  Assert (next_free_entry == nonzero_elements,
448  ExcInternalError());
449 
450  // set iterator-past-the-end
451  rowstart[rows] = next_row_start;
452 
453  // set colnums to the newly allocated array and delete previous content
454  // in the process
455  colnums = std::move(new_colnums);
456 
457  // store the size
458  max_vec_len = nonzero_elements;
459 
460  compressed = true;
461 }
462 
463 
464 
465 void
467 {
468  // first determine row lengths for each row. if the matrix is quadratic,
469  // then we might have to add an additional entry for the diagonal, if that
470  // is not yet present. as we have to call compress anyway later on, don't
471  // bother to check whether that diagonal entry is in a certain row or not
472  const bool do_diag_optimize = (sp.n_rows() == sp.n_cols());
473  std::vector<unsigned int> row_lengths (sp.n_rows());
474  for (size_type i=0; i<sp.n_rows(); ++i)
475  {
476  row_lengths[i] = sp.row_length(i);
477  if (do_diag_optimize && !sp.exists(i,i))
478  ++row_lengths[i];
479  }
480  reinit (sp.n_rows(), sp.n_cols(), row_lengths);
481 
482  // now enter all the elements into the matrix, if there are any. note that
483  // if the matrix is quadratic, then we already have the diagonal element
484  // preallocated
485  if (n_rows() != 0 && n_cols() != 0)
486  for (size_type row = 0; row<sp.n_rows(); ++row)
487  {
488  size_type *cols = &colnums[rowstart[row]] + (do_diag_optimize ? 1 : 0);
489  typename SparsityPattern::iterator col_num = sp.begin (row),
490  end_row = sp.end (row);
491 
492  for (; col_num != end_row; ++col_num)
493  {
494  const size_type col = col_num->column();
495  if ((col!=row) || !do_diag_optimize)
496  *cols++ = col;
497  }
498  }
499 
500  // do not need to compress the sparsity pattern since we already have
501  // allocated the right amount of data, and the SparsityPattern data is
502  // sorted, too.
503  compressed = true;
504 }
505 
506 
507 
508 // Use a special implementation for DynamicSparsityPattern where we can use
509 // the column_number method to gain faster access to the
510 // entries. DynamicSparsityPattern::iterator can show quadratic complexity in
511 // case many rows are empty and the begin() method needs to jump to the next
512 // free row. Otherwise, the code is exactly the same as above.
513 void
515 {
516  const bool do_diag_optimize = (dsp.n_rows() == dsp.n_cols());
517  std::vector<unsigned int> row_lengths (dsp.n_rows());
518  for (size_type i=0; i<dsp.n_rows(); ++i)
519  {
520  row_lengths[i] = dsp.row_length(i);
521  if (do_diag_optimize && !dsp.exists(i,i))
522  ++row_lengths[i];
523  }
524  reinit (dsp.n_rows(), dsp.n_cols(), row_lengths);
525 
526  if (n_rows() != 0 && n_cols() != 0)
527  for (size_type row = 0; row<dsp.n_rows(); ++row)
528  {
529  size_type *cols = &colnums[rowstart[row]] + (do_diag_optimize ? 1 : 0);
530  const unsigned int row_length = dsp.row_length(row);
531  for (unsigned int index=0; index < row_length; ++index)
532  {
533  const size_type col = dsp.column_number(row, index);
534  if ((col!=row) || !do_diag_optimize)
535  *cols++ = col;
536  }
537  }
538 
539  // do not need to compress the sparsity pattern since we already have
540  // allocated the right amount of data, and the SparsityPatternType data is sorted,
541  // too.
542  compressed = true;
543 }
544 
545 
546 
547 template <typename number>
549 {
550  // first init with the number of entries per row. if this matrix is square
551  // then we also have to allocate memory for the diagonal entry, unless we
552  // have already counted it
553  const bool matrix_is_square = (matrix.m() == matrix.n());
554 
555  std::vector<unsigned int> entries_per_row (matrix.m(), 0);
556  for (size_type row=0; row<matrix.m(); ++row)
557  {
558  for (size_type col=0; col<matrix.n(); ++col)
559  if (matrix(row,col) != 0)
560  ++entries_per_row[row];
561  if (matrix_is_square
562  &&
563  (matrix(row,row) == 0))
564  ++entries_per_row[row];
565  }
566 
567  reinit (matrix.m(), matrix.n(), entries_per_row);
568 
569  // now set entries. if we enter entries row by row, then we'll get
570  // quadratic complexity in the number of entries per row. this is
571  // not usually a problem (we don't usually create dense matrices),
572  // but there are cases where it matters -- so we may as well be
573  // gentler and hand over a whole row of entries at a time
574  std::vector<size_type> column_indices;
575  column_indices.reserve (entries_per_row.size() > 0
576  ?
577  *std::max_element (entries_per_row.begin(),
578  entries_per_row.end())
579  :
580  0);
581  for (size_type row=0; row<matrix.m(); ++row)
582  {
583  column_indices.resize(entries_per_row[row]);
584 
585  size_type current_index = 0;
586  for (size_type col=0; col<matrix.n(); ++col)
587  if (matrix(row,col) != 0)
588  {
589  column_indices[current_index] = col;
590  ++current_index;
591  }
592  else
593  // the (row,col) entry is zero; check if we need to add it
594  // anyway because it's the diagonal entry of a square
595  // matrix
596  if (matrix_is_square
597  &&
598  (col == row))
599  {
600  column_indices[current_index] = row;
601  ++current_index;
602  }
603 
604  // check that we really added the correct number of indices
605  Assert (current_index == entries_per_row[row], ExcInternalError());
606 
607  // now bulk add all of these entries
608  add_entries(row, column_indices.begin(), column_indices.end(), true);
609  }
610 
611  // finally compress
612  compress ();
613 }
614 
615 
616 void
618  const size_type n,
619  const std::vector<unsigned int> &row_lengths)
620 {
621  reinit(m, n, make_slice(row_lengths));
622 }
623 
624 
625 
626 
627 bool
629 {
630  // let's try to be on the safe side of life by using multiple possibilities
631  // in the check for emptiness... (sorry for this kludge -- emptying matrices
632  // and freeing memory was not present in the original implementation and I
633  // don't know at how many places I missed something in adding it, so I try
634  // to be cautious. wb)
635  if ((rowstart==nullptr) || (rows==0) || (cols==0))
636  {
637  Assert (rowstart==nullptr, ExcInternalError());
638  Assert (rows==0, ExcInternalError());
639  Assert (cols==0, ExcInternalError());
640  Assert (colnums==nullptr, ExcInternalError());
642 
643  return true;
644  };
645  return false;
646 }
647 
648 
649 
652 {
653  // if compress() has not yet been called, we can get the maximum number of
654  // elements per row using the stored value
655  if (!compressed)
656  return max_row_length;
657 
658  // if compress() was called, we use a better algorithm which gives us a
659  // sharp bound
660  size_type m = 0;
661  for (size_type i=1; i<=rows; ++i)
662  m = std::max (m, static_cast<size_type>(rowstart[i]-rowstart[i-1]));
663 
664  return m;
665 }
666 
667 
668 
671  const size_type j) const
672 {
673  Assert ((rowstart!=nullptr) && (colnums!=nullptr), ExcEmptyObject());
674  Assert (i<rows, ExcIndexRange(i,0,rows));
675  Assert (j<cols, ExcIndexRange(j,0,cols));
677 
678  // let's see whether there is something in this line
679  if (rowstart[i] == rowstart[i+1])
680  return invalid_entry;
681 
682  // If special storage of diagonals was requested, we can get the diagonal
683  // element faster by this query.
684  if (store_diagonal_first_in_row && (i==j))
685  return rowstart[i];
686 
687  // all other entries are sorted, so we can use a binary search algorithm
688  //
689  // note that the entries are only sorted upon compression, so this would
690  // fail for non-compressed sparsity patterns; however, that is why the
691  // Assertion is at the top of this function, so it may not be called for
692  // noncompressed structures.
693  const size_type *sorted_region_start = (store_diagonal_first_in_row ?
694  &colnums[rowstart[i]+1] :
695  &colnums[rowstart[i]]);
696  const size_type *const p
697  = Utilities::lower_bound<const size_type *> (sorted_region_start,
698  &colnums[rowstart[i+1]],
699  j);
700  if ((p != &colnums[rowstart[i+1]]) && (*p == j))
701  return (p - colnums.get());
702  else
703  return invalid_entry;
704 }
705 
706 
707 
708 void
710  const size_type j)
711 {
712  Assert ((rowstart!=nullptr) && (colnums!=nullptr), ExcEmptyObject());
713  Assert (i<rows, ExcIndexRange(i,0,rows));
714  Assert (j<cols, ExcIndexRange(j,0,cols));
716 
717  for (std::size_t k=rowstart[i]; k<rowstart[i+1]; k++)
718  {
719  // entry already exists
720  if (colnums[k] == j) return;
721  // empty entry found, put new entry here
722  if (colnums[k] == invalid_entry)
723  {
724  colnums[k] = j;
725  return;
726  };
727  };
728 
729  // if we came thus far, something went wrong: there was not enough space in
730  // this line
731  Assert (false, ExcNotEnoughSpace(i, rowstart[i+1]-rowstart[i]));
732 }
733 
734 
735 
736 template <typename ForwardIterator>
737 void
739  ForwardIterator begin,
740  ForwardIterator end,
741  const bool indices_are_sorted)
742 {
743  if (indices_are_sorted == true)
744  {
745  if (begin != end)
746  {
747  ForwardIterator it = begin;
748  bool has_larger_entries = false;
749  // skip diagonal
750  std::size_t k=rowstart[row]+store_diagonal_first_in_row;
751  for ( ; k<rowstart[row+1]; k++)
752  if (colnums[k] == invalid_entry)
753  break;
754  else if (colnums[k] >= *it)
755  {
756  has_larger_entries = true;
757  break;
758  }
759  if (has_larger_entries == false)
760  for ( ; it != end; ++it)
761  {
762  if (store_diagonal_first_in_row && *it == row)
763  continue;
764  Assert (k <= rowstart[row+1],
765  ExcNotEnoughSpace(row, rowstart[row+1]-rowstart[row]));
766  colnums[k++] = *it;
767  }
768  else
769  // cannot just append the new range at the end, forward to the
770  // other function
771  for (ForwardIterator p = begin; p != end; ++p)
772  add (row, *p);
773  }
774  }
775  else
776  {
777  // forward to the other function.
778  for (ForwardIterator it = begin; it != end; ++it)
779  add (row, *it);
780  }
781 }
782 
783 
784 bool
786 {
787  Assert ((rowstart!=nullptr) && (colnums!=nullptr), ExcEmptyObject());
788  Assert (i<rows, ExcIndexRange(i,0,rows));
789  Assert (j<cols, ExcIndexRange(j,0,cols));
790 
791  for (size_type k=rowstart[i]; k<rowstart[i+1]; k++)
792  {
793  // entry already exists
794  if (colnums[k] == j) return true;
795  }
796  return false;
797 }
798 
799 
800 
803 {
804  Assert ((rowstart!=nullptr) && (colnums!=nullptr), ExcEmptyObject());
805  Assert (i<rows, ExcIndexRange(i,0,rows));
806  Assert (j<cols, ExcIndexRange(j,0,cols));
807 
808  for (size_type k=rowstart[i]; k<rowstart[i+1]; k++)
809  {
810  // entry exists
811  if (colnums[k] == j) return k-rowstart[i];
812  }
814 }
815 
816 
817 
818 std::pair<SparsityPattern::size_type, SparsityPattern::size_type>
819 SparsityPattern::matrix_position (const std::size_t global_index) const
820 {
821  Assert (compressed == true, ExcNotCompressed());
822  Assert (global_index < n_nonzero_elements(),
823  ExcIndexRange (global_index, 0, n_nonzero_elements()));
824 
825  // first find the row in which the entry is located. for this note that the
826  // rowstart array indexes the global indices at which each row starts. since
827  // it is sorted, and since there is an element for the one-past-last row, we
828  // can simply use a bisection search on it
829  const size_type row
830  = (std::upper_bound (rowstart.get(), rowstart.get() + rows, global_index)
831  - rowstart.get() - 1);
832 
833  // now, the column index is simple since that is what the colnums array
834  // stores:
835  const size_type col = colnums[global_index];
836 
837  // so return the respective pair
838  return std::make_pair (row,col);
839 }
840 
841 
842 
843 void
845 {
846  Assert ((rowstart!=nullptr) && (colnums!=nullptr), ExcEmptyObject());
848  // Note that we only require a quadratic matrix here, no special treatment
849  // of diagonals
851 
852  // loop over all elements presently in the sparsity pattern and add the
853  // transpose element. note:
854  //
855  // 1. that the sparsity pattern changes which we work on, but not the
856  // present row
857  //
858  // 2. that the @p{add} function can be called on elements that already exist
859  // without any harm
860  for (size_type row=0; row<rows; ++row)
861  for (size_type k=rowstart[row]; k<rowstart[row+1]; k++)
862  {
863  // check whether we are at the end of the entries of this row. if so,
864  // go to next row
865  if (colnums[k] == invalid_entry)
866  break;
867 
868  // otherwise add the transpose entry if this is not the diagonal (that
869  // would not harm, only take time to check up)
870  if (colnums[k] != row)
871  add (colnums[k], row);
872  };
873 }
874 
875 
876 
877 void
878 SparsityPattern::print (std::ostream &out) const
879 {
880  Assert ((rowstart!=nullptr) && (colnums!=nullptr), ExcEmptyObject());
881 
882  AssertThrow (out, ExcIO());
883 
884  for (size_type i=0; i<rows; ++i)
885  {
886  out << '[' << i;
887  for (size_type j=rowstart[i]; j<rowstart[i+1]; ++j)
888  if (colnums[j] != invalid_entry)
889  out << ',' << colnums[j];
890  out << ']' << std::endl;
891  }
892 
893  AssertThrow (out, ExcIO());
894 }
895 
896 
897 
898 void
899 SparsityPattern::print_gnuplot (std::ostream &out) const
900 {
901  Assert ((rowstart!=nullptr) && (colnums!=nullptr), ExcEmptyObject());
902 
903  AssertThrow (out, ExcIO());
904 
905  for (size_type i=0; i<rows; ++i)
906  for (size_type j=rowstart[i]; j<rowstart[i+1]; ++j)
907  if (colnums[j] != invalid_entry)
908  // while matrix entries are usually written (i,j), with i vertical and
909  // j horizontal, gnuplot output is x-y, that is we have to exchange
910  // the order of output
911  out << colnums[j] << " " << -static_cast<signed int>(i) << std::endl;
912 
913  AssertThrow (out, ExcIO());
914 }
915 
916 void
917 SparsityPattern::print_svg (std::ostream &out) const
918 {
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"
924  " <![CDATA[\n"
925  " rect.pixel {\n"
926  " fill: #ff0000;\n"
927  " }\n"
928  " ]]>\n"
929  " </style>\n\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";
933 
935  it = this->begin(),
936  end = this->end();
937  for (; it!=end; ++it)
938  {
939  out << " <rect class=\"pixel\" x=\"" << it->column()+1
940  << "\" y=\"" << it->row()+1
941  << "\" width=\".9\" height=\".9\"/>\n";
942  }
943  out << "</svg>" << std::endl;
944 
945 }
946 
947 
948 
949 
952 {
953  Assert ((rowstart!=nullptr) && (colnums!=nullptr), ExcEmptyObject());
954  size_type b=0;
955  for (size_type i=0; i<rows; ++i)
956  for (size_type j=rowstart[i]; j<rowstart[i+1]; ++j)
957  if (colnums[j] != invalid_entry)
958  {
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]));
961  }
962  else
963  // leave if at the end of the entries of this line
964  break;
965  return b;
966 }
967 
968 
969 void
970 SparsityPattern::block_write (std::ostream &out) const
971 {
972  AssertThrow (out, ExcIO());
973 
974  // first the simple objects, bracketed in [...]
975  out << '[' << max_dim << ' '
976  << rows << ' '
977  << cols << ' '
978  << max_vec_len << ' '
979  << max_row_length << ' '
980  << compressed << ' '
981  << store_diagonal_first_in_row << "][";
982  // then write out real data
983  out.write (reinterpret_cast<const char *>(rowstart.get()),
984  reinterpret_cast<const char *>(rowstart.get() + max_dim + 1)
985  - reinterpret_cast<const char *>(rowstart.get()));
986  out << "][";
987  out.write (reinterpret_cast<const char *>(colnums.get()),
988  reinterpret_cast<const char *>(colnums.get() + max_vec_len)
989  - reinterpret_cast<const char *>(colnums.get()));
990  out << ']';
991 
992  AssertThrow (out, ExcIO());
993 }
994 
995 
996 
997 void
998 SparsityPattern::block_read (std::istream &in)
999 {
1000  AssertThrow (in, ExcIO());
1001 
1002  char c;
1003 
1004  // first read in simple data
1005  in >> c;
1006  AssertThrow (c == '[', ExcIO());
1007  in >> max_dim
1008  >> rows
1009  >> cols
1010  >> max_vec_len
1011  >> max_row_length
1012  >> compressed
1014 
1015  in >> c;
1016  AssertThrow (c == ']', ExcIO());
1017  in >> c;
1018  AssertThrow (c == '[', ExcIO());
1019 
1020  // reallocate space
1021  rowstart = std_cxx14::make_unique<std::size_t[]> (max_dim+1);
1022  colnums = std_cxx14::make_unique<size_type[]> (max_vec_len);
1023 
1024  // then read data
1025  in.read (reinterpret_cast<char *>(rowstart.get()),
1026  reinterpret_cast<char *>(rowstart.get() + max_dim + 1)
1027  - reinterpret_cast<char *>(rowstart.get()));
1028  in >> c;
1029  AssertThrow (c == ']', ExcIO());
1030  in >> c;
1031  AssertThrow (c == '[', ExcIO());
1032  in.read (reinterpret_cast<char *>(colnums.get()),
1033  reinterpret_cast<char *>(colnums.get() + max_vec_len)
1034  - reinterpret_cast<char *>(colnums.get()));
1035  in >> c;
1036  AssertThrow (c == ']', ExcIO());
1037 }
1038 
1039 
1040 
1041 std::size_t
1043 {
1044  return (max_dim * sizeof(size_type) +
1045  sizeof(*this) +
1046  max_vec_len * sizeof(size_type));
1047 }
1048 
1049 
1050 
1051 // explicit instantiations
1052 template void SparsityPattern::copy_from<float> (const FullMatrix<float> &);
1053 template void SparsityPattern::copy_from<double> (const FullMatrix<double> &);
1054 
1055 template void SparsityPattern::add_entries<const SparsityPattern::size_type *> (const size_type,
1056  const size_type *,
1057  const size_type *,
1058  const bool);
1059 #ifndef DEAL_II_VECTOR_ITERATOR_IS_POINTER
1060 template void SparsityPattern::add_entries<std::vector<SparsityPattern::size_type>::const_iterator>
1061 (const size_type,
1062  std::vector<size_type>::const_iterator,
1063  std::vector<size_type>::const_iterator,
1064  const bool);
1065 #endif
1066 template void SparsityPattern::add_entries<std::vector<SparsityPattern::size_type>::iterator>
1067 (const size_type,
1068  std::vector<size_type>::iterator,
1069  std::vector<size_type>::iterator,
1070  const bool);
1071 
1072 DEAL_II_NAMESPACE_CLOSE
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:891
void block_write(std::ostream &out) const
const types::global_dof_index invalid_size_type
Definition: types.h:182
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false)
static ::ExceptionBase & ExcMatrixIsCompressed()
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
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
iterator end() const
bool exists(const size_type i, const size_type j) const
bool empty() const
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
size_type n_cols() const
bool is_compressed() const
std::size_t max_vec_len
size_type n_rows() 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)
Definition: exceptions.h:1142
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
std::unique_ptr< std::size_t[]> rowstart
static ::ExceptionBase & ExcEmptyObject()
void block_read(std::istream &in)
iterator begin() const
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()