]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove the deprecated RefinementListener approach in the Triangulation class.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 4 Jan 2015 21:27:34 +0000 (15:27 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 9 Jan 2015 06:51:53 +0000 (00:51 -0600)
This method has been supplanted by the far more versatile signals approach.

doc/news/changes.h
include/deal.II/grid/tria.h
source/grid/tria.cc
tests/deal.II/refinement_listener_01.cc [deleted file]
tests/deal.II/refinement_listener_01.output [deleted file]

index 9a962faa5677ab999fdb54c26394195b3fa013d9..cc4cbc701c42843474d805e4b5abf7adc3b55742 100644 (file)
@@ -92,6 +92,10 @@ inconvenience this causes.
   </li>
 </ol>
 
+- The refinement listener concept of the Triangulation class. This
+  approach to getting notified about what happens to triangulations
+  has been superseded by the signals defined by the triangulation
+  class.
 
 <!-- ----------- GENERAL IMPROVEMENTS ----------------- -->
 
index ad344e9d734c0b8e7c9c81dc955b0c7e14621332..22d249d7a1bc3b4e84135787661862fd4ebd41cf 100644 (file)
@@ -1431,99 +1431,6 @@ public:
 
   typedef typename IteratorSelector::hex_iterator         hex_iterator;
   typedef typename IteratorSelector::active_hex_iterator  active_hex_iterator;
-
-  /**
-   * Base class for refinement listeners. Other classes, which need to be
-   * informed about refinements of the Triangulation, can be derived from
-   * RefinementListener.
-   *
-   * @note The use of this class has been superseded by the signals mechanism.
-   * See the general documentation of the Triangulation class for more
-   * information.
-   *
-   * @deprecated
-   */
-  class RefinementListener
-  {
-  public:
-    /**
-     * Destructor. Does nothing, but is declared virtual because this class
-     * also has virtual functions.
-     *
-     * @note The use of this class has been superseded by the signals
-     * mechanism. See the general documentation of the Triangulation class for
-     * more information.
-     *
-     * @deprecated
-     */
-    virtual ~RefinementListener ();
-
-    /**
-     * Before refinement is actually performed, the triangulation class calls
-     * this method on all objects derived from this class and registered with
-     * the triangulation.
-     *
-     * @note The use of this class has been superseded by the signals
-     * mechanism. See the general documentation of the Triangulation class for
-     * more information.
-     *
-     * @deprecated
-     */
-    virtual
-    void
-    pre_refinement_notification (const Triangulation<dim, spacedim> &tria);
-
-    /**
-     * After refinement is actually performed, the triangulation class calls
-     * this method on all objects derived from this class and registered with
-     * the triangulation.
-     *
-     * @note The use of this class has been superseded by the signals
-     * mechanism. See the general documentation of the Triangulation class for
-     * more information.
-     *
-     * @deprecated
-     */
-    virtual
-    void
-    post_refinement_notification (const Triangulation<dim, spacedim> &tria);
-
-    /**
-     * At the end of a call to copy_triangulation() the Triangulation class
-     * calls this method on all objects derived from this class and registered
-     * with the original Triangulation @p old_tria so that they might
-     * subscribe to the copied one @p new_tria as well, if that is desired. By
-     * default this method does nothing, a different behavior has to be
-     * implemented in derived classes.
-     *
-     * @note The use of this class has been superseded by the signals
-     * mechanism. See the general documentation of the Triangulation class for
-     * more information.
-     *
-     * @deprecated
-     */
-    virtual
-    void
-    copy_notification (const Triangulation<dim, spacedim> &old_tria,
-                       const Triangulation<dim, spacedim> &new_tria);
-
-    /**
-     * At the end of a call to create_triangulation() the Triangulation class
-     * calls this method on all objects derived from this class and registered
-     * with the current Triangulation object. By default this method does
-     * nothing, a different behavior has to be implemented in derived classes.
-     *
-     * @note The use of this class has been superseded by the signals
-     * mechanism. See the general documentation of the Triangulation class for
-     * more information.
-     *
-     * @deprecated
-     */
-    virtual
-    void
-    create_notification (const Triangulation<dim, spacedim> &tria);
-  };
-
   /**
    * A structure that is used as an exception object by the
    * create_triangulation() function to indicate which cells among the coarse
@@ -1982,31 +1889,6 @@ public:
    * @{
    */
 
-  /**
-   * Add a RefinementListener. Adding listeners to the Triangulation allows
-   * other classes to be informed when the Triangulation is refined.
-   *
-   * @note The use of this function has been superseded by the signals
-   * mechanism.  See the general documentation of the Triangulation class for
-   * more information.
-   *
-   * @deprecated
-   */
-  void add_refinement_listener (RefinementListener &listener) const DEAL_II_DEPRECATED;
-
-  /**
-   * Remove a RefinementListener. When some class needs no longer to be
-   * informed about refinements, the listener should be removed from the
-   * Triangulation.
-   *
-   * @note The use of this function has been superseded by the signals
-   * mechanism.  See the general documentation of the Triangulation class for
-   * more information.
-   *
-   * @deprecated
-   */
-  void remove_refinement_listener (RefinementListener &listener) const DEAL_II_DEPRECATED;
-
   /**
    * A structure that has boost::signal objects for a number of actions that a
    * triangulation can do to itself. See the general documentation of the
@@ -3252,24 +3134,7 @@ private:
    */
   std::map<unsigned int, types::manifold_id> *vertex_to_manifold_id_map_1d;
 
-
-  /**
-   * A map that correlates each refinement listener that has been added
-   * through the outdated RefinementListener interface via
-   * add_refinement_listener(), with the new-style boost::signal connections
-   * for each of the member function. We need to keep this list around so that
-   * we can later terminate the connection again when someone calls
-   * remove_refinement_listener().
-   *
-   * The data type is a multimap since, although this would be weird, the same
-   * object may add itself multiple times as a listener.
-   */
-  mutable
-  std::multimap<const RefinementListener *, std::vector<boost::signals2::connection> >
-  refinement_listener_map;
-
-  // make a couple of classes
-  // friends
+  // make a couple of classes friends
   template <int,int,int> friend class TriaAccessorBase;
   template <int,int,int> friend class TriaAccessor;
   friend class TriaAccessor<0, 1, spacedim>;
index 107b65eab21723ef439e313a5bf368bbe916b41a..de79df8f7aaa57fa75a358e9e592f0f620ed3580 100644 (file)
@@ -12997,12 +12997,6 @@ Triangulation<dim, spacedim>::memory_consumption () const
 
 
 
-template<int dim, int spacedim>
-Triangulation<dim, spacedim>::RefinementListener::~RefinementListener ()
-{}
-
-
-
 template<int dim, int spacedim>
 Triangulation<dim, spacedim>::DistortedCellList::~DistortedCellList () throw ()
 {
@@ -13012,93 +13006,6 @@ Triangulation<dim, spacedim>::DistortedCellList::~DistortedCellList () throw ()
 }
 
 
-
-
-template<int dim, int spacedim>
-void Triangulation<dim, spacedim>::
-RefinementListener::pre_refinement_notification (const Triangulation<dim, spacedim> &)
-{}
-
-
-
-template<int dim, int spacedim>
-void Triangulation<dim, spacedim>::
-RefinementListener::post_refinement_notification (const Triangulation<dim, spacedim> &)
-{}
-
-
-
-template<int dim, int spacedim>
-void Triangulation<dim, spacedim>::
-RefinementListener::copy_notification (const Triangulation<dim, spacedim> &,
-                                       const Triangulation<dim, spacedim> &)
-{}
-
-
-
-template<int dim, int spacedim>
-void Triangulation<dim, spacedim>::
-RefinementListener::create_notification (const Triangulation<dim, spacedim> &)
-{}
-
-
-
-template<int dim, int spacedim>
-void
-Triangulation<dim, spacedim>::add_refinement_listener (RefinementListener &listener) const
-{
-  // in this compatibility mode with the old-style refinement
-  // listeners, an external class presents itself as one that may or
-  // may not have overloaded all of the functions that the
-  // RefinementListener class has. consequently, we need to connect
-  // each of its functions to the relevant signals. for those
-  // functions that haven't been overloaded, that means that
-  // triggering the signal yields a call to the function in the
-  // RefinementListener base class which simply does nothing
-  std::vector<boost::signals2::connection> connections;
-
-  connections.push_back
-  (signals.create.connect (std_cxx11::bind (&RefinementListener::create_notification,
-                                            std_cxx11::ref(listener),
-                                            std_cxx11::cref(*this))));
-  connections.push_back
-  (signals.copy.connect (std_cxx11::bind (&RefinementListener::copy_notification,
-                                          std_cxx11::ref(listener),
-                                          std_cxx11::cref(*this),
-                                          std_cxx11::_1)));
-  connections.push_back
-  (signals.pre_refinement.connect (std_cxx11::bind (&RefinementListener::pre_refinement_notification,
-                                                    std_cxx11::ref(listener),
-                                                    std_cxx11::cref(*this))));
-  connections.push_back
-  (signals.post_refinement.connect (std_cxx11::bind (&RefinementListener::post_refinement_notification,
-                                                     std_cxx11::ref(listener),
-                                                     std_cxx11::cref(*this))));
-
-  // now push the set of connections into the map
-  refinement_listener_map.insert (std::make_pair(&listener, connections));
-}
-
-
-
-template<int dim, int spacedim>
-void
-Triangulation<dim, spacedim>::remove_refinement_listener (RefinementListener &listener) const
-{
-  Assert (refinement_listener_map.find (&listener) != refinement_listener_map.end(),
-          ExcMessage("You try to remove a refinement listener that does "
-                     "not appear to have been added previously."));
-
-  // get the element of the map, and terminate these connections. then
-  // erase the element from the list
-  std::vector<boost::signals2::connection> connections
-    = refinement_listener_map.find(&listener)->second;
-  for (unsigned int i=0; i<connections.size(); ++i)
-    connections[i].disconnect ();
-
-  refinement_listener_map.erase (refinement_listener_map.find (&listener));
-}
-
 template <>
 const Manifold<2,1> &Triangulation<2, 1>::get_manifold(const types::manifold_id) const
 {
diff --git a/tests/deal.II/refinement_listener_01.cc b/tests/deal.II/refinement_listener_01.cc
deleted file mode 100644 (file)
index 7d268c7..0000000
+++ /dev/null
@@ -1,123 +0,0 @@
-// ---------------------------------------------------------------------
-//
-// Copyright (C) 2006 - 2013 by the deal.II authors
-//
-// This file is part of the deal.II library.
-//
-// The deal.II library is free software; you can use it, redistribute
-// it, and/or modify it under the terms of the GNU Lesser General
-// Public License as published by the Free Software Foundation; either
-// version 2.1 of the License, or (at your option) any later version.
-// The full text of the license can be found in the file LICENSE at
-// the top level of the deal.II distribution.
-//
-// ---------------------------------------------------------------------
-
-
-
-// test the various refinement listener functions
-
-
-#include "../tests.h"
-#include <deal.II/grid/tria.h>
-#include <deal.II/grid/grid_generator.h>
-#include <deal.II/base/logstream.h>
-
-#include <fstream>
-#include <iomanip>
-#include <cstdio>
-
-std::ofstream logfile("output");
-
-
-template <int dim, int spacedim = dim>
-class RefinementListener :
-  public Triangulation<dim,spacedim>::RefinementListener
-{
-public:
-  RefinementListener (const std::string &prefix)
-    :
-    prefix(prefix)
-  {}
-
-  virtual
-  void
-  pre_refinement_notification (const Triangulation<dim, spacedim> &tria)
-  {
-    deallog << prefix << ' ' << "Pre-refinement: " << tria.n_active_cells() << std::endl;
-  }
-
-  virtual
-  void
-  post_refinement_notification (const Triangulation<dim, spacedim> &tria)
-  {
-    deallog << prefix << ' ' << "Post-refinement: " << tria.n_active_cells() << std::endl;
-  }
-
-  virtual
-  void
-  copy_notification (const Triangulation<dim, spacedim> &old_tria,
-                     const Triangulation<dim, spacedim> &new_tria)
-  {
-    deallog << prefix << ' ' << "Copy: "
-            << old_tria.n_active_cells() << ' '
-            << new_tria.n_active_cells() << std::endl;
-  }
-
-  virtual
-  void
-  create_notification (const Triangulation<dim, spacedim> &tria)
-  {
-    deallog << prefix << ' ' << "Create: " << tria.n_active_cells() << std::endl;
-  }
-
-private:
-  std::string prefix;
-};
-
-
-template <int dim>
-void test ()
-{
-  deallog << dim << "D" << std::endl;
-
-  Triangulation<dim> tria_1, tria_2;
-
-  GridGenerator::hyper_cube(tria_2);
-
-  RefinementListener<dim> rl_1 ("tria_1");
-  RefinementListener<dim> rl_2 ("tria_2");
-  tria_1.add_refinement_listener (rl_1);
-  tria_2.add_refinement_listener (rl_2);
-
-  // this should print the create note
-  GridGenerator::hyper_cube(tria_1);
-
-  // this should print the pre- and
-  // post-refinement note
-  tria_1.refine_global (1);
-
-  // this should print the copy note
-  tria_1.clear ();
-  tria_1.copy_triangulation (tria_2);
-
-  // no longer print anything
-  tria_1.remove_refinement_listener (rl_1);
-  tria_1.refine_global (2);
-}
-
-
-int main ()
-{
-  deallog << std::setprecision(2);
-  logfile << std::setprecision(2);
-  deallog.attach(logfile);
-  deallog.depth_console(0);
-  deallog.threshold_double(1.e-10);
-
-  test<1> ();
-  test<2> ();
-  test<3> ();
-
-  return 0;
-}
diff --git a/tests/deal.II/refinement_listener_01.output b/tests/deal.II/refinement_listener_01.output
deleted file mode 100644 (file)
index 4a4c71d..0000000
+++ /dev/null
@@ -1,19 +0,0 @@
-
-DEAL::1D
-DEAL::tria_1 Create: 1
-DEAL::tria_1 Pre-refinement: 1
-DEAL::tria_1 Post-refinement: 2
-DEAL::tria_2 Copy: 1 1
-DEAL::tria_1 Create: 1
-DEAL::2D
-DEAL::tria_1 Create: 1
-DEAL::tria_1 Pre-refinement: 1
-DEAL::tria_1 Post-refinement: 4
-DEAL::tria_2 Copy: 1 1
-DEAL::tria_1 Create: 1
-DEAL::3D
-DEAL::tria_1 Create: 1
-DEAL::tria_1 Pre-refinement: 1
-DEAL::tria_1 Post-refinement: 8
-DEAL::tria_2 Copy: 1 1
-DEAL::tria_1 Create: 1

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.