Reference documentation for deal.II version 9.3.3
\(\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
39#ifdef signals
40# error \
41 "The name 'signals' is already defined. You are most likely using the QT library \
42and 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
48
49namespace mg
50{
60 struct Signals
61 {
67 boost::signals2::signal<void(const bool before)> transfer_to_mg;
68
74 boost::signals2::signal<void(const bool before)> transfer_to_global;
75
85 boost::signals2::signal<void(const bool before, const unsigned int level)>
87
96 boost::signals2::signal<void(const bool before, const unsigned int level)>
98
104 boost::signals2::signal<void(const bool before, const unsigned int level)>
106
115 boost::signals2::signal<void(const bool before, const unsigned int level)>
117
124 boost::signals2::signal<void(const bool before, const unsigned int level)>
126 };
127} // namespace mg
128
151template <typename VectorType>
152class Multigrid : public Subscriptor
153{
154public:
158 enum Cycle
159 {
165 f_cycle
166 };
167
168 using vector_type = VectorType;
169 using const_vector_type = const VectorType;
170
185 const unsigned int minlevel = 0,
186 const unsigned int maxlevel = numbers::invalid_unsigned_int,
188
192 void
193 reinit(const unsigned int minlevel, const unsigned int maxlevel);
194
199 void
201
212 void
214
227 void
230
243 void
246
250 unsigned int
252
256 unsigned int
258
264 void
265 set_maxlevel(const unsigned int);
266
282 void
283 set_minlevel(const unsigned int level, bool relative = false);
284
289
293 boost::signals2::connection
295 const std::function<void(const bool, const unsigned int)> &slot);
296
300 boost::signals2::connection
302 const std::function<void(const bool, const unsigned int)> &slot);
303
307 boost::signals2::connection
309 const std::function<void(const bool, const unsigned int)> &slot);
310
314 boost::signals2::connection
316 const std::function<void(const bool, const unsigned int)> &slot);
317
321 boost::signals2::connection
323 const std::function<void(const bool, const unsigned int)> &slot);
324
325private:
330
337 void
338 level_v_step(const unsigned int level);
339
347 void
348 level_step(const unsigned int level, Cycle cycle);
349
354
358 unsigned int minlevel;
359
363 unsigned int maxlevel;
364
365public:
371
376
377private:
382
387
388
393
399
405
411
417
424
432
439
446
447 template <int dim, class OtherVectorType, class TRANSFER>
448 friend class PreconditionMG;
449};
450
451
466template <int dim, typename VectorType, class TRANSFER>
468{
469public:
474 PreconditionMG(const DoFHandler<dim> &dof_handler,
476 const TRANSFER & transfer);
477
482 PreconditionMG(const std::vector<const DoFHandler<dim> *> &dof_handler,
484 const TRANSFER & transfer);
485
489 bool
490 empty() const;
491
498 template <class OtherVectorType>
499 void
500 vmult(OtherVectorType &dst, const OtherVectorType &src) const;
501
506 template <class OtherVectorType>
507 void
508 vmult_add(OtherVectorType &dst, const OtherVectorType &src) const;
509
515 template <class OtherVectorType>
516 void
517 Tvmult(OtherVectorType &dst, const OtherVectorType &src) const;
518
524 template <class OtherVectorType>
525 void
526 Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const;
527
535 locally_owned_range_indices(const unsigned int block = 0) const;
536
544 locally_owned_domain_indices(const unsigned int block = 0) const;
545
551
555 boost::signals2::connection
556 connect_transfer_to_mg(const std::function<void(bool)> &slot);
557
561 boost::signals2::connection
562 connect_transfer_to_global(const std::function<void(bool)> &slot);
563
564private:
568 std::vector<SmartPointer<const DoFHandler<dim>,
571
576 std::vector<const DoFHandler<dim> *> dof_handler_vector_raw;
577
583
589
595
600};
601
604#ifndef DOXYGEN
605/* --------------------------- inline functions --------------------- */
606
607
608template <typename VectorType>
610 const MGCoarseGridBase<VectorType> &coarse,
611 const MGTransferBase<VectorType> & transfer,
612 const MGSmootherBase<VectorType> & pre_smooth,
613 const MGSmootherBase<VectorType> &post_smooth,
614 const unsigned int min_level,
615 const unsigned int max_level,
616 Cycle cycle)
617 : cycle_type(cycle)
618 , matrix(&matrix, typeid(*this).name())
619 , coarse(&coarse, typeid(*this).name())
620 , transfer(&transfer, typeid(*this).name())
621 , pre_smooth(&pre_smooth, typeid(*this).name())
622 , post_smooth(&post_smooth, typeid(*this).name())
623 , edge_out(nullptr, typeid(*this).name())
624 , edge_in(nullptr, typeid(*this).name())
625 , edge_down(nullptr, typeid(*this).name())
626 , edge_up(nullptr, typeid(*this).name())
627{
628 if (max_level == numbers::invalid_unsigned_int)
629 maxlevel = matrix.get_maxlevel();
630 else
631 maxlevel = max_level;
632 reinit(min_level, maxlevel);
633}
634
635
636
637template <typename VectorType>
638inline unsigned int
640{
641 return maxlevel;
642}
643
644
645
646template <typename VectorType>
647inline unsigned int
649{
650 return minlevel;
651}
652
653
654/* --------------------------- inline functions --------------------- */
655
656
657namespace internal
658{
659 namespace PreconditionMGImplementation
660 {
661 template <int dim,
662 typename VectorType,
663 class TRANSFER,
664 typename OtherVectorType>
665 typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
666 vmult(
667 const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
668 ::Multigrid<VectorType> & multigrid,
669 const TRANSFER & transfer,
670 OtherVectorType & dst,
671 const OtherVectorType & src,
672 const bool uses_dof_handler_vector,
673 const typename ::mg::Signals &signals,
674 int)
675 {
676 signals.transfer_to_mg(true);
677 if (uses_dof_handler_vector)
678 transfer.copy_to_mg(dof_handler_vector, multigrid.defect, src);
679 else
680 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
681 signals.transfer_to_mg(false);
682
683 multigrid.cycle();
684
685 signals.transfer_to_global(true);
686 if (uses_dof_handler_vector)
687 transfer.copy_from_mg(dof_handler_vector, dst, multigrid.solution);
688 else
689 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.solution);
690 signals.transfer_to_global(false);
691 }
692
693 template <int dim,
694 typename VectorType,
695 class TRANSFER,
696 typename OtherVectorType>
697 void
698 vmult(
699 const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
700 ::Multigrid<VectorType> & multigrid,
701 const TRANSFER & transfer,
702 OtherVectorType & dst,
703 const OtherVectorType & src,
704 const bool uses_dof_handler_vector,
705 const typename ::mg::Signals &signals,
706 ...)
707 {
708 (void)uses_dof_handler_vector;
709 Assert(!uses_dof_handler_vector, ExcInternalError());
710
711 signals.transfer_to_mg(true);
712 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
713 signals.transfer_to_mg(false);
714
715 multigrid.cycle();
716
717 signals.transfer_to_global(true);
718 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.solution);
719 signals.transfer_to_global(false);
720 }
721
722 template <int dim,
723 typename VectorType,
724 class TRANSFER,
725 typename OtherVectorType>
726 typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
727 vmult_add(
728 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 int)
736 {
737 signals.transfer_to_mg(true);
738 if (uses_dof_handler_vector)
739 transfer.copy_to_mg(dof_handler_vector, multigrid.defect, src);
740 else
741 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
742 signals.transfer_to_mg(false);
743
744 multigrid.cycle();
745
746 signals.transfer_to_global(true);
747 if (uses_dof_handler_vector)
748 transfer.copy_from_mg_add(dof_handler_vector, dst, multigrid.solution);
749 else
750 transfer.copy_from_mg_add(*dof_handler_vector[0],
751 dst,
752 multigrid.solution);
753 signals.transfer_to_global(false);
754 }
755
756 template <int dim,
757 typename VectorType,
758 class TRANSFER,
759 typename OtherVectorType>
760 void
761 vmult_add(
762 const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
763 ::Multigrid<VectorType> & multigrid,
764 const TRANSFER & transfer,
765 OtherVectorType & dst,
766 const OtherVectorType & src,
767 const bool uses_dof_handler_vector,
768 const typename ::mg::Signals &signals,
769 ...)
770 {
771 (void)uses_dof_handler_vector;
772 Assert(!uses_dof_handler_vector, ExcInternalError());
773
774 signals.transfer_to_mg(true);
775 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
776 signals.transfer_to_mg(false);
777
778 multigrid.cycle();
779
780 signals.transfer_to_global(true);
781 transfer.copy_from_mg_add(*dof_handler_vector[0],
782 dst,
783 multigrid.solution);
784 signals.transfer_to_global(false);
785 }
786 } // namespace PreconditionMGImplementation
787} // namespace internal
788
789template <int dim, typename VectorType, class TRANSFER>
791 const DoFHandler<dim> &dof_handler,
793 const TRANSFER & transfer)
794 : dof_handler_vector(1, &dof_handler)
795 , dof_handler_vector_raw(1, &dof_handler)
796 , multigrid(&mg)
797 , transfer(&transfer)
798 , uses_dof_handler_vector(false)
799{}
800
801template <int dim, typename VectorType, class TRANSFER>
803 const std::vector<const DoFHandler<dim> *> &dof_handler,
805 const TRANSFER & transfer)
806 : dof_handler_vector(dof_handler.size())
807 , dof_handler_vector_raw(dof_handler.size())
808 , multigrid(&mg)
809 , transfer(&transfer)
810 , uses_dof_handler_vector(true)
811{
812 for (unsigned int i = 0; i < dof_handler.size(); ++i)
813 {
814 dof_handler_vector[i] = dof_handler[i];
815 dof_handler_vector_raw[i] = dof_handler[i];
816 }
817}
818
819template <int dim, typename VectorType, class TRANSFER>
820inline bool
822{
823 return false;
824}
825
826template <int dim, typename VectorType, class TRANSFER>
827template <class OtherVectorType>
828void
830 OtherVectorType & dst,
831 const OtherVectorType &src) const
832{
833 internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,
834 *multigrid,
835 *transfer,
836 dst,
837 src,
838 uses_dof_handler_vector,
839 this->signals,
840 0);
841}
842
843
844template <int dim, typename VectorType, class TRANSFER>
847 const unsigned int block) const
848{
849 AssertIndexRange(block, dof_handler_vector.size());
850 return dof_handler_vector[block]->locally_owned_dofs();
851}
852
853
854template <int dim, typename VectorType, class TRANSFER>
857 const unsigned int block) const
858{
859 AssertIndexRange(block, dof_handler_vector.size());
860 return dof_handler_vector[block]->locally_owned_dofs();
861}
862
863
864
865template <int dim, typename VectorType, class TRANSFER>
868{
869 // currently parallel GMG works with parallel triangulations only,
870 // so it should be a safe bet to use it to query MPI communicator:
871 const Triangulation<dim> &tria = dof_handler_vector[0]->get_triangulation();
873 dynamic_cast<const parallel::TriangulationBase<dim> *>(&tria);
874 Assert(ptria != nullptr, ExcInternalError());
875 return ptria->get_communicator();
876}
877
878
879
880template <int dim, typename VectorType, class TRANSFER>
881boost::signals2::connection
883 const std::function<void(bool)> &slot)
884{
885 return this->signals.transfer_to_mg.connect(slot);
886}
887
888
889
890template <int dim, typename VectorType, class TRANSFER>
891boost::signals2::connection
893 const std::function<void(bool)> &slot)
894{
895 return this->signals.transfer_to_global.connect(slot);
896}
897
898
899
900template <int dim, typename VectorType, class TRANSFER>
901template <class OtherVectorType>
902void
904 OtherVectorType & dst,
905 const OtherVectorType &src) const
906{
907 internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,
908 *multigrid,
909 *transfer,
910 dst,
911 src,
912 uses_dof_handler_vector,
913 this->signals,
914 0);
915}
916
917
918template <int dim, typename VectorType, class TRANSFER>
919template <class OtherVectorType>
920void
922 const OtherVectorType &) const
923{
924 Assert(false, ExcNotImplemented());
925}
926
927
928template <int dim, typename VectorType, class TRANSFER>
929template <class OtherVectorType>
930void
932 OtherVectorType &,
933 const OtherVectorType &) const
934{
935 Assert(false, ExcNotImplemented());
936}
937
938#endif // DOXYGEN
939
941
942#endif
@ v_cycle
The V-cycle.
Definition: multigrid.h:161
@ w_cycle
The W-cycle.
Definition: multigrid.h:163
@ f_cycle
The F-cycle.
Definition: multigrid.h:165
MGLevelObject< VectorType > defect2
Definition: multigrid.h:386
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:404
unsigned int minlevel
Definition: multigrid.h:358
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_up
Definition: multigrid.h:445
mg::Signals signals
Definition: multigrid.h:329
SmartPointer< const MGMatrixBase< VectorType > > edge_in
Definition: multigrid.h:431
SmartPointer< const MGMatrixBase< VectorType > > edge_out
Definition: multigrid.h:423
MGLevelObject< VectorType > solution
Definition: multigrid.h:375
unsigned int maxlevel
Definition: multigrid.h:363
VectorType vector_type
Definition: multigrid.h:168
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > matrix
Definition: multigrid.h:392
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)
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:169
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > pre_smooth
Definition: multigrid.h:410
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > post_smooth
Definition: multigrid.h:416
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:398
Cycle cycle_type
Definition: multigrid.h:353
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:438
MGLevelObject< VectorType > t
Definition: multigrid.h:381
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:370
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:582
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:588
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:570
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)
PreconditionMG(const DoFHandler< dim > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer)
mg::Signals signals
Definition: multigrid.h:599
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:594
MPI_Comm get_mpi_communicator() const
std::vector< const DoFHandler< dim > * > dof_handler_vector_raw
Definition: multigrid.h:576
Triangulation< dim, spacedim > & get_triangulation()
virtual MPI_Comm get_communicator() const override
Definition: tria_base.cc:140
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
unsigned int level
Definition: grid_out.cc:4590
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
static ::ExceptionBase & ExcInternalError()
@ 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:196
boost::signals2::signal< void(const bool before, const unsigned int level)> prolongation
Definition: multigrid.h:105
boost::signals2::signal< void(const bool before, const unsigned int level)> pre_smoother_step
Definition: multigrid.h:116
boost::signals2::signal< void(const bool before, const unsigned int level)> coarse_solve
Definition: multigrid.h:86
boost::signals2::signal< void(const bool before, const unsigned int level)> post_smoother_step
Definition: multigrid.h:125
boost::signals2::signal< void(const bool before)> transfer_to_global
Definition: multigrid.h:74
boost::signals2::signal< void(const bool before, const unsigned int level)> restriction
Definition: multigrid.h:97
boost::signals2::signal< void(const bool before)> transfer_to_mg
Definition: multigrid.h:67