ExcMessage("The Triangulation you are using is empty!"));
Assert(ff.size() > 0, ExcMessage("The given hp::FECollection is empty!"));
+#ifdef DEBUG
+ // make sure that the provided FE collection is large enough to
+ // cover all FE indices presently in use on the mesh
+ if ((hp_cell_active_fe_indices.size() > 0) &&
+ (hp_cell_future_fe_indices.size() > 0))
+ {
+ Assert(hp_capability_enabled, ExcInternalError());
+
+ for (const auto &cell :
+ this->active_cell_iterators() | IteratorFilters::LocallyOwnedCell())
+ {
+ Assert(cell->active_fe_index() < ff.size(),
+ ExcInvalidFEIndex(cell->active_fe_index(), ff.size()));
+ Assert(cell->future_fe_index() < ff.size(),
+ ExcInvalidFEIndex(cell->future_fe_index(), ff.size()));
+ }
+ }
+#endif
+
//
// register the new finite element collection
//
// on both its own cells and all ghost cells
internal::hp::DoFHandlerImplementation::Implementation::
communicate_active_fe_indices(*this);
-
-#ifdef DEBUG
- // make sure that the FE collection is large enough to
- // cover all FE indices presently in use on the mesh
- for (const auto &cell : this->active_cell_iterators())
- {
- if (!cell->is_artificial())
- Assert(cell->active_fe_index() < this->fe_collection.size(),
- ExcInvalidFEIndex(cell->active_fe_index(),
- this->fe_collection.size()));
- if (cell->is_locally_owned())
- Assert(cell->future_fe_index() < this->fe_collection.size(),
- ExcInvalidFEIndex(cell->future_fe_index(),
- this->fe_collection.size()));
- }
-#endif
}
{
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2023 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// warn users that they should have nonzero active FE indices in non-hp mode
+
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+test()
+{
+ Triangulation<dim> triangulation;
+ GridGenerator::hyper_cube(triangulation);
+
+ FE_Q<dim> fe(1);
+
+ DoFHandler<dim> dof_handler(triangulation);
+
+ // Choose index out of range here.
+ dof_handler.begin_active()->set_active_fe_index(1);
+
+ try
+ {
+ dof_handler.distribute_dofs(fe);
+ }
+ catch (ExceptionBase &e)
+ {
+ deallog << e.what() << std::endl;
+ }
+}
+
+
+int
+main()
+{
+ deal_II_exceptions::disable_abort_on_exception();
+
+ initlog();
+
+ test<2>();
+
+ deallog << "OK" << std::endl;
+}
--- /dev/null
+
+DEAL::
+--------------------------------------------------------
+An error occurred in file <dof_handler.cc> in function
+ void dealii::DoFHandler<dim, spacedim>::distribute_dofs(const dealii::hp::FECollection<dim, spacedim>&) [with int dim = 2; int spacedim = 2]
+The violated condition was:
+ cell->active_fe_index() < ff.size()
+Additional information:
+ The mesh contains a cell with an active FE index of 1, but the finite
+ element collection only has 1 elements
+--------------------------------------------------------
+
+DEAL::OK