]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Remove memory leak and generally improve infrastructure.
authorwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 29 Nov 2001 15:07:08 +0000 (15:07 +0000)
committerwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 29 Nov 2001 15:07:08 +0000 (15:07 +0000)
git-svn-id: https://svn.dealii.org/trunk@5310 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/multigrid/mg_dof_handler.h
deal.II/deal.II/source/multigrid/mg_dof_handler.cc

index 56442cfb6571910aa53a4d7e691a7c649af27fc1..70a82b66f647891c4c7a6294d1a5ba3803eb5535 100644 (file)
@@ -197,6 +197,13 @@ class MGDoFHandler : public DoFHandler<dim>
     virtual void distribute_dofs (const FiniteElement<dim> &, 
                                  const unsigned int offset = 0);
 
+                                    /**
+                                     * Clear all data of this object
+                                     * and call the respective
+                                     * function of the base class.
+                                     */
+    virtual void clear ();
+
                                     /**
                                      * Actually do the renumbering based on
                                      * a list of new dof numbers for all the
@@ -907,24 +914,27 @@ class MGDoFHandler : public DoFHandler<dim>
     };
 
 
-/**
- * Distribute dofs on the given cell,
- * with new dofs starting with index
- * @p{next_free_dof}. Return the next
- * unused index number. The finite
- * element used is the one given to
- * @p{distribute_dofs}, which is copied
- * to @p{selected_fe}.
- *
- * This function is excluded from the
- * @p{distribute_dofs} function since
- * it can not be implemented dimension
- * independent.
- *
- * Note that unlike for the usual dofs,
- * here all cells and not only active
- * ones are allowed.
- */
+                                    /**
+                                     * Distribute dofs on the given
+                                     * cell, with new dofs starting
+                                     * with index
+                                     * @p{next_free_dof}. Return the
+                                     * next unused index number. The
+                                     * finite element used is the one
+                                     * given to @p{distribute_dofs},
+                                     * which is copied to
+                                     * @p{selected_fe}.
+                                     *
+                                     * This function is excluded from
+                                     * the @p{distribute_dofs}
+                                     * function since it can not be
+                                     * implemented dimension
+                                     * independent.
+                                     *
+                                     * Note that unlike for the usual
+                                     * dofs, here all cells and not
+                                     * only active ones are allowed.
+                                     */
     unsigned int distribute_dofs_on_cell (cell_iterator &cell,
                                          unsigned int   next_free_dof);
     
@@ -934,6 +944,11 @@ class MGDoFHandler : public DoFHandler<dim>
                                      */
     void reserve_space ();
 
+                                    /**
+                                     * Free all used memory.
+                                     */
+    void clear_space ();
+    
                                     /**
                                      * Space to store the DoF numbers for the
                                      * different levels. Unlike the @p{levels}
index d6fc4e1dc94f145687bdb4bee2c49cd577555933..ac0b253e75a84edb7899e367bb0c8d5f492864a1 100644 (file)
@@ -85,7 +85,10 @@ MGDoFHandler<dim>::MGDoFHandler (Triangulation<dim> &tria) :
 
 
 template <int dim>
-MGDoFHandler<dim>::~MGDoFHandler () {};
+MGDoFHandler<dim>::~MGDoFHandler ()
+{
+  clear ();
+};
 
 
 #if deal_II_dimension == 1
@@ -1181,7 +1184,7 @@ void MGDoFHandler<dim>::distribute_dofs (const FiniteElement<dim> &fe,
   DoFHandler<dim>::distribute_dofs (fe, offset);
 
 
-// reserve space for the MG dof numbers
+                                  // reserve space for the MG dof numbers
   reserve_space ();
   mg_used_dofs.resize (tria->n_levels(), 0);
 
@@ -1369,6 +1372,19 @@ MGDoFHandler<3>::distribute_dofs_on_cell (cell_iterator &cell,
 #endif
 
 
+template <int dim>
+void
+MGDoFHandler<dim>::clear ()
+{
+                                  // release own memory
+  clear_space ();
+
+                                  // let base class release its mem
+                                  // as well
+  DoFHandler<dim>::clear ();  
+};
+
+
 template <int dim>
 unsigned int MGDoFHandler<dim>::n_dofs (const unsigned int level) const {
   Assert (level < mg_used_dofs.size(), ExcInvalidLevel(level));
@@ -1503,19 +1519,7 @@ void MGDoFHandler<1>::reserve_space () {
 
                                   //////////////////////////
                                   // DESTRUCTION
-  
-                                   // delete all levels and set them up
-                                   // newly, since vectors are
-                                   // troublesome if you want to change
-                                   // their size
-  for (unsigned int i=0; i<mg_levels.size(); ++i)
-    delete mg_levels[i];
-  mg_levels.resize (0);
-
-                                  // also delete vector of vertex indices
-                                  // this calls the destructor which
-                                  // must free the space
-  mg_vertex_dofs.resize (0);
+  clear_space ();
 
                                   ////////////////////////////
                                   // CONSTRUCTION
@@ -1588,20 +1592,7 @@ void MGDoFHandler<2>::reserve_space () {
   
                                   ////////////////////////////
                                   // DESTRUCTION
-
-                                   // delete all levels and set them up
-                                   // newly, since vectors are
-                                   // troublesome if you want to change
-                                   // their size
-  for (unsigned int i=0; i<mg_levels.size(); ++i)
-    delete mg_levels[i];
-  mg_levels.clear ();
-
-                                  // also delete vector of vertex indices
-                                  // this calls the destructor which
-                                  // must free the space
-  mg_vertex_dofs.clear ();
-
+  clear_space ();
 
                                   ////////////////////////////
                                   // CONSTRUCTION
@@ -1679,20 +1670,7 @@ void MGDoFHandler<3>::reserve_space () {
   
                                   ////////////////////////////
                                   // DESTRUCTION
-
-                                   // delete all levels and set them up
-                                   // newly, since vectors are
-                                   // troublesome if you want to change
-                                   // their size
-  for (unsigned int i=0; i<mg_levels.size(); ++i)
-    delete mg_levels[i];
-  mg_levels.clear ();
-
-                                  // also delete vector of vertex indices
-                                  // this calls the destructor which
-                                  // must free the space
-  mg_vertex_dofs.clear ();
-
+  clear_space ();
 
                                   ////////////////////////////
                                   // CONSTRUCTION
@@ -1762,5 +1740,25 @@ void MGDoFHandler<3>::reserve_space () {
 #endif
 
 
+template <int dim>
+void MGDoFHandler<dim>::clear_space ()
+{
+                                   // delete all levels and set them up
+                                   // newly, since vectors are
+                                   // troublesome if you want to change
+                                   // their size
+  for (unsigned int i=0; i<mg_levels.size(); ++i)
+    delete mg_levels[i];
+  mg_levels.clear ();
+
+                                  // also delete vector of vertex
+                                  // indices
+  typename std::vector<MGVertexDoFs> tmp;
+  std::swap (mg_vertex_dofs, tmp);
+};
+
+
+
 // explicit instantiations
 template class MGDoFHandler<deal_II_dimension>;
+

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.