Reference documentation for deal.II version 9.0.0
block_sparsity_pattern.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2017 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 #ifndef dealii_block_sparsity_pattern_h
17 #define dealii_block_sparsity_pattern_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/table.h>
23 #include <deal.II/base/subscriptor.h>
24 #include <deal.II/base/smartpointer.h>
25 #include <deal.II/lac/sparsity_pattern.h>
26 #include <deal.II/lac/trilinos_sparsity_pattern.h>
27 #include <deal.II/lac/dynamic_sparsity_pattern.h>
28 #include <deal.II/lac/block_indices.h>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
33 template <typename number> class BlockSparseMatrix;
35 
73 template <typename SparsityPatternType>
75 {
76 public:
81 
91 
98 
105  const size_type n_block_columns);
106 
114 
119 
133  void reinit (const size_type n_block_rows,
134  const size_type n_block_columns);
135 
142 
149  void collect_sizes ();
150 
154  SparsityPatternType &
155  block (const size_type row,
156  const size_type column);
157 
158 
163  const SparsityPatternType &
164  block (const size_type row,
165  const size_type column) const;
166 
171  const BlockIndices &
172  get_row_indices () const;
173 
178  const BlockIndices &
179  get_column_indices () const;
180 
185  void compress ();
186 
190  size_type n_block_rows () const;
191 
195  size_type n_block_cols () const;
196 
203  bool empty () const;
204 
211 
221  void add (const size_type i, const size_type j);
222 
233  template <typename ForwardIterator>
234  void add_entries (const size_type row,
235  ForwardIterator begin,
236  ForwardIterator end,
237  const bool indices_are_sorted = false);
238 
243  size_type n_rows () const;
244 
250  size_type n_cols () const;
251 
255  bool exists (const size_type i, const size_type j) const;
256 
261  unsigned int row_length (const size_type row) const;
262 
274  size_type n_nonzero_elements () const;
275 
281  void print (std::ostream &out) const;
282 
290  void print_gnuplot (std::ostream &out) const;
291 
301  int, int, int, int,
302  << "The blocks [" << arg1 << ',' << arg2 << "] and ["
303  << arg3 << ',' << arg4 << "] have differing row numbers.");
308  int, int, int, int,
309  << "The blocks [" << arg1 << ',' << arg2 << "] and ["
310  << arg3 << ',' << arg4 << "] have differing column numbers.");
312 
313 protected:
314 
319 
324 
329 
335 
341 
342 private:
347  std::vector<size_type > counter_within_block;
348 
353  std::vector<std::vector<size_type > > block_column_indices;
354 
359  template <typename number>
360  friend class BlockSparseMatrix;
361 };
362 
363 
364 
376 class BlockSparsityPattern : public BlockSparsityPatternBase<SparsityPattern>
377 {
378 public:
379 
385  BlockSparsityPattern () = default;
386 
393  const size_type n_columns);
394 
398  void reinit (const size_type n_block_rows,
399  const size_type n_block_columns);
400 
413  void reinit (const BlockIndices &row_indices,
414  const BlockIndices &col_indices,
415  const std::vector<std::vector<unsigned int> > &row_lengths);
416 
417 
422  bool is_compressed () const;
423 
428  std::size_t memory_consumption () const;
429 
435  void copy_from (const BlockDynamicSparsityPattern &dsp);
436 };
437 
438 
439 
482 class BlockDynamicSparsityPattern : public BlockSparsityPatternBase<DynamicSparsityPattern>
483 {
484 public:
485 
491  BlockDynamicSparsityPattern () = default;
492 
499  const size_type n_columns);
500 
507  BlockDynamicSparsityPattern (const std::vector<size_type> &row_block_sizes,
508  const std::vector<size_type> &col_block_sizes);
509 
518  BlockDynamicSparsityPattern (const std::vector<IndexSet> &partitioning);
519 
525  const BlockIndices &col_indices);
526 
527 
536  void reinit (const std::vector<size_type> &row_block_sizes,
537  const std::vector<size_type> &col_block_sizes);
538 
543  void reinit(const std::vector<IndexSet> &partitioning);
544 
550  void reinit (const BlockIndices &row_indices, const BlockIndices &col_indices);
551 
556  size_type column_number (const size_type row,
557  const unsigned int index) const;
558 
563 };
564 
568 #ifdef DEAL_II_WITH_TRILINOS
569 
570 
571 namespace TrilinosWrappers
572 {
573 
598  public ::BlockSparsityPatternBase<SparsityPattern>
599  {
600  public:
601 
607  BlockSparsityPattern () = default;
608 
615  const size_type n_columns);
616 
623  BlockSparsityPattern (const std::vector<size_type> &row_block_sizes,
624  const std::vector<size_type> &col_block_sizes);
625 
636  DEAL_II_DEPRECATED
637  BlockSparsityPattern (const std::vector<Epetra_Map> &parallel_partitioning);
638 
646  BlockSparsityPattern (const std::vector<IndexSet> &parallel_partitioning,
647  const MPI_Comm &communicator = MPI_COMM_WORLD);
648 
660  BlockSparsityPattern (const std::vector<IndexSet> &row_parallel_partitioning,
661  const std::vector<IndexSet> &column_parallel_partitioning,
662  const std::vector<IndexSet> &writeable_rows,
663  const MPI_Comm &communicator = MPI_COMM_WORLD);
664 
674  void reinit (const std::vector<size_type> &row_block_sizes,
675  const std::vector<size_type> &col_block_sizes);
676 
684  DEAL_II_DEPRECATED
685  void reinit (const std::vector<Epetra_Map> &parallel_partitioning);
686 
691  void reinit (const std::vector<IndexSet> &parallel_partitioning,
692  const MPI_Comm &communicator = MPI_COMM_WORLD);
693 
699  void reinit (const std::vector<IndexSet> &row_parallel_partitioning,
700  const std::vector<IndexSet> &column_parallel_partitioning,
701  const MPI_Comm &communicator = MPI_COMM_WORLD);
702 
711  void reinit (const std::vector<IndexSet> &row_parallel_partitioning,
712  const std::vector<IndexSet> &column_parallel_partitioning,
713  const std::vector<IndexSet> &writeable_rows,
714  const MPI_Comm &communicator = MPI_COMM_WORLD);
715 
720  };
721 
724 } /* namespace TrilinosWrappers */
725 
726 #endif
727 
728 /*---------------------- Template functions -----------------------------------*/
729 
730 
731 
732 template <typename SparsityPatternType>
733 inline
734 SparsityPatternType &
736  const size_type column)
737 {
738  Assert (row<rows, ExcIndexRange(row,0,rows));
739  Assert (column<columns, ExcIndexRange(column,0,columns));
740  return *sub_objects[row][column];
741 }
742 
743 
744 
745 template <typename SparsityPatternType>
746 inline
747 const SparsityPatternType &
749  const size_type column) const
750 {
751  Assert (row<rows, ExcIndexRange(row,0,rows));
752  Assert (column<columns, ExcIndexRange(column,0,columns));
753  return *sub_objects[row][column];
754 }
755 
756 
757 
758 template <typename SparsityPatternType>
759 inline
760 const BlockIndices &
762 {
763  return row_indices;
764 }
765 
766 
767 
768 template <typename SparsityPatternType>
769 inline
770 const BlockIndices &
772 {
773  return column_indices;
774 }
775 
776 
777 
778 template <typename SparsityPatternType>
779 inline
780 void
782  const size_type j)
783 {
784  // if you get an error here, are
785  // you sure you called
786  // <tt>collect_sizes()</tt> before?
787  const std::pair<size_type,size_type>
788  row_index = row_indices.global_to_local (i),
789  col_index = column_indices.global_to_local (j);
790  sub_objects[row_index.first][col_index.first]->add (row_index.second,
791  col_index.second);
792 }
793 
794 
795 
796 template <typename SparsityPatternType>
797 template <typename ForwardIterator>
798 void
800  ForwardIterator begin,
801  ForwardIterator end,
802  const bool indices_are_sorted)
803 {
804  // Resize scratch arrays
805  if (block_column_indices.size() < this->n_block_cols())
806  {
807  block_column_indices.resize (this->n_block_cols());
808  counter_within_block.resize (this->n_block_cols());
809  }
810 
811  const size_type n_cols = static_cast<size_type>(end - begin);
812 
813  // Resize sub-arrays to n_cols. This
814  // is a bit wasteful, but we resize
815  // only a few times (then the maximum
816  // row length won't increase that
817  // much any more). At least we know
818  // that all arrays are going to be of
819  // the same size, so we can check
820  // whether the size of one is large
821  // enough before actually going
822  // through all of them.
823  if (block_column_indices[0].size() < n_cols)
824  for (size_type i=0; i<this->n_block_cols(); ++i)
825  block_column_indices[i].resize(n_cols);
826 
827  // Reset the number of added elements
828  // in each block to zero.
829  for (size_type i=0; i<this->n_block_cols(); ++i)
830  counter_within_block[i] = 0;
831 
832  // Go through the column indices to
833  // find out which portions of the
834  // values should be set in which
835  // block of the matrix. We need to
836  // touch all the data, since we can't
837  // be sure that the data of one block
838  // is stored contiguously (in fact,
839  // indices will be intermixed when it
840  // comes from an element matrix).
841  for (ForwardIterator it = begin; it != end; ++it)
842  {
843  const size_type col = *it;
844 
845  const std::pair<size_type, size_type>
846  col_index = this->column_indices.global_to_local(col);
847 
848  const size_type local_index = counter_within_block[col_index.first]++;
849 
850  block_column_indices[col_index.first][local_index] = col_index.second;
851  }
852 
853 #ifdef DEBUG
854  // If in debug mode, do a check whether
855  // the right length has been obtained.
856  size_type length = 0;
857  for (size_type i=0; i<this->n_block_cols(); ++i)
858  length += counter_within_block[i];
859  Assert (length == n_cols, ExcInternalError());
860 #endif
861 
862  // Now we found out about where the
863  // individual columns should start and
864  // where we should start reading out
865  // data. Now let's write the data into
866  // the individual blocks!
867  const std::pair<size_type, size_type>
868  row_index = this->row_indices.global_to_local (row);
869  for (size_type block_col=0; block_col<n_block_cols(); ++block_col)
870  {
871  if (counter_within_block[block_col] == 0)
872  continue;
873  sub_objects[row_index.first][block_col]->
874  add_entries (row_index.second,
875  block_column_indices[block_col].begin(),
876  block_column_indices[block_col].begin()+counter_within_block[block_col],
877  indices_are_sorted);
878  }
879 }
880 
881 
882 
883 template <typename SparsityPatternType>
884 inline
885 bool
887  const size_type j) const
888 {
889  // if you get an error here, are
890  // you sure you called
891  // <tt>collect_sizes()</tt> before?
892  const std::pair<size_type, size_type>
893  row_index = row_indices.global_to_local (i),
894  col_index = column_indices.global_to_local (j);
895  return sub_objects[row_index.first][col_index.first]->exists (row_index.second,
896  col_index.second);
897 }
898 
899 
900 
901 template <typename SparsityPatternType>
902 inline
903 unsigned int
905 row_length (const size_type row) const
906 {
907  const std::pair<size_type, size_type>
908  row_index = row_indices.global_to_local (row);
909 
910  unsigned int c = 0;
911 
912  for (size_type b=0; b<rows; ++b)
913  c += sub_objects[row_index.first][b]->row_length (row_index.second);
914 
915  return c;
916 }
917 
918 
919 
920 template <typename SparsityPatternType>
921 inline
924 {
925  return columns;
926 }
927 
928 
929 
930 template <typename SparsityPatternType>
931 inline
934 {
935  return rows;
936 }
937 
938 
939 inline
942  const unsigned int index) const
943 {
944  // .first= ith block, .second = jth row in that block
945  const std::pair<size_type,size_type >
946  row_index = row_indices.global_to_local (row);
947 
948  Assert(index<row_length(row), ExcIndexRange(index, 0, row_length(row)));
949 
950  size_type c = 0;
951  size_type block_columns = 0; //sum of n_cols for all blocks to the left
952  for (unsigned int b=0; b<columns; ++b)
953  {
954  unsigned int rowlen = sub_objects[row_index.first][b]->row_length (row_index.second);
955  if (index<c+rowlen)
956  return block_columns+sub_objects[row_index.first][b]->column_number(row_index.second, index-c);
957  c += rowlen;
958  block_columns += sub_objects[row_index.first][b]->n_cols();
959  }
960 
961  Assert(false, ExcInternalError());
962  return 0;
963 }
964 
965 
966 inline
967 void
969  const size_type n_block_rows,
970  const size_type n_block_columns)
971 {
973  n_block_rows, n_block_columns);
974 }
975 
976 
977 DEAL_II_NAMESPACE_CLOSE
978 
979 #endif
static ::ExceptionBase & ExcIncompatibleRowNumbers(int arg1, int arg2, int arg3, int arg4)
void reinit(const std::vector< size_type > &row_block_sizes, const std::vector< size_type > &col_block_sizes)
size_type column_number(const size_type row, const unsigned int index) const
void reinit(const std::vector< size_type > &row_block_sizes, const std::vector< size_type > &col_block_sizes)
static const size_type invalid_entry
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
SparsityPatternType & block(const size_type row, const size_type column)
bool exists(const size_type i, const size_type j) const
std::vector< size_type > counter_within_block
void reinit(const size_type n_block_rows, const size_type n_block_columns)
types::global_dof_index size_type
void reinit(const size_type n_block_rows, const size_type n_block_columns)
static const size_type invalid_entry
std::size_t memory_consumption() const
void add(const size_type i, const size_type j)
unsigned int global_dof_index
Definition: types.h:88
void print(std::ostream &out) const
#define Assert(cond, exc)
Definition: exceptions.h:1142
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_sorted=false)
const BlockIndices & get_row_indices() const
BlockSparsityPatternBase & operator=(const BlockSparsityPatternBase &)
static ::ExceptionBase & ExcIncompatibleColNumbers(int arg1, int arg2, int arg3, int arg4)
std::vector< std::vector< size_type > > block_column_indices
Table< 2, SmartPointer< SparsityPatternType, BlockSparsityPatternBase< SparsityPatternType > > > sub_objects
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition: exceptions.h:382
unsigned int row_length(const size_type row) const
void print_gnuplot(std::ostream &out) const
Definition: table.h:34
void copy_from(const BlockDynamicSparsityPattern &dsp)
BlockSparsityPattern()=default
std::pair< unsigned int, size_type > global_to_local(const size_type i) const
const BlockIndices & get_column_indices() const
static ::ExceptionBase & ExcInternalError()