]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use make_unique/shared in dofs
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sat, 10 Feb 2018 14:46:59 +0000 (15:46 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 15 Feb 2018 21:00:31 +0000 (22:00 +0100)
source/dofs/dof_handler.cc
source/dofs/dof_handler_policy.cc

index 502af9ab8fc9c1dee8be2267aab1f5f20b912ed4..a519a22583e97e578e243c1c76221016ea8bb80d 100644 (file)
@@ -319,7 +319,7 @@ namespace internal
                      numbers::invalid_dof_index);
           }
 
-        dof_handler.faces.reset (new internal::DoFHandler::DoFFaces<2>);
+        dof_handler.faces = std_cxx14::make_unique<internal::DoFHandler::DoFFaces<2>> ();
         // avoid access to n_raw_lines when there are no cells
         if (dof_handler.tria->n_cells() > 0)
           {
@@ -354,7 +354,7 @@ namespace internal
                      dof_handler.get_fe().dofs_per_cell,
                      numbers::invalid_dof_index);
           }
-        dof_handler.faces.reset (new internal::DoFHandler::DoFFaces<3>);
+        dof_handler.faces = std_cxx14::make_unique<internal::DoFHandler::DoFFaces<3>> ();
 
         // avoid access to n_raw_lines when there are no cells
         if (dof_handler.tria->n_cells() > 0)
@@ -439,11 +439,11 @@ namespace internal
 
         for (unsigned int i = 0; i < n_levels; ++i)
           {
-            dof_handler.mg_levels.emplace_back (new internal::DoFHandler::DoFLevel<2>);
+            dof_handler.mg_levels.emplace_back (std_cxx14::make_unique<internal::DoFHandler::DoFLevel<2>> ());
             dof_handler.mg_levels.back ()->dof_object.dofs = std::vector<types::global_dof_index> (tria.n_raw_quads (i) * fe.dofs_per_quad, numbers::invalid_dof_index);
           }
 
-        dof_handler.mg_faces.reset (new internal::DoFHandler::DoFFaces<2>);
+        dof_handler.mg_faces = std_cxx14::make_unique<internal::DoFHandler::DoFFaces<2>> ();
         dof_handler.mg_faces->lines.dofs = std::vector<types::global_dof_index> (tria.n_raw_lines () * fe.dofs_per_line, numbers::invalid_dof_index);
 
         const unsigned int &n_vertices = tria.n_vertices ();
@@ -498,11 +498,11 @@ namespace internal
 
         for (unsigned int i = 0; i < n_levels; ++i)
           {
-            dof_handler.mg_levels.emplace_back (new internal::DoFHandler::DoFLevel<3>);
+            dof_handler.mg_levels.emplace_back (std_cxx14::make_unique<internal::DoFHandler::DoFLevel<3>> ());
             dof_handler.mg_levels.back ()->dof_object.dofs = std::vector<types::global_dof_index> (tria.n_raw_hexs (i) * fe.dofs_per_hex, numbers::invalid_dof_index);
           }
 
-        dof_handler.mg_faces.reset (new internal::DoFHandler::DoFFaces<3>);
+        dof_handler.mg_faces = std_cxx14::make_unique<internal::DoFHandler::DoFFaces<3>> ();
         dof_handler.mg_faces->lines.dofs = std::vector<types::global_dof_index> (tria.n_raw_lines () * fe.dofs_per_line, numbers::invalid_dof_index);
         dof_handler.mg_faces->quads.dofs = std::vector<types::global_dof_index> (tria.n_raw_quads () * fe.dofs_per_quad, numbers::invalid_dof_index);
 
@@ -661,13 +661,13 @@ DoFHandler<dim,spacedim>::DoFHandler (const Triangulation<dim,spacedim> &tria)
   if (dynamic_cast<const parallel::shared::Triangulation< dim, spacedim>*>
       (&tria)
       != nullptr)
-    policy.reset (new internal::DoFHandler::Policy::ParallelShared<DoFHandler<dim,spacedim> >(*this));
+    policy = std_cxx14::make_unique<internal::DoFHandler::Policy::ParallelShared<DoFHandler<dim,spacedim> >> (*this);
   else if (dynamic_cast<const parallel::distributed::Triangulation< dim, spacedim >*>
            (&tria)
            == nullptr)
-    policy.reset (new internal::DoFHandler::Policy::Sequential<DoFHandler<dim,spacedim> >(*this));
+    policy = std_cxx14::make_unique<internal::DoFHandler::Policy::Sequential<DoFHandler<dim,spacedim> >> (*this);
   else
-    policy.reset (new internal::DoFHandler::Policy::ParallelDistributed<DoFHandler<dim,spacedim> >(*this));
+    policy = std_cxx14::make_unique<internal::DoFHandler::Policy::ParallelDistributed<DoFHandler<dim,spacedim> >> (*this);
 }
 
 
@@ -707,11 +707,11 @@ initialize(const Triangulation<dim,spacedim> &t,
 
   // decide whether we need a sequential or a parallel distributed policy
   if (dynamic_cast<const parallel::shared::Triangulation< dim, spacedim>*> (&t) != nullptr)
-    policy.reset (new internal::DoFHandler::Policy::ParallelShared<DoFHandler<dim,spacedim> >(*this));
+    policy = std_cxx14::make_unique<internal::DoFHandler::Policy::ParallelShared<DoFHandler<dim,spacedim> >> (*this);
   else if (dynamic_cast<const parallel::distributed::Triangulation< dim, spacedim >*> (&t) != nullptr)
-    policy.reset (new internal::DoFHandler::Policy::ParallelDistributed<DoFHandler<dim,spacedim> >(*this));
+    policy = std_cxx14::make_unique<internal::DoFHandler::Policy::ParallelDistributed<DoFHandler<dim,spacedim> >> (*this);
   else
-    policy.reset (new internal::DoFHandler::Policy::Sequential<DoFHandler<dim,spacedim> >(*this));
+    policy = std_cxx14::make_unique<internal::DoFHandler::Policy::Sequential<DoFHandler<dim,spacedim> >> (*this);
 
   distribute_dofs(fe);
 }
@@ -1308,7 +1308,7 @@ void DoFHandler<dim, spacedim>::MGVertexDoFs::init (const unsigned int cl,
       const unsigned int n_levels = finest_level - coarsest_level + 1;
       const unsigned int n_indices = n_levels * dofs_per_vertex;
 
-      indices.reset (new types::global_dof_index[n_indices]);
+      indices = std_cxx14::make_unique<types::global_dof_index[]> (n_indices);
       std::fill (indices.get(), indices.get()+n_indices,
                  numbers::invalid_dof_index);
     }
index 35132174b87dcf60666b3c9ab21cc506ad959850..8b20aa899be69d71dd132b43ea81233f2c62bc61 100644 (file)
@@ -40,6 +40,7 @@
 #include <boost/serialization/array.hpp>
 #endif
 
+#include <memory>
 #include <set>
 #include <algorithm>
 #include <numeric>
@@ -160,7 +161,7 @@ namespace internal
         void
         ensure_existence_of_dof_identities (const FiniteElement<dim,spacedim> &fe1,
                                             const FiniteElement<dim,spacedim> &fe2,
-                                            std::shared_ptr<DoFIdentities> &identities)
+                                            std::unique_ptr<DoFIdentities> &identities)
         {
           // see if we need to fill this entry, or whether it already
           // exists
@@ -170,25 +171,22 @@ namespace internal
                 {
                 case 0:
                 {
-                  identities =
-                    std::shared_ptr<DoFIdentities>
-                    (new DoFIdentities(fe1.hp_vertex_dof_identities(fe2)));
+                  identities = std_cxx14::make_unique<DoFIdentities>
+                               (fe1.hp_vertex_dof_identities(fe2));
                   break;
                 }
 
                 case 1:
                 {
-                  identities =
-                    std::shared_ptr<DoFIdentities>
-                    (new DoFIdentities(fe1.hp_line_dof_identities(fe2)));
+                  identities = std_cxx14::make_unique<DoFIdentities>
+                               (fe1.hp_line_dof_identities(fe2));
                   break;
                 }
 
                 case 2:
                 {
-                  identities =
-                    std::shared_ptr<DoFIdentities>
-                    (new DoFIdentities(fe1.hp_quad_dof_identities(fe2)));
+                  identities = std_cxx14::make_unique<DoFIdentities>
+                               (fe1.hp_quad_dof_identities(fe2));
                   break;
                 }
 
@@ -606,7 +604,7 @@ namespace internal
           // it is not clear whether this is actually necessary for
           // vertices at all, I can't think of a finite element that
           // would make that necessary...
-          dealii::Table<2,std::shared_ptr<DoFIdentities> >
+          dealii::Table<2,std::unique_ptr<DoFIdentities> >
           vertex_dof_identities (dof_handler.get_fe_collection().size(),
                                  dof_handler.get_fe_collection().size());
 
@@ -773,7 +771,7 @@ namespace internal
           // is to first treat all pairs of finite elements that have *identical* dofs,
           // and then only deal with those that are not identical of which we can
           // handle at most 2
-          dealii::Table<2,std::shared_ptr<DoFIdentities> >
+          dealii::Table<2,std::unique_ptr<DoFIdentities> >
           line_dof_identities (dof_handler.finite_elements->size(),
                                dof_handler.finite_elements->size());
 
@@ -1005,7 +1003,7 @@ namespace internal
           // trouble. note that this only happens for lines in 3d and
           // higher, and for quads only in 4d and higher, so this
           // isn't a particularly frequent case
-          dealii::Table<2,std::shared_ptr<DoFIdentities> >
+          dealii::Table<2,std::unique_ptr<DoFIdentities> >
           quad_dof_identities (dof_handler.finite_elements->size(),
                                dof_handler.finite_elements->size());
 

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.