16 #ifndef dealii_multigrid_h 17 #define dealii_multigrid_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/subscriptor.h> 22 #include <deal.II/base/smartpointer.h> 23 #include <deal.II/dofs/dof_handler.h> 24 #include <deal.II/lac/sparse_matrix.h> 25 #include <deal.II/lac/vector.h> 26 #include <deal.II/multigrid/mg_base.h> 27 #include <deal.II/base/mg_level_object.h> 28 #include <deal.II/distributed/tria.h> 32 DEAL_II_NAMESPACE_OPEN
67 boost::signals2::signal<void (const bool before, const unsigned int level)>
coarse_solve;
74 boost::signals2::signal<void (const bool before, const unsigned int level)>
restriction;
81 boost::signals2::signal<void (const bool before, const unsigned int level)>
prolongation;
88 boost::signals2::signal<void (const bool before, const unsigned int level)>
pre_smoother_step;
123 template <
typename VectorType>
140 typedef VectorType vector_type;
141 typedef const VectorType const_vector_type;
277 bool relative =
false);
296 boost::signals2::connection
302 boost::signals2::connection
308 boost::signals2::connection
314 boost::signals2::connection
320 boost::signals2::connection
445 template <
int dim,
class OtherVectorType,
class TRANSFER>
friend class PreconditionMG;
465 template <
int dim,
typename VectorType,
class TRANSFER>
496 template <
class OtherVectorType>
497 void vmult (OtherVectorType &dst,
498 const OtherVectorType &src)
const;
504 template <
class OtherVectorType>
506 const OtherVectorType &src)
const;
513 template <
class OtherVectorType>
514 void Tvmult (OtherVectorType &dst,
515 const OtherVectorType &src)
const;
522 template <
class OtherVectorType>
524 const OtherVectorType &src)
const;
550 boost::signals2::connection
556 boost::signals2::connection
599 template <
typename VectorType>
607 const unsigned int min_level,
608 const unsigned int max_level,
613 matrix(&matrix, typeid(*this).name()),
614 coarse(&coarse, typeid(*this).name()),
615 transfer(&transfer, typeid(*this).name()),
616 pre_smooth(&pre_smooth, typeid(*this).name()),
617 post_smooth(&post_smooth, typeid(*this).name()),
618 edge_down(nullptr, typeid(*this).name()),
619 edge_up(nullptr, typeid(*this).name()),
622 const unsigned int dof_handler_max_level
625 maxlevel = dof_handler_max_level;
627 maxlevel = max_level;
629 reinit(minlevel, maxlevel);
634 template <
typename VectorType>
640 const unsigned int min_level,
641 const unsigned int max_level,
646 coarse(&coarse, typeid(*this).name()),
647 transfer(&transfer, typeid(*this).name()),
648 pre_smooth(&pre_smooth, typeid(*this).name()),
649 post_smooth(&post_smooth, typeid(*this).name()),
650 edge_out(nullptr, typeid(*this).name()),
651 edge_in(nullptr, typeid(*this).name()),
652 edge_down(nullptr, typeid(*this).name()),
653 edge_up(nullptr, typeid(*this).name()),
657 maxlevel =
matrix.get_maxlevel();
659 maxlevel = max_level;
660 reinit(min_level, maxlevel);
665 template <
typename VectorType>
675 template <
typename VectorType>
689 namespace PreconditionMGImplementation
691 template <
int dim,
typename VectorType,
class TRANSFER,
typename OtherVectorType>
692 typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
693 vmult(
const std::vector<const ::DoFHandler<dim>*> &dof_handler_vector,
695 const TRANSFER &transfer,
696 OtherVectorType &dst,
697 const OtherVectorType &src,
698 const bool uses_dof_handler_vector,
699 const typename ::mg::Signals &signals,
int)
701 signals.transfer_to_mg(
true);
702 if (uses_dof_handler_vector)
703 transfer.copy_to_mg(dof_handler_vector,
707 transfer.copy_to_mg(*dof_handler_vector[0],
710 signals.transfer_to_mg(
false);
714 signals.transfer_to_global(
true);
715 if (uses_dof_handler_vector)
716 transfer.copy_from_mg(dof_handler_vector,
720 transfer.copy_from_mg(*dof_handler_vector[0],
723 signals.transfer_to_global(
false);
726 template <
int dim,
typename VectorType,
class TRANSFER,
typename OtherVectorType>
728 vmult(
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, ...)
736 (void) uses_dof_handler_vector;
739 signals.transfer_to_mg(
true);
740 transfer.copy_to_mg(*dof_handler_vector[0],
743 signals.transfer_to_mg(
false);
747 signals.transfer_to_global(
true);
748 transfer.copy_from_mg(*dof_handler_vector[0],
751 signals.transfer_to_global(
false);
754 template <
int dim,
typename VectorType,
class TRANSFER,
typename OtherVectorType>
755 typename std::enable_if<TRANSFER::supports_dof_handler_vector>::type
756 vmult_add(
const std::vector<const ::DoFHandler<dim>*> &dof_handler_vector,
758 const TRANSFER &transfer,
759 OtherVectorType &dst,
760 const OtherVectorType &src,
761 const bool uses_dof_handler_vector,
762 const typename ::mg::Signals &signals,
int)
764 signals.transfer_to_mg(
true);
765 if (uses_dof_handler_vector)
766 transfer.copy_to_mg(dof_handler_vector,
770 transfer.copy_to_mg(*dof_handler_vector[0],
773 signals.transfer_to_mg(
false);
777 signals.transfer_to_global(
true);
778 if (uses_dof_handler_vector)
779 transfer.copy_from_mg_add(dof_handler_vector,
783 transfer.copy_from_mg_add(*dof_handler_vector[0],
786 signals.transfer_to_global(
false);
789 template <
int dim,
typename VectorType,
class TRANSFER,
typename OtherVectorType>
791 vmult_add(
const std::vector<const ::DoFHandler<dim>*> &dof_handler_vector,
793 const TRANSFER &transfer,
794 OtherVectorType &dst,
795 const OtherVectorType &src,
796 const bool uses_dof_handler_vector,
797 const typename ::mg::Signals &signals, ...)
799 (void) uses_dof_handler_vector;
802 signals.transfer_to_mg(
true);
803 transfer.copy_to_mg(*dof_handler_vector[0],
806 signals.transfer_to_mg(
false);
810 signals.transfer_to_global(
true);
811 transfer.copy_from_mg_add(*dof_handler_vector[0],
814 signals.transfer_to_global(
false);
819 template <
int dim,
typename VectorType,
class TRANSFER>
823 const TRANSFER &transfer)
825 dof_handler_vector(1,&dof_handler),
826 dof_handler_vector_raw(1,&dof_handler),
829 uses_dof_handler_vector(false)
832 template <
int dim,
typename VectorType,
class TRANSFER>
836 const TRANSFER &transfer)
838 dof_handler_vector(dof_handler.size()),
839 dof_handler_vector_raw(dof_handler.size()),
842 uses_dof_handler_vector(true)
844 for (
unsigned int i = 0; i< dof_handler.size() ; ++i)
846 dof_handler_vector[i] = dof_handler[i];
847 dof_handler_vector_raw[i] = dof_handler[i];
851 template <
int dim,
typename VectorType,
class TRANSFER>
858 template <
int dim,
typename VectorType,
class TRANSFER>
859 template <
class OtherVectorType>
862 (OtherVectorType &dst,
863 const OtherVectorType &src)
const 865 internal::PreconditionMGImplementation::vmult(dof_handler_vector_raw,*multigrid,*transfer,
866 dst,src,uses_dof_handler_vector,
871 template <
int dim,
typename VectorType,
class TRANSFER>
876 return dof_handler_vector[block]->locally_owned_dofs();
880 template <
int dim,
typename VectorType,
class TRANSFER>
885 return dof_handler_vector[block]->locally_owned_dofs();
890 template <
int dim,
typename VectorType,
class TRANSFER>
904 template <
int dim,
typename VectorType,
class TRANSFER>
905 boost::signals2::connection
909 return this->signals.transfer_to_mg.connect(slot);
914 template <
int dim,
typename VectorType,
class TRANSFER>
915 boost::signals2::connection
918 const std::function<
void (
bool)> &slot)
920 return this->signals.transfer_to_global.connect(slot);
925 template <
int dim,
typename VectorType,
class TRANSFER>
926 template <
class OtherVectorType>
929 (OtherVectorType &dst,
930 const OtherVectorType &src)
const 932 internal::PreconditionMGImplementation::vmult_add(dof_handler_vector_raw,*multigrid,*transfer,
933 dst,src,uses_dof_handler_vector,
938 template <
int dim,
typename VectorType,
class TRANSFER>
939 template <
class OtherVectorType>
943 const OtherVectorType &)
const 949 template <
int dim,
typename VectorType,
class TRANSFER>
950 template <
class OtherVectorType>
954 const OtherVectorType &)
const 961 DEAL_II_NAMESPACE_CLOSE
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > matrix
boost::signals2::connection connect_transfer_to_global(const std::function< void(bool)> &slot)
boost::signals2::signal< void(const bool before)> transfer_to_mg
SmartPointer< const TRANSFER, PreconditionMG< dim, VectorType, TRANSFER > > transfer
static const unsigned int invalid_unsigned_int
SmartPointer< const MGMatrixBase< VectorType > > edge_in
void vmult(OtherVectorType &dst, const OtherVectorType &src) const
MPI_Comm get_mpi_communicator() const
Contents is actually a matrix.
void set_edge_matrices(const MGMatrixBase< VectorType > &edge_out, const MGMatrixBase< VectorType > &edge_in)
unsigned int get_maxlevel() const
boost::signals2::signal< void(const bool before, const unsigned int level)> restriction
#define AssertIndexRange(index, range)
void set_maxlevel(const unsigned int)
MGLevelObject< VectorType > t
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::signal< void(const bool before, const unsigned int level)> coarse_solve
const bool uses_dof_handler_vector
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_up
void set_debug(const unsigned int)
void level_v_step(const unsigned int level)
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > post_smooth
IndexSet locally_owned_range_indices(const unsigned int block=0) const
#define Assert(cond, exc)
std::vector< SmartPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TRANSFER > > > dof_handler_vector
boost::signals2::signal< void(const bool before)> transfer_to_global
virtual MPI_Comm get_communicator() const
void Tvmult(OtherVectorType &dst, const OtherVectorType &src) const
SmartPointer< Multigrid< VectorType >, PreconditionMG< dim, VectorType, TRANSFER > > multigrid
boost::signals2::connection connect_restriction(const std::function< void(const bool, const unsigned int)> &slot)
unsigned int get_minlevel() const
SmartPointer< const MGMatrixBase< VectorType > > edge_out
void vmult_add(OtherVectorType &dst, const OtherVectorType &src) const
SmartPointer< const MGTransferBase< VectorType >, Multigrid< VectorType > > transfer
void set_minlevel(const unsigned int level, bool relative=false)
void set_edge_flux_matrices(const MGMatrixBase< VectorType > &edge_down, const MGMatrixBase< VectorType > &edge_up)
boost::signals2::connection connect_transfer_to_mg(const std::function< void(bool)> &slot)
void Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const
MGLevelObject< VectorType > solution
boost::signals2::connection connect_pre_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
boost::signals2::signal< void(const bool before, const unsigned int level)> prolongation
boost::signals2::connection connect_prolongation(const std::function< void(const bool, const unsigned int)> &slot)
std::vector< const DoFHandler< dim > * > dof_handler_vector_raw
boost::signals2::signal< void(const bool before, const unsigned int level)> post_smoother_step
IndexSet locally_owned_domain_indices(const unsigned int block=0) const
void level_step(const unsigned int level, Cycle cycle)
virtual unsigned int n_global_levels() const
boost::signals2::connection connect_post_smoother_step(const std::function< void(const bool, const unsigned int)> &slot)
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > pre_smooth
static ::ExceptionBase & ExcNotImplemented()
const Triangulation< dim, spacedim > & get_triangulation() const
Multigrid(const DoFHandler< dim > &mg_dof_handler, 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)
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_down
MGLevelObject< VectorType > defect2
MGLevelObject< VectorType > defect
PreconditionMG(const DoFHandler< dim > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer)
boost::signals2::signal< void(const bool before, const unsigned int level)> pre_smoother_step
static ::ExceptionBase & ExcInternalError()
Triangulation< dim, spacedim > & get_triangulation()
void reinit(const unsigned int minlevel, const unsigned int maxlevel)