From: Guido Kanschat Date: Wed, 4 Apr 2018 15:04:12 +0000 (+0200) Subject: Move Signals from Multigrid class to mg namespace and use more common names for transfer X-Git-Tag: v9.0.0-rc1~235 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F6129%2Fhead;p=dealii.git Move Signals from Multigrid class to mg namespace and use more common names for transfer --- diff --git a/include/deal.II/multigrid/multigrid.h b/include/deal.II/multigrid/multigrid.h index dfc811ae2d..57e0f329c1 100644 --- a/include/deal.II/multigrid/multigrid.h +++ b/include/deal.II/multigrid/multigrid.h @@ -34,102 +34,106 @@ DEAL_II_NAMESPACE_OPEN /*!@addtogroup mg */ /*@{*/ -/** - * Implementation of the multigrid method. The implementation supports both - * continuous and discontinuous elements and follows the procedure described in - * the @ref mg_paper "multigrid paper by Janssen and Kanschat". - * - * The function which starts a multigrid cycle on the finest level is cycle(). - * Depending on the cycle type chosen with the constructor (see enum Cycle), - * this function triggers one of the cycles level_v_step() or level_step(), - * where the latter one can do different types of cycles. - * - * Using this class, it is expected that the right hand side has been - * converted from a vector living on the locally finest level to a multilevel - * vector. This is a nontrivial operation, usually initiated automatically by - * the class PreconditionMG and performed by the classes derived from - * MGTransferBase. - * - * @note The interface of this class is still very clumsy. In particular, you - * will have to set up quite a few auxiliary objects before you can use it. - * Unfortunately, it seems that this can be avoided only be restricting the - * flexibility of this class in an unacceptable way. - * - * @author Guido Kanschat, 1999 - 2005 - */ -template -class Multigrid : public Subscriptor +namespace mg { -public: /** - * List of implemented cycle types. - */ - enum Cycle - { - /// The V-cycle - v_cycle, - /// The W-cycle - w_cycle, - /// The F-cycle - f_cycle - }; - - /** - * A structure that has boost::signal objects for a number of actions that - * happen for the typical usage of this class. + * A structure containing boost::signal objects for optional processing in multigrid solvers. * - * For documentation on signals, see - * http://www.boost.org/doc/libs/release/libs/signals2 . + * Each of these signals is called twice, once before and once after + * the action is performed. The two function calls differ in the + * boolean argument @p before, which is true the first time and + * false the second. */ struct Signals { /** - * This signal is triggered before (@p start is true) and after (@p start is + * This signal is triggered before (@p before is true) and after (@p before is * false) the call to MGTransfer::copy_to_mg which transfers the vector * given to it to a multi-level vector. */ - boost::signals2::signal transfer_to_mg; + boost::signals2::signal transfer_to_mg; /** - * This signal is triggered before (@p start is true) and after (@p start is + * This signal is triggered before (@p before is true) and after (@p before is * false) the call to MGTransfer::copy_from_mg which transfers the * multi-level vector given to it to a normal vector. */ - boost::signals2::signal transfer_to_global; + boost::signals2::signal transfer_to_global; /** - * This signal is triggered before (@p start is true) and after (@p start is + * This signal is triggered before (@p before is true) and after (@p before is * false) the call to the coarse solver on @p level. */ - boost::signals2::signal coarse_solve; + boost::signals2::signal coarse_solve; /** - * This signal is triggered before (@p start is true) and after (@p start is + * This signal is triggered before (@p before is true) and after (@p before is * false) the call to MGTransfer::restrict_and_add() which restricts a * vector from @p level to the next coarser one (@p level - 1). */ - boost::signals2::signal transfer_down_to; + boost::signals2::signal restriction; /** - * This signal is triggered before (@p start is true) and after (@p start is + * This signal is triggered before (@p before is true) and after (@p before is * false) the call to MGTransfer::prolongate() which prolongs a vector to * @p level from the next coarser one (@p level - 1). */ - boost::signals2::signal transfer_up_to; + boost::signals2::signal prolongation; /** - * This signal is triggered before (@p start is true) and after (@p start is + * This signal is triggered before (@p before is true) and after (@p before is * false) the call to a pre-smoothing step via MGPreSmoother::apply() on * @p level. */ - boost::signals2::signal pre_smoother_step; + boost::signals2::signal pre_smoother_step; /** - * This signal is triggered before (@p start is true) and after (@p start is + * This signal is triggered before (@p before is true) and after (@p before is * false) the call to a post-smoothing step via MGPostSmoother::apply() on * @p level. */ - boost::signals2::signal post_smoother_step; + boost::signals2::signal post_smoother_step; + }; +} + +/** + * Implementation of the multigrid method. The implementation supports both + * continuous and discontinuous elements and follows the procedure described in + * the @ref mg_paper "multigrid paper by Janssen and Kanschat". + * + * The function which starts a multigrid cycle on the finest level is cycle(). + * Depending on the cycle type chosen with the constructor (see enum Cycle), + * this function triggers one of the cycles level_v_step() or level_step(), + * where the latter one can do different types of cycles. + * + * Using this class, it is expected that the right hand side has been + * converted from a vector living on the locally finest level to a multilevel + * vector. This is a nontrivial operation, usually initiated automatically by + * the class PreconditionMG and performed by the classes derived from + * MGTransferBase. + * + * @note The interface of this class is still very clumsy. In particular, you + * will have to set up quite a few auxiliary objects before you can use it. + * Unfortunately, it seems that this can be avoided only be restricting the + * flexibility of this class in an unacceptable way. + * + * @author Guido Kanschat, 1999 - 2005 + */ +template +class Multigrid : public Subscriptor +{ +public: + /** + * List of implemented cycle types. + */ + enum Cycle + { + /// The V-cycle + v_cycle, + /// The W-cycle + w_cycle, + /// The F-cycle + f_cycle }; typedef VectorType vector_type; @@ -277,37 +281,40 @@ public: void set_cycle(Cycle); /** + * @deprecated Debug output will go away. Use signals instead. + * * Set the debug level. Higher values will create more debugging output * during the multigrid cycles. */ + DEAL_II_DEPRECATED void set_debug (const unsigned int); /** - * Connect a function to Signals::coarse_solve. + * Connect a function to mg::Signals::coarse_solve. */ boost::signals2::connection connect_coarse_solve(const std::function &slot); /** - * Connect a function to Signals::transfer_down_to. + * Connect a function to mg::Signals::restriction. */ boost::signals2::connection - connect_transfer_down_to(const std::function &slot); + connect_restriction(const std::function &slot); /** - * Connect a function to Signals::transfer_up_to. + * Connect a function to mg::Signals::prolongation. */ boost::signals2::connection - connect_transfer_up_to(const std::function &slot); + connect_prolongation(const std::function &slot); /** - * Connect a function to Signals::pre_smoother_step. + * Connect a function to mg::Signals::pre_smoother_step. */ boost::signals2::connection connect_pre_smoother_step(const std::function &slot); /** - * Connect a function to Signals::post_smoother_step. + * Connect a function to mg::Signals::post_smoother_step. */ boost::signals2::connection connect_post_smoother_step(const std::function &slot); @@ -317,7 +324,7 @@ private: /** * Signals for the various actions that the Multigrid algorithm uses. */ - Signals signals; + mg::Signals signals; /** * The V-cycle multigrid method. level is the level the function @@ -537,13 +544,13 @@ public: MPI_Comm get_mpi_communicator() const; /** - * Connect a function to Signals::transfer_to_mg. + * Connect a function to mg::Signals::transfer_to_mg. */ boost::signals2::connection connect_transfer_to_mg(const std::function &slot); /** - * Connect a function to Signals::transfer_to_global. + * Connect a function to mg::Signals::transfer_to_global. */ boost::signals2::connection connect_transfer_to_global(const std::function &slot); @@ -579,7 +586,7 @@ private: /** * Signals used by this object */ - typename Multigrid::Signals signals; + mg::Signals signals; }; /*@}*/ @@ -688,7 +695,7 @@ namespace internal OtherVectorType &dst, const OtherVectorType &src, const bool uses_dof_handler_vector, - const typename dealii::Multigrid::Signals &signals, int) + const typename dealii::mg::Signals &signals, int) { signals.transfer_to_mg(true); if (uses_dof_handler_vector) @@ -723,7 +730,7 @@ namespace internal OtherVectorType &dst, const OtherVectorType &src, const bool uses_dof_handler_vector, - const typename dealii::Multigrid::Signals &signals, ...) + const typename dealii::mg::Signals &signals, ...) { (void) uses_dof_handler_vector; Assert (!uses_dof_handler_vector, ExcInternalError()); @@ -751,9 +758,9 @@ namespace internal OtherVectorType &dst, const OtherVectorType &src, const bool uses_dof_handler_vector, - const typename dealii::Multigrid::Signals &signals, int) + const typename dealii::mg::Signals &signals, int) { - signals.transfer_to_mg(0, true); + signals.transfer_to_mg(true); if (uses_dof_handler_vector) transfer.copy_to_mg(dof_handler_vector, multigrid.defect, @@ -762,11 +769,11 @@ namespace internal transfer.copy_to_mg(*dof_handler_vector[0], multigrid.defect, src); - signals.transfer_to_mg(0, false); + signals.transfer_to_mg(false); multigrid.cycle(); - signals.transfer_to_global(0, true); + signals.transfer_to_global(true); if (uses_dof_handler_vector) transfer.copy_from_mg_add(dof_handler_vector, dst, @@ -775,7 +782,7 @@ namespace internal transfer.copy_from_mg_add(*dof_handler_vector[0], dst, multigrid.solution); - signals.transfer_to_global(0, false); + signals.transfer_to_global(false); } template @@ -786,7 +793,7 @@ namespace internal OtherVectorType &dst, const OtherVectorType &src, const bool uses_dof_handler_vector, - const typename dealii::Multigrid::Signals &signals, ...) + const typename dealii::mg::Signals &signals, ...) { (void) uses_dof_handler_vector; Assert (!uses_dof_handler_vector, ExcInternalError()); diff --git a/include/deal.II/multigrid/multigrid.templates.h b/include/deal.II/multigrid/multigrid.templates.h index 3d384e5fd4..d3f3fb418d 100644 --- a/include/deal.II/multigrid/multigrid.templates.h +++ b/include/deal.II/multigrid/multigrid.templates.h @@ -160,17 +160,17 @@ Multigrid::level_v_step (const unsigned int level) defect[level-1] -= t[level-1]; } - this->signals.transfer_down_to(true, level); + this->signals.restriction(true, level); transfer->restrict_and_add(level, defect[level-1], t[level]); - this->signals.transfer_down_to(false, level); + this->signals.restriction(false, level); // do recursion level_v_step(level-1); // do coarse grid correction - this->signals.transfer_up_to(true, level); + this->signals.prolongation(true, level); transfer->prolongate(level, t[level], solution[level-1]); - this->signals.transfer_up_to(false, level); + this->signals.prolongation(false, level); if (debug>2) deallog << "Prolongate norm " << t[level].l2_norm() << std::endl; @@ -394,9 +394,9 @@ connect_coarse_solve(const std::function template boost::signals2::connection Multigrid:: -connect_transfer_down_to(const std::function &slot) +connect_restriction(const std::function &slot) { - return this->signals.transfer_down_to.connect(slot); + return this->signals.restriction.connect(slot); } @@ -404,9 +404,9 @@ connect_transfer_down_to(const std::function boost::signals2::connection Multigrid:: -connect_transfer_up_to(const std::function &slot) +connect_prolongation(const std::function &slot) { - return this->signals.transfer_up_to.connect(slot); + return this->signals.prolongation.connect(slot); } diff --git a/tests/multigrid/events_01.cc b/tests/multigrid/events_01.cc index 3f00a66802..37937feae3 100644 --- a/tests/multigrid/events_01.cc +++ b/tests/multigrid/events_01.cc @@ -402,16 +402,16 @@ namespace Step50 }; - auto print_transfer_down_to = [](const bool start, const unsigned int level) + auto print_restriction = [](const bool start, const unsigned int level) { - deallog << "transfer_down_to " << (start?"started":"finished") + deallog << "restriction " << (start?"started":"finished") << " on level " << level << std::endl; }; - auto print_transfer_up_to = [](const bool start, const unsigned int level) + auto print_prolongation = [](const bool start, const unsigned int level) { - deallog << "transfer_up_to " << (start?"started":"finished") + deallog << "prolongation " << (start?"started":"finished") << " on level " << level << std::endl; }; @@ -436,8 +436,8 @@ namespace Step50 preconditioner.connect_transfer_to_global(print_transfer_to_global); mg.connect_coarse_solve(print_coarse_solve); - mg.connect_transfer_down_to(print_transfer_down_to); - mg.connect_transfer_up_to(print_transfer_up_to); + mg.connect_restriction(print_restriction); + mg.connect_prolongation(print_prolongation); mg.connect_pre_smoother_step(print_pre_smoother_step); mg.connect_post_smoother_step(print_post_smoother_step); diff --git a/tests/multigrid/events_01.with_mpi=true.with_p4est=true.with_trilinos=true.mpirun=1.output b/tests/multigrid/events_01.with_mpi=true.with_p4est=true.with_trilinos=true.mpirun=1.output index bf4004d436..b44455da11 100644 --- a/tests/multigrid/events_01.with_mpi=true.with_p4est=true.with_trilinos=true.mpirun=1.output +++ b/tests/multigrid/events_01.with_mpi=true.with_p4est=true.with_trilinos=true.mpirun=1.output @@ -6,20 +6,20 @@ DEAL::transfer_to_mg started DEAL::transfer_to_mg finished DEAL::pre_smoother_step started on level 2 DEAL::pre_smoother_step finished on level 2 -DEAL::transfer_down_to started on level 2 -DEAL::transfer_down_to finished on level 2 +DEAL::restriction started on level 2 +DEAL::restriction finished on level 2 DEAL::pre_smoother_step started on level 1 DEAL::pre_smoother_step finished on level 1 -DEAL::transfer_down_to started on level 1 -DEAL::transfer_down_to finished on level 1 +DEAL::restriction started on level 1 +DEAL::restriction finished on level 1 DEAL::coarse_solve started on level 0 DEAL::coarse_solve finished on level 0 -DEAL::transfer_up_to started on level 1 -DEAL::transfer_up_to finished on level 1 +DEAL::prolongation started on level 1 +DEAL::prolongation finished on level 1 DEAL::post_smoother_step started on level 1 DEAL::post_smoother_step finished on level 1 -DEAL::transfer_up_to started on level 2 -DEAL::transfer_up_to finished on level 2 +DEAL::prolongation started on level 2 +DEAL::prolongation finished on level 2 DEAL::post_smoother_step started on level 2 DEAL::post_smoother_step finished on level 2 DEAL::transfer_to_global started @@ -33,20 +33,20 @@ DEAL::transfer_to_mg started DEAL::transfer_to_mg finished DEAL::pre_smoother_step started on level 2 DEAL::pre_smoother_step finished on level 2 -DEAL::transfer_down_to started on level 2 -DEAL::transfer_down_to finished on level 2 +DEAL::restriction started on level 2 +DEAL::restriction finished on level 2 DEAL::pre_smoother_step started on level 1 DEAL::pre_smoother_step finished on level 1 -DEAL::transfer_down_to started on level 1 -DEAL::transfer_down_to finished on level 1 +DEAL::restriction started on level 1 +DEAL::restriction finished on level 1 DEAL::coarse_solve started on level 0 DEAL::coarse_solve finished on level 0 -DEAL::transfer_up_to started on level 1 -DEAL::transfer_up_to finished on level 1 +DEAL::prolongation started on level 1 +DEAL::prolongation finished on level 1 DEAL::post_smoother_step started on level 1 DEAL::post_smoother_step finished on level 1 -DEAL::transfer_up_to started on level 2 -DEAL::transfer_up_to finished on level 2 +DEAL::prolongation started on level 2 +DEAL::prolongation finished on level 2 DEAL::post_smoother_step started on level 2 DEAL::post_smoother_step finished on level 2 DEAL::transfer_to_global started @@ -60,20 +60,20 @@ DEAL::transfer_to_mg started DEAL::transfer_to_mg finished DEAL::pre_smoother_step started on level 2 DEAL::pre_smoother_step finished on level 2 -DEAL::transfer_down_to started on level 2 -DEAL::transfer_down_to finished on level 2 +DEAL::restriction started on level 2 +DEAL::restriction finished on level 2 DEAL::pre_smoother_step started on level 1 DEAL::pre_smoother_step finished on level 1 -DEAL::transfer_down_to started on level 1 -DEAL::transfer_down_to finished on level 1 +DEAL::restriction started on level 1 +DEAL::restriction finished on level 1 DEAL::coarse_solve started on level 0 DEAL::coarse_solve finished on level 0 -DEAL::transfer_up_to started on level 1 -DEAL::transfer_up_to finished on level 1 +DEAL::prolongation started on level 1 +DEAL::prolongation finished on level 1 DEAL::post_smoother_step started on level 1 DEAL::post_smoother_step finished on level 1 -DEAL::transfer_up_to started on level 2 -DEAL::transfer_up_to finished on level 2 +DEAL::prolongation started on level 2 +DEAL::prolongation finished on level 2 DEAL::post_smoother_step started on level 2 DEAL::post_smoother_step finished on level 2 DEAL::transfer_to_global started