Reference documentation for deal.II version 9.2.0
\(\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 - 2020 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 
41 
42 namespace mg
43 {
53  struct Signals
54  {
60  boost::signals2::signal<void(const bool before)> transfer_to_mg;
61 
67  boost::signals2::signal<void(const bool before)> transfer_to_global;
68 
74  boost::signals2::signal<void(const bool before, const unsigned int level)>
76 
82  boost::signals2::signal<void(const bool before, const unsigned int level)>
84 
90  boost::signals2::signal<void(const bool before, const unsigned int level)>
92 
98  boost::signals2::signal<void(const bool before, const unsigned int level)>
100 
107  boost::signals2::signal<void(const bool before, const unsigned int level)>
109  };
110 } // namespace mg
111 
136 template <typename VectorType>
137 class Multigrid : public Subscriptor
138 {
139 public:
143  enum Cycle
144  {
151  };
152 
155 
170  const unsigned int minlevel = 0,
171  const unsigned int maxlevel = numbers::invalid_unsigned_int,
172  Cycle cycle = v_cycle);
173 
177  void
178  reinit(const unsigned int minlevel, const unsigned int maxlevel);
179 
184  void
185  cycle();
186 
197  void
198  vcycle();
199 
212  void
215 
228  void
231 
235  unsigned int
236  get_maxlevel() const;
237 
241  unsigned int
242  get_minlevel() const;
243 
249  void
250  set_maxlevel(const unsigned int);
251 
267  void
268  set_minlevel(const unsigned int level, bool relative = false);
269 
273  void set_cycle(Cycle);
274 
278  boost::signals2::connection
280  const std::function<void(const bool, const unsigned int)> &slot);
281 
285  boost::signals2::connection
287  const std::function<void(const bool, const unsigned int)> &slot);
288 
292  boost::signals2::connection
294  const std::function<void(const bool, const unsigned int)> &slot);
295 
299  boost::signals2::connection
301  const std::function<void(const bool, const unsigned int)> &slot);
302 
306  boost::signals2::connection
308  const std::function<void(const bool, const unsigned int)> &slot);
309 
310 private:
315 
322  void
323  level_v_step(const unsigned int level);
324 
332  void
333  level_step(const unsigned int level, Cycle cycle);
334 
339 
343  unsigned int minlevel;
344 
348  unsigned int maxlevel;
349 
350 public:
356 
361 
362 private:
367 
372 
373 
378 
384 
390 
396 
402 
409 
417 
424 
431 
432  template <int dim, class OtherVectorType, class TRANSFER>
433  friend class PreconditionMG;
434 };
435 
436 
453 template <int dim, typename VectorType, class TRANSFER>
455 {
456 public:
461  PreconditionMG(const DoFHandler<dim> &dof_handler,
463  const TRANSFER & transfer);
464 
469  PreconditionMG(const std::vector<const DoFHandler<dim> *> &dof_handler,
471  const TRANSFER & transfer);
472 
476  bool
477  empty() const;
478 
485  template <class OtherVectorType>
486  void
487  vmult(OtherVectorType &dst, const OtherVectorType &src) const;
488 
493  template <class OtherVectorType>
494  void
495  vmult_add(OtherVectorType &dst, const OtherVectorType &src) const;
496 
502  template <class OtherVectorType>
503  void
504  Tvmult(OtherVectorType &dst, const OtherVectorType &src) const;
505 
511  template <class OtherVectorType>
512  void
513  Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const;
514 
521  IndexSet
522  locally_owned_range_indices(const unsigned int block = 0) const;
523 
530  IndexSet
531  locally_owned_domain_indices(const unsigned int block = 0) const;
532 
536  const MPI_Comm &
537  get_mpi_communicator() const;
538 
542  boost::signals2::connection
543  connect_transfer_to_mg(const std::function<void(bool)> &slot);
544 
548  boost::signals2::connection
549  connect_transfer_to_global(const std::function<void(bool)> &slot);
550 
551 private:
555  std::vector<SmartPointer<const DoFHandler<dim>,
558 
563  std::vector<const DoFHandler<dim> *> dof_handler_vector_raw;
564 
570 
576 
582 
587 };
588 
591 #ifndef DOXYGEN
592 /* --------------------------- inline functions --------------------- */
593 
594 
595 template <typename VectorType>
597  const MGCoarseGridBase<VectorType> &coarse,
598  const MGTransferBase<VectorType> & transfer,
599  const MGSmootherBase<VectorType> & pre_smooth,
600  const MGSmootherBase<VectorType> &post_smooth,
601  const unsigned int min_level,
602  const unsigned int max_level,
603  Cycle cycle)
604  : cycle_type(cycle)
605  , matrix(&matrix, typeid(*this).name())
606  , coarse(&coarse, typeid(*this).name())
607  , transfer(&transfer, typeid(*this).name())
608  , pre_smooth(&pre_smooth, typeid(*this).name())
609  , post_smooth(&post_smooth, typeid(*this).name())
610  , edge_out(nullptr, typeid(*this).name())
611  , edge_in(nullptr, typeid(*this).name())
612  , edge_down(nullptr, typeid(*this).name())
613  , edge_up(nullptr, typeid(*this).name())
614 {
615  if (max_level == numbers::invalid_unsigned_int)
616  maxlevel = matrix.get_maxlevel();
617  else
618  maxlevel = max_level;
619  reinit(min_level, maxlevel);
620 }
621 
622 
623 
624 template <typename VectorType>
625 inline unsigned int
627 {
628  return maxlevel;
629 }
630 
631 
632 
633 template <typename VectorType>
634 inline unsigned int
636 {
637  return minlevel;
638 }
639 
640 
641 /* --------------------------- inline functions --------------------- */
642 
643 
644 namespace internal
645 {
646  namespace PreconditionMGImplementation
647  {
648  template <int dim,
649  typename VectorType,
650  class TRANSFER,
651  typename OtherVectorType>
652  typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
653  vmult(
654  const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
655  ::Multigrid<VectorType> & multigrid,
656  const TRANSFER & transfer,
657  OtherVectorType & dst,
658  const OtherVectorType & src,
659  const bool uses_dof_handler_vector,
660  const typename ::mg::Signals &signals,
661  int)
662  {
663  signals.transfer_to_mg(true);
664  if (uses_dof_handler_vector)
665  transfer.copy_to_mg(dof_handler_vector, multigrid.defect, src);
666  else
667  transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
668  signals.transfer_to_mg(false);
669 
670  multigrid.cycle();
671 
672  signals.transfer_to_global(true);
673  if (uses_dof_handler_vector)
674  transfer.copy_from_mg(dof_handler_vector, dst, multigrid.solution);
675  else
676  transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.solution);
677  signals.transfer_to_global(false);
678  }
679 
680  template <int dim,
681  typename VectorType,
682  class TRANSFER,
683  typename OtherVectorType>
684  void
685  vmult(
686  const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
687  ::Multigrid<VectorType> & multigrid,
688  const TRANSFER & transfer,
689  OtherVectorType & dst,
690  const OtherVectorType & src,
691  const bool uses_dof_handler_vector,
692  const typename ::mg::Signals &signals,
693  ...)
694  {
695  (void)uses_dof_handler_vector;
696  Assert(!uses_dof_handler_vector, ExcInternalError());
697 
698  signals.transfer_to_mg(true);
699  transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
700  signals.transfer_to_mg(false);
701 
702  multigrid.cycle();
703 
704  signals.transfer_to_global(true);
705  transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.solution);
706  signals.transfer_to_global(false);
707  }
708 
709  template <int dim,
710  typename VectorType,
711  class TRANSFER,
712  typename OtherVectorType>
713  typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
714  vmult_add(
715  const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
716  ::Multigrid<VectorType> & multigrid,
717  const TRANSFER & transfer,
718  OtherVectorType & dst,
719  const OtherVectorType & src,
720  const bool uses_dof_handler_vector,
721  const typename ::mg::Signals &signals,
722  int)
723  {
724  signals.transfer_to_mg(true);
725  if (uses_dof_handler_vector)
726  transfer.copy_to_mg(dof_handler_vector, multigrid.defect, src);
727  else
728  transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
729  signals.transfer_to_mg(false);
730 
731  multigrid.cycle();
732 
733  signals.transfer_to_global(true);
734  if (uses_dof_handler_vector)
735  transfer.copy_from_mg_add(dof_handler_vector, dst, multigrid.solution);
736  else
737  transfer.copy_from_mg_add(*dof_handler_vector[0],
738  dst,
739  multigrid.solution);
740  signals.transfer_to_global(false);
741  }
742 
743  template <int dim,
744  typename VectorType,
745  class TRANSFER,
746  typename OtherVectorType>
747  void
748  vmult_add(
749  const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
750  ::Multigrid<VectorType> & multigrid,
751  const TRANSFER & transfer,
752  OtherVectorType & dst,
753  const OtherVectorType & src,
754  const bool uses_dof_handler_vector,
755  const typename ::mg::Signals &signals,
756  ...)
757  {
758  (void)uses_dof_handler_vector;
759  Assert(!uses_dof_handler_vector, ExcInternalError());
760 
761  signals.transfer_to_mg(true);
762  transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
763  signals.transfer_to_mg(false);
764 
765  multigrid.cycle();
766 
767  signals.transfer_to_global(true);
768  transfer.copy_from_mg_add(*dof_handler_vector[0],
769  dst,
770  multigrid.solution);
771  signals.transfer_to_global(false);
772  }
773  } // namespace PreconditionMGImplementation
774 } // namespace internal
775 
776 template <int dim, typename VectorType, class TRANSFER>
778  const DoFHandler<dim> &dof_handler,
780  const TRANSFER & transfer)
781  : dof_handler_vector(1, &dof_handler)
782  , dof_handler_vector_raw(1, &dof_handler)
783  , multigrid(&mg)
784  , transfer(&transfer)
785  , uses_dof_handler_vector(false)
786 {}
787 
788 template <int dim, typename VectorType, class TRANSFER>
790  const std::vector<const DoFHandler<dim> *> &dof_handler,
792  const TRANSFER & transfer)
793  : dof_handler_vector(dof_handler.size())
794  , dof_handler_vector_raw(dof_handler.size())
795  , multigrid(&mg)
796  , transfer(&transfer)
797  , uses_dof_handler_vector(true)
798 {
799  for (unsigned int i = 0; i < dof_handler.size(); ++i)
800  {
801  dof_handler_vector[i] = dof_handler[i];
802  dof_handler_vector_raw[i] = dof_handler[i];
803  }
804 }
805 
806 template <int dim, typename VectorType, class TRANSFER>
807 inline bool
809 {
810  return false;
811 }
812 
813 template <int dim, typename VectorType, class TRANSFER>
814 template <class OtherVectorType>
815 void
817  OtherVectorType & dst,
818  const OtherVectorType &src) const
819 {
820  internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,
821  *multigrid,
822  *transfer,
823  dst,
824  src,
825  uses_dof_handler_vector,
826  this->signals,
827  0);
828 }
829 
830 
831 template <int dim, typename VectorType, class TRANSFER>
832 IndexSet
834  const unsigned int block) const
835 {
836  AssertIndexRange(block, dof_handler_vector.size());
837  return dof_handler_vector[block]->locally_owned_dofs();
838 }
839 
840 
841 template <int dim, typename VectorType, class TRANSFER>
842 IndexSet
844  const unsigned int block) const
845 {
846  AssertIndexRange(block, dof_handler_vector.size());
847  return dof_handler_vector[block]->locally_owned_dofs();
848 }
849 
850 
851 
852 template <int dim, typename VectorType, class TRANSFER>
853 const MPI_Comm &
855 {
856  // currently parallel GMG works with parallel triangulations only,
857  // so it should be a safe bet to use it to query MPI communicator:
858  const Triangulation<dim> &tria = dof_handler_vector[0]->get_triangulation();
859  const parallel::TriangulationBase<dim> *ptria =
860  dynamic_cast<const parallel::TriangulationBase<dim> *>(&tria);
861  Assert(ptria != nullptr, ExcInternalError());
862  return ptria->get_communicator();
863 }
864 
865 
866 
867 template <int dim, typename VectorType, class TRANSFER>
868 boost::signals2::connection
870  const std::function<void(bool)> &slot)
871 {
872  return this->signals.transfer_to_mg.connect(slot);
873 }
874 
875 
876 
877 template <int dim, typename VectorType, class TRANSFER>
878 boost::signals2::connection
880  const std::function<void(bool)> &slot)
881 {
882  return this->signals.transfer_to_global.connect(slot);
883 }
884 
885 
886 
887 template <int dim, typename VectorType, class TRANSFER>
888 template <class OtherVectorType>
889 void
891  OtherVectorType & dst,
892  const OtherVectorType &src) const
893 {
894  internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,
895  *multigrid,
896  *transfer,
897  dst,
898  src,
899  uses_dof_handler_vector,
900  this->signals,
901  0);
902 }
903 
904 
905 template <int dim, typename VectorType, class TRANSFER>
906 template <class OtherVectorType>
907 void
909  const OtherVectorType &) const
910 {
911  Assert(false, ExcNotImplemented());
912 }
913 
914 
915 template <int dim, typename VectorType, class TRANSFER>
916 template <class OtherVectorType>
917 void
919  OtherVectorType &,
920  const OtherVectorType &) const
921 {
922  Assert(false, ExcNotImplemented());
923 }
924 
925 #endif // DOXYGEN
926 
928 
929 #endif
IndexSet
Definition: index_set.h:74
mg_base.h
sparse_matrix.h
PreconditionMG::vmult
void vmult(OtherVectorType &dst, const OtherVectorType &src) const
mg::Signals::post_smoother_step
boost::signals2::signal< void(const bool before, const unsigned int level)> post_smoother_step
Definition: multigrid.h:108
Multigrid::connect_pre_smoother_step
boost::signals2::connection connect_pre_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
PreconditionMG::signals
mg::Signals signals
Definition: multigrid.h:586
Multigrid::t
MGLevelObject< VectorType > t
Definition: multigrid.h:366
Multigrid::pre_smooth
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > pre_smooth
Definition: multigrid.h:395
MGTransferBase
Definition: mg_base.h:177
Multigrid::Multigrid
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)
Multigrid::reinit
void reinit(const unsigned int minlevel, const unsigned int maxlevel)
MGMatrixBase
Definition: mg_base.h:49
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Triangulation::get_triangulation
Triangulation< dim, spacedim > & get_triangulation()
Definition: tria.cc:13338
Triangulation< dim >
Multigrid::edge_in
SmartPointer< const MGMatrixBase< VectorType > > edge_in
Definition: multigrid.h:416
Multigrid::Cycle
Cycle
Definition: multigrid.h:143
VectorType
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
Multigrid::cycle_type
Cycle cycle_type
Definition: multigrid.h:338
MGCoarseGridBase
Definition: mg_base.h:111
MPI_Comm
PreconditionMG::connect_transfer_to_global
boost::signals2::connection connect_transfer_to_global(const std::function< void(bool)> &slot)
mg
Definition: mg.h:81
PreconditionMG::locally_owned_range_indices
IndexSet locally_owned_range_indices(const unsigned int block=0) const
Multigrid::transfer
SmartPointer< const MGTransferBase< VectorType >, Multigrid< VectorType > > transfer
Definition: multigrid.h:389
PreconditionMG::uses_dof_handler_vector
const bool uses_dof_handler_vector
Definition: multigrid.h:581
Multigrid::solution
MGLevelObject< VectorType > solution
Definition: multigrid.h:360
Multigrid
Definition: multigrid.h:137
DoFHandler< dim >
Multigrid::connect_coarse_solve
boost::signals2::connection connect_coarse_solve(const std::function< void(const bool, const unsigned int)> &slot)
PreconditionMG::Tvmult
void Tvmult(OtherVectorType &dst, const OtherVectorType &src) const
Subscriptor
Definition: subscriptor.h:62
Multigrid::signals
mg::Signals signals
Definition: multigrid.h:314
level
unsigned int level
Definition: grid_out.cc:4355
PreconditionMG::multigrid
SmartPointer< Multigrid< VectorType >, PreconditionMG< dim, VectorType, TRANSFER > > multigrid
Definition: multigrid.h:569
PreconditionMG::get_mpi_communicator
const MPI_Comm & get_mpi_communicator() const
mg::Signals::pre_smoother_step
boost::signals2::signal< void(const bool before, const unsigned int level)> pre_smoother_step
Definition: multigrid.h:99
Multigrid::w_cycle
@ w_cycle
The W-cycle.
Definition: multigrid.h:148
Multigrid::edge_down
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_down
Definition: multigrid.h:423
PreconditionMG::empty
bool empty() const
Multigrid::edge_out
SmartPointer< const MGMatrixBase< VectorType > > edge_out
Definition: multigrid.h:408
Multigrid::defect
MGLevelObject< VectorType > defect
Definition: multigrid.h:355
PreconditionMG::vmult_add
void vmult_add(OtherVectorType &dst, const OtherVectorType &src) const
Multigrid::level_v_step
void level_v_step(const unsigned int level)
internal::reinit
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:621
subscriptor.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
PreconditionMG::dof_handler_vector_raw
std::vector< const DoFHandler< dim > * > dof_handler_vector_raw
Definition: multigrid.h:563
Multigrid::f_cycle
@ f_cycle
The F-cycle.
Definition: multigrid.h:150
Multigrid::matrix
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > matrix
Definition: multigrid.h:377
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
Multigrid::coarse
SmartPointer< const MGCoarseGridBase< VectorType >, Multigrid< VectorType > > coarse
Definition: multigrid.h:383
MGSmootherBase
Definition: mg_base.h:245
smartpointer.h
Multigrid::set_maxlevel
void set_maxlevel(const unsigned int)
mg_level_object.h
Multigrid::connect_prolongation
boost::signals2::connection connect_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
mg::Signals::transfer_to_global
boost::signals2::signal< void(const bool before)> transfer_to_global
Definition: multigrid.h:67
mg::Signals::prolongation
boost::signals2::signal< void(const bool before, const unsigned int level)> prolongation
Definition: multigrid.h:91
mg::Signals::transfer_to_mg
boost::signals2::signal< void(const bool before)> transfer_to_mg
Definition: multigrid.h:60
PreconditionMG::connect_transfer_to_mg
boost::signals2::connection connect_transfer_to_mg(const std::function< void(bool)> &slot)
Multigrid::set_edge_matrices
void set_edge_matrices(const MGMatrixBase< VectorType > &edge_out, const MGMatrixBase< VectorType > &edge_in)
mg::Signals::coarse_solve
boost::signals2::signal< void(const bool before, const unsigned int level)> coarse_solve
Definition: multigrid.h:75
Multigrid::connect_restriction
boost::signals2::connection connect_restriction(const std::function< void(const bool, const unsigned int)> &slot)
Multigrid::post_smooth
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > post_smooth
Definition: multigrid.h:401
SmartPointer
Definition: smartpointer.h:68
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
tria.h
Multigrid::level_step
void level_step(const unsigned int level, Cycle cycle)
Multigrid::set_edge_flux_matrices
void set_edge_flux_matrices(const MGMatrixBase< VectorType > &edge_down, const MGMatrixBase< VectorType > &edge_up)
parallel::TriangulationBase::get_communicator
virtual const MPI_Comm & get_communicator() const
Definition: tria_base.cc:138
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
PreconditionMG
Definition: multigrid.h:454
MGLevelObject< VectorType >
Multigrid::connect_post_smoother_step
boost::signals2::connection connect_post_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
vector.h
dof_handler.h
mg::Signals
Definition: multigrid.h:53
PreconditionMG::locally_owned_domain_indices
IndexSet locally_owned_domain_indices(const unsigned int block=0) const
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
PreconditionMG::transfer
SmartPointer< const TRANSFER, PreconditionMG< dim, VectorType, TRANSFER > > transfer
Definition: multigrid.h:575
Multigrid::v_cycle
@ v_cycle
The V-cycle.
Definition: multigrid.h:146
Multigrid::get_minlevel
unsigned int get_minlevel() const
Multigrid::maxlevel
unsigned int maxlevel
Definition: multigrid.h:348
parallel::TriangulationBase
Definition: tria_base.h:78
PreconditionMG::Tvmult_add
void Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const
PreconditionMG::PreconditionMG
PreconditionMG(const DoFHandler< dim > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer)
config.h
Multigrid::set_minlevel
void set_minlevel(const unsigned int level, bool relative=false)
internal
Definition: aligned_vector.h:369
Multigrid::vcycle
void vcycle()
Multigrid::edge_up
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_up
Definition: multigrid.h:430
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Multigrid::get_maxlevel
unsigned int get_maxlevel() const
Multigrid::minlevel
unsigned int minlevel
Definition: multigrid.h:343
mg::Signals::restriction
boost::signals2::signal< void(const bool before, const unsigned int level)> restriction
Definition: multigrid.h:83
PreconditionMG::dof_handler_vector
std::vector< SmartPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TRANSFER > > > dof_handler_vector
Definition: multigrid.h:557
Multigrid::cycle
void cycle()
Multigrid::set_cycle
void set_cycle(Cycle)
Multigrid::defect2
MGLevelObject< VectorType > defect2
Definition: multigrid.h:371