Reference documentation for deal.II version 8.5.1
block_sparsity_pattern.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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 #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 
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 
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 
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  BlockSparsityPattern (const std::vector<Epetra_Map> &parallel_partitioning) DEAL_II_DEPRECATED;
637 
645  BlockSparsityPattern (const std::vector<IndexSet> &parallel_partitioning,
646  const MPI_Comm &communicator = MPI_COMM_WORLD);
647 
659  BlockSparsityPattern (const std::vector<IndexSet> &row_parallel_partitioning,
660  const std::vector<IndexSet> &column_parallel_partitioning,
661  const std::vector<IndexSet> &writeable_rows,
662  const MPI_Comm &communicator = MPI_COMM_WORLD);
663 
673  void reinit (const std::vector<size_type> &row_block_sizes,
674  const std::vector<size_type> &col_block_sizes);
675 
683  void reinit (const std::vector<Epetra_Map> &parallel_partitioning) DEAL_II_DEPRECATED;
684 
689  void reinit (const std::vector<IndexSet> &parallel_partitioning,
690  const MPI_Comm &communicator = MPI_COMM_WORLD);
691 
697  void reinit (const std::vector<IndexSet> &row_parallel_partitioning,
698  const std::vector<IndexSet> &column_parallel_partitioning,
699  const MPI_Comm &communicator = MPI_COMM_WORLD);
700 
709  void reinit (const std::vector<IndexSet> &row_parallel_partitioning,
710  const std::vector<IndexSet> &column_parallel_partitioning,
711  const std::vector<IndexSet> &writeable_rows,
712  const MPI_Comm &communicator = MPI_COMM_WORLD);
713 
718  };
719 
722 } /* namespace TrilinosWrappers */
723 
724 #endif
725 
726 /*---------------------- Template functions -----------------------------------*/
727 
728 
729 
730 template <typename SparsityPatternType>
731 inline
732 SparsityPatternType &
734  const size_type column)
735 {
736  Assert (row<rows, ExcIndexRange(row,0,rows));
737  Assert (column<columns, ExcIndexRange(column,0,columns));
738  return *sub_objects[row][column];
739 }
740 
741 
742 
743 template <typename SparsityPatternType>
744 inline
745 const SparsityPatternType &
747  const size_type column) const
748 {
749  Assert (row<rows, ExcIndexRange(row,0,rows));
750  Assert (column<columns, ExcIndexRange(column,0,columns));
751  return *sub_objects[row][column];
752 }
753 
754 
755 
756 template <typename SparsityPatternType>
757 inline
758 const BlockIndices &
760 {
761  return row_indices;
762 }
763 
764 
765 
766 template <typename SparsityPatternType>
767 inline
768 const BlockIndices &
770 {
771  return column_indices;
772 }
773 
774 
775 
776 template <typename SparsityPatternType>
777 inline
778 void
780  const size_type j)
781 {
782  // if you get an error here, are
783  // you sure you called
784  // <tt>collect_sizes()</tt> before?
785  const std::pair<size_type,size_type>
786  row_index = row_indices.global_to_local (i),
787  col_index = column_indices.global_to_local (j);
788  sub_objects[row_index.first][col_index.first]->add (row_index.second,
789  col_index.second);
790 }
791 
792 
793 
794 template <typename SparsityPatternType>
795 template <typename ForwardIterator>
796 void
798  ForwardIterator begin,
799  ForwardIterator end,
800  const bool indices_are_sorted)
801 {
802  // Resize scratch arrays
803  if (block_column_indices.size() < this->n_block_cols())
804  {
805  block_column_indices.resize (this->n_block_cols());
806  counter_within_block.resize (this->n_block_cols());
807  }
808 
809  const size_type n_cols = static_cast<size_type>(end - begin);
810 
811  // Resize sub-arrays to n_cols. This
812  // is a bit wasteful, but we resize
813  // only a few times (then the maximum
814  // row length won't increase that
815  // much any more). At least we know
816  // that all arrays are going to be of
817  // the same size, so we can check
818  // whether the size of one is large
819  // enough before actually going
820  // through all of them.
821  if (block_column_indices[0].size() < n_cols)
822  for (size_type i=0; i<this->n_block_cols(); ++i)
823  block_column_indices[i].resize(n_cols);
824 
825  // Reset the number of added elements
826  // in each block to zero.
827  for (size_type i=0; i<this->n_block_cols(); ++i)
828  counter_within_block[i] = 0;
829 
830  // Go through the column indices to
831  // find out which portions of the
832  // values should be set in which
833  // block of the matrix. We need to
834  // touch all the data, since we can't
835  // be sure that the data of one block
836  // is stored contiguously (in fact,
837  // indices will be intermixed when it
838  // comes from an element matrix).
839  for (ForwardIterator it = begin; it != end; ++it)
840  {
841  const size_type col = *it;
842 
843  const std::pair<size_type , size_type>
844  col_index = this->column_indices.global_to_local(col);
845 
846  const size_type local_index = counter_within_block[col_index.first]++;
847 
848  block_column_indices[col_index.first][local_index] = col_index.second;
849  }
850 
851 #ifdef DEBUG
852  // If in debug mode, do a check whether
853  // the right length has been obtained.
854  size_type length = 0;
855  for (size_type i=0; i<this->n_block_cols(); ++i)
856  length += counter_within_block[i];
857  Assert (length == n_cols, ExcInternalError());
858 #endif
859 
860  // Now we found out about where the
861  // individual columns should start and
862  // where we should start reading out
863  // data. Now let's write the data into
864  // the individual blocks!
865  const std::pair<size_type , size_type>
866  row_index = this->row_indices.global_to_local (row);
867  for (size_type block_col=0; block_col<n_block_cols(); ++block_col)
868  {
869  if (counter_within_block[block_col] == 0)
870  continue;
871  sub_objects[row_index.first][block_col]->
872  add_entries (row_index.second,
873  block_column_indices[block_col].begin(),
874  block_column_indices[block_col].begin()+counter_within_block[block_col],
875  indices_are_sorted);
876  }
877 }
878 
879 
880 
881 template <typename SparsityPatternType>
882 inline
883 bool
885  const size_type j) const
886 {
887  // if you get an error here, are
888  // you sure you called
889  // <tt>collect_sizes()</tt> before?
890  const std::pair<size_type , size_type>
891  row_index = row_indices.global_to_local (i),
892  col_index = column_indices.global_to_local (j);
893  return sub_objects[row_index.first][col_index.first]->exists (row_index.second,
894  col_index.second);
895 }
896 
897 
898 
899 template <typename SparsityPatternType>
900 inline
901 unsigned int
903 row_length (const size_type row) const
904 {
905  const std::pair<size_type , size_type>
906  row_index = row_indices.global_to_local (row);
907 
908  unsigned int c = 0;
909 
910  for (size_type b=0; b<rows; ++b)
911  c += sub_objects[row_index.first][b]->row_length (row_index.second);
912 
913  return c;
914 }
915 
916 
917 
918 template <typename SparsityPatternType>
919 inline
922 {
923  return columns;
924 }
925 
926 
927 
928 template <typename SparsityPatternType>
929 inline
932 {
933  return rows;
934 }
935 
936 
937 inline
940  const unsigned int index) const
941 {
942  // .first= ith block, .second = jth row in that block
943  const std::pair<size_type ,size_type >
944  row_index = row_indices.global_to_local (row);
945 
946  Assert(index<row_length(row), ExcIndexRange(index, 0, row_length(row)));
947 
948  size_type c = 0;
949  size_type block_columns = 0; //sum of n_cols for all blocks to the left
950  for (unsigned int b=0; b<columns; ++b)
951  {
952  unsigned int rowlen = sub_objects[row_index.first][b]->row_length (row_index.second);
953  if (index<c+rowlen)
954  return block_columns+sub_objects[row_index.first][b]->column_number(row_index.second, index-c);
955  c += rowlen;
956  block_columns += sub_objects[row_index.first][b]->n_cols();
957  }
958 
959  Assert(false, ExcInternalError());
960  return 0;
961 }
962 
963 
964 inline
965 void
967  const size_type n_block_rows,
968  const size_type n_block_columns)
969 {
971  n_block_rows, n_block_columns);
972 }
973 
974 
975 DEAL_II_NAMESPACE_CLOSE
976 
977 #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:313
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:600
unsigned int row_length(const size_type row) const
void print_gnuplot(std::ostream &out) const
Definition: table.h:33
void copy_from(const BlockDynamicSparsityPattern &dsp)
std::pair< unsigned int, size_type > global_to_local(const size_type i) const
const BlockIndices & get_column_indices() const
static ::ExceptionBase & ExcInternalError()