Reference documentation for deal.II version 9.0.0
filtered_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 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 #ifndef dealii_filtered_matrix_h
17 #define dealii_filtered_matrix_h
18 
19 
20 
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/smartpointer.h>
23 #include <deal.II/base/thread_management.h>
24 #include <deal.II/base/memory_consumption.h>
25 #include <deal.II/lac/pointer_matrix.h>
26 #include <deal.II/lac/vector_memory.h>
27 #include <vector>
28 #include <algorithm>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 template <class VectorType> class FilteredMatrixBlock;
33 
194 template <typename VectorType>
195 class DEAL_II_DEPRECATED FilteredMatrix : public Subscriptor
196 {
197 public:
199 
204 
208  class Accessor
209  {
214  Accessor (const FilteredMatrix<VectorType> *matrix,
215  const size_type index);
216 
217  public:
221  size_type row() const;
222 
226  size_type column() const;
227 
231  double value() const;
232 
233  private:
237  void advance ();
238 
243 
248  /*
249  * Make enclosing class a
250  * friend.
251  */
252  friend class const_iterator;
253  };
254 
259  {
260  public:
265  const size_type index);
266 
270  const_iterator &operator++ ();
271 
275  const_iterator &operator++ (int);
276 
280  const Accessor &operator* () const;
281 
285  const Accessor *operator-> () const;
286 
290  bool operator == (const const_iterator &) const;
294  bool operator != (const const_iterator &) const;
295 
300  bool operator < (const const_iterator &) const;
301 
306  bool operator > (const const_iterator &) const;
307 
308  private:
313  };
314 
319  typedef std::pair<size_type, double> IndexValuePair;
320 
329  FilteredMatrix ();
330 
335  FilteredMatrix (const FilteredMatrix &fm);
336 
345  template <typename MatrixType>
346  FilteredMatrix (const MatrixType &matrix,
347  const bool expect_constrained_source = false);
348 
353 
363  template <typename MatrixType>
364  void initialize (const MatrixType &m,
365  const bool expect_constrained_source = false);
366 
370  void clear ();
372 
380  void add_constraint (const size_type i, const double v);
381 
397  template <class ConstraintList>
398  void add_constraints (const ConstraintList &new_constraints);
399 
403  void clear_constraints ();
405 
416  DEAL_II_DEPRECATED
417  void apply_constraints (VectorType &v,
418  const bool matrix_is_symmetric) const;
423  void apply_constraints (VectorType &v) const;
424 
429  void vmult (VectorType &dst,
430  const VectorType &src) const;
431 
437  void Tvmult (VectorType &dst,
438  const VectorType &src) const;
439 
447  void vmult_add (VectorType &dst,
448  const VectorType &src) const;
449 
457  void Tvmult_add (VectorType &dst,
458  const VectorType &src) const;
460 
468  const_iterator begin () const;
472  const_iterator end () const;
474 
480  std::size_t memory_consumption () const;
481 
482 private:
495 
501  typedef typename std::vector<IndexValuePair>::const_iterator const_index_value_iterator;
502 
508  {
512  bool operator () (const IndexValuePair &i1,
513  const IndexValuePair &i2) const;
514  };
515 
519  std::shared_ptr<PointerMatrixBase<VectorType> > matrix;
520 
525  std::vector<IndexValuePair> constraints;
526 
531  void pre_filter (VectorType &v) const;
532 
538  void post_filter (const VectorType &in,
539  VectorType &out) const;
540 
541  friend class Accessor;
545  friend class FilteredMatrixBlock<VectorType>;
546 };
547 
549 /*---------------------- Inline functions -----------------------------------*/
550 
551 
552 //--------------------------------Iterators--------------------------------------//
553 
554 template <typename VectorType>
555 inline
558  const size_type index)
559  :
560  matrix(matrix),
561  index(index)
562 {
563  Assert (index <= matrix->constraints.size(),
564  ExcIndexRange(index, 0, matrix->constraints.size()));
565 }
566 
567 
568 
569 template <typename VectorType>
570 inline
573 {
574  return matrix->constraints[index].first;
575 }
576 
577 
578 
579 template <typename VectorType>
580 inline
583 {
584  return matrix->constraints[index].first;
585 }
586 
587 
588 
589 template <typename VectorType>
590 inline
591 double
593 {
594  return matrix->constraints[index].second;
595 }
596 
597 
598 
599 template <typename VectorType>
600 inline
601 void
603 {
604  Assert (index < matrix->constraints.size(), ExcIteratorPastEnd());
605  ++index;
606 }
607 
608 
609 
610 
611 template <typename VectorType>
612 inline
615  const size_type index)
616  :
617  accessor(matrix, index)
618 {}
619 
620 
621 
622 template <typename VectorType>
623 inline
626 {
627  accessor.advance();
628  return *this;
629 }
630 
631 
632 template <typename number>
633 inline
634 const typename FilteredMatrix<number>::Accessor &
636 {
637  return accessor;
638 }
639 
640 
641 template <typename number>
642 inline
643 const typename FilteredMatrix<number>::Accessor *
645 {
646  return &accessor;
647 }
648 
649 
650 template <typename number>
651 inline
652 bool
654 operator == (const const_iterator &other) const
655 {
656  return (accessor.index == other.accessor.index
657  && accessor.matrix == other.accessor.matrix);
658 }
659 
660 
661 template <typename number>
662 inline
663 bool
665 operator != (const const_iterator &other) const
666 {
667  return ! (*this == other);
668 }
669 
670 
671 
672 //------------------------------- FilteredMatrix ---------------------------------------//
673 
674 template <typename number>
675 inline
678 {
679  return const_iterator(this, 0);
680 }
681 
682 
683 template <typename number>
684 inline
687 {
688  return const_iterator(this, constraints.size());
689 }
690 
691 
692 template <typename VectorType>
693 inline
694 bool
697  const IndexValuePair &i2) const
698 {
699  return (i1.first < i2.first);
700 }
701 
702 
703 
704 template <typename VectorType>
705 template <typename MatrixType>
706 inline
707 void
708 FilteredMatrix<VectorType>::initialize (const MatrixType &m, bool ecs)
709 {
710  matrix.reset (new_pointer_matrix_base(m, VectorType()));
711 
713 }
714 
715 
716 
717 template <typename VectorType>
718 inline
720  :
722 {}
723 
724 
725 
726 template <typename VectorType>
727 inline
729  :
730  Subscriptor(),
731  expect_constrained_source(fm.expect_constrained_source),
732  matrix(fm.matrix),
733  constraints (fm.constraints)
734 {}
735 
736 
737 
738 template <typename VectorType>
739 template <typename MatrixType>
740 inline
742 FilteredMatrix (const MatrixType &m,
743  const bool ecs)
744  :
745  expect_constrained_source (false)
746 {
747  initialize (m, ecs);
748 }
749 
750 
751 
752 template <typename VectorType>
753 inline
756 {
757  matrix = fm.matrix;
758  expect_constrained_source = fm.expect_constrained_source;
759  constraints = fm.constraints;
760  return *this;
761 }
762 
763 
764 
765 template <typename VectorType>
766 inline
767 void
768 FilteredMatrix<VectorType>::add_constraint (const size_type index, const double value)
769 {
770  // add new constraint to end
771  constraints.push_back(IndexValuePair(index, value));
772 }
773 
774 
775 
776 template <typename VectorType>
777 template <class ConstraintList>
778 inline
779 void
780 FilteredMatrix<VectorType>::add_constraints (const ConstraintList &new_constraints)
781 {
782  // add new constraints to end
783  const size_type old_size = constraints.size();
784  constraints.reserve (old_size + new_constraints.size());
785  constraints.insert (constraints.end(),
786  new_constraints.begin(),
787  new_constraints.end());
788  // then merge the two arrays to
789  // form one sorted one
790  std::inplace_merge (constraints.begin(),
791  constraints.begin()+old_size,
792  constraints.end(),
793  PairComparison());
794 }
795 
796 
797 
798 template <typename VectorType>
799 inline
800 void
802 {
803  // swap vectors to release memory
804  std::vector<IndexValuePair> empty;
805  constraints.swap (empty);
806 }
807 
808 
809 
810 template <typename VectorType>
811 inline
812 void
814 {
815  clear_constraints();
816  matrix.reset();
817 }
818 
819 
820 
821 template <typename VectorType>
822 inline
823 void
825 (VectorType &v,
826  const bool /* matrix_is_symmetric */) const
827 {
828  apply_constraints(v);
829 }
830 
831 
832 template <typename VectorType>
833 inline
834 void
836 {
838  typename VectorMemory<VectorType>::Pointer tmp_vector(mem);
839  tmp_vector->reinit(v);
840  const_index_value_iterator i = constraints.begin();
841  const const_index_value_iterator e = constraints.end();
842  for (; i!=e; ++i)
843  {
844  AssertIsFinite(i->second);
845  (*tmp_vector)(i->first) = -i->second;
846  }
847 
848  // This vmult is without bc, to get
849  // the rhs correction in a correct
850  // way.
851  matrix->vmult_add(v, *tmp_vector);
852  // finally set constrained
853  // entries themselves
854  for (i=constraints.begin(); i!=e; ++i)
855  {
856  AssertIsFinite(i->second);
857  v(i->first) = i->second;
858  }
859 }
860 
861 
862 template <typename VectorType>
863 inline
864 void
866 {
867  // iterate over all constraints and
868  // zero out value
869  const_index_value_iterator i = constraints.begin();
870  const const_index_value_iterator e = constraints.end();
871  for (; i!=e; ++i)
872  v(i->first) = 0;
873 }
874 
875 
876 
877 template <typename VectorType>
878 inline
879 void
881  VectorType &out) const
882 {
883  // iterate over all constraints and
884  // set value correctly
885  const_index_value_iterator i = constraints.begin();
886  const const_index_value_iterator e = constraints.end();
887  for (; i!=e; ++i)
888  {
889  AssertIsFinite(in(i->first));
890  out(i->first) = in(i->first);
891  }
892 }
893 
894 
895 
896 template <typename VectorType>
897 inline
898 void
899 FilteredMatrix<VectorType>::vmult (VectorType &dst, const VectorType &src) const
900 {
901  if (!expect_constrained_source)
902  {
904  typename VectorMemory<VectorType>::Pointer tmp_vector (mem);
905  // first copy over src vector and
906  // pre-filter
907  tmp_vector->reinit(src, true);
908  *tmp_vector = src;
909  pre_filter (*tmp_vector);
910  // then let matrix do its work
911  matrix->vmult (dst, *tmp_vector);
912  }
913  else
914  {
915  matrix->vmult (dst, src);
916  }
917 
918  // finally do post-filtering
919  post_filter (src, dst);
920 }
921 
922 
923 
924 template <typename VectorType>
925 inline
926 void
927 FilteredMatrix<VectorType>::Tvmult (VectorType &dst, const VectorType &src) const
928 {
929  if (!expect_constrained_source)
930  {
932  typename VectorMemory<VectorType>::Pointer tmp_vector (mem);
933  // first copy over src vector and
934  // pre-filter
935  tmp_vector->reinit(src, true);
936  *tmp_vector = src;
937  pre_filter (*tmp_vector);
938  // then let matrix do its work
939  matrix->Tvmult (dst, *tmp_vector);
940  }
941  else
942  {
943  matrix->Tvmult (dst, src);
944  }
945 
946  // finally do post-filtering
947  post_filter (src, dst);
948 }
949 
950 
951 
952 template <typename VectorType>
953 inline
954 void
955 FilteredMatrix<VectorType>::vmult_add (VectorType &dst, const VectorType &src) const
956 {
957  if (!expect_constrained_source)
958  {
960  typename VectorMemory<VectorType>::Pointer tmp_vector (mem);
961  // first copy over src vector and
962  // pre-filter
963  tmp_vector->reinit(src, true);
964  *tmp_vector = src;
965  pre_filter (*tmp_vector);
966  // then let matrix do its work
967  matrix->vmult_add (dst, *tmp_vector);
968  }
969  else
970  {
971  matrix->vmult_add (dst, src);
972  }
973 
974  // finally do post-filtering
975  post_filter (src, dst);
976 }
977 
978 
979 
980 template <typename VectorType>
981 inline
982 void
983 FilteredMatrix<VectorType>::Tvmult_add (VectorType &dst, const VectorType &src) const
984 {
985  if (!expect_constrained_source)
986  {
988  typename VectorMemory<VectorType>::Pointer tmp_vector (mem);
989  // first copy over src vector and
990  // pre-filter
991  tmp_vector->reinit(src, true);
992  *tmp_vector = src;
993  pre_filter (*tmp_vector);
994  // then let matrix do its work
995  matrix->Tvmult_add (dst, *tmp_vector);
996  }
997  else
998  {
999  matrix->Tvmult_add (dst, src);
1000  }
1001 
1002  // finally do post-filtering
1003  post_filter (src, dst);
1004 }
1005 
1006 
1007 
1008 template <typename VectorType>
1009 inline
1010 std::size_t
1012 {
1013  return (MemoryConsumption::memory_consumption (matrix) +
1015 }
1016 
1017 
1018 
1019 DEAL_II_NAMESPACE_CLOSE
1020 
1021 #endif
1022 /*---------------------------- filtered_matrix.h ---------------------------*/
void initialize(const MatrixType &m, const bool expect_constrained_source=false)
const FilteredMatrix< VectorType > * matrix
const Accessor * operator->() const
const_iterator begin() const
types::global_dof_index size_type
Subscriptor & operator=(const Subscriptor &)
Definition: subscriptor.cc:129
void vmult_add(VectorType &dst, const VectorType &src) const
size_type column() const
void pre_filter(VectorType &v) const
std::vector< IndexValuePair >::const_iterator const_index_value_iterator
std::shared_ptr< PointerMatrixBase< VectorType > > matrix
void apply_constraints(VectorType &v, const bool matrix_is_symmetric) const
void post_filter(const VectorType &in, VectorType &out) const
std::vector< IndexValuePair > constraints
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
SymmetricTensor< rank_, dim, Number > operator*(const SymmetricTensor< rank_, dim, Number > &t, const Number &factor)
FilteredMatrix & operator=(const FilteredMatrix &fm)
Accessor(const FilteredMatrix< VectorType > *matrix, const size_type index)
void add_constraint(const size_type i, const double v)
PointerMatrixBase< VectorType > * new_pointer_matrix_base(MatrixType &matrix, const VectorType &, const char *name="PointerMatrixAux")
const Accessor & operator*() const
void vmult(VectorType &dst, const VectorType &src) const
unsigned int global_dof_index
Definition: types.h:88
void Tvmult(VectorType &dst, const VectorType &src) const
#define Assert(cond, exc)
Definition: exceptions.h:1142
bool expect_constrained_source
bool operator()(const IndexValuePair &i1, const IndexValuePair &i2) const
bool operator!=(const const_iterator &) const
bool operator==(const const_iterator &) const
void Tvmult_add(VectorType &dst, const VectorType &src) const
void add_constraints(const ConstraintList &new_constraints)
void clear_constraints()
static ::ExceptionBase & ExcIteratorPastEnd()
std::pair< size_type, double > IndexValuePair
const_iterator(const FilteredMatrix< VectorType > *matrix, const size_type index)
#define AssertIsFinite(number)
Definition: exceptions.h:1299
std::size_t memory_consumption() const
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
const_iterator end() const