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."
85 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
96 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
104 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
115 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
124 boost::signals2::signal<void(
const bool before,
const unsigned int level)>
151template <
typename VectorType>
293 boost::signals2::connection
295 const std::function<
void(
const bool,
const unsigned int)> &slot);
300 boost::signals2::connection
302 const std::function<
void(
const bool,
const unsigned int)> &slot);
307 boost::signals2::connection
309 const std::function<
void(
const bool,
const unsigned int)> &slot);
314 boost::signals2::connection
316 const std::function<
void(
const bool,
const unsigned int)> &slot);
321 boost::signals2::connection
323 const std::function<
void(
const bool,
const unsigned int)> &slot);
447 template <
int dim,
class OtherVectorType,
class TRANSFER>
466template <
int dim,
typename VectorType,
class TRANSFER>
498 template <
class OtherVectorType>
500 vmult(OtherVectorType &dst,
const OtherVectorType &src)
const;
506 template <
class OtherVectorType>
508 vmult_add(OtherVectorType &dst,
const OtherVectorType &src)
const;
515 template <
class OtherVectorType>
517 Tvmult(OtherVectorType &dst,
const OtherVectorType &src)
const;
524 template <
class OtherVectorType>
526 Tvmult_add(OtherVectorType &dst,
const OtherVectorType &src)
const;
555 boost::signals2::connection
561 boost::signals2::connection
568 std::vector<SmartPointer<const DoFHandler<dim>,
608template <
typename VectorType>
614 const unsigned int min_level,
615 const unsigned int max_level,
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())
629 maxlevel =
matrix.get_maxlevel();
631 maxlevel = max_level;
632 reinit(min_level, maxlevel);
637template <
typename VectorType>
646template <
typename VectorType>
659 namespace PreconditionMGImplementation
664 typename OtherVectorType>
665 typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
667 const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
669 const TRANSFER & transfer,
670 OtherVectorType & dst,
671 const OtherVectorType & src,
672 const bool uses_dof_handler_vector,
673 const typename ::mg::Signals &signals,
676 signals.transfer_to_mg(
true);
677 if (uses_dof_handler_vector)
678 transfer.copy_to_mg(dof_handler_vector, multigrid.
defect, src);
680 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
681 signals.transfer_to_mg(
false);
685 signals.transfer_to_global(
true);
686 if (uses_dof_handler_vector)
687 transfer.copy_from_mg(dof_handler_vector, dst, multigrid.
solution);
689 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.
solution);
690 signals.transfer_to_global(
false);
696 typename OtherVectorType>
699 const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
701 const TRANSFER & transfer,
702 OtherVectorType & dst,
703 const OtherVectorType & src,
704 const bool uses_dof_handler_vector,
705 const typename ::mg::Signals &signals,
708 (void)uses_dof_handler_vector;
711 signals.transfer_to_mg(
true);
712 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
713 signals.transfer_to_mg(
false);
717 signals.transfer_to_global(
true);
718 transfer.copy_from_mg(*dof_handler_vector[0], dst, multigrid.
solution);
719 signals.transfer_to_global(
false);
725 typename OtherVectorType>
726 typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
728 const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
730 const TRANSFER & transfer,
731 OtherVectorType & dst,
732 const OtherVectorType & src,
733 const bool uses_dof_handler_vector,
734 const typename ::mg::Signals &signals,
737 signals.transfer_to_mg(
true);
738 if (uses_dof_handler_vector)
739 transfer.copy_to_mg(dof_handler_vector, multigrid.
defect, src);
741 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
742 signals.transfer_to_mg(
false);
746 signals.transfer_to_global(
true);
747 if (uses_dof_handler_vector)
748 transfer.copy_from_mg_add(dof_handler_vector, dst, multigrid.
solution);
750 transfer.copy_from_mg_add(*dof_handler_vector[0],
753 signals.transfer_to_global(
false);
759 typename OtherVectorType>
762 const std::vector<const ::DoFHandler<dim> *> &dof_handler_vector,
764 const TRANSFER & transfer,
765 OtherVectorType & dst,
766 const OtherVectorType & src,
767 const bool uses_dof_handler_vector,
768 const typename ::mg::Signals &signals,
771 (void)uses_dof_handler_vector;
774 signals.transfer_to_mg(
true);
775 transfer.copy_to_mg(*dof_handler_vector[0], multigrid.
defect, src);
776 signals.transfer_to_mg(
false);
780 signals.transfer_to_global(
true);
781 transfer.copy_from_mg_add(*dof_handler_vector[0],
784 signals.transfer_to_global(
false);
789template <
int dim,
typename VectorType,
class TRANSFER>
793 const TRANSFER & transfer)
794 : dof_handler_vector(1, &dof_handler)
795 , dof_handler_vector_raw(1, &dof_handler)
797 , transfer(&transfer)
798 , uses_dof_handler_vector(false)
801template <
int dim,
typename VectorType,
class TRANSFER>
805 const TRANSFER & transfer)
806 : dof_handler_vector(dof_handler.size())
807 , dof_handler_vector_raw(dof_handler.size())
809 , transfer(&transfer)
810 , uses_dof_handler_vector(true)
812 for (
unsigned int i = 0; i < dof_handler.size(); ++i)
819template <
int dim,
typename VectorType,
class TRANSFER>
826template <
int dim,
typename VectorType,
class TRANSFER>
827template <
class OtherVectorType>
830 OtherVectorType & dst,
831 const OtherVectorType &src)
const
833 internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,
838 uses_dof_handler_vector,
844template <
int dim,
typename VectorType,
class TRANSFER>
847 const unsigned int block)
const
850 return dof_handler_vector[block]->locally_owned_dofs();
854template <
int dim,
typename VectorType,
class TRANSFER>
857 const unsigned int block)
const
860 return dof_handler_vector[block]->locally_owned_dofs();
865template <
int dim,
typename VectorType,
class TRANSFER>
880template <
int dim,
typename VectorType,
class TRANSFER>
881boost::signals2::connection
883 const std::function<
void(
bool)> &slot)
885 return this->signals.transfer_to_mg.connect(slot);
890template <
int dim,
typename VectorType,
class TRANSFER>
891boost::signals2::connection
893 const std::function<
void(
bool)> &slot)
895 return this->signals.transfer_to_global.connect(slot);
900template <
int dim,
typename VectorType,
class TRANSFER>
901template <
class OtherVectorType>
904 OtherVectorType & dst,
905 const OtherVectorType &src)
const
907 internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,
912 uses_dof_handler_vector,
918template <
int dim,
typename VectorType,
class TRANSFER>
919template <
class OtherVectorType>
922 const OtherVectorType &)
const
928template <
int dim,
typename VectorType,
class TRANSFER>
929template <
class OtherVectorType>
933 const OtherVectorType &)
const
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
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)
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
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
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)
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()
@ matrix
Contents is actually a matrix.
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
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)> post_smoother_step
boost::signals2::signal< void(const bool before)> transfer_to_global
boost::signals2::signal< void(const bool before, const unsigned int level)> restriction
boost::signals2::signal< void(const bool before)> transfer_to_mg