]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use make_unique/shared in fe
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sat, 10 Feb 2018 14:22:44 +0000 (15:22 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 15 Feb 2018 21:00:31 +0000 (22:00 +0100)
include/deal.II/fe/fe_tools_interpolate.templates.h
source/fe/fe_q_bubbles.cc
source/fe/fe_values.cc
source/fe/mapping_c1.cc
source/fe/mapping_q.cc
source/fe/mapping_q_eulerian.cc

index 744b39d12b0ed704ddb5def50cbcff867af7ca7d..afd472b2e7c03b95ee02e6f2f6da47b1c0a9206e 100644 (file)
@@ -19,6 +19,7 @@
 
 #include <deal.II/base/quadrature_lib.h>
 #include <deal.II/base/qprojector.h>
+#include <deal.II/base/std_cxx14/memory.h>
 #include <deal.II/base/thread_management.h>
 #include <deal.II/base/utilities.h>
 #include <deal.II/lac/vector.h>
@@ -110,12 +111,12 @@ namespace FETools
     Vector<typename OutVector::value_type> u1_local(DoFTools::max_dofs_per_cell(dof1));
     Vector<typename OutVector::value_type> u2_local(DoFTools::max_dofs_per_cell(dof2));
 
-    // have a map for interpolation
-    // matrices. shared_ptr make sure
-    // that memory is released again
+    // have a map for interpolation matrices.
+    // Using a unique_ptr makes sure that the
+    // memory is released again automatically.
     std::map<const FiniteElement<dim,spacedim> *,
         std::map<const FiniteElement<dim,spacedim> *,
-        std::shared_ptr<FullMatrix<double> > > >
+        std::unique_ptr<FullMatrix<double> > > >
         interpolation_matrices;
 
     typename DoFHandlerType1<dim,spacedim>::active_cell_iterator cell1 = dof1.begin_active(),
@@ -176,15 +177,15 @@ namespace FETools
           // there
           if (interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()].get() == nullptr)
             {
-              std::shared_ptr<FullMatrix<double> >
-              interpolation_matrix (new FullMatrix<double> (dofs_per_cell2,
-                                                            dofs_per_cell1));
-              interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()]
-                = interpolation_matrix;
+              auto interpolation_matrix = std_cxx14::make_unique<FullMatrix<double> >
+                                          (dofs_per_cell2, dofs_per_cell1);
 
               get_interpolation_matrix(cell1->get_fe(),
                                        cell2->get_fe(),
                                        *interpolation_matrix);
+
+              interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()]
+                = std::move(interpolation_matrix);
             }
 
           cell1->get_dof_values(u1, u1_local);
@@ -291,7 +292,7 @@ namespace FETools
     // dof1 to the back_interpolation
     // matrices
     std::map<const FiniteElement<dim> *,
-        std::shared_ptr<FullMatrix<double> > > interpolation_matrices;
+        std::unique_ptr<FullMatrix<double> > > interpolation_matrices;
 
     for (; cell!=endc; ++cell)
       if ((cell->subdomain_id() == subdomain_id)
@@ -320,9 +321,8 @@ namespace FETools
           // matrix is available
           if (interpolation_matrices[&cell->get_fe()] == nullptr)
             {
-              interpolation_matrices[&cell->get_fe()] =
-                std::shared_ptr<FullMatrix<double> >
-                (new FullMatrix<double>(dofs_per_cell1, dofs_per_cell1));
+              interpolation_matrices[&cell->get_fe()] = std_cxx14::make_unique<FullMatrix<double> >
+                                                        (dofs_per_cell1, dofs_per_cell1);
               get_back_interpolation_matrix(cell->get_fe(), fe2,
                                             *interpolation_matrices[&cell->get_fe()]);
             }
index e57b3731a54e4161f62924957ab144aef50a25c2..e070b7cbecadbdf57adf892c2185c4abc4dfd6d7 100644 (file)
@@ -61,20 +61,20 @@ namespace internal
           {
           case 1:
             if (spacedim==1)
-              q_fine.reset(new QGauss<dim> (degree+1));
+              q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
             else if (spacedim==2)
-              q_fine.reset(new QAnisotropic<dim>(QGauss<1>(degree+1), q_dummy));
+              q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), q_dummy);
             else
-              q_fine.reset(new QAnisotropic<dim>(QGauss<1>(degree+1), q_dummy, q_dummy));
+              q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), q_dummy, q_dummy);
             break;
           case 2:
             if (spacedim==2)
-              q_fine.reset(new QGauss<dim> (degree+1));
+              q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
             else
-              q_fine.reset(new QAnisotropic<dim>(QGauss<1>(degree+1), QGauss<1>(degree+1), q_dummy));
+              q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), QGauss<1>(degree+1), q_dummy);
             break;
           case 3:
-            q_fine.reset(new QGauss<dim> (degree+1));
+            q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
             break;
           default:
             Assert(false, ExcInternalError());
index 130916c3338c658147eb8c8bc716e60c76e9bcaa..0810a05d1c40277e0a8398f8684ef1ae176edcd9 100644 (file)
@@ -19,6 +19,7 @@
 #include <deal.II/base/numbers.h>
 #include <deal.II/base/quadrature.h>
 #include <deal.II/base/signaling_nan.h>
+#include <deal.II/base/std_cxx14/memory.h>
 #include <deal.II/differentiation/ad.h>
 #include <deal.II/lac/vector.h>
 #include <deal.II/lac/block_vector.h>
@@ -39,7 +40,6 @@
 #include <deal.II/fe/fe.h>
 
 #include <iomanip>
-#include <memory>
 #include <type_traits>
 
 
@@ -3979,7 +3979,7 @@ FEValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
   if (flags & update_mapping)
     this->mapping_data.reset (mapping_get_data.return_value());
   else
-    this->mapping_data.reset (new typename Mapping<dim,spacedim>::InternalDataBase());
+    this->mapping_data = std_cxx14::make_unique<typename Mapping<dim,spacedim>::InternalDataBase> ();
 }
 
 
@@ -4010,7 +4010,7 @@ namespace
       }
     else
       // if the types don't match, there is nothing we can do here
-      present_cell.reset (new Type(new_cell));
+      present_cell = std_cxx14::make_unique<Type> (new_cell);
   }
 }
 
@@ -4229,7 +4229,7 @@ FEFaceValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
   if (flags & update_mapping)
     this->mapping_data.reset (mapping_get_data.return_value());
   else
-    this->mapping_data.reset (new typename Mapping<dim,spacedim>::InternalDataBase());
+    this->mapping_data = std_cxx14::make_unique<typename Mapping<dim,spacedim>::InternalDataBase> ();
 }
 
 
@@ -4397,7 +4397,7 @@ FESubfaceValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
   if (flags & update_mapping)
     this->mapping_data.reset (mapping_get_data.return_value());
   else
-    this->mapping_data.reset (new typename Mapping<dim,spacedim>::InternalDataBase());
+    this->mapping_data = std_cxx14::make_unique<typename Mapping<dim,spacedim>::InternalDataBase> ();
 }
 
 
index 0532e863afa380325df8d3089592754b0d64120a..a2a6266328b54a3fd05a97e157fb15aed8faf90f 100644 (file)
@@ -47,7 +47,7 @@ MappingC1<dim,spacedim>::MappingC1 ()
   //
   // we only need to replace the Qp mapping because that's the one that's
   // used on boundary cells where it matters
-  this->qp_mapping.reset (new MappingC1<dim,spacedim>::MappingC1Generic());
+  this->qp_mapping = std::make_shared<MappingC1<dim,spacedim>::MappingC1Generic>();
 }
 
 
index 25f641e662b684bd2c3c3bb7bc8f7b2dd6b75c85..122916c5d2d6dccf9bbd2f2adc2ef91a04db8fab 100644 (file)
@@ -72,13 +72,13 @@ MappingQ<dim,spacedim>::MappingQ (const unsigned int degree,
                               (dim != spacedim)),
   // create a Q1 mapping for use on interior cells (if necessary)
   // or to create a good initial guess in transform_real_to_unit_cell()
-  q1_mapping (new MappingQGeneric<dim,spacedim>(1)),
+  q1_mapping (std::make_shared<MappingQGeneric<dim,spacedim> >(1)),
 
   // create a Q_p mapping; if p=1, simply share the Q_1 mapping already
   // created via the shared_ptr objects
   qp_mapping (this->polynomial_degree>1
               ?
-              std::make_shared<const MappingQGeneric<dim,spacedim>>(degree)
+              std::make_shared<const MappingQGeneric<dim,spacedim> >(degree)
               :
               q1_mapping)
 {}
index f7ccc78f4a87c396aef0d9fc0eef96b8ce90d6c6..a1b315f17f70526d16cb7c53188bbe10d1432c8e 100644 (file)
@@ -70,11 +70,11 @@ MappingQEulerian (const unsigned int              degree,
   // reset the q1 mapping we use for interior cells (and previously
   // set by the MappingQ constructor) to a MappingQ1Eulerian with the
   // current vector
-  this->q1_mapping.reset (new MappingQ1Eulerian<dim,VectorType,spacedim>(euler_dof_handler,
-                          euler_vector));
+  this->q1_mapping = std::make_shared<MappingQ1Eulerian<dim,VectorType,spacedim>>
+                     (euler_dof_handler, euler_vector);
 
   // also reset the qp mapping pointer with our own class
-  this->qp_mapping.reset (new MappingQEulerianGeneric(degree,*this));
+  this->qp_mapping = std::make_shared<MappingQEulerianGeneric>(degree,*this);
 }
 
 

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.