]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use a std_cxx11::unique_ptr in FEValues.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sat, 8 Aug 2015 02:58:21 +0000 (21:58 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sat, 8 Aug 2015 19:31:17 +0000 (14:31 -0500)
Specifically, use it for the internal data objects of mapping and
finite element. These pointers are internally created and used in the
FEValues class. They are not shared by any external interface. Using
a unique_ptr ensures no memory leaks.

include/deal.II/fe/fe_values.h
source/fe/fe_values.cc

index 62d0f27bcf30555a068bc5b6e7d00eecb60b0a5c..179e007c5c1616b9e82747713b868a55b44a5605 100644 (file)
@@ -2334,12 +2334,12 @@ protected:
   /**
    * Internal data of mapping.
    */
-  SmartPointer<typename Mapping<dim,spacedim>::InternalDataBase,FEValuesBase<dim,spacedim> > mapping_data;
+  std_cxx11::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase> mapping_data;
 
   /**
    * Internal data of finite element.
    */
-  SmartPointer<typename FiniteElement<dim,spacedim>::InternalDataBase,FEValuesBase<dim,spacedim> > fe_data;
+  std_cxx11::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase> fe_data;
 
   /**
    * Original update flags handed to the constructor of FEValues.
index 40548fbae9a02ea58689a8b7a5630e030df38234..6f18dc0fbec5962485741b7e1ca66fdc3e4b30cb 100644 (file)
@@ -2161,8 +2161,6 @@ FEValuesBase<dim,spacedim>::FEValuesBase (const unsigned int n_q_points,
   dofs_per_cell (dofs_per_cell),
   mapping(&mapping, typeid(*this).name()),
   fe(&fe, typeid(*this).name()),
-  mapping_data(0, typeid(*this).name()),
-  fe_data(0, typeid(*this).name()),
   fe_values_views_cache (*this)
 {
   Assert (n_q_points > 0,
@@ -2177,25 +2175,6 @@ FEValuesBase<dim,spacedim>::FEValuesBase (const unsigned int n_q_points,
 template <int dim, int spacedim>
 FEValuesBase<dim,spacedim>::~FEValuesBase ()
 {
-  // delete those fields that were
-  // created by the mapping and
-  // finite element objects,
-  // respectively, but of which we
-  // have assumed ownership
-  if (fe_data != 0)
-    {
-      typename FiniteElement<dim,spacedim>::InternalDataBase *tmp1=fe_data;
-      fe_data=0;
-      delete tmp1;
-    }
-
-  if (mapping_data != 0)
-    {
-      typename Mapping<dim,spacedim>::InternalDataBase *tmp1=mapping_data;
-      mapping_data=0;
-      delete tmp1;
-    }
-
   tria_listener.disconnect ();
 }
 
@@ -3427,8 +3406,8 @@ FEValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
   // FE and the Mapping can store
   // intermediate data used across
   // calls to reinit
-  this->mapping_data = this->mapping->get_data(flags, quadrature);
-  this->fe_data      = this->fe->get_data(flags, *this->mapping, quadrature);
+  this->mapping_data.reset (this->mapping->get_data(flags, quadrature));
+  this->fe_data.reset (this->fe->get_data(flags, *this->mapping, quadrature));
 
   // initialize the base classes
   internal::FEValues::MappingRelatedData<dim,spacedim>::initialize(this->n_quadrature_points, flags);
@@ -3657,8 +3636,8 @@ FEFaceValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
   // FE and the Mapping can store
   // intermediate data used across
   // calls to reinit
-  this->mapping_data = this->mapping->get_face_data(flags, this->quadrature);
-  this->fe_data      = this->fe->get_face_data(flags, *this->mapping, this->quadrature);
+  this->mapping_data.reset (this->mapping->get_face_data(flags, this->quadrature));
+  this->fe_data.reset (this->fe->get_face_data(flags, *this->mapping, this->quadrature));
 
   // initialize the base classes
   internal::FEValues::MappingRelatedData<dim,spacedim>::initialize(this->n_quadrature_points, flags);
@@ -3806,10 +3785,10 @@ FESubfaceValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
   // FE and the Mapping can store
   // intermediate data used across
   // calls to reinit
-  this->mapping_data = this->mapping->get_subface_data(flags, this->quadrature);
-  this->fe_data      = this->fe->get_subface_data(flags,
+  this->mapping_data.reset (this->mapping->get_subface_data(flags, this->quadrature));
+  this->fe_data.reset (this->fe->get_subface_data(flags,
                                                   *this->mapping,
-                                                  this->quadrature);
+                                                  this->quadrature));
 
   // initialize the base classes
   internal::FEValues::MappingRelatedData<dim,spacedim>::initialize(this->n_quadrature_points, flags);

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.