--- /dev/null
+// ---------------------------------------------------------------------\r
+// $Id: no_flux_11.cc 31349 2013-10-20 19:07:06Z maier $\r
+//\r
+// Copyright (C) 2012 - 2013 by the deal.II authors\r
+//\r
+// This file is part of the deal.II library.\r
+//\r
+// The deal.II library is free software; you can use it, redistribute\r
+// it, and/or modify it under the terms of the GNU Lesser General\r
+// Public License as published by the Free Software Foundation; either\r
+// version 2.1 of the License, or (at your option) any later version.\r
+// The full text of the license can be found in the file LICENSE at\r
+// the top level of the deal.II distribution.\r
+//\r
+// ---------------------------------------------------------------------\r
+\r
+// bug report from mailing list from 11/15/2013 (simplified). no_normal_flux throws\r
+// an exception:\r
+/*\r
+--------------------------------------------------------\r
+An error occurred in line <4314> of file </scratch/deal-trunk/deal.II/include/deal.II/numerics/vector_tools.templates.h> in function\r
+ void dealii::VectorTools::compute_no_normal_flux_constraints(const DH<dim, spacedim>&, unsigned int, const std::set<unsigned char>&, dealii::ConstraintMatrix&, const dealii::Mapping<dim, spacedim>&) [with int dim = 3; DH = dealii::DoFHandler; int spacedim = 3]\r
+The violated condition was: \r
+ vector_dofs.dof_indices[d] < dof_handler.n_dofs()\r
+The name and call sequence of the exception was:\r
+ ExcInternalError()\r
+Additional Information: \r
+(none)\r
+*/ \r
+\r
+#include "../tests.h"\r
+\r
+#include <deal.II/base/function.h>\r
+#include <deal.II/base/logstream.h>\r
+#include <deal.II/lac/constraint_matrix.h>\r
+#include <deal.II/lac/compressed_sparsity_pattern.h>\r
+#include <deal.II/grid/tria.h>\r
+#include <deal.II/grid/grid_generator.h>\r
+#include <deal.II/dofs/dof_handler.h>\r
+#include <deal.II/fe/fe_q.h>\r
+#include <deal.II/numerics/vector_tools.h>\r
+\r
+using namespace dealii;\r
+\r
+template <int dim>\r
+void test()\r
+{\r
+ Triangulation<dim> triangulation;\r
+ GridGenerator::hyper_cube (triangulation,-1.0,1.0);\r
+ triangulation.begin_active()->face(1)->set_all_boundary_indicators(1);\r
+\r
+ FE_Q<dim> fe(1);\r
+ DoFHandler<dim> dof_handler(triangulation);\r
+ dof_handler.distribute_dofs(fe);\r
+ \r
+ ConstraintMatrix constraints;\r
+ std::set<types::boundary_id> no_normal_flux_boundaries;\r
+ no_normal_flux_boundaries.insert (1);\r
+ VectorTools::compute_no_normal_flux_constraints (dof_handler,\r
+ 0,\r
+ no_normal_flux_boundaries,\r
+ constraints);\r
+\r
+ constraints.close();\r
+ constraints.print(deallog.get_file_stream ());\r
+}\r
+\r
+\r
+int main ()\r
+{\r
+ initlog();\r
+ deallog.depth_console (0);\r
+\r
+ test<3> ();\r
+}\r