Reference documentation for deal.II version 8.5.1
matrix_block.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2007 - 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__matrix_block_h
17 #define dealii__matrix_block_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/smartpointer.h>
21 #include <deal.II/base/std_cxx11/shared_ptr.h>
22 #include <deal.II/base/memory_consumption.h>
23 #include <deal.II/base/mg_level_object.h>
24 #include <deal.II/lac/block_indices.h>
25 #include <deal.II/lac/block_sparsity_pattern.h>
26 #include <deal.II/lac/sparse_matrix.h>
27 #include <deal.II/lac/full_matrix.h>
28 #include <deal.II/algorithms/any_data.h>
29 
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 template <typename MatrixType> class MatrixBlock;
34 
35 namespace internal
36 {
37  template <typename MatrixType>
38  void
39  reinit(MatrixBlock<MatrixType> &v,
40  const BlockSparsityPattern &p);
41 
42  template <typename number>
43  void
44  reinit(MatrixBlock<::SparseMatrix<number> > &v,
45  const BlockSparsityPattern &p);
46 }
47 
103 template <typename MatrixType>
104 class MatrixBlock
105  : public Subscriptor
106 {
107 public:
112 
116  typedef typename MatrixType::value_type value_type;
117 
121  MatrixBlock();
122 
127 
133 
141  void reinit(const BlockSparsityPattern &sparsity);
142 
143  operator MatrixType &();
144  operator const MatrixType &() const;
145 
150  void add (const size_type i,
151  const size_type j,
152  const typename MatrixType::value_type value);
153 
169  template <typename number>
170  void add (const std::vector<size_type> &indices,
171  const FullMatrix<number> &full_matrix,
172  const bool elide_zero_values = true);
173 
188  template <typename number>
189  void add (const std::vector<size_type> &row_indices,
190  const std::vector<size_type> &col_indices,
191  const FullMatrix<number> &full_matrix,
192  const bool elide_zero_values = true);
193 
210  template <typename number>
211  void add (const size_type row_index,
212  const std::vector<size_type> &col_indices,
213  const std::vector<number> &values,
214  const bool elide_zero_values = true);
215 
225  template <typename number>
226  void add (const size_type row,
227  const size_type n_cols,
228  const size_type *col_indices,
229  const number *values,
230  const bool elide_zero_values = true,
231  const bool col_indices_are_sorted = false);
232 
238  template<class VectorType>
239  void vmult (VectorType &w, const VectorType &v) const;
240 
246  template<class VectorType>
247  void vmult_add (VectorType &w, const VectorType &v) const;
248 
254  template<class VectorType>
255  void Tvmult (VectorType &w, const VectorType &v) const;
256 
262  template<class VectorType>
263  void Tvmult_add (VectorType &w, const VectorType &v) const;
264 
268  std::size_t memory_consumption () const;
269 
275  << "Block index " << arg1 << " does not match " << arg2);
276 
287 
291  MatrixType matrix;
292 
293 private:
305 
306  template <class OTHER_MatrixType>
307  friend
308  void ::internal::reinit(MatrixBlock<OTHER_MatrixType> &,
309  const BlockSparsityPattern &);
310 
311  template <typename number>
312  friend
313  void internal::reinit(MatrixBlock<::SparseMatrix<number> > &v,
314  const BlockSparsityPattern &p);
315 };
316 
317 
327 template <typename MatrixType>
329  :
330  private AnyData
331 {
332 public:
337 
342 
347  typedef std_cxx11::shared_ptr<value_type> ptr_type;
348 
353  void add(size_type row, size_type column, const std::string &name);
354 
359  void reinit(const BlockSparsityPattern &sparsity);
360 
371  void clear (bool really_clean = false);
372 
376  std::size_t memory_consumption () const;
377 
381  const value_type &block(size_type i) const;
382 
387 
391  MatrixType &matrix(size_type i);
392 
396  using AnyData::subscribe;
397  using AnyData::unsubscribe;
398  using AnyData::size;
399  using AnyData::name;
400 };
401 
402 
412 template <typename MatrixType>
414  : public Subscriptor
415 {
416 public:
421 
435  MGMatrixBlockVector(const bool edge_matrices = false,
436  const bool edge_flux_matrices = false);
437 
441  unsigned int size () const;
442 
448  void add(size_type row, size_type column, const std::string &name);
449 
464  void reinit_edge(const MGLevelObject<BlockSparsityPattern> &sparsity);
472 
483  void clear (bool really_clean = false);
484 
488  const value_type &block(size_type i) const;
489 
494 
499  const value_type &block_in(size_type i) const;
500 
505 
510  const value_type &block_out(size_type i) const;
511 
516 
521  const value_type &block_up(size_type i) const;
522 
527 
532  const value_type &block_down(size_type i) const;
533 
538 
542  std::size_t memory_consumption () const;
543 private:
545  void clear_object(AnyData &);
546 
548  const bool edge_matrices;
549 
551  const bool edge_flux_matrices;
552 
563 };
564 
565 
566 //----------------------------------------------------------------------//
567 
568 namespace internal
569 {
570  template <typename MatrixType>
571  void
572  reinit(MatrixBlock<MatrixType> &v,
573  const BlockSparsityPattern &p)
574  {
577  }
578 
579 
580  template <typename number>
581  void
582  reinit(MatrixBlock<::SparseMatrix<number> > &v,
583  const BlockSparsityPattern &p)
584  {
587  v.matrix.reinit(p.block(v.row, v.column));
588  }
589 }
590 
591 
592 template <typename MatrixType>
593 inline
595  :
596  row(numbers::invalid_size_type),
597  column(numbers::invalid_size_type)
598 {}
599 
600 
601 template <typename MatrixType>
602 inline
604  :
605  Subscriptor(),
606  row(M.row),
607  column(M.column),
608  matrix(M.matrix),
609  row_indices(M.row_indices),
610  column_indices(M.column_indices)
611 {}
612 
613 
614 template <typename MatrixType>
615 inline
617  :
618  row(i), column(j)
619 {}
620 
621 
622 template <typename MatrixType>
623 inline
624 void
626 {
627  internal::reinit(*this, sparsity);
628 }
629 
630 
631 template <typename MatrixType>
632 inline
634 {
635  return matrix;
636 }
637 
638 
639 template <typename MatrixType>
640 inline
641 MatrixBlock<MatrixType>::operator const MatrixType &() const
642 {
643  return matrix;
644 }
645 
646 
647 template <typename MatrixType>
648 inline void
650  const size_type gj,
651  const typename MatrixType::value_type value)
652 {
653  Assert(row_indices.size() != 0, ExcNotInitialized());
654  Assert(column_indices.size() != 0, ExcNotInitialized());
655 
656  const std::pair<unsigned int, size_type> bi
657  = row_indices.global_to_local(gi);
658  const std::pair<unsigned int, size_type> bj
659  = column_indices.global_to_local(gj);
660 
661  Assert (bi.first == row, ExcBlockIndexMismatch(bi.first, row));
662  Assert (bj.first == column, ExcBlockIndexMismatch(bj.first, column));
663 
664  matrix.add(bi.second, bj.second, value);
665 }
666 
667 
668 template <typename MatrixType>
669 template <typename number>
670 inline
671 void
672 MatrixBlock<MatrixType>::add (const std::vector<size_type> &r_indices,
673  const std::vector<size_type> &c_indices,
674  const FullMatrix<number> &values,
675  const bool elide_zero_values)
676 {
677  Assert(row_indices.size() != 0, ExcNotInitialized());
678  Assert(column_indices.size() != 0, ExcNotInitialized());
679 
680  AssertDimension (r_indices.size(), values.m());
681  AssertDimension (c_indices.size(), values.n());
682 
683  for (size_type i=0; i<row_indices.size(); ++i)
684  add (r_indices[i], c_indices.size(), &c_indices[0], &values(i,0),
685  elide_zero_values);
686 }
687 
688 
689 template <typename MatrixType>
690 template <typename number>
691 inline
692 void
694  const size_type n_cols,
695  const size_type *col_indices,
696  const number *values,
697  const bool,
698  const bool)
699 {
700  Assert(row_indices.size() != 0, ExcNotInitialized());
701  Assert(column_indices.size() != 0, ExcNotInitialized());
702 
703  const std::pair<unsigned int, size_type> bi
704  = row_indices.global_to_local(b_row);
705 
706  // In debug mode, we check whether
707  // all indices are in the correct
708  // block.
709 
710  // Actually, for the time being, we
711  // leave it at this. While it may
712  // not be the most efficient way,
713  // it is at least thread safe.
714 //#ifdef DEBUG
715  Assert(bi.first == row, ExcBlockIndexMismatch(bi.first, row));
716 
717  for (size_type j=0; j<n_cols; ++j)
718  {
719  const std::pair<unsigned int, size_type> bj
720  = column_indices.global_to_local(col_indices[j]);
721  Assert(bj.first == column, ExcBlockIndexMismatch(bj.first, column));
722 
723  matrix.add(bi.second, bj.second, values[j]);
724  }
725 //#endif
726 }
727 
728 
729 template <typename MatrixType>
730 template <typename number>
731 inline
732 void
733 MatrixBlock<MatrixType>::add (const std::vector<size_type> &indices,
734  const FullMatrix<number> &values,
735  const bool elide_zero_values)
736 {
737  Assert(row_indices.size() != 0, ExcNotInitialized());
738  Assert(column_indices.size() != 0, ExcNotInitialized());
739 
740  AssertDimension (indices.size(), values.m());
741  Assert (values.n() == values.m(), ExcNotQuadratic());
742 
743  for (size_type i=0; i<indices.size(); ++i)
744  add (indices[i], indices.size(), &indices[0], &values(i,0),
745  elide_zero_values);
746 }
747 
748 
749 
750 template <typename MatrixType>
751 template <typename number>
752 inline
753 void
755  const std::vector<size_type> &col_indices,
756  const std::vector<number> &values,
757  const bool elide_zero_values)
758 {
759  Assert(row_indices.size() != 0, ExcNotInitialized());
760  Assert(column_indices.size() != 0, ExcNotInitialized());
761 
762  AssertDimension (col_indices.size(), values.size());
763  add (row, col_indices.size(), &col_indices[0], &values[0],
764  elide_zero_values);
765 }
766 
767 
768 template <typename MatrixType>
769 template <class VectorType>
770 inline
771 void
772 MatrixBlock<MatrixType>::vmult (VectorType &w, const VectorType &v) const
773 {
774  matrix.vmult(w,v);
775 }
776 
777 
778 template <typename MatrixType>
779 template <class VectorType>
780 inline
781 void
782 MatrixBlock<MatrixType>::vmult_add (VectorType &w, const VectorType &v) const
783 {
784  matrix.vmult_add(w,v);
785 }
786 
787 
788 template <typename MatrixType>
789 template <class VectorType>
790 inline
791 void
792 MatrixBlock<MatrixType>::Tvmult (VectorType &w, const VectorType &v) const
793 {
794  matrix.Tvmult(w,v);
795 }
796 
797 
798 template <typename MatrixType>
799 template <class VectorType>
800 inline
801 void
802 MatrixBlock<MatrixType>::Tvmult_add (VectorType &w, const VectorType &v) const
803 {
804  matrix.Tvmult_add(w,v);
805 }
806 
807 
808 template <typename MatrixType>
809 inline
810 std::size_t
812 {
813  return (sizeof(*this)
815  - sizeof(matrix));
816 }
817 
818 //----------------------------------------------------------------------//
819 
820 template <typename MatrixType>
821 inline void
823  size_type column,
824  const std::string &name)
825 {
826  ptr_type p(new value_type(row, column));
827  AnyData::add(p, name);
828 }
829 
830 
831 template <typename MatrixType>
832 inline void
834 {
835  for (size_type i=0; i<this->size(); ++i)
836  {
837  block(i).reinit(sparsity);
838  }
839 }
840 
841 
842 template <typename MatrixType>
843 inline void
845 {
846  if (really_clean)
847  {
848  Assert(false, ExcNotImplemented());
849  }
850  else
851  {
852  for (size_type i=0; i<this->size(); ++i)
853  matrix(i).clear();
854  }
855 }
856 
857 
858 
859 template <typename MatrixType>
860 inline const MatrixBlock<MatrixType> &
862 {
863  return *this->read<ptr_type>(i);
864 }
865 
866 
867 template <typename MatrixType>
870 {
871  return *this->entry<ptr_type>(i);
872 }
873 
874 
875 template <typename MatrixType>
876 inline MatrixType &
878 {
879  return this->entry<ptr_type>(i)->matrix;
880 }
881 
882 
883 
884 //----------------------------------------------------------------------//
885 
886 template <typename MatrixType>
887 inline
889  :
890  edge_matrices(e),
891  edge_flux_matrices(f)
892 {}
893 
894 
895 template <typename MatrixType>
896 inline
897 unsigned int
899 {
900  return matrices.size();
901 }
902 
903 
904 template <typename MatrixType>
905 inline void
907  size_type row, size_type column,
908  const std::string &name)
909 {
911  p[0].row = row;
912  p[0].column = column;
913 
914  matrices.add(p, name);
915  if (edge_matrices)
916  {
917  matrices_in.add(p, name);
918  matrices_out.add(p, name);
919  }
920  if (edge_flux_matrices)
921  {
922  flux_matrices_up.add(p, name);
923  flux_matrices_down.add(p, name);
924  }
925 }
926 
927 
928 template <typename MatrixType>
931 {
932  return *matrices.read<const MGLevelObject<MatrixType>* >(i);
933 }
934 
935 
936 template <typename MatrixType>
939 {
940  return *matrices.entry<MGLevelObject<MatrixType>* >(i);
941 }
942 
943 
944 template <typename MatrixType>
947 {
948  return *matrices_in.read<const MGLevelObject<MatrixType>* >(i);
949 }
950 
951 
952 template <typename MatrixType>
955 {
956  return *matrices_in.entry<MGLevelObject<MatrixType>* >(i);
957 }
958 
959 
960 template <typename MatrixType>
963 {
964  return *matrices_out.read<const MGLevelObject<MatrixType>* >(i);
965 }
966 
967 
968 template <typename MatrixType>
971 {
972  return *matrices_out.entry<MGLevelObject<MatrixType>* >(i);
973 }
974 
975 
976 template <typename MatrixType>
979 {
980  return *flux_matrices_up.read<const MGLevelObject<MatrixType>* >(i);
981 }
982 
983 
984 template <typename MatrixType>
987 {
988  return *flux_matrices_up.entry<MGLevelObject<MatrixType>* >(i);
989 }
990 
991 
992 template <typename MatrixType>
995 {
996  return *flux_matrices_down.read<const MGLevelObject<MatrixType>* >(i);
997 }
998 
999 
1000 template <typename MatrixType>
1003 {
1004  return *flux_matrices_down.entry<MGLevelObject<MatrixType>* >(i);
1005 }
1006 
1007 
1008 template <typename MatrixType>
1009 inline void
1011 {
1012  for (size_type i=0; i<this->size(); ++i)
1013  {
1014  MGLevelObject<MatrixBlock<MatrixType> > &o = block(i);
1015  const size_type row = o[o.min_level()].row;
1016  const size_type col = o[o.min_level()].column;
1017 
1018  o.resize(sparsity.min_level(), sparsity.max_level());
1019  for (size_type level = o.min_level(); level <= o.max_level(); ++level)
1020  {
1021  o[level].row = row;
1022  o[level].column = col;
1023  internal::reinit(o[level], sparsity[level]);
1024  }
1025  }
1026 }
1027 
1028 
1029 template <typename MatrixType>
1030 inline void
1032 {
1033  for (size_type i=0; i<this->size(); ++i)
1034  {
1035  MGLevelObject<MatrixBlock<MatrixType> > &o = block(i);
1036  const size_type row = o[o.min_level()].row;
1037  const size_type col = o[o.min_level()].column;
1038 
1039  block_in(i).resize(sparsity.min_level(), sparsity.max_level());
1040  block_out(i).resize(sparsity.min_level(), sparsity.max_level());
1041  for (size_type level = o.min_level(); level <= o.max_level(); ++level)
1042  {
1043  block_in(i)[level].row = row;
1044  block_in(i)[level].column = col;
1045  internal::reinit(block_in(i)[level], sparsity[level]);
1046  block_out(i)[level].row = row;
1047  block_out(i)[level].column = col;
1048  internal::reinit(block_out(i)[level], sparsity[level]);
1049  }
1050  }
1051 }
1052 
1053 
1054 template <typename MatrixType>
1055 inline void
1058 {
1059  for (size_type i=0; i<this->size(); ++i)
1060  {
1061  MGLevelObject<MatrixBlock<MatrixType> > &o = block(i);
1062  const size_type row = o[o.min_level()].row;
1063  const size_type col = o[o.min_level()].column;
1064 
1065  block_up(i).resize(sparsity.min_level(), sparsity.max_level());
1066  block_down(i).resize(sparsity.min_level(), sparsity.max_level());
1067  for (size_type level = o.min_level(); level <= o.max_level(); ++level)
1068  {
1069  block_up(i)[level].row = row;
1070  block_up(i)[level].column = col;
1071  internal::reinit(block_up(i)[level], sparsity[level]);
1072  block_down(i)[level].row = row;
1073  block_down(i)[level].column = col;
1074  internal::reinit(block_down(i)[level], sparsity[level]);
1075  }
1076 
1077  }
1078 }
1079 
1080 
1081 template <typename MatrixType>
1082 inline void
1084 {
1085  for (size_type i=0; i<mo.size(); ++i)
1086  {
1088  for (size_type level = o.min_level(); level <= o.max_level(); ++level)
1089  o[level].matrix.clear();
1090  }
1091 }
1092 
1093 
1094 template <typename MatrixType>
1095 inline void
1097 {
1098  if (really_clean)
1099  {
1100  Assert(false, ExcNotImplemented());
1101  }
1102  else
1103  {
1104  clear_object(matrices);
1105  clear_object(matrices_in);
1106  clear_object(matrices_out);
1107  clear_object(flux_matrices_up);
1108  clear_object(flux_matrices_down);
1109  }
1110 }
1111 
1112 
1113 
1114 
1115 
1116 DEAL_II_NAMESPACE_CLOSE
1117 
1118 #endif
void reinit_matrix(const MGLevelObject< BlockSparsityPattern > &sparsity)
unsigned int max_level() const
size_type m() const
const bool edge_matrices
Flag for storing matrices_in and matrices_out.
Definition: matrix_block.h:548
type entry(const std::string &name)
Access to stored data object by name.
Definition: any_data.h:343
AnyData matrices
The level matrices.
Definition: matrix_block.h:554
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
types::global_dof_index size_type
Definition: matrix_block.h:111
MatrixBlock< MatrixType > value_type
Definition: matrix_block.h:341
AnyData matrices_out
The matrix from the refinement edge to the interior of a level.
Definition: matrix_block.h:558
MatrixType & matrix(size_type i)
Definition: matrix_block.h:877
const value_type & block(size_type i) const
Definition: matrix_block.h:861
void Tvmult(VectorType &w, const VectorType &v) const
Definition: matrix_block.h:792
void vmult(VectorType &w, const VectorType &v) const
Definition: matrix_block.h:772
const value_type & block(size_type i) const
Definition: matrix_block.h:930
void subscribe(const char *identifier=0) const
Definition: subscriptor.cc:173
AnyData matrices_in
The matrix from the interior of a level to the refinement edge.
Definition: matrix_block.h:556
types::global_dof_index size_type
Definition: matrix_block.h:336
static ::ExceptionBase & ExcNotInitialized()
void add(const size_type i, const size_type j, const typename MatrixType::value_type value)
Definition: matrix_block.h:649
size_type column
Definition: matrix_block.h:286
MatrixType::value_type value_type
Definition: matrix_block.h:116
unsigned int min_level() const
SparsityPatternType & block(const size_type row, const size_type column)
const bool edge_flux_matrices
Flag for storing flux_matrices_up and flux_matrices_down.
Definition: matrix_block.h:551
void clear(bool really_clean=false)
Definition: matrix_block.h:844
const value_type & block_out(size_type i) const
Definition: matrix_block.h:962
void add(size_type row, size_type column, const std::string &name)
Definition: matrix_block.h:822
void reinit_edge(const MGLevelObject< BlockSparsityPattern > &sparsity)
MGMatrixBlockVector(const bool edge_matrices=false, const bool edge_flux_matrices=false)
Definition: matrix_block.h:888
size_type n() const
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
AnyData flux_matrices_down
The DG flux from a level to the lower level.
Definition: matrix_block.h:560
const BlockIndices & get_row_indices() const
std::size_t memory_consumption() const
Definition: matrix_block.h:811
types::global_dof_index size_type
Definition: matrix_block.h:420
void unsubscribe(const char *identifier=0) const
Definition: subscriptor.cc:200
BlockIndices column_indices
Definition: matrix_block.h:304
const value_type & block_in(size_type i) const
Definition: matrix_block.h:946
void vmult_add(VectorType &w, const VectorType &v) const
Definition: matrix_block.h:782
BlockIndices row_indices
Definition: matrix_block.h:298
MatrixType matrix
Definition: matrix_block.h:291
static ::ExceptionBase & ExcNotQuadratic()
void add(size_type row, size_type column, const std::string &name)
Definition: matrix_block.h:906
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcBlockIndexMismatch(size_type arg1, size_type arg2)
void clear(bool really_clean=false)
unsigned int size() const
Definition: matrix_block.h:898
void reinit_edge_flux(const MGLevelObject< BlockSparsityPattern > &sparsity)
void reinit(const BlockSparsityPattern &sparsity)
Definition: matrix_block.h:833
const value_type & block_up(size_type i) const
Definition: matrix_block.h:978
std::size_t memory_consumption() const
MGLevelObject< MatrixBlock< MatrixType > > value_type
Definition: matrix_block.h:425
void add(type entry, const std::string &name)
Add a new data object.
Definition: any_data.h:432
unsigned int size() const
Number of stored data objects.
Definition: any_data.h:204
void reinit(const BlockSparsityPattern &sparsity)
Definition: matrix_block.h:625
AnyData flux_matrices_up
The DG flux from the lower level to a level.
Definition: matrix_block.h:562
const value_type & block_down(size_type i) const
Definition: matrix_block.h:994
static ::ExceptionBase & ExcNotImplemented()
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel)
void clear_object(AnyData &)
Clear one of the matrix objects.
std_cxx11::shared_ptr< value_type > ptr_type
Definition: matrix_block.h:347
const BlockIndices & get_column_indices() const
void Tvmult_add(VectorType &w, const VectorType &v) const
Definition: matrix_block.h:802
const std::string & name(const unsigned int i) const
Name of object at index.
Definition: any_data.h:299
std::size_t memory_consumption() const
size_type row
Definition: matrix_block.h:281