deal.II version GIT relicensing-2206-gaa53ff9447 2024-12-02 09:10:00+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\}}\)
Loading...
Searching...
No Matches
multigrid.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1999 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_multigrid_h
16#define dealii_multigrid_h
17
18
19#include <deal.II/base/config.h>
20
24
26
28
30#include <deal.II/lac/vector.h>
31
33
34#include <vector>
35
37
38#ifdef signals
39# error \
40 "The name 'signals' is already defined. You are most likely using the QT library \
41and using the 'signals' keyword. You can either #include the Qt headers (or any conflicting headers) \
42*after* the deal.II headers or you can define the 'QT_NO_KEYWORDS' macro and use the 'Q_SIGNALS' macro."
43#endif
44
50namespace mg
51{
61 struct Signals
62 {
68 boost::signals2::signal<void(const bool before)> transfer_to_mg;
69
75 boost::signals2::signal<void(const bool before)> transfer_to_global;
76
86 boost::signals2::signal<void(const bool before, const unsigned int level)>
88
97 boost::signals2::signal<void(const bool before, const unsigned int level)>
99
105 boost::signals2::signal<void(const bool before, const unsigned int level)>
107
116 boost::signals2::signal<void(const bool before, const unsigned int level)>
118
124 boost::signals2::signal<void(const bool before, const unsigned int level)>
126
132 boost::signals2::signal<void(const bool before, const unsigned int level)>
134
139 boost::signals2::signal<void(const bool before, const unsigned int level)>
141 };
142} // namespace mg
143
161template <typename VectorType>
163{
164public:
177
178 using vector_type = VectorType;
179 using const_vector_type = const VectorType;
180
195 const unsigned int minlevel = 0,
196 const unsigned int maxlevel = numbers::invalid_unsigned_int,
198
202 void
203 reinit(const unsigned int minlevel, const unsigned int maxlevel);
204
209 void
211
222 void
224
237 void
240
247 void
249
262 void
265
269 unsigned int
271
275 unsigned int
277
283 void
284 set_maxlevel(const unsigned int);
285
301 void
302 set_minlevel(const unsigned int level, bool relative = false);
303
308
312 boost::signals2::connection
314 const std::function<void(const bool, const unsigned int)> &slot);
315
319 boost::signals2::connection
321 const std::function<void(const bool, const unsigned int)> &slot);
322
326 boost::signals2::connection
328 const std::function<void(const bool, const unsigned int)> &slot);
329
333 boost::signals2::connection
335 const std::function<void(const bool, const unsigned int)> &slot);
336
340 boost::signals2::connection
342 const std::function<void(const bool, const unsigned int)> &slot);
343
347 boost::signals2::connection
349 const std::function<void(const bool, const unsigned int)> &slot);
350
354 boost::signals2::connection
356 const std::function<void(const bool, const unsigned int)> &slot);
357
358private:
363
370 void
371 level_v_step(const unsigned int level);
372
380 void
381 level_step(const unsigned int level, Cycle cycle);
382
387
391 unsigned int minlevel;
392
396 unsigned int maxlevel;
397
398public:
404
409
410private:
415
420
421
426
432
438
444
450
457
465
473
481
482 template <int dim, typename OtherVectorType, typename TransferType>
483 friend class PreconditionMG;
484};
485
486
501template <int dim, typename VectorType, typename TransferType>
503{
504public:
509 PreconditionMG(const DoFHandler<dim> &dof_handler,
511 const TransferType &transfer);
512
517 PreconditionMG(const std::vector<const DoFHandler<dim> *> &dof_handler,
519 const TransferType &transfer);
520
524 bool
525 empty() const;
526
533 template <typename OtherVectorType>
534 void
535 vmult(OtherVectorType &dst, const OtherVectorType &src) const;
536
541 template <typename OtherVectorType>
542 void
543 vmult_add(OtherVectorType &dst, const OtherVectorType &src) const;
544
550 template <typename OtherVectorType>
551 void
552 Tvmult(OtherVectorType &dst, const OtherVectorType &src) const;
553
559 template <typename OtherVectorType>
560 void
561 Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const;
562
570 locally_owned_range_indices(const unsigned int block = 0) const;
571
579 locally_owned_domain_indices(const unsigned int block = 0) const;
580
586
590 boost::signals2::connection
591 connect_transfer_to_mg(const std::function<void(bool)> &slot);
592
596 boost::signals2::connection
597 connect_transfer_to_global(const std::function<void(bool)> &slot);
598
604
610
611private:
615 std::vector<ObserverPointer<const DoFHandler<dim>,
618
623 std::vector<const DoFHandler<dim> *> dof_handler_vector_raw;
624
631
635 ObserverPointer<const TransferType,
638
644
649};
650
653#ifndef DOXYGEN
654/* --------------------------- inline functions --------------------- */
655
656
657template <typename VectorType>
659 const MGCoarseGridBase<VectorType> &coarse,
660 const MGTransferBase<VectorType> &transfer,
661 const MGSmootherBase<VectorType> &pre_smooth,
662 const MGSmootherBase<VectorType> &post_smooth,
663 const unsigned int min_level,
664 const unsigned int max_level,
665 Cycle cycle)
666 : cycle_type(cycle)
667 , matrix(&matrix, typeid(*this).name())
668 , coarse(&coarse, typeid(*this).name())
669 , transfer(&transfer, typeid(*this).name())
670 , pre_smooth(&pre_smooth, typeid(*this).name())
671 , post_smooth(&post_smooth, typeid(*this).name())
672 , edge_out(nullptr, typeid(*this).name())
673 , edge_in(nullptr, typeid(*this).name())
674 , edge_down(nullptr, typeid(*this).name())
675 , edge_up(nullptr, typeid(*this).name())
676{
677 if (max_level == numbers::invalid_unsigned_int)
678 maxlevel = matrix.get_maxlevel();
679 else
680 maxlevel = max_level;
681 reinit(min_level, maxlevel);
682}
683
684
685
686template <typename VectorType>
687inline unsigned int
689{
690 return maxlevel;
691}
692
693
694
695template <typename VectorType>
696inline unsigned int
698{
699 return minlevel;
700}
701
702
703/* --------------------------- inline functions --------------------- */
704
705
706namespace internal
707{
708 namespace PreconditionMGImplementation
709 {
710 template <int dim,
711 typename VectorType,
712 typename TransferType,
713 typename OtherVectorType>
714 std::enable_if_t<TransferType::supports_dof_handler_vector>
715 vmult(const std::vector<const DoFHandler<dim> *> &dof_handler_vector,
716 Multigrid<VectorType> &multigrid,
717 const TransferType &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(dof_handler_vector, dst, multigrid.solution);
736 else
737 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.solution);
738 signals.transfer_to_global(false);
739 }
740
741 template <int dim,
742 typename VectorType,
743 typename TransferType,
744 typename OtherVectorType>
745 void
746 vmult(const std::vector<const DoFHandler<dim> *> &dof_handler_vector,
747 Multigrid<VectorType> &multigrid,
748 const TransferType &transfer,
749 OtherVectorType &dst,
750 const OtherVectorType &src,
751 const bool uses_dof_handler_vector,
752 const typename ::mg::Signals &signals,
753 ...)
754 {
755 (void)uses_dof_handler_vector;
756 Assert(!uses_dof_handler_vector, ExcInternalError());
757
758 signals.transfer_to_mg(true);
759 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
760 signals.transfer_to_mg(false);
761
762 multigrid.cycle();
763
764 signals.transfer_to_global(true);
765 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.solution);
766 signals.transfer_to_global(false);
767 }
768
769 template <int dim,
770 typename VectorType,
771 typename TransferType,
772 typename OtherVectorType>
773 std::enable_if_t<TransferType::supports_dof_handler_vector>
774 vmult_add(const std::vector<const DoFHandler<dim> *> &dof_handler_vector,
775 Multigrid<VectorType> &multigrid,
776 const TransferType &transfer,
777 OtherVectorType &dst,
778 const OtherVectorType &src,
779 const bool uses_dof_handler_vector,
780 const typename ::mg::Signals &signals,
781 int)
782 {
783 signals.transfer_to_mg(true);
784 if (uses_dof_handler_vector)
785 transfer.copy_to_mg(dof_handler_vector, multigrid.defect, src);
786 else
787 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
788 signals.transfer_to_mg(false);
789
790 multigrid.cycle();
791
792 signals.transfer_to_global(true);
793 if (uses_dof_handler_vector)
794 transfer.copy_from_mg_add(dof_handler_vector, dst, multigrid.solution);
795 else
796 transfer.copy_from_mg_add(*dof_handler_vector[0],
797 dst,
798 multigrid.solution);
799 signals.transfer_to_global(false);
800 }
801
802 template <int dim,
803 typename VectorType,
804 typename TransferType,
805 typename OtherVectorType>
806 void
807 vmult_add(const std::vector<const DoFHandler<dim> *> &dof_handler_vector,
808 Multigrid<VectorType> &multigrid,
809 const TransferType &transfer,
810 OtherVectorType &dst,
811 const OtherVectorType &src,
812 const bool uses_dof_handler_vector,
813 const typename ::mg::Signals &signals,
814 ...)
815 {
816 (void)uses_dof_handler_vector;
817 Assert(!uses_dof_handler_vector, ExcInternalError());
818
819 signals.transfer_to_mg(true);
820 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src);
821 signals.transfer_to_mg(false);
822
823 multigrid.cycle();
824
825 signals.transfer_to_global(true);
826 transfer.copy_from_mg_add(*dof_handler_vector[0],
827 dst,
828 multigrid.solution);
829 signals.transfer_to_global(false);
830 }
831 } // namespace PreconditionMGImplementation
832} // namespace internal
833
834template <int dim, typename VectorType, typename TransferType>
836 const DoFHandler<dim> &dof_handler,
838 const TransferType &transfer)
839 : dof_handler_vector(1, &dof_handler)
840 , dof_handler_vector_raw(1, &dof_handler)
841 , multigrid(&mg)
842 , transfer(&transfer)
843 , uses_dof_handler_vector(false)
844{}
845
846template <int dim, typename VectorType, typename TransferType>
848 const std::vector<const DoFHandler<dim> *> &dof_handler,
850 const TransferType &transfer)
851 : dof_handler_vector(dof_handler.size())
852 , dof_handler_vector_raw(dof_handler.size())
853 , multigrid(&mg)
854 , transfer(&transfer)
855 , uses_dof_handler_vector(true)
856{
857 for (unsigned int i = 0; i < dof_handler.size(); ++i)
858 {
859 dof_handler_vector[i] = dof_handler[i];
860 dof_handler_vector_raw[i] = dof_handler[i];
861 }
862}
863
864template <int dim, typename VectorType, typename TransferType>
865inline bool
867{
868 return false;
869}
870
871template <int dim, typename VectorType, typename TransferType>
872template <typename OtherVectorType>
873void
875 OtherVectorType &dst,
876 const OtherVectorType &src) const
877{
878 internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,
879 *multigrid,
880 *transfer,
881 dst,
882 src,
883 uses_dof_handler_vector,
884 this->signals,
885 0);
886}
887
888
889template <int dim, typename VectorType, typename TransferType>
892 const unsigned int block) const
893{
894 AssertIndexRange(block, dof_handler_vector.size());
895 return dof_handler_vector[block]->locally_owned_dofs();
896}
897
898
899template <int dim, typename VectorType, typename TransferType>
902 const unsigned int block) const
903{
904 AssertIndexRange(block, dof_handler_vector.size());
905 return dof_handler_vector[block]->locally_owned_dofs();
906}
907
908
909
910template <int dim, typename VectorType, typename TransferType>
913{
914 // currently parallel GMG works with parallel triangulations only,
915 // so it should be a safe bet to use it to query MPI communicator:
916 const Triangulation<dim> &tria = dof_handler_vector[0]->get_triangulation();
918 dynamic_cast<const parallel::TriangulationBase<dim> *>(&tria);
919 Assert(ptria != nullptr, ExcInternalError());
920 return ptria->get_mpi_communicator();
921}
922
923
924
925template <int dim, typename VectorType, typename TransferType>
926boost::signals2::connection
928 const std::function<void(bool)> &slot)
929{
930 return this->signals.transfer_to_mg.connect(slot);
931}
932
933
934
935template <int dim, typename VectorType, typename TransferType>
936boost::signals2::connection
938 const std::function<void(bool)> &slot)
939{
940 return this->signals.transfer_to_global.connect(slot);
941}
942
943
944
945template <int dim, typename VectorType, typename TransferType>
946template <typename OtherVectorType>
947void
949 OtherVectorType &dst,
950 const OtherVectorType &src) const
951{
952 internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,
953 *multigrid,
954 *transfer,
955 dst,
956 src,
957 uses_dof_handler_vector,
958 this->signals,
959 0);
960}
961
962
963template <int dim, typename VectorType, typename TransferType>
964template <typename OtherVectorType>
965void
967 OtherVectorType &,
968 const OtherVectorType &) const
969{
971}
972
973
974template <int dim, typename VectorType, typename TransferType>
975template <typename OtherVectorType>
976void
978 OtherVectorType &,
979 const OtherVectorType &) const
980{
982}
983
984
985template <int dim, typename VectorType, typename TransferType>
988{
989 return *this->multigrid;
990}
991
992
993template <int dim, typename VectorType, typename TransferType>
996{
997 return *this->multigrid;
998}
999
1000#endif // DOXYGEN
1001
1003
1004#endif
@ v_cycle
The V-cycle.
Definition multigrid.h:171
@ w_cycle
The W-cycle.
Definition multigrid.h:173
@ f_cycle
The F-cycle.
Definition multigrid.h:175
ObserverPointer< const MGMatrixBase< VectorType > > edge_out
Definition multigrid.h:456
MGLevelObject< VectorType > defect2
Definition multigrid.h:419
void cycle()
void set_edge_flux_matrices(const MGMatrixBase< VectorType > &edge_down, const MGMatrixBase< VectorType > &edge_up)
unsigned int minlevel
Definition multigrid.h:391
mg::Signals signals
Definition multigrid.h:362
MGLevelObject< VectorType > solution
Definition multigrid.h:408
ObserverPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > post_smooth
Definition multigrid.h:449
boost::signals2::connection connect_edge_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
ObserverPointer< const MGMatrixBase< VectorType > > edge_in
Definition multigrid.h:464
unsigned int maxlevel
Definition multigrid.h:396
VectorType vector_type
Definition multigrid.h:178
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
ObserverPointer< const MGCoarseGridBase< VectorType >, Multigrid< VectorType > > coarse
Definition multigrid.h:431
ObserverPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > pre_smooth
Definition multigrid.h:443
const VectorType const_vector_type
Definition multigrid.h:179
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)
ObserverPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > matrix
Definition multigrid.h:425
Cycle cycle_type
Definition multigrid.h:386
boost::signals2::connection connect_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
void set_edge_in_matrix(const MGMatrixBase< VectorType > &edge_in)
MGLevelObject< VectorType > t
Definition multigrid.h:414
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:403
void vcycle()
ObserverPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_down
Definition multigrid.h:472
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)
ObserverPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_up
Definition multigrid.h:480
boost::signals2::connection connect_pre_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
ObserverPointer< const MGTransferBase< VectorType >, Multigrid< VectorType > > transfer
Definition multigrid.h:437
std::vector< const DoFHandler< dim > * > dof_handler_vector_raw
Definition multigrid.h:623
void Tvmult(OtherVectorType &dst, const OtherVectorType &src) const
IndexSet locally_owned_domain_indices(const unsigned int block=0) const
ObserverPointer< Multigrid< VectorType >, PreconditionMG< dim, VectorType, TransferType > > multigrid
Definition multigrid.h:630
const Multigrid< VectorType > & get_multigrid() const
boost::signals2::connection connect_transfer_to_global(const std::function< void(bool)> &slot)
mg::Signals signals
Definition multigrid.h:648
boost::signals2::connection connect_transfer_to_mg(const std::function< void(bool)> &slot)
PreconditionMG(const DoFHandler< dim > &dof_handler, Multigrid< VectorType > &mg, const TransferType &transfer)
ObserverPointer< const TransferType, PreconditionMG< dim, VectorType, TransferType > > transfer
Definition multigrid.h:637
MPI_Comm get_mpi_communicator() const
const bool uses_dof_handler_vector
Definition multigrid.h:643
void vmult_add(OtherVectorType &dst, const OtherVectorType &src) const
bool empty() const
IndexSet locally_owned_range_indices(const unsigned int block=0) const
std::vector< ObserverPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TransferType > > > dof_handler_vector
Definition multigrid.h:617
PreconditionMG(const std::vector< const DoFHandler< dim > * > &dof_handler, Multigrid< VectorType > &mg, const TransferType &transfer)
Multigrid< VectorType > & get_multigrid()
void Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const
void vmult(OtherVectorType &dst, const OtherVectorType &src) const
Triangulation< dim, spacedim > & get_triangulation()
virtual MPI_Comm get_mpi_communicator() const override
Definition tria_base.cc:160
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
unsigned int level
Definition grid_out.cc:4626
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
#define DEAL_II_NOT_IMPLEMENTED()
std::size_t size
Definition mpi.cc:734
Tpetra::Vector< Number, LO, GO, NodeType< MemorySpace > > VectorType
Definition mg.h:81
static const unsigned int invalid_unsigned_int
Definition types.h:220
boost::signals2::signal< void(const bool before, const unsigned int level)> prolongation
Definition multigrid.h:106
boost::signals2::signal< void(const bool before, const unsigned int level)> pre_smoother_step
Definition multigrid.h:117
boost::signals2::signal< void(const bool before, const unsigned int level)> coarse_solve
Definition multigrid.h:87
boost::signals2::signal< void(const bool before, const unsigned int level)> edge_prolongation
Definition multigrid.h:140
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:75
boost::signals2::signal< void(const bool before, const unsigned int level)> residual_step
Definition multigrid.h:133
boost::signals2::signal< void(const bool before, const unsigned int level)> restriction
Definition multigrid.h:98
boost::signals2::signal< void(const bool before)> transfer_to_mg
Definition multigrid.h:68