]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use std::unique_ptr in hp::DoFHandler. 4551/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Mon, 26 Jun 2017 23:43:47 +0000 (17:43 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Mon, 26 Jun 2017 23:43:47 +0000 (17:43 -0600)
This replaces raw pointers. A similar change was made last week to the ::DoFHandler
class already.

include/deal.II/hp/dof_handler.h
source/hp/dof_handler.cc

index 93925d423c46b95dc58b03c9ac4219250e4d1889..81e4da9484c75e6bc9283149d3b0d2b3df06fab4 100644 (file)
@@ -249,11 +249,27 @@ namespace hp
      */
     DoFHandler (const Triangulation<dim,spacedim> &tria);
 
+    /**
+     * Copy constructor. DoFHandler objects are large and expensive.
+     * They should not be copied, in particular not by accident, but
+     * rather deliberately constructed. As a consequence, this constructor
+     * is explicitly removed from the interface of this class.
+     */
+    DoFHandler (const DoFHandler &) = delete;
+
     /**
      * Destructor.
      */
     virtual ~DoFHandler ();
 
+    /**
+     * Copy operator. DoFHandler objects are large and expensive.
+     * They should not be copied, in particular not by accident, but
+     * rather deliberately constructed. As a consequence, this operator
+     * is explicitly removed from the interface of this class.
+     */
+    DoFHandler &operator = (const DoFHandler &) = delete;
+
     /**
      * Go through the triangulation and "distribute" the degrees of freedoms
      * needed for the given finite element. "Distributing" degrees of freedom
@@ -720,32 +736,23 @@ namespace hp
 
   private:
 
-    /**
-     * Copy constructor. I can see no reason why someone might want to use it,
-     * so I don't provide it. Since this class has pointer members, making it
-     * private prevents the compiler to provide it's own, incorrect one if
-     * anyone chose to copy such an object.
-     */
-    DoFHandler (const DoFHandler &);
-
-    /**
-     * Copy operator. I can see no reason why someone might want to use it, so
-     * I don't provide it. Since this class has pointer members, making it
-     * private prevents the compiler to provide it's own, incorrect one if
-     * anyone chose to copy such an object.
-     */
-    DoFHandler &operator = (const DoFHandler &);
-
     /**
      * Free all used memory.
      */
     void clear_space ();
 
-    template<int structdim>
-    types::global_dof_index get_dof_index (const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index) const;
+    template <int structdim>
+    types::global_dof_index get_dof_index (const unsigned int obj_level,
+                                           const unsigned int obj_index,
+                                           const unsigned int fe_index,
+                                           const unsigned int local_index) const;
 
-    template<int structdim>
-    void set_dof_index (const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index) const;
+    template <int structdim>
+    void set_dof_index (const unsigned int obj_level,
+                        const unsigned int obj_index,
+                        const unsigned int fe_index,
+                        const unsigned int local_index,
+                        const types::global_dof_index global_index) const;
 
     /**
      * Create default tables for the active_fe_indices in the
@@ -814,13 +821,13 @@ namespace hp
      * Space to store the DoF numbers for the different levels. Analogous to
      * the <tt>levels[]</tt> tree of the Triangulation objects.
      */
-    std::vector<dealii::internal::hp::DoFLevel *> levels;
+    std::vector<std::unique_ptr<dealii::internal::hp::DoFLevel> > levels;
 
     /**
      * Space to store the DoF numbers for the faces. Analogous to the
      * <tt>faces</tt> pointer of the Triangulation objects.
      */
-    dealii::internal::hp::DoFIndicesOnFaces<dim> *faces;
+    std::unique_ptr<dealii::internal::hp::DoFIndicesOnFaces<dim> > faces;
 
     /**
      * A structure that contains all sorts of numbers that characterize the
@@ -865,7 +872,7 @@ namespace hp
      * refinement, i.e. from between when pre_refinement_action is called and
      * when post_refinement_action runs.
      */
-    std::vector<std::vector<bool> *> has_children;
+    std::vector<std::unique_ptr<std::vector<bool> > > has_children;
 
     /**
      * A list of connections with which this object connects to the
index 3bc530f43581d5c4741320d9b5e36f503c8486b6..e9fa9c081ba219da3be8f6b5cb5ab8b85959e6c7 100644 (file)
@@ -571,7 +571,7 @@ namespace internal
 
             for (unsigned int level=0; level<dof_handler.tria->n_levels(); ++level)
               {
-                dof_handler.levels.push_back (new internal::hp::DoFLevel);
+                dof_handler.levels.emplace_back (new internal::hp::DoFLevel);
                 std::swap (active_fe_backup[level],
                            dof_handler.levels[level]->active_fe_indices);
               }
@@ -699,11 +699,11 @@ namespace internal
 
             for (unsigned int level=0; level<dof_handler.tria->n_levels(); ++level)
               {
-                dof_handler.levels.push_back (new internal::hp::DoFLevel);
+                dof_handler.levels.emplace_back (new internal::hp::DoFLevel);
                 std::swap (active_fe_backup[level],
                            dof_handler.levels[level]->active_fe_indices);
               }
-            dof_handler.faces = new internal::hp::DoFIndicesOnFaces<2>;
+            dof_handler.faces.reset (new internal::hp::DoFIndicesOnFaces<2>);
           }
 
           // QUAD (CELL) DOFs
@@ -1076,11 +1076,11 @@ namespace internal
 
             for (unsigned int level=0; level<dof_handler.tria->n_levels(); ++level)
               {
-                dof_handler.levels.push_back (new internal::hp::DoFLevel);
+                dof_handler.levels.emplace_back (new internal::hp::DoFLevel);
                 std::swap (active_fe_backup[level],
                            dof_handler.levels[level]->active_fe_indices);
               }
-            dof_handler.faces = new internal::hp::DoFIndicesOnFaces<3>;
+            dof_handler.faces.reset (new internal::hp::DoFIndicesOnFaces<3>);
           }
 
           // HEX (CELL) DOFs
@@ -3148,7 +3148,7 @@ namespace hp
     // Create sufficiently many
     // hp::DoFLevels.
     while (levels.size () < tria->n_levels ())
-      levels.push_back (new dealii::internal::hp::DoFLevel);
+      levels.emplace_back (new dealii::internal::hp::DoFLevel);
 
     // then make sure that on each
     // level we have the appropriate
@@ -3199,8 +3199,7 @@ namespace hp
     for (unsigned int i=0; i<levels.size(); ++i)
       {
         const unsigned int cells_on_level = tria->n_raw_cells(i);
-        std::vector<bool> *has_children_level =
-          new std::vector<bool> (cells_on_level);
+        std::unique_ptr<std::vector<bool> > has_children_level(new std::vector<bool> (cells_on_level));
 
         // Check for each cell, if it has children. in 1d,
         // we don't store refinement cases, so use the 'children'
@@ -3220,7 +3219,7 @@ namespace hp
                                      std::placeholders::_1,
                                      static_cast<unsigned char>(RefinementCase<dim>::no_refinement)));
 
-        has_children.push_back (has_children_level);
+        has_children.emplace_back (std::move(has_children_level));
       }
   }
 
@@ -3235,13 +3234,13 @@ namespace hp
     // Normally only one level is added, but if this Triangulation
     // is created by copy_triangulation, it can be more than one level.
     while (levels.size () < tria->n_levels ())
-      levels.push_back (new dealii::internal::hp::DoFLevel);
+      levels.emplace_back (new dealii::internal::hp::DoFLevel);
 
     // Coarsening can lead to the loss
     // of levels. Hence remove them.
     while (levels.size () > tria->n_levels ())
       {
-        delete levels[levels.size ()-1];
+        // drop the last element. that also releases the memory pointed to
         levels.pop_back ();
       }
 
@@ -3306,11 +3305,6 @@ namespace hp
       }
 
     // Free buffer objects
-    std::vector<std::vector<bool> *>::iterator children_level;
-    for (children_level = has_children.begin ();
-         children_level != has_children.end ();
-         ++children_level)
-      delete (*children_level);
     has_children.clear ();
   }
 
@@ -3344,11 +3338,8 @@ namespace hp
   template<int dim, int spacedim>
   void DoFHandler<dim,spacedim>::clear_space ()
   {
-    for (unsigned int i=0; i<levels.size(); ++i)
-      delete levels[i];
-    levels.resize (0);
-    delete faces;
-    faces = nullptr;
+    levels.clear ();
+    faces.reset ();
 
     {
       std::vector<types::global_dof_index> tmp;

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.