Reference documentation for deal.II version 9.0.0
multigrid.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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_multigrid_h
17 #define dealii_multigrid_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/subscriptor.h>
22 #include <deal.II/base/smartpointer.h>
23 #include <deal.II/dofs/dof_handler.h>
24 #include <deal.II/lac/sparse_matrix.h>
25 #include <deal.II/lac/vector.h>
26 #include <deal.II/multigrid/mg_base.h>
27 #include <deal.II/base/mg_level_object.h>
28 #include <deal.II/distributed/tria.h>
29 
30 #include <vector>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
36 
37 namespace mg
38 {
47  struct Signals
48  {
54  boost::signals2::signal<void (const bool before)> transfer_to_mg;
55 
61  boost::signals2::signal<void (const bool before)> transfer_to_global;
62 
67  boost::signals2::signal<void (const bool before, const unsigned int level)> coarse_solve;
68 
74  boost::signals2::signal<void (const bool before, const unsigned int level)> restriction;
75 
81  boost::signals2::signal<void (const bool before, const unsigned int level)> prolongation;
82 
88  boost::signals2::signal<void (const bool before, const unsigned int level)> pre_smoother_step;
89 
95  boost::signals2::signal<void (const bool before, const unsigned int level)> post_smoother_step;
96  };
97 }
98 
123 template <typename VectorType>
124 class Multigrid : public Subscriptor
125 {
126 public:
130  enum Cycle
131  {
138  };
139 
140  typedef VectorType vector_type;
141  typedef const VectorType const_vector_type;
142 
158  template <int dim>
159 #ifndef _MSC_VER
160  DEAL_II_DEPRECATED
161 #endif
162  Multigrid(const DoFHandler<dim> &mg_dof_handler,
168  const unsigned int minlevel = 0,
169  const unsigned int maxlevel = numbers::invalid_unsigned_int,
170  Cycle cycle = v_cycle);
171 
186  const unsigned int minlevel = 0,
187  const unsigned int maxlevel = numbers::invalid_unsigned_int,
188  Cycle cycle = v_cycle);
189 
193  void reinit (const unsigned int minlevel,
194  const unsigned int maxlevel);
195 
200  void cycle ();
201 
212  void vcycle ();
213 
228 
243 
247  unsigned int get_maxlevel() const;
248 
252  unsigned int get_minlevel() const;
253 
259  void set_maxlevel (const unsigned int);
260 
276  void set_minlevel (const unsigned int level,
277  bool relative = false);
278 
282  void set_cycle(Cycle);
283 
290  DEAL_II_DEPRECATED
291  void set_debug (const unsigned int);
292 
296  boost::signals2::connection
297  connect_coarse_solve(const std::function<void (const bool, const unsigned int)> &slot);
298 
302  boost::signals2::connection
303  connect_restriction(const std::function<void (const bool, const unsigned int)> &slot);
304 
308  boost::signals2::connection
309  connect_prolongation(const std::function<void (const bool, const unsigned int)> &slot);
310 
314  boost::signals2::connection
315  connect_pre_smoother_step(const std::function<void (const bool, const unsigned int)> &slot);
316 
320  boost::signals2::connection
321  connect_post_smoother_step(const std::function<void (const bool, const unsigned int)> &slot);
322 
323 private:
324 
329 
336  void level_v_step (const unsigned int level);
337 
345  void level_step (const unsigned int level, Cycle cycle);
346 
351 
355  unsigned int minlevel;
356 
360  unsigned int maxlevel;
361 
362 public:
368 
373 
374 private:
379 
384 
385 
390 
395 
400 
405 
410 
417 
425 
432 
439 
443  unsigned int debug;
444 
445  template <int dim, class OtherVectorType, class TRANSFER> friend class PreconditionMG;
446 };
447 
448 
465 template <int dim, typename VectorType, class TRANSFER>
467 {
468 public:
473  PreconditionMG(const DoFHandler<dim> &dof_handler,
475  const TRANSFER &transfer);
476 
481  PreconditionMG(const std::vector<const DoFHandler<dim>*> &dof_handler,
483  const TRANSFER &transfer);
484 
488  bool empty () const;
489 
496  template <class OtherVectorType>
497  void vmult (OtherVectorType &dst,
498  const OtherVectorType &src) const;
499 
504  template <class OtherVectorType>
505  void vmult_add (OtherVectorType &dst,
506  const OtherVectorType &src) const;
507 
513  template <class OtherVectorType>
514  void Tvmult (OtherVectorType &dst,
515  const OtherVectorType &src) const;
516 
522  template <class OtherVectorType>
523  void Tvmult_add (OtherVectorType &dst,
524  const OtherVectorType &src) const;
525 
532  IndexSet locally_owned_range_indices(const unsigned int block=0) const;
533 
540  IndexSet locally_owned_domain_indices(const unsigned int block=0) const;
541 
545  MPI_Comm get_mpi_communicator() const;
546 
550  boost::signals2::connection
551  connect_transfer_to_mg(const std::function<void (bool)> &slot);
552 
556  boost::signals2::connection
557  connect_transfer_to_global(const std::function<void (bool)> &slot);
558 
559 private:
563  std::vector<SmartPointer<const DoFHandler<dim>,PreconditionMG<dim,VectorType,TRANSFER> > > dof_handler_vector;
564 
569  std::vector<const DoFHandler<dim>*> dof_handler_vector_raw;
570 
575 
580 
586 
591 };
592 
595 #ifndef DOXYGEN
596 /* --------------------------- inline functions --------------------- */
597 
598 
599 template <typename VectorType>
600 template <int dim>
602  const MGMatrixBase<VectorType> &matrix,
603  const MGCoarseGridBase<VectorType> &coarse,
604  const MGTransferBase<VectorType> &transfer,
605  const MGSmootherBase<VectorType> &pre_smooth,
606  const MGSmootherBase<VectorType> &post_smooth,
607  const unsigned int min_level,
608  const unsigned int max_level,
609  Cycle cycle)
610  :
611  cycle_type(cycle),
612  minlevel(min_level),
613  matrix(&matrix, typeid(*this).name()),
614  coarse(&coarse, typeid(*this).name()),
615  transfer(&transfer, typeid(*this).name()),
616  pre_smooth(&pre_smooth, typeid(*this).name()),
617  post_smooth(&post_smooth, typeid(*this).name()),
618  edge_down(nullptr, typeid(*this).name()),
619  edge_up(nullptr, typeid(*this).name()),
620  debug(0)
621 {
622  const unsigned int dof_handler_max_level
623  = mg_dof_handler.get_triangulation().n_global_levels()-1;
624  if (max_level == numbers::invalid_unsigned_int)
625  maxlevel = dof_handler_max_level;
626  else
627  maxlevel = max_level;
628 
629  reinit(minlevel, maxlevel);
630 }
631 
632 
633 
634 template <typename VectorType>
636  const MGCoarseGridBase<VectorType> &coarse,
637  const MGTransferBase<VectorType> &transfer,
638  const MGSmootherBase<VectorType> &pre_smooth,
639  const MGSmootherBase<VectorType> &post_smooth,
640  const unsigned int min_level,
641  const unsigned int max_level,
642  Cycle cycle)
643  :
644  cycle_type(cycle),
645  matrix(&matrix, typeid(*this).name()),
646  coarse(&coarse, typeid(*this).name()),
647  transfer(&transfer, typeid(*this).name()),
648  pre_smooth(&pre_smooth, typeid(*this).name()),
649  post_smooth(&post_smooth, typeid(*this).name()),
650  edge_out(nullptr, typeid(*this).name()),
651  edge_in(nullptr, typeid(*this).name()),
652  edge_down(nullptr, typeid(*this).name()),
653  edge_up(nullptr, typeid(*this).name()),
654  debug(0)
655 {
656  if (max_level == numbers::invalid_unsigned_int)
657  maxlevel = matrix.get_maxlevel();
658  else
659  maxlevel = max_level;
660  reinit(min_level, maxlevel);
661 }
662 
663 
664 
665 template <typename VectorType>
666 inline
667 unsigned int
669 {
670  return maxlevel;
671 }
672 
673 
674 
675 template <typename VectorType>
676 inline
677 unsigned int
679 {
680  return minlevel;
681 }
682 
683 
684 /* --------------------------- inline functions --------------------- */
685 
686 
687 namespace internal
688 {
689  namespace PreconditionMGImplementation
690  {
691  template <int dim, typename VectorType, class TRANSFER, typename OtherVectorType>
692  typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
693  vmult(const std::vector<const ::DoFHandler<dim>*> &dof_handler_vector,
694  ::Multigrid<VectorType> &multigrid,
695  const TRANSFER &transfer,
696  OtherVectorType &dst,
697  const OtherVectorType &src,
698  const bool uses_dof_handler_vector,
699  const typename ::mg::Signals &signals, int)
700  {
701  signals.transfer_to_mg(true);
702  if (uses_dof_handler_vector)
703  transfer.copy_to_mg(dof_handler_vector,
704  multigrid.defect,
705  src);
706  else
707  transfer.copy_to_mg(*dof_handler_vector[0],
708  multigrid.defect,
709  src);
710  signals.transfer_to_mg(false);
711 
712  multigrid.cycle();
713 
714  signals.transfer_to_global(true);
715  if (uses_dof_handler_vector)
716  transfer.copy_from_mg(dof_handler_vector,
717  dst,
718  multigrid.solution);
719  else
720  transfer.copy_from_mg(*dof_handler_vector[0],
721  dst,
722  multigrid.solution);
723  signals.transfer_to_global(false);
724  }
725 
726  template <int dim, typename VectorType, class TRANSFER, typename OtherVectorType>
727  void
728  vmult(const std::vector<const ::DoFHandler<dim>*> &dof_handler_vector,
729  ::Multigrid<VectorType> &multigrid,
730  const TRANSFER &transfer,
731  OtherVectorType &dst,
732  const OtherVectorType &src,
733  const bool uses_dof_handler_vector,
734  const typename ::mg::Signals &signals, ...)
735  {
736  (void) uses_dof_handler_vector;
737  Assert (!uses_dof_handler_vector, ExcInternalError());
738 
739  signals.transfer_to_mg(true);
740  transfer.copy_to_mg(*dof_handler_vector[0],
741  multigrid.defect,
742  src);
743  signals.transfer_to_mg(false);
744 
745  multigrid.cycle();
746 
747  signals.transfer_to_global(true);
748  transfer.copy_from_mg(*dof_handler_vector[0],
749  dst,
750  multigrid.solution);
751  signals.transfer_to_global(false);
752  }
753 
754  template <int dim, typename VectorType, class TRANSFER, typename OtherVectorType>
755  typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
756  vmult_add(const std::vector<const ::DoFHandler<dim>*> &dof_handler_vector,
757  ::Multigrid<VectorType> &multigrid,
758  const TRANSFER &transfer,
759  OtherVectorType &dst,
760  const OtherVectorType &src,
761  const bool uses_dof_handler_vector,
762  const typename ::mg::Signals &signals, int)
763  {
764  signals.transfer_to_mg(true);
765  if (uses_dof_handler_vector)
766  transfer.copy_to_mg(dof_handler_vector,
767  multigrid.defect,
768  src);
769  else
770  transfer.copy_to_mg(*dof_handler_vector[0],
771  multigrid.defect,
772  src);
773  signals.transfer_to_mg(false);
774 
775  multigrid.cycle();
776 
777  signals.transfer_to_global(true);
778  if (uses_dof_handler_vector)
779  transfer.copy_from_mg_add(dof_handler_vector,
780  dst,
781  multigrid.solution);
782  else
783  transfer.copy_from_mg_add(*dof_handler_vector[0],
784  dst,
785  multigrid.solution);
786  signals.transfer_to_global(false);
787  }
788 
789  template <int dim, typename VectorType, class TRANSFER, typename OtherVectorType>
790  void
791  vmult_add(const std::vector<const ::DoFHandler<dim>*> &dof_handler_vector,
792  ::Multigrid<VectorType> &multigrid,
793  const TRANSFER &transfer,
794  OtherVectorType &dst,
795  const OtherVectorType &src,
796  const bool uses_dof_handler_vector,
797  const typename ::mg::Signals &signals, ...)
798  {
799  (void) uses_dof_handler_vector;
800  Assert (!uses_dof_handler_vector, ExcInternalError());
801 
802  signals.transfer_to_mg(true);
803  transfer.copy_to_mg(*dof_handler_vector[0],
804  multigrid.defect,
805  src);
806  signals.transfer_to_mg(false);
807 
808  multigrid.cycle();
809 
810  signals.transfer_to_global(true);
811  transfer.copy_from_mg_add(*dof_handler_vector[0],
812  dst,
813  multigrid.solution);
814  signals.transfer_to_global(false);
815  }
816  }
817 }
818 
819 template <int dim, typename VectorType, class TRANSFER>
821 ::PreconditionMG(const DoFHandler<dim> &dof_handler,
823  const TRANSFER &transfer)
824  :
825  dof_handler_vector(1,&dof_handler),
826  dof_handler_vector_raw(1,&dof_handler),
827  multigrid(&mg),
828  transfer(&transfer),
829  uses_dof_handler_vector(false)
830 {}
831 
832 template <int dim, typename VectorType, class TRANSFER>
834 ::PreconditionMG(const std::vector<const DoFHandler<dim>*> &dof_handler,
836  const TRANSFER &transfer)
837  :
838  dof_handler_vector(dof_handler.size()),
839  dof_handler_vector_raw(dof_handler.size()),
840  multigrid(&mg),
841  transfer(&transfer),
842  uses_dof_handler_vector(true)
843 {
844  for (unsigned int i = 0; i< dof_handler.size() ; ++i)
845  {
846  dof_handler_vector[i] = dof_handler[i];
847  dof_handler_vector_raw[i] = dof_handler[i];
848  }
849 }
850 
851 template <int dim, typename VectorType, class TRANSFER>
852 inline bool
854 {
855  return false;
856 }
857 
858 template <int dim, typename VectorType, class TRANSFER>
859 template <class OtherVectorType>
860 void
862 (OtherVectorType &dst,
863  const OtherVectorType &src) const
864 {
865  internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,*multigrid,*transfer,
866  dst,src,uses_dof_handler_vector,
867  this->signals, 0);
868 }
869 
870 
871 template <int dim, typename VectorType, class TRANSFER>
872 IndexSet
874 {
875  AssertIndexRange(block,dof_handler_vector.size());
876  return dof_handler_vector[block]->locally_owned_dofs();
877 }
878 
879 
880 template <int dim, typename VectorType, class TRANSFER>
881 IndexSet
883 {
884  AssertIndexRange(block,dof_handler_vector.size());
885  return dof_handler_vector[block]->locally_owned_dofs();
886 }
887 
888 
889 
890 template <int dim, typename VectorType, class TRANSFER>
891 MPI_Comm
893 {
894  // currently parallel GMG works with distributed Triangulation only,
895  // so it should be a safe bet to use it to query MPI communicator:
896  const Triangulation<dim> &tria = dof_handler_vector[0]->get_triangulation();
898  Assert (ptria != nullptr, ExcInternalError());
899  return ptria->get_communicator ();
900 }
901 
902 
903 
904 template <int dim, typename VectorType, class TRANSFER>
905 boost::signals2::connection
907 connect_transfer_to_mg(const std::function<void (bool)> &slot)
908 {
909  return this->signals.transfer_to_mg.connect(slot);
910 }
911 
912 
913 
914 template <int dim, typename VectorType, class TRANSFER>
915 boost::signals2::connection
918  const std::function<void (bool)> &slot)
919 {
920  return this->signals.transfer_to_global.connect(slot);
921 }
922 
923 
924 
925 template <int dim, typename VectorType, class TRANSFER>
926 template <class OtherVectorType>
927 void
929 (OtherVectorType &dst,
930  const OtherVectorType &src) const
931 {
932  internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,*multigrid,*transfer,
933  dst,src,uses_dof_handler_vector,
934  this->signals, 0);
935 }
936 
937 
938 template <int dim, typename VectorType, class TRANSFER>
939 template <class OtherVectorType>
940 void
942 (OtherVectorType &,
943  const OtherVectorType &) const
944 {
945  Assert(false, ExcNotImplemented());
946 }
947 
948 
949 template <int dim, typename VectorType, class TRANSFER>
950 template <class OtherVectorType>
951 void
953 (OtherVectorType &,
954  const OtherVectorType &) const
955 {
956  Assert(false, ExcNotImplemented());
957 }
958 
959 #endif // DOXYGEN
960 
961 DEAL_II_NAMESPACE_CLOSE
962 
963 #endif
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > matrix
Definition: multigrid.h:389
boost::signals2::connection connect_transfer_to_global(const std::function< void(bool)> &slot)
boost::signals2::signal< void(const bool before)> transfer_to_mg
Definition: multigrid.h:54
void vcycle()
void cycle()
SmartPointer< const TRANSFER, PreconditionMG< dim, VectorType, TRANSFER > > transfer
Definition: multigrid.h:579
static const unsigned int invalid_unsigned_int
Definition: types.h:173
SmartPointer< const MGMatrixBase< VectorType > > edge_in
Definition: multigrid.h:424
void vmult(OtherVectorType &dst, const OtherVectorType &src) const
MPI_Comm get_mpi_communicator() const
Contents is actually a matrix.
void set_edge_matrices(const MGMatrixBase< VectorType > &edge_out, const MGMatrixBase< VectorType > &edge_in)
Definition: mg.h:81
unsigned int get_maxlevel() const
boost::signals2::signal< void(const bool before, const unsigned int level)> restriction
Definition: multigrid.h:74
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
void set_maxlevel(const unsigned int)
MGLevelObject< VectorType > t
Definition: multigrid.h:378
boost::signals2::connection connect_coarse_solve(const std::function< void(const bool, const unsigned int)> &slot)
bool empty() const
SmartPointer< const MGCoarseGridBase< VectorType >, Multigrid< VectorType > > coarse
Definition: multigrid.h:394
boost::signals2::signal< void(const bool before, const unsigned int level)> coarse_solve
Definition: multigrid.h:67
const bool uses_dof_handler_vector
Definition: multigrid.h:585
unsigned int debug
Definition: multigrid.h:443
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_up
Definition: multigrid.h:438
void set_debug(const unsigned int)
void level_v_step(const unsigned int level)
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > post_smooth
Definition: multigrid.h:409
mg::Signals signals
Definition: multigrid.h:328
unsigned int maxlevel
Definition: multigrid.h:360
IndexSet locally_owned_range_indices(const unsigned int block=0) const
#define Assert(cond, exc)
Definition: exceptions.h:1142
std::vector< SmartPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TRANSFER > > > dof_handler_vector
Definition: multigrid.h:563
boost::signals2::signal< void(const bool before)> transfer_to_global
Definition: multigrid.h:61
The F-cycle.
Definition: multigrid.h:137
virtual MPI_Comm get_communicator() const
Definition: tria_base.cc:146
void Tvmult(OtherVectorType &dst, const OtherVectorType &src) const
SmartPointer< Multigrid< VectorType >, PreconditionMG< dim, VectorType, TRANSFER > > multigrid
Definition: multigrid.h:574
boost::signals2::connection connect_restriction(const std::function< void(const bool, const unsigned int)> &slot)
unsigned int get_minlevel() const
SmartPointer< const MGMatrixBase< VectorType > > edge_out
Definition: multigrid.h:416
void vmult_add(OtherVectorType &dst, const OtherVectorType &src) const
The W-cycle.
Definition: multigrid.h:135
SmartPointer< const MGTransferBase< VectorType >, Multigrid< VectorType > > transfer
Definition: multigrid.h:399
void set_minlevel(const unsigned int level, bool relative=false)
The V-cycle.
Definition: multigrid.h:133
void set_edge_flux_matrices(const MGMatrixBase< VectorType > &edge_down, const MGMatrixBase< VectorType > &edge_up)
boost::signals2::connection connect_transfer_to_mg(const std::function< void(bool)> &slot)
void Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const
unsigned int minlevel
Definition: multigrid.h:355
MGLevelObject< VectorType > solution
Definition: multigrid.h:372
void set_cycle(Cycle)
boost::signals2::connection connect_pre_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
boost::signals2::signal< void(const bool before, const unsigned int level)> prolongation
Definition: multigrid.h:81
boost::signals2::connection connect_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
std::vector< const DoFHandler< dim > * > dof_handler_vector_raw
Definition: multigrid.h:569
boost::signals2::signal< void(const bool before, const unsigned int level)> post_smoother_step
Definition: multigrid.h:95
IndexSet locally_owned_domain_indices(const unsigned int block=0) const
void level_step(const unsigned int level, Cycle cycle)
virtual unsigned int n_global_levels() const
boost::signals2::connection connect_post_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
mg::Signals signals
Definition: multigrid.h:590
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > pre_smooth
Definition: multigrid.h:404
static ::ExceptionBase & ExcNotImplemented()
const Triangulation< dim, spacedim > & get_triangulation() const
Multigrid(const DoFHandler< dim > &mg_dof_handler, const MGMatrixBase< VectorType > &matrix, const MGCoarseGridBase< VectorType > &coarse, const MGTransferBase< VectorType > &transfer, const MGSmootherBase< VectorType > &pre_smooth, const MGSmootherBase< VectorType > &post_smooth, const unsigned int minlevel=0, const unsigned int maxlevel=numbers::invalid_unsigned_int, Cycle cycle=v_cycle)
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_down
Definition: multigrid.h:431
MGLevelObject< VectorType > defect2
Definition: multigrid.h:383
MGLevelObject< VectorType > defect
Definition: multigrid.h:367
PreconditionMG(const DoFHandler< dim > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer)
boost::signals2::signal< void(const bool before, const unsigned int level)> pre_smoother_step
Definition: multigrid.h:88
Cycle cycle_type
Definition: multigrid.h:350
static ::ExceptionBase & ExcInternalError()
Triangulation< dim, spacedim > & get_triangulation()
Definition: tria.cc:11695
void reinit(const unsigned int minlevel, const unsigned int maxlevel)