Assert (tria->n_levels() > 0,
ExcMessage("The Triangulation you are using is empty!"));
- fe_collection = std_cxx14::make_unique<hp::FECollection<dim, spacedim>>(ff);
+ // Only recreate the FECollection if we don't already store
+ // the exact same FiniteElement object.
+ if (fe_collection == nullptr || &((*fe_collection)[0]) != &ff)
+ fe_collection = std_cxx14::make_unique<hp::FECollection<dim, spacedim>>(ff);
// delete all levels and set them
// up newly. note that we still
--- /dev/null
+//
+// Copyright (C) 2018 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// Test that calling hp::DoFhandler::distribute_dofs again with the same
+// hp::FECollection doesn't recreate the copy.
+
+
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+
+#include "../tests.h"
+
+int main()
+{
+ initlog();
+
+ constexpr int dim = 2;
+ constexpr int spacedim = 2;
+
+ Triangulation<dim,spacedim> tria;
+ GridGenerator::hyper_cube (tria);
+ FE_Q<dim,spacedim> fe(1);
+
+ DoFHandler<dim,spacedim> dh(tria);
+ dh.distribute_dofs(fe);
+
+ SmartPointer<const FiniteElement<dim,spacedim>> fe_p(&dh.get_fe());
+
+ dh.distribute_dofs(*fe_p);
+
+ deallog << "OK" << std::endl;
+}
--- /dev/null
+
+DEAL::OK
--- /dev/null
+//
+// Copyright (C) 2018 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// Test that calling hp::DoFhandler::distribute_dofs again with the same
+// hp::FECollection doesn't recreate the copy.
+
+
+#include <deal.II/hp/dof_handler.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+
+#include "../tests.h"
+
+int main()
+{
+ initlog();
+
+ constexpr int dim = 2;
+ constexpr int spacedim = 2;
+
+ Triangulation<dim,spacedim> tria;
+ GridGenerator::hyper_cube (tria);
+ FE_Q<dim,spacedim> fe(1);
+ hp::FECollection<dim, spacedim> fe_collection(fe);
+
+ hp::DoFHandler<dim,spacedim> dh(tria);
+ dh.distribute_dofs(fe_collection);
+
+ SmartPointer<const hp::FECollection<dim,spacedim>> fe_p(&dh.get_fe_collection());
+
+ dh.distribute_dofs(*fe_p);
+
+ deallog << "OK" << std::endl;
+}
--- /dev/null
+
+DEAL::OK