15#ifndef dealii_multigrid_h
16#define dealii_multigrid_h
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."
86 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
97 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
105 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
116 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
124 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
132 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
139 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
161template <
typename VectorType>
312 boost::signals2::connection
314 const std::function<
void(
const bool,
const unsigned int)> &slot);
319 boost::signals2::connection
321 const std::function<
void(
const bool,
const unsigned int)> &slot);
326 boost::signals2::connection
328 const std::function<
void(
const bool,
const unsigned int)> &slot);
333 boost::signals2::connection
335 const std::function<
void(
const bool,
const unsigned int)> &slot);
340 boost::signals2::connection
342 const std::function<
void(
const bool,
const unsigned int)> &slot);
347 boost::signals2::connection
349 const std::function<
void(
const bool,
const unsigned int)> &slot);
354 boost::signals2::connection
356 const std::function<
void(
const bool,
const unsigned int)> &slot);
480 template <
int dim,
typename OtherVectorType,
typename TransferType>
499template <
int dim,
typename VectorType,
typename TransferType>
531 template <
typename OtherVectorType>
533 vmult(OtherVectorType &dst,
const OtherVectorType &src)
const;
539 template <
typename OtherVectorType>
541 vmult_add(OtherVectorType &dst,
const OtherVectorType &src)
const;
548 template <
typename OtherVectorType>
550 Tvmult(OtherVectorType &dst,
const OtherVectorType &src)
const;
557 template <
typename OtherVectorType>
559 Tvmult_add(OtherVectorType &dst,
const OtherVectorType &src)
const;
588 boost::signals2::connection
594 boost::signals2::connection
613 std::vector<SmartPointer<const DoFHandler<dim>,
655template <
typename VectorType>
661 const unsigned int min_level,
662 const unsigned int max_level,
665 , matrix(&matrix, typeid(*this).name())
666 , coarse(&coarse, typeid(*this).name())
667 , transfer(&transfer, typeid(*this).name())
668 , pre_smooth(&pre_smooth, typeid(*this).name())
669 , post_smooth(&post_smooth, typeid(*this).name())
670 , edge_out(nullptr, typeid(*this).name())
671 , edge_in(nullptr, typeid(*this).name())
672 , edge_down(nullptr, typeid(*this).name())
673 , edge_up(nullptr, typeid(*this).name())
676 maxlevel = matrix.get_maxlevel();
678 maxlevel = max_level;
679 reinit(min_level, maxlevel);
684template <
typename VectorType>
693template <
typename VectorType>
706 namespace PreconditionMGImplementation
710 typename TransferType,
711 typename OtherVectorType>
712 std::enable_if_t<TransferType::supports_dof_handler_vector>
715 const TransferType &transfer,
716 OtherVectorType &dst,
717 const OtherVectorType &src,
718 const bool uses_dof_handler_vector,
719 const typename ::mg::Signals &signals,
722 signals.transfer_to_mg(
true);
723 if (uses_dof_handler_vector)
724 transfer.copy_to_mg(dof_handler_vector, multigrid.
defect, src);
726 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
727 signals.transfer_to_mg(
false);
731 signals.transfer_to_global(
true);
732 if (uses_dof_handler_vector)
733 transfer.copy_from_mg(dof_handler_vector, dst, multigrid.
solution);
735 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.
solution);
736 signals.transfer_to_global(
false);
741 typename TransferType,
742 typename OtherVectorType>
746 const TransferType &transfer,
747 OtherVectorType &dst,
748 const OtherVectorType &src,
749 const bool uses_dof_handler_vector,
750 const typename ::mg::Signals &signals,
753 (void)uses_dof_handler_vector;
756 signals.transfer_to_mg(
true);
757 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
758 signals.transfer_to_mg(
false);
762 signals.transfer_to_global(
true);
763 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.
solution);
764 signals.transfer_to_global(
false);
769 typename TransferType,
770 typename OtherVectorType>
771 std::enable_if_t<TransferType::supports_dof_handler_vector>
772 vmult_add(
const std::vector<
const DoFHandler<dim> *> &dof_handler_vector,
774 const TransferType &transfer,
775 OtherVectorType &dst,
776 const OtherVectorType &src,
777 const bool uses_dof_handler_vector,
778 const typename ::mg::Signals &signals,
781 signals.transfer_to_mg(
true);
782 if (uses_dof_handler_vector)
783 transfer.copy_to_mg(dof_handler_vector, multigrid.
defect, src);
785 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
786 signals.transfer_to_mg(
false);
790 signals.transfer_to_global(
true);
791 if (uses_dof_handler_vector)
792 transfer.copy_from_mg_add(dof_handler_vector, dst, multigrid.
solution);
794 transfer.copy_from_mg_add(*dof_handler_vector[0],
797 signals.transfer_to_global(
false);
802 typename TransferType,
803 typename OtherVectorType>
805 vmult_add(
const std::vector<
const DoFHandler<dim> *> &dof_handler_vector,
807 const TransferType &transfer,
808 OtherVectorType &dst,
809 const OtherVectorType &src,
810 const bool uses_dof_handler_vector,
811 const typename ::mg::Signals &signals,
814 (void)uses_dof_handler_vector;
817 signals.transfer_to_mg(
true);
818 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
819 signals.transfer_to_mg(
false);
823 signals.transfer_to_global(
true);
824 transfer.copy_from_mg_add(*dof_handler_vector[0],
827 signals.transfer_to_global(
false);
832template <
int dim,
typename VectorType,
typename TransferType>
836 const TransferType &transfer)
837 : dof_handler_vector(1, &dof_handler)
838 , dof_handler_vector_raw(1, &dof_handler)
840 , transfer(&transfer)
841 , uses_dof_handler_vector(false)
844template <
int dim,
typename VectorType,
typename TransferType>
848 const TransferType &transfer)
849 : dof_handler_vector(dof_handler.size())
850 , dof_handler_vector_raw(dof_handler.size())
852 , transfer(&transfer)
853 , uses_dof_handler_vector(true)
855 for (
unsigned int i = 0; i < dof_handler.size(); ++i)
862template <
int dim,
typename VectorType,
typename TransferType>
869template <
int dim,
typename VectorType,
typename TransferType>
870template <
typename OtherVectorType>
873 OtherVectorType &dst,
874 const OtherVectorType &src)
const
876 internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,
881 uses_dof_handler_vector,
887template <
int dim,
typename VectorType,
typename TransferType>
890 const unsigned int block)
const
893 return dof_handler_vector[block]->locally_owned_dofs();
897template <
int dim,
typename VectorType,
typename TransferType>
900 const unsigned int block)
const
903 return dof_handler_vector[block]->locally_owned_dofs();
908template <
int dim,
typename VectorType,
typename TransferType>
923template <
int dim,
typename VectorType,
typename TransferType>
924boost::signals2::connection
926 const std::function<
void(
bool)> &slot)
928 return this->signals.transfer_to_mg.connect(slot);
933template <
int dim,
typename VectorType,
typename TransferType>
934boost::signals2::connection
936 const std::function<
void(
bool)> &slot)
938 return this->signals.transfer_to_global.connect(slot);
943template <
int dim,
typename VectorType,
typename TransferType>
944template <
typename OtherVectorType>
947 OtherVectorType &dst,
948 const OtherVectorType &src)
const
950 internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,
955 uses_dof_handler_vector,
961template <
int dim,
typename VectorType,
typename TransferType>
962template <
typename OtherVectorType>
966 const OtherVectorType &)
const
972template <
int dim,
typename VectorType,
typename TransferType>
973template <
typename OtherVectorType>
977 const OtherVectorType &)
const
983template <
int dim,
typename VectorType,
typename TransferType>
987 return *this->multigrid;
991template <
int dim,
typename VectorType,
typename TransferType>
995 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)
std::vector< const DoFHandler< dim > * > dof_handler_vector_raw
void Tvmult(OtherVectorType &dst, const OtherVectorType &src) const
SmartPointer< const TransferType, PreconditionMG< dim, VectorType, TransferType > > transfer
IndexSet locally_owned_domain_indices(const unsigned int block=0) const
const Multigrid< VectorType > & get_multigrid() const
boost::signals2::connection connect_transfer_to_global(const std::function< void(bool)> &slot)
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)
MPI_Comm get_mpi_communicator() const
const bool uses_dof_handler_vector
void vmult_add(OtherVectorType &dst, const OtherVectorType &src) const
IndexSet locally_owned_range_indices(const unsigned int block=0) const
SmartPointer< Multigrid< VectorType >, PreconditionMG< dim, VectorType, TransferType > > multigrid
PreconditionMG(const std::vector< const DoFHandler< dim > * > &dof_handler, Multigrid< VectorType > &mg, const TransferType &transfer)
Multigrid< VectorType > & get_multigrid()
std::vector< SmartPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TransferType > > > dof_handler_vector
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_communicator() const override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
#define DEAL_II_NOT_IMPLEMENTED()
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