Reference documentation for deal.II version 8.5.1
chunk_sparsity_pattern.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2016 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/lac/chunk_sparsity_pattern.h>
18 #include <deal.II/lac/dynamic_sparsity_pattern.h>
19 #include <deal.II/lac/full_matrix.h>
20 
21 
22 DEAL_II_NAMESPACE_OPEN
23 
24 
26 {
27  reinit (0,0,0,1);
28 }
29 
30 
31 
33  :
34  Subscriptor(),
35  chunk_size (s.chunk_size),
36  sparsity_pattern(s.sparsity_pattern)
37 {
38  Assert (s.rows==0 && s.cols==0,
39  ExcMessage("This constructor can only be called if the provided argument "
40  "is the sparsity pattern for an empty matrix. This constructor can "
41  "not be used to copy-construct a non-empty sparsity pattern."));
42 
43  reinit (0,0,0,0);
44 }
45 
46 
47 
49  const size_type n,
50  const size_type max_per_row,
51  const size_type chunk_size)
52 {
54 
55  reinit (m,n,max_per_row, chunk_size);
56 }
57 
58 
59 
61  const size_type m,
62  const size_type n,
63  const std::vector<size_type> &row_lengths,
64  const size_type chunk_size)
65 {
67 
68  reinit (m, n, row_lengths, chunk_size);
69 }
70 
71 
72 
74  const size_type max_per_row,
75  const size_type chunk_size)
76 {
77  reinit (n, n, max_per_row, chunk_size);
78 }
79 
80 
81 
83  const size_type m,
84  const std::vector<size_type > &row_lengths,
85  const size_type chunk_size)
86 {
88 
89  reinit (m, m, row_lengths, chunk_size);
90 }
91 
92 
93 
95 {}
96 
97 
98 
101 {
102  Assert (s.rows==0 && s.cols==0,
103  ExcMessage("This operator can only be called if the provided argument "
104  "is the sparsity pattern for an empty matrix. This operator can "
105  "not be used to copy a non-empty sparsity pattern."));
106 
107  Assert (rows==0 && cols==0,
108  ExcMessage("This operator can only be called if the current object is "
109  "empty."));
110 
111  // perform the checks in the underlying object as well
113 
114  return *this;
115 }
116 
117 
118 
119 void
121  const size_type n,
122  const size_type max_per_row,
123  const size_type chunk_size)
124 {
126 
127  // simply map this function to the other @p{reinit} function
128  const std::vector<size_type> row_lengths (m, max_per_row);
129  reinit (m, n, row_lengths, chunk_size);
130 }
131 
132 
133 
134 void
136  const size_type m,
137  const size_type n,
138  const VectorSlice<const std::vector<size_type> > &row_lengths,
139  const size_type chunk_size)
140 {
141  Assert (row_lengths.size() == m, ExcInvalidNumber (m));
143 
144  rows = m;
145  cols = n;
146 
147  this->chunk_size = chunk_size;
148 
149  // pass down to the necessary information to the underlying object. we need
150  // to calculate how many chunks we need: we need to round up (m/chunk_size)
151  // and (n/chunk_size). rounding up in integer arithmetic equals
152  // ((m+chunk_size-1)/chunk_size):
153  const size_type m_chunks = (m+chunk_size-1) / chunk_size,
154  n_chunks = (n+chunk_size-1) / chunk_size;
155 
156  // compute the maximum number of chunks in each row. the passed array
157  // denotes the number of entries in each row of the big matrix -- in the
158  // worst case, these are all in independent chunks, so we have to calculate
159  // it as follows (as an example: let chunk_size==2, row_lengths={2,2,...},
160  // and entries in row zero at columns {0,2} and for row one at {4,6} -->
161  // we'll need 4 chunks for the first chunk row!) :
162  std::vector<unsigned int> chunk_row_lengths (m_chunks, 0);
163  for (size_type i=0; i<m; ++i)
164  chunk_row_lengths[i/chunk_size] += row_lengths[i];
165 
166  // for the case that the reduced sparsity pattern optimizes the diagonal but
167  // the actual sparsity pattern does not, need to take one more entry in the
168  // row to fit the user-required entry
169  if (m != n && m_chunks == n_chunks)
170  for (unsigned int i=0; i<m_chunks; ++i)
171  ++chunk_row_lengths[i];
172 
173  sparsity_pattern.reinit (m_chunks,
174  n_chunks,
175  chunk_row_lengths);
176 }
177 
178 
179 
180 void
182 {
184 }
185 
186 
187 
188 template <typename SparsityPatternType>
189 void
190 ChunkSparsityPattern::copy_from (const SparsityPatternType &dsp,
191  const size_type chunk_size)
192 {
194  this->chunk_size = chunk_size;
195  rows = dsp.n_rows();
196  cols = dsp.n_cols();
197 
198  // simple case: just use the given sparsity pattern
199  if (chunk_size == 1)
200  {
202  return;
203  }
204 
205  // create a temporary compressed sparsity pattern that collects all entries
206  // from the input sparsity pattern and then initialize the underlying small
207  // sparsity pattern
208  const size_type m_chunks = (dsp.n_rows()+chunk_size-1) / chunk_size,
209  n_chunks = (dsp.n_cols()+chunk_size-1) / chunk_size;
210  DynamicSparsityPattern temporary_sp(m_chunks, n_chunks);
211 
212  for (size_type row = 0; row<dsp.n_rows(); ++row)
213  {
214  const size_type reduced_row = row/chunk_size;
215 
216  // TODO: This could be made more efficient if we cached the
217  // previous column and only called add() if the previous and the
218  // current column lead to different chunk columns
219  for (typename SparsityPatternType::iterator col_num = dsp.begin(row);
220  col_num != dsp.end(row); ++col_num)
221  temporary_sp.add (reduced_row, col_num->column()/chunk_size);
222  }
223 
224  sparsity_pattern.copy_from (temporary_sp);
225 }
226 
227 
228 
229 
230 template <typename number>
232  const size_type chunk_size)
233 {
235 
236  // count number of entries per row, then initialize the underlying sparsity
237  // pattern. remember to also allocate space for the diagonal entry (if that
238  // hasn't happened yet) if m==n since we always allocate that for diagonal
239  // matrices
240  std::vector<size_type> entries_per_row (matrix.m(), 0);
241  for (size_type row=0; row<matrix.m(); ++row)
242  {
243  for (size_type col=0; col<matrix.n(); ++col)
244  if (matrix(row,col) != 0)
245  ++entries_per_row[row];
246 
247  if ((matrix.m() == matrix.n())
248  &&
249  (matrix(row,row) == 0))
250  ++entries_per_row[row];
251  }
252 
253  reinit (matrix.m(), matrix.n(),
254  entries_per_row,
255  chunk_size);
256 
257  // then actually fill it
258  for (size_type row=0; row<matrix.m(); ++row)
259  for (size_type col=0; col<matrix.n(); ++col)
260  if (matrix(row,col) != 0)
261  add (row,col);
262 
263  // finally compress
264  compress ();
265 }
266 
267 
268 
269 void
271  const size_type m,
272  const size_type n,
273  const std::vector<size_type> &row_lengths,
274  const size_type chunk_size)
275 {
277 
278  reinit(m, n, make_slice(row_lengths), chunk_size);
279 }
280 
281 
282 
283 namespace internal
284 {
285  namespace
286  {
287  template <typename SparsityPatternType>
288  void copy_sparsity (const SparsityPatternType &src,
289  SparsityPattern &dst)
290  {
291  dst.copy_from(src);
292  }
293 
294  void copy_sparsity (const SparsityPattern &src,
295  SparsityPattern &dst)
296  {
297  dst = src;
298  }
299  }
300 }
301 
302 
303 
304 template <typename Sparsity>
305 void
307 (const unsigned int m,
308  const unsigned int n,
309  const Sparsity &sparsity_pattern_for_chunks,
310  const unsigned int chunk_size_in,
311  const bool)
312 {
313  Assert (m > (sparsity_pattern_for_chunks.n_rows()-1) * chunk_size_in &&
314  m <= sparsity_pattern_for_chunks.n_rows() * chunk_size_in,
315  ExcMessage("Number of rows m is not compatible with chunk size "
316  "and number of rows in sparsity pattern for the chunks."));
317  Assert (n > (sparsity_pattern_for_chunks.n_cols()-1) * chunk_size_in &&
318  n <= sparsity_pattern_for_chunks.n_cols() * chunk_size_in,
319  ExcMessage("Number of columns m is not compatible with chunk size "
320  "and number of columns in sparsity pattern for the chunks."));
321 
322  internal::copy_sparsity(sparsity_pattern_for_chunks, sparsity_pattern);
323  chunk_size = chunk_size_in;
324  rows = m;
325  cols = n;
326 }
327 
328 
329 
330 bool
332 {
333  return sparsity_pattern.empty();
334 }
335 
336 
337 
340 {
342 }
343 
344 
345 
346 void
348  const size_type j)
349 {
352 
354 }
355 
356 
357 bool
359  const size_type j) const
360 {
361  Assert (i<rows, ExcIndexRange(i,0,rows));
362  Assert (j<cols, ExcIndexRange(j,0,cols));
363 
365  j/chunk_size);
366 }
367 
368 
369 
370 void
372 {
373  // matrix must be square. note that the for some matrix sizes, the current
374  // sparsity pattern may not be square even if the underlying sparsity
375  // pattern is (e.g. a 10x11 matrix with chunk_size 4)
377 
379 }
380 
381 
382 
385 {
386  Assert (i<rows, ExcIndexRange(i,0,rows));
387 
388  // find out if we did padding and if this row is affected by it
389  if (n_cols() % chunk_size == 0)
391  else
392  // if columns don't align, then just iterate over all chunks and see
393  // what this leads to
394  {
397  unsigned int n = 0;
398  for ( ; p != end; ++p)
399  if (p->column() != sparsity_pattern.n_cols() - 1)
400  n += chunk_size;
401  else
402  n += (n_cols() % chunk_size);
403  return n;
404  }
405 }
406 
407 
408 
411 {
412  if ((n_rows() % chunk_size == 0)
413  &&
414  (n_cols() % chunk_size == 0))
416  chunk_size *
417  chunk_size);
418  else
419  // some of the chunks reach beyond the extent of this matrix. this
420  // requires a somewhat more complicated computations, in particular if the
421  // columns don't align
422  {
423  if ((n_rows() % chunk_size != 0)
424  &&
425  (n_cols() % chunk_size == 0))
426  {
427  // columns align with chunks, but not rows
429  chunk_size *
430  chunk_size;
431  n -= (sparsity_pattern.n_rows() * chunk_size - n_rows()) *
433  chunk_size;
434  return n;
435  }
436 
437  else
438  {
439  // if columns don't align, then just iterate over all chunks and see
440  // what this leads to. follow the advice in the documentation of the
441  // sparsity pattern iterators to do the loop over individual rows,
442  // rather than all elements
443  size_type n = 0;
444 
445  for (size_type row = 0; row < sparsity_pattern.n_rows(); ++row)
446  {
448  for (; p!=sparsity_pattern.end(row); ++p)
449  if ((row != sparsity_pattern.n_rows() - 1)
450  &&
451  (p->column() != sparsity_pattern.n_cols() - 1))
452  n += chunk_size * chunk_size;
453  else if ((row == sparsity_pattern.n_rows() - 1)
454  &&
455  (p->column() != sparsity_pattern.n_cols() - 1))
456  // last chunk row, but not last chunk column. only a smaller
457  // number (n_rows % chunk_size) of rows actually exist
458  n += (n_rows() % chunk_size) * chunk_size;
459  else if ((row != sparsity_pattern.n_rows() - 1)
460  &&
461  (p->column() == sparsity_pattern.n_cols() - 1))
462  // last chunk column, but not row
463  n += (n_cols() % chunk_size) * chunk_size;
464  else
465  // bottom right chunk
466  n += (n_cols() % chunk_size) *
467  (n_rows() % chunk_size);
468  }
469 
470  return n;
471  }
472  }
473 }
474 
475 
476 
477 void
478 ChunkSparsityPattern::print (std::ostream &out) const
479 {
481  ExcEmptyObject());
482 
483  AssertThrow (out, ExcIO());
484 
485  for (size_type i=0; i<sparsity_pattern.rows; ++i)
486  for (size_type d=0;
487  (d<chunk_size) && (i*chunk_size + d < n_rows());
488  ++d)
489  {
490  out << '[' << i *chunk_size+d;
492  j<sparsity_pattern.rowstart[i+1]; ++j)
494  for (size_type e=0;
495  ((e<chunk_size) &&
497  ++e)
498  out << ',' << sparsity_pattern.colnums[j]*chunk_size+e;
499  out << ']' << std::endl;
500  }
501 
502  AssertThrow (out, ExcIO());
503 }
504 
505 
506 
507 void
508 ChunkSparsityPattern::print_gnuplot (std::ostream &out) const
509 {
512 
513  AssertThrow (out, ExcIO());
514 
515  // for each entry in the underlying sparsity pattern, repeat everything
516  // chunk_size x chunk_size times
517  for (size_type i=0; i<sparsity_pattern.rows; ++i)
519  j<sparsity_pattern.rowstart[i+1]; ++j)
521  for (size_type d=0;
522  ((d<chunk_size) &&
524  ++d)
525  for (size_type e=0;
526  (e<chunk_size) && (i*chunk_size + e < n_rows());
527  ++e)
528  // while matrix entries are usually written (i,j), with i vertical
529  // and j horizontal, gnuplot output is x-y, that is we have to
530  // exchange the order of output
531  out << sparsity_pattern.colnums[j]*chunk_size+d << " "
532  << -static_cast<signed int>(i*chunk_size+e)
533  << std::endl;
534 
535  AssertThrow (out, ExcIO());
536 }
537 
538 
539 
542 {
543  // calculate the bandwidth from that of the underlying sparsity
544  // pattern. note that even if the bandwidth of that is zero, then the
545  // bandwidth of the chunky pattern is chunk_size-1, if it is 1 then the
546  // chunky pattern has chunk_size+(chunk_size-1), etc
547  //
548  // we'll cut it off at max(n(),m())
549  return std::min (sparsity_pattern.bandwidth()*chunk_size
550  + (chunk_size-1),
551  std::max(n_rows(), n_cols()));
552 }
553 
554 
555 
556 bool
558 {
559  if (chunk_size == 1)
561  else
562  return false;
563 }
564 
565 
566 
567 void
568 ChunkSparsityPattern::block_write (std::ostream &out) const
569 {
570  AssertThrow (out, ExcIO());
571 
572  // first the simple objects, bracketed in [...]
573  out << '['
574  << rows << ' '
575  << cols << ' '
576  << chunk_size << ' '
577  << "][";
578  // then the underlying sparsity pattern
580  out << ']';
581 
582  AssertThrow (out, ExcIO());
583 }
584 
585 
586 
587 void
589 {
590  AssertThrow (in, ExcIO());
591 
592  char c;
593 
594  // first read in simple data
595  in >> c;
596  AssertThrow (c == '[', ExcIO());
597  in >> rows
598  >> cols
599  >> chunk_size;
600 
601  in >> c;
602  AssertThrow (c == ']', ExcIO());
603  in >> c;
604  AssertThrow (c == '[', ExcIO());
605 
606  // then read the underlying sparsity pattern
608 
609  in >> c;
610  AssertThrow (c == ']', ExcIO());
611 }
612 
613 
614 
615 std::size_t
617 {
618  return (sizeof(*this) +
620 }
621 
622 
623 
624 // explicit instantiations
625 template
626 void ChunkSparsityPattern::copy_from<DynamicSparsityPattern> (const DynamicSparsityPattern &,
627  const size_type);
628 template
629 void ChunkSparsityPattern::create_from<SparsityPattern>
630 (const unsigned int,
631  const unsigned int,
632  const SparsityPattern &,
633  const unsigned int,
634  const bool);
635 template
636 void ChunkSparsityPattern::create_from<DynamicSparsityPattern>
637 (const unsigned int,
638  const unsigned int,
639  const DynamicSparsityPattern &,
640  const unsigned int,
641  const bool);
642 template
643 void ChunkSparsityPattern::copy_from<float> (const FullMatrix<float> &,
644  const size_type);
645 template
646 void ChunkSparsityPattern::copy_from<double> (const FullMatrix<double> &,
647  const size_type);
648 
649 DEAL_II_NAMESPACE_CLOSE
void block_write(std::ostream &out) const
void block_write(std::ostream &out) const
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end, const size_type chunk_size)
std::size_t memory_consumption() const
static ::ExceptionBase & ExcEmptyObject()
std::size_t n_nonzero_elements() const
static ::ExceptionBase & ExcInvalidNumber(int arg1)
static ::ExceptionBase & ExcInvalidIndex(int arg1, int arg2)
static ::ExceptionBase & ExcIO()
static const size_type invalid_entry
void add(const size_type i, const size_type j)
iterator end() const
std::size_t memory_consumption() const
void block_read(std::istream &in)
iterator end() const
bool exists(const size_type i, const size_type j) const
bool exists(const size_type i, const size_type j) const
bool empty() const
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
void add(const size_type i, const size_type j)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
size_type n_cols() const
SparsityPattern sparsity_pattern
size_type n_rows() const
size_type bandwidth() const
size_type n_cols() const
static ::ExceptionBase & ExcMessage(std::string arg1)
void add(const size_type i, const size_type j)
#define Assert(cond, exc)
Definition: exceptions.h:313
size_type row_length(const size_type row) const
size_type n_rows() const
size_type max_entries_per_row() const
size_type n_nonzero_elements() const
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
ChunkSparsityPattern & operator=(const ChunkSparsityPattern &)
void reinit(const size_type m, const size_type n, const unsigned int max_per_row)
std::size_t * rowstart
types::global_dof_index size_type
static ::ExceptionBase & ExcNotQuadratic()
bool stores_only_added_elements() const
void block_read(std::istream &in)
iterator begin() const
void print(std::ostream &out) const
void create_from(const unsigned int m, const unsigned int n, const Sparsity &sparsity_pattern_for_chunks, const unsigned int chunk_size, const bool optimize_diagonal=true)
void print_gnuplot(std::ostream &out) const
void reinit(const size_type m, const size_type n, const size_type max_per_row, const size_type chunk_size)
unsigned int row_length(const size_type row) const
size_type max_entries_per_row() const