]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Avoid dangling pointer in initializer_list
authorDaniel Arndt <arndtd@ornl.gov>
Thu, 19 Sep 2019 02:07:14 +0000 (20:07 -0600)
committerDaniel Arndt <arndtd@ornl.gov>
Thu, 19 Sep 2019 02:07:14 +0000 (20:07 -0600)
include/deal.II/base/patterns.h
include/deal.II/hp/fe_collection.h
include/deal.II/hp/q_collection.h
tests/mpi/hp_unify_dof_indices_07.cc
tests/vector_tools/integrate_difference_04_hp.cc

index e369cacde0e1568abb9b8fa04472bcd7be6ffad4..413781ec151623b92ca8820819ed8e27fff1b737 100644 (file)
@@ -1439,7 +1439,7 @@ namespace Patterns
                   "are derived from PatternBase");
     static_assert(sizeof...(ps) > 0,
                   "The number of PatternTypes must be greater than zero!");
-    auto pattern_pointers = {(static_cast<const PatternBase *>(&ps))...};
+    const auto &pattern_pointers = {(static_cast<const PatternBase *>(&ps))...};
     for (const auto p : pattern_pointers)
       patterns.push_back(p->clone());
   }
index a422679dd4974d351445b80644858c7435edcaeb..613d9392ee459f63071915f3d2ce8532ff858e5e 100644 (file)
@@ -834,7 +834,8 @@ namespace hp
     // loop over all of the given arguments and add the finite elements to
     // this collection. Inlining the definition of fe_pointers causes internal
     // compiler errors on GCC 7.1.1 so we define it separately:
-    const auto fe_pointers = {&fes...};
+    const auto &fe_pointers =
+      std::initializer_list<const FiniteElement<dim, spacedim> *>{&fes...};
     for (const auto p : fe_pointers)
       push_back(*p);
   }
index 1bba0eff849205480fb12bf4974a2d553497d1fb..ecdc5d8384f8f2cde2cbf3cf2faeb5b489d89c11 100644 (file)
@@ -156,7 +156,8 @@ namespace hp
     // loop over all of the given arguments and add the quadrature objects to
     // this collection. Inlining the definition of q_pointers causes internal
     // compiler errors on GCC 7.1.1 so we define it separately:
-    const auto q_pointers = {&quadrature_objects...};
+    const auto &q_pointers =
+      std::initializer_list<const Quadrature<dim> *>{&quadrature_objects...};
     for (const auto p : q_pointers)
       push_back(*p);
   }
index 8bd1d7b6fdcc4ddbc53a78838a2529a2147ef745..4991bf1781d014c0117d859dff5a1417edf92331 100644 (file)
@@ -70,9 +70,7 @@ test()
   Assert(triangulation.n_global_active_cells() == 2, ExcInternalError());
   Assert(triangulation.n_active_cells() == 2, ExcInternalError());
 
-  hp::FECollection<dim> fe;
-  fe.push_back(FE_Q<dim>(2));
-  fe.push_back(FE_Q<dim>(2));
+  hp::FECollection<dim> fe(FE_Q<dim>(2), FE_Q<dim>(2));
 
   hp::DoFHandler<dim> dof_handler(triangulation);
   if (dof_handler.begin_active()->is_locally_owned())
index 531c56c5bec97c99378e8dadbd832cfdebd07b80..db4aecb441d9cc7e72a7445f2b3ee949f105f3cb 100644 (file)
@@ -101,8 +101,7 @@ test(VectorTools::NormType norm, double value, double exp = 2.0)
   solution = interpolated;
 
   Vector<double>       cellwise_errors(tria.n_active_cells());
-  hp::QCollection<dim> quadrature;
-  quadrature.push_back(QIterated<dim>(QTrapez<1>(), 5));
+  hp::QCollection<dim> quadrature(QIterated<dim>(QTrapez<1>(), 5));
 
   const dealii::Function<dim, double> *w = nullptr;
   VectorTools::integrate_difference(dofh,

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.