Reference documentation for deal.II version GIT b8135fa6eb 2023-03-25 15:55:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
multigrid.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2022 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_multigrid_h
17 #define dealii_multigrid_h
18 
19 
20 #include <deal.II/base/config.h>
21 
25 
27 
29 
31 #include <deal.II/lac/vector.h>
32 
34 
35 #include <vector>
36 
38 
39 #ifdef signals
40 # error \
41  "The name 'signals' is already defined. You are most likely using the QT library \
42 and using the 'signals' keyword. You can either #include the Qt headers (or any conflicting headers) \
43 *after* the deal.II headers or you can define the 'QT_NO_KEYWORDS' macro and use the 'Q_SIGNALS' macro."
44 #endif
45 
51 namespace mg
52 {
62  struct Signals
63  {
69  boost::signals2::signal<void(const bool before)> transfer_to_mg;
70 
76  boost::signals2::signal<void(const bool before)> transfer_to_global;
77 
87  boost::signals2::signal<void(const bool before, const unsigned int level)>
89 
98  boost::signals2::signal<void(const bool before, const unsigned int level)>
100 
106  boost::signals2::signal<void(const bool before, const unsigned int level)>
108 
117  boost::signals2::signal<void(const bool before, const unsigned int level)>
119 
125  boost::signals2::signal<void(const bool before, const unsigned int level)>
127 
133  boost::signals2::signal<void(const bool before, const unsigned int level)>
135 
140  boost::signals2::signal<void(const bool before, const unsigned int level)>
142  };
143 } // namespace mg
144 
162 template <typename VectorType>
163 class Multigrid : public Subscriptor
164 {
165 public:
169  enum Cycle
170  {
176  f_cycle
177  };
178 
179  using vector_type = VectorType;
180  using const_vector_type = const VectorType;
181 
196  const unsigned int minlevel = 0,
197  const unsigned int maxlevel = numbers::invalid_unsigned_int,
198  Cycle cycle = v_cycle);
199 
203  void
204  reinit(const unsigned int minlevel, const unsigned int maxlevel);
205 
210  void
211  cycle();
212 
223  void
225 
238  void
241 
248  void
250 
263  void
266 
270  unsigned int
271  get_maxlevel() const;
272 
276  unsigned int
277  get_minlevel() const;
278 
284  void
285  set_maxlevel(const unsigned int);
286 
302  void
303  set_minlevel(const unsigned int level, bool relative = false);
304 
309 
313  boost::signals2::connection
315  const std::function<void(const bool, const unsigned int)> &slot);
316 
320  boost::signals2::connection
322  const std::function<void(const bool, const unsigned int)> &slot);
323 
327  boost::signals2::connection
329  const std::function<void(const bool, const unsigned int)> &slot);
330 
334  boost::signals2::connection
336  const std::function<void(const bool, const unsigned int)> &slot);
337 
341  boost::signals2::connection
343  const std::function<void(const bool, const unsigned int)> &slot);
344 
348  boost::signals2::connection
350  const std::function<void(const bool, const unsigned int)> &slot);
351 
355  boost::signals2::connection
357  const std::function<void(const bool, const unsigned int)> &slot);
358 
359 private:
364 
371  void
372  level_v_step(const unsigned int level);
373 
381  void
382  level_step(const unsigned int level, Cycle cycle);
383 
388 
392  unsigned int minlevel;
393 
397  unsigned int maxlevel;
398 
399 public:
405 
410 
411 private:
416 
421 
422 
427 
433 
439 
445 
451 
458 
466 
473 
480 
481  template <int dim, class OtherVectorType, class TRANSFER>
482  friend class PreconditionMG;
483 };
484 
485 
500 template <int dim, typename VectorType, class TRANSFER>
502 {
503 public:
508  PreconditionMG(const DoFHandler<dim> &dof_handler,
510  const TRANSFER & transfer);
511 
516  PreconditionMG(const std::vector<const DoFHandler<dim> *> &dof_handler,
518  const TRANSFER & transfer);
519 
523  bool
524  empty() const;
525 
532  template <class OtherVectorType>
533  void
534  vmult(OtherVectorType &dst, const OtherVectorType &src) const;
535 
540  template <class OtherVectorType>
541  void
542  vmult_add(OtherVectorType &dst, const OtherVectorType &src) const;
543 
549  template <class OtherVectorType>
550  void
551  Tvmult(OtherVectorType &dst, const OtherVectorType &src) const;
552 
558  template <class OtherVectorType>
559  void
560  Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const;
561 
568  IndexSet
569  locally_owned_range_indices(const unsigned int block = 0) const;
570 
577  IndexSet
578  locally_owned_domain_indices(const unsigned int block = 0) const;
579 
583  MPI_Comm
585 
589  boost::signals2::connection
590  connect_transfer_to_mg(const std::function<void(bool)> &slot);
591 
595  boost::signals2::connection
596  connect_transfer_to_global(const std::function<void(bool)> &slot);
597 
603 
607  const Multigrid<VectorType> &
608  get_multigrid() const;
609 
610 private:
614  std::vector<SmartPointer<const DoFHandler<dim>,
617 
622  std::vector<const DoFHandler<dim> *> dof_handler_vector_raw;
623 
629 
635 
641 
646 };
647 
650 #ifndef DOXYGEN
651 /* --------------------------- inline functions --------------------- */
652 
653 
654 template <typename VectorType>
656  const MGCoarseGridBase<VectorType> &coarse,
657  const MGTransferBase<VectorType> & transfer,
658  const MGSmootherBase<VectorType> & pre_smooth,
659  const MGSmootherBase<VectorType> &post_smooth,
660  const unsigned int min_level,
661  const unsigned int max_level,
662  Cycle cycle)
663  : cycle_type(cycle)
664  , matrix(&matrix, typeid(*this).name())
665  , coarse(&coarse, typeid(*this).name())
666  , transfer(&transfer, typeid(*this).name())
667  , pre_smooth(&pre_smooth, typeid(*this).name())
668  , post_smooth(&post_smooth, typeid(*this).name())
669  , edge_out(nullptr, typeid(*this).name())
670  , edge_in(nullptr, typeid(*this).name())
671  , edge_down(nullptr, typeid(*this).name())
672  , edge_up(nullptr, typeid(*this).name())
673 {
674  if (max_level == numbers::invalid_unsigned_int)
675  maxlevel = matrix.get_maxlevel();
676  else
677  maxlevel = max_level;
678  reinit(min_level, maxlevel);
679 }
680 
681 
682 
683 template <typename VectorType>
684 inline unsigned int
686 {
687  return maxlevel;
688 }
689 
690 
691 
692 template <typename VectorType>
693 inline unsigned int
695 {
696  return minlevel;
697 }
698 
699 
700 /* --------------------------- inline functions --------------------- */
701 
702 
703 namespace internal
704 {
705  namespace PreconditionMGImplementation
706  {
707  template <int dim,
708  typename VectorType,
709  class TRANSFER,
710  typename OtherVectorType>
711  std::enable_if_t<TRANSFER::supports_dof_handler_vector>
712  vmult(const std::vector<const DoFHandler<dim> *> &dof_handler_vector,
713  Multigrid<VectorType> & multigrid,
714  const TRANSFER & transfer,
715  OtherVectorType & dst,
716  const OtherVectorType & src,
717  const bool uses_dof_handler_vector,
718  const typename ::mg::Signals & signals,
719  int)
720  {
721  signals.transfer_to_mg(true);
722  if (uses_dof_handler_vector)
723  transfer.copy_to_mg(dof_handler_vector, multigrid.defect, src);
724  else
725  transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
726  signals.transfer_to_mg(false);
727 
728  multigrid.cycle();
729 
730  signals.transfer_to_global(true);
731  if (uses_dof_handler_vector)
732  transfer.copy_from_mg(dof_handler_vector, dst, multigrid.solution);
733  else
734  transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.solution);
735  signals.transfer_to_global(false);
736  }
737 
738  template <int dim,
739  typename VectorType,
740  class TRANSFER,
741  typename OtherVectorType>
742  void
743  vmult(const std::vector<const DoFHandler<dim> *> &dof_handler_vector,
744  Multigrid<VectorType> & multigrid,
745  const TRANSFER & transfer,
746  OtherVectorType & dst,
747  const OtherVectorType & src,
748  const bool uses_dof_handler_vector,
749  const typename ::mg::Signals & signals,
750  ...)
751  {
752  (void)uses_dof_handler_vector;
753  Assert(!uses_dof_handler_vector, ExcInternalError());
754 
755  signals.transfer_to_mg(true);
756  transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
757  signals.transfer_to_mg(false);
758 
759  multigrid.cycle();
760 
761  signals.transfer_to_global(true);
762  transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.solution);
763  signals.transfer_to_global(false);
764  }
765 
766  template <int dim,
767  typename VectorType,
768  class TRANSFER,
769  typename OtherVectorType>
770  std::enable_if_t<TRANSFER::supports_dof_handler_vector>
771  vmult_add(const std::vector<const DoFHandler<dim> *> &dof_handler_vector,
772  Multigrid<VectorType> & multigrid,
773  const TRANSFER & transfer,
774  OtherVectorType & dst,
775  const OtherVectorType & src,
776  const bool uses_dof_handler_vector,
777  const typename ::mg::Signals &signals,
778  int)
779  {
780  signals.transfer_to_mg(true);
781  if (uses_dof_handler_vector)
782  transfer.copy_to_mg(dof_handler_vector, multigrid.defect, src);
783  else
784  transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
785  signals.transfer_to_mg(false);
786 
787  multigrid.cycle();
788 
789  signals.transfer_to_global(true);
790  if (uses_dof_handler_vector)
791  transfer.copy_from_mg_add(dof_handler_vector, dst, multigrid.solution);
792  else
793  transfer.copy_from_mg_add(*dof_handler_vector[0],
794  dst,
795  multigrid.solution);
796  signals.transfer_to_global(false);
797  }
798 
799  template <int dim,
800  typename VectorType,
801  class TRANSFER,
802  typename OtherVectorType>
803  void
804  vmult_add(const std::vector<const DoFHandler<dim> *> &dof_handler_vector,
805  Multigrid<VectorType> & multigrid,
806  const TRANSFER & transfer,
807  OtherVectorType & dst,
808  const OtherVectorType & src,
809  const bool uses_dof_handler_vector,
810  const typename ::mg::Signals &signals,
811  ...)
812  {
813  (void)uses_dof_handler_vector;
814  Assert(!uses_dof_handler_vector, ExcInternalError());
815 
816  signals.transfer_to_mg(true);
817  transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
818  signals.transfer_to_mg(false);
819 
820  multigrid.cycle();
821 
822  signals.transfer_to_global(true);
823  transfer.copy_from_mg_add(*dof_handler_vector[0],
824  dst,
825  multigrid.solution);
826  signals.transfer_to_global(false);
827  }
828  } // namespace PreconditionMGImplementation
829 } // namespace internal
830 
831 template <int dim, typename VectorType, class TRANSFER>
833  const DoFHandler<dim> &dof_handler,
835  const TRANSFER & transfer)
836  : dof_handler_vector(1, &dof_handler)
837  , dof_handler_vector_raw(1, &dof_handler)
838  , multigrid(&mg)
839  , transfer(&transfer)
840  , uses_dof_handler_vector(false)
841 {}
842 
843 template <int dim, typename VectorType, class TRANSFER>
845  const std::vector<const DoFHandler<dim> *> &dof_handler,
847  const TRANSFER & transfer)
848  : dof_handler_vector(dof_handler.size())
849  , dof_handler_vector_raw(dof_handler.size())
850  , multigrid(&mg)
851  , transfer(&transfer)
852  , uses_dof_handler_vector(true)
853 {
854  for (unsigned int i = 0; i < dof_handler.size(); ++i)
855  {
856  dof_handler_vector[i] = dof_handler[i];
857  dof_handler_vector_raw[i] = dof_handler[i];
858  }
859 }
860 
861 template <int dim, typename VectorType, class TRANSFER>
862 inline bool
864 {
865  return false;
866 }
867 
868 template <int dim, typename VectorType, class TRANSFER>
869 template <class OtherVectorType>
870 void
872  OtherVectorType & dst,
873  const OtherVectorType &src) const
874 {
875  internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,
876  *multigrid,
877  *transfer,
878  dst,
879  src,
880  uses_dof_handler_vector,
881  this->signals,
882  0);
883 }
884 
885 
886 template <int dim, typename VectorType, class TRANSFER>
887 IndexSet
889  const unsigned int block) const
890 {
891  AssertIndexRange(block, dof_handler_vector.size());
892  return dof_handler_vector[block]->locally_owned_dofs();
893 }
894 
895 
896 template <int dim, typename VectorType, class TRANSFER>
897 IndexSet
899  const unsigned int block) const
900 {
901  AssertIndexRange(block, dof_handler_vector.size());
902  return dof_handler_vector[block]->locally_owned_dofs();
903 }
904 
905 
906 
907 template <int dim, typename VectorType, class TRANSFER>
908 MPI_Comm
910 {
911  // currently parallel GMG works with parallel triangulations only,
912  // so it should be a safe bet to use it to query MPI communicator:
913  const Triangulation<dim> &tria = dof_handler_vector[0]->get_triangulation();
914  const parallel::TriangulationBase<dim> *ptria =
915  dynamic_cast<const parallel::TriangulationBase<dim> *>(&tria);
916  Assert(ptria != nullptr, ExcInternalError());
917  return ptria->get_communicator();
918 }
919 
920 
921 
922 template <int dim, typename VectorType, class TRANSFER>
923 boost::signals2::connection
925  const std::function<void(bool)> &slot)
926 {
927  return this->signals.transfer_to_mg.connect(slot);
928 }
929 
930 
931 
932 template <int dim, typename VectorType, class TRANSFER>
933 boost::signals2::connection
935  const std::function<void(bool)> &slot)
936 {
937  return this->signals.transfer_to_global.connect(slot);
938 }
939 
940 
941 
942 template <int dim, typename VectorType, class TRANSFER>
943 template <class OtherVectorType>
944 void
946  OtherVectorType & dst,
947  const OtherVectorType &src) const
948 {
949  internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,
950  *multigrid,
951  *transfer,
952  dst,
953  src,
954  uses_dof_handler_vector,
955  this->signals,
956  0);
957 }
958 
959 
960 template <int dim, typename VectorType, class TRANSFER>
961 template <class OtherVectorType>
962 void
964  const OtherVectorType &) const
965 {
966  Assert(false, ExcNotImplemented());
967 }
968 
969 
970 template <int dim, typename VectorType, class TRANSFER>
971 template <class OtherVectorType>
972 void
974  OtherVectorType &,
975  const OtherVectorType &) const
976 {
977  Assert(false, ExcNotImplemented());
978 }
979 
980 
981 template <int dim, typename VectorType, class TRANSFER>
984 {
985  return *this->multigrid;
986 }
987 
988 
989 template <int dim, typename VectorType, class TRANSFER>
990 const Multigrid<VectorType> &
992 {
993  return *this->multigrid;
994 }
995 
996 #endif // DOXYGEN
997 
999 
1000 #endif
@ v_cycle
The V-cycle.
Definition: multigrid.h:172
@ w_cycle
The W-cycle.
Definition: multigrid.h:174
@ f_cycle
The F-cycle.
Definition: multigrid.h:176
MGLevelObject< VectorType > defect2
Definition: multigrid.h:420
void cycle()
void set_edge_flux_matrices(const MGMatrixBase< VectorType > &edge_down, const MGMatrixBase< VectorType > &edge_up)
SmartPointer< const MGTransferBase< VectorType >, Multigrid< VectorType > > transfer
Definition: multigrid.h:438
unsigned int minlevel
Definition: multigrid.h:392
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_up
Definition: multigrid.h:479
mg::Signals signals
Definition: multigrid.h:363
SmartPointer< const MGMatrixBase< VectorType > > edge_in
Definition: multigrid.h:465
SmartPointer< const MGMatrixBase< VectorType > > edge_out
Definition: multigrid.h:457
MGLevelObject< VectorType > solution
Definition: multigrid.h:409
boost::signals2::connection connect_edge_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
unsigned int maxlevel
Definition: multigrid.h:397
VectorType vector_type
Definition: multigrid.h:179
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > matrix
Definition: multigrid.h:426
void set_maxlevel(const unsigned int)
unsigned int get_minlevel() const
void set_edge_matrices(const MGMatrixBase< VectorType > &edge_out, const MGMatrixBase< VectorType > &edge_in)
boost::signals2::connection connect_residual_step(const std::function< void(const bool, const unsigned int)> &slot)
void level_step(const unsigned int level, Cycle cycle)
void set_cycle(Cycle)
unsigned int get_maxlevel() const
const VectorType const_vector_type
Definition: multigrid.h:180
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > pre_smooth
Definition: multigrid.h:444
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > post_smooth
Definition: multigrid.h:450
boost::signals2::connection connect_post_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
boost::signals2::connection connect_coarse_solve(const std::function< void(const bool, const unsigned int)> &slot)
SmartPointer< const MGCoarseGridBase< VectorType >, Multigrid< VectorType > > coarse
Definition: multigrid.h:432
Cycle cycle_type
Definition: multigrid.h:387
boost::signals2::connection connect_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_down
Definition: multigrid.h:472
void set_edge_in_matrix(const MGMatrixBase< VectorType > &edge_in)
MGLevelObject< VectorType > t
Definition: multigrid.h:415
void set_minlevel(const unsigned int level, bool relative=false)
void level_v_step(const unsigned int level)
Multigrid(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)
MGLevelObject< VectorType > defect
Definition: multigrid.h:404
void vcycle()
boost::signals2::connection connect_restriction(const std::function< void(const bool, const unsigned int)> &slot)
void reinit(const unsigned int minlevel, const unsigned int maxlevel)
boost::signals2::connection connect_pre_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
SmartPointer< Multigrid< VectorType >, PreconditionMG< dim, VectorType, TRANSFER > > multigrid
Definition: multigrid.h:628
void vmult_add(OtherVectorType &dst, const OtherVectorType &src) const
void vmult(OtherVectorType &dst, const OtherVectorType &src) const
IndexSet locally_owned_range_indices(const unsigned int block=0) const
SmartPointer< const TRANSFER, PreconditionMG< dim, VectorType, TRANSFER > > transfer
Definition: multigrid.h:634
PreconditionMG(const std::vector< const DoFHandler< dim > * > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer)
std::vector< SmartPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TRANSFER > > > dof_handler_vector
Definition: multigrid.h:616
bool empty() const
void Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const
void Tvmult(OtherVectorType &dst, const OtherVectorType &src) const
boost::signals2::connection connect_transfer_to_mg(const std::function< void(bool)> &slot)
const Multigrid< VectorType > & get_multigrid() const
PreconditionMG(const DoFHandler< dim > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer)
mg::Signals signals
Definition: multigrid.h:645
IndexSet locally_owned_domain_indices(const unsigned int block=0) const
boost::signals2::connection connect_transfer_to_global(const std::function< void(bool)> &slot)
const bool uses_dof_handler_vector
Definition: multigrid.h:640
MPI_Comm get_mpi_communicator() const
Multigrid< VectorType > & get_multigrid()
std::vector< const DoFHandler< dim > * > dof_handler_vector_raw
Definition: multigrid.h:622
Triangulation< dim, spacedim > & get_triangulation()
virtual MPI_Comm get_communicator() const override
Definition: tria_base.cc:143
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
unsigned int level
Definition: grid_out.cc:4618
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
Definition: exceptions.h:1586
static ::ExceptionBase & ExcNotImplemented()
#define AssertIndexRange(index, range)
Definition: exceptions.h:1827
@ matrix
Contents is actually a matrix.
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618
Definition: mg.h:82
static const unsigned int invalid_unsigned_int
Definition: types.h:213
boost::signals2::signal< void(const bool before, const unsigned int level)> prolongation
Definition: multigrid.h:107
boost::signals2::signal< void(const bool before, const unsigned int level)> pre_smoother_step
Definition: multigrid.h:118
boost::signals2::signal< void(const bool before, const unsigned int level)> coarse_solve
Definition: multigrid.h:88
boost::signals2::signal< void(const bool before, const unsigned int level)> edge_prolongation
Definition: multigrid.h:141
boost::signals2::signal< void(const bool before, const unsigned int level)> post_smoother_step
Definition: multigrid.h:126
boost::signals2::signal< void(const bool before)> transfer_to_global
Definition: multigrid.h:76
boost::signals2::signal< void(const bool before, const unsigned int level)> residual_step
Definition: multigrid.h:134
boost::signals2::signal< void(const bool before, const unsigned int level)> restriction
Definition: multigrid.h:99
boost::signals2::signal< void(const bool before)> transfer_to_mg
Definition: multigrid.h:69
const ::Triangulation< dim, spacedim > & tria