16#ifndef dealii_multigrid_h
17#define dealii_multigrid_h
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."
87 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
98 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
106 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
117 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
125 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
133 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
140 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
162template <
typename VectorType>
313 boost::signals2::connection
315 const std::function<
void(
const bool,
const unsigned int)> &slot);
320 boost::signals2::connection
322 const std::function<
void(
const bool,
const unsigned int)> &slot);
327 boost::signals2::connection
329 const std::function<
void(
const bool,
const unsigned int)> &slot);
334 boost::signals2::connection
336 const std::function<
void(
const bool,
const unsigned int)> &slot);
341 boost::signals2::connection
343 const std::function<
void(
const bool,
const unsigned int)> &slot);
348 boost::signals2::connection
350 const std::function<
void(
const bool,
const unsigned int)> &slot);
355 boost::signals2::connection
357 const std::function<
void(
const bool,
const unsigned int)> &slot);
481 template <
int dim,
class OtherVectorType,
class TRANSFER>
500template <
int dim,
typename VectorType,
class TRANSFER>
532 template <
class OtherVectorType>
534 vmult(OtherVectorType &dst,
const OtherVectorType &src)
const;
540 template <
class OtherVectorType>
542 vmult_add(OtherVectorType &dst,
const OtherVectorType &src)
const;
549 template <
class OtherVectorType>
551 Tvmult(OtherVectorType &dst,
const OtherVectorType &src)
const;
558 template <
class OtherVectorType>
560 Tvmult_add(OtherVectorType &dst,
const OtherVectorType &src)
const;
589 boost::signals2::connection
595 boost::signals2::connection
614 std::vector<SmartPointer<const DoFHandler<dim>,
654template <
typename VectorType>
660 const unsigned int min_level,
661 const unsigned int max_level,
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())
675 maxlevel = matrix.get_maxlevel();
677 maxlevel = max_level;
678 reinit(min_level, maxlevel);
683template <
typename VectorType>
692template <
typename VectorType>
705 namespace PreconditionMGImplementation
710 typename OtherVectorType>
711 std::enable_if_t<TRANSFER::supports_dof_handler_vector>
714 const TRANSFER & transfer,
715 OtherVectorType & dst,
716 const OtherVectorType & src,
717 const bool uses_dof_handler_vector,
718 const typename ::mg::Signals & signals,
721 signals.transfer_to_mg(
true);
722 if (uses_dof_handler_vector)
723 transfer.copy_to_mg(dof_handler_vector, multigrid.
defect, src);
725 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
726 signals.transfer_to_mg(
false);
730 signals.transfer_to_global(
true);
731 if (uses_dof_handler_vector)
732 transfer.copy_from_mg(dof_handler_vector, dst, multigrid.
solution);
734 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.
solution);
735 signals.transfer_to_global(
false);
741 typename OtherVectorType>
745 const TRANSFER & transfer,
746 OtherVectorType & dst,
747 const OtherVectorType & src,
748 const bool uses_dof_handler_vector,
749 const typename ::mg::Signals & signals,
752 (void)uses_dof_handler_vector;
755 signals.transfer_to_mg(
true);
756 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
757 signals.transfer_to_mg(
false);
761 signals.transfer_to_global(
true);
762 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.
solution);
763 signals.transfer_to_global(
false);
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,
773 const TRANSFER & transfer,
774 OtherVectorType & dst,
775 const OtherVectorType & src,
776 const bool uses_dof_handler_vector,
777 const typename ::mg::Signals &signals,
780 signals.transfer_to_mg(
true);
781 if (uses_dof_handler_vector)
782 transfer.copy_to_mg(dof_handler_vector, multigrid.
defect, src);
784 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
785 signals.transfer_to_mg(
false);
789 signals.transfer_to_global(
true);
790 if (uses_dof_handler_vector)
791 transfer.copy_from_mg_add(dof_handler_vector, dst, multigrid.
solution);
793 transfer.copy_from_mg_add(*dof_handler_vector[0],
796 signals.transfer_to_global(
false);
802 typename OtherVectorType>
804 vmult_add(
const std::vector<
const DoFHandler<dim> *> &dof_handler_vector,
806 const TRANSFER & transfer,
807 OtherVectorType & dst,
808 const OtherVectorType & src,
809 const bool uses_dof_handler_vector,
810 const typename ::mg::Signals &signals,
813 (void)uses_dof_handler_vector;
816 signals.transfer_to_mg(
true);
817 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
818 signals.transfer_to_mg(
false);
822 signals.transfer_to_global(
true);
823 transfer.copy_from_mg_add(*dof_handler_vector[0],
826 signals.transfer_to_global(
false);
831template <
int dim,
typename VectorType,
class TRANSFER>
835 const TRANSFER & transfer)
836 : dof_handler_vector(1, &dof_handler)
837 , dof_handler_vector_raw(1, &dof_handler)
839 , transfer(&transfer)
840 , uses_dof_handler_vector(false)
843template <
int dim,
typename VectorType,
class TRANSFER>
847 const TRANSFER & transfer)
848 : dof_handler_vector(dof_handler.size())
849 , dof_handler_vector_raw(dof_handler.size())
851 , transfer(&transfer)
852 , uses_dof_handler_vector(true)
854 for (
unsigned int i = 0; i < dof_handler.size(); ++i)
861template <
int dim,
typename VectorType,
class TRANSFER>
868template <
int dim,
typename VectorType,
class TRANSFER>
869template <
class OtherVectorType>
872 OtherVectorType & dst,
873 const OtherVectorType &src)
const
875 internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,
880 uses_dof_handler_vector,
886template <
int dim,
typename VectorType,
class TRANSFER>
889 const unsigned int block)
const
892 return dof_handler_vector[block]->locally_owned_dofs();
896template <
int dim,
typename VectorType,
class TRANSFER>
899 const unsigned int block)
const
902 return dof_handler_vector[block]->locally_owned_dofs();
907template <
int dim,
typename VectorType,
class TRANSFER>
922template <
int dim,
typename VectorType,
class TRANSFER>
923boost::signals2::connection
925 const std::function<
void(
bool)> &slot)
927 return this->signals.transfer_to_mg.connect(slot);
932template <
int dim,
typename VectorType,
class TRANSFER>
933boost::signals2::connection
935 const std::function<
void(
bool)> &slot)
937 return this->signals.transfer_to_global.connect(slot);
942template <
int dim,
typename VectorType,
class TRANSFER>
943template <
class OtherVectorType>
946 OtherVectorType & dst,
947 const OtherVectorType &src)
const
949 internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,
954 uses_dof_handler_vector,
960template <
int dim,
typename VectorType,
class TRANSFER>
961template <
class OtherVectorType>
964 const OtherVectorType &)
const
970template <
int dim,
typename VectorType,
class TRANSFER>
971template <
class OtherVectorType>
975 const OtherVectorType &)
const
981template <
int dim,
typename VectorType,
class TRANSFER>
985 return *this->multigrid;
989template <
int dim,
typename VectorType,
class TRANSFER>
993 return *this->multigrid;
MGLevelObject< VectorType > defect2
void set_edge_flux_matrices(const MGMatrixBase< VectorType > &edge_down, const MGMatrixBase< VectorType > &edge_up)
SmartPointer< const MGTransferBase< VectorType >, Multigrid< VectorType > > transfer
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_up
SmartPointer< const MGMatrixBase< VectorType > > edge_in
SmartPointer< const MGMatrixBase< VectorType > > edge_out
MGLevelObject< VectorType > solution
boost::signals2::connection connect_edge_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > matrix
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)
unsigned int get_maxlevel() const
const VectorType const_vector_type
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > pre_smooth
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > post_smooth
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
boost::signals2::connection connect_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_down
void set_edge_in_matrix(const MGMatrixBase< VectorType > &edge_in)
MGLevelObject< VectorType > t
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
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
Multigrid< VectorType > & get_multigrid()
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
PreconditionMG(const std::vector< const DoFHandler< dim > * > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer)
const Multigrid< VectorType > & get_multigrid() const
std::vector< SmartPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TRANSFER > > > dof_handler_vector
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)
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
MPI_Comm get_mpi_communicator() const
std::vector< const DoFHandler< dim > * > dof_handler_vector_raw
Triangulation< dim, spacedim > & get_triangulation()
virtual MPI_Comm get_communicator() const override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static const unsigned int invalid_unsigned_int
boost::signals2::signal< void(const bool before, const unsigned int level)> prolongation
boost::signals2::signal< void(const bool before, const unsigned int level)> pre_smoother_step
boost::signals2::signal< void(const bool before, const unsigned int level)> coarse_solve
boost::signals2::signal< void(const bool before, const unsigned int level)> edge_prolongation
boost::signals2::signal< void(const bool before, const unsigned int level)> post_smoother_step
boost::signals2::signal< void(const bool before)> transfer_to_global
boost::signals2::signal< void(const bool before, const unsigned int level)> residual_step
boost::signals2::signal< void(const bool before, const unsigned int level)> restriction
boost::signals2::signal< void(const bool before)> transfer_to_mg
const ::Triangulation< dim, spacedim > & tria