]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add test for project_boundary_values_curl_conforming
authorAlexander Grayver <agrayver@erdw.ethz.ch>
Wed, 13 Aug 2014 16:05:23 +0000 (18:05 +0200)
committerAlexander Grayver <agrayver@erdw.ethz.ch>
Wed, 13 Aug 2014 16:05:23 +0000 (18:05 +0200)
tests/deal.II/project_bv_curl_conf_02.cc [new file with mode: 0755]

diff --git a/tests/deal.II/project_bv_curl_conf_02.cc b/tests/deal.II/project_bv_curl_conf_02.cc
new file mode 100755 (executable)
index 0000000..b150d06
--- /dev/null
@@ -0,0 +1,131 @@
+// ---------------------------------------------------------------------\r
+// $Id$\r
+//\r
+// Copyright (C) 2014 by the Alexander Grayver & 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
+// The test checks that project_boundary_values_curl_conforming\r
+// works correctly for high-order FE_Nedelec elements used via\r
+// FESystem. This requires the produced constraints to be the same\r
+// for FE_Nedelec and FESystem(FE_Nedelec, 1).\r
+\r
+#include "../tests.h"\r
+#include <deal.II/base/function.h>\r
+#include <deal.II/base/logstream.h>\r
+#include <deal.II/dofs/dof_handler.h>\r
+#include <deal.II/grid/grid_generator.h>\r
+#include <deal.II/grid/tria.h>\r
+#include <deal.II/fe/fe_nedelec.h>\r
+#include <deal.II/fe/fe_system.h>\r
+#include <deal.II/lac/constraint_matrix.h>\r
+#include <deal.II/numerics/vector_tools.h>\r
+\r
+std::ofstream logfile ("output");\r
+\r
+template <int dim>\r
+class BoundaryFunction: public Function<dim>\r
+{\r
+public:\r
+  BoundaryFunction ();\r
+  virtual void vector_value (const Point<dim> &p, Vector<double> &values) const;\r
+};\r
+\r
+template <int dim>\r
+BoundaryFunction<dim>::BoundaryFunction (): Function<dim> (dim)\r
+{\r
+}\r
+\r
+template <int dim>\r
+void BoundaryFunction<dim>::vector_value (const Point<dim> &,\r
+                                          Vector<double> &values) const\r
+{\r
+  for (unsigned int d = 0; d < dim; ++d)\r
+    values (d) = d + 1.0;\r
+}\r
+\r
+template <int dim>\r
+void test_boundary_values (const FiniteElement<dim> &fe, ConstraintMatrix &constraints)\r
+{\r
+  Triangulation<dim> triangulation;\r
+  GridGenerator::subdivided_hyper_cube (triangulation, 2);\r
+  DoFHandler<dim> dof_handler (triangulation);\r
+  dof_handler.distribute_dofs (fe);\r
+  BoundaryFunction<dim> boundary_function;\r
+  constraints.clear ();\r
+  VectorTools::project_boundary_values_curl_conforming (dof_handler, 0, boundary_function, 0, constraints);\r
+  constraints.close ();\r
+}\r
+\r
+template <int dim>\r
+void test(unsigned order)\r
+{\r
+  deallog << "dim:" << dim << " order:" << order << "\t";\r
+\r
+  ConstraintMatrix constraints_fe, constraints_fes;\r
+\r
+  {\r
+    FE_Nedelec<3> fe (order);\r
+    test_boundary_values (fe, constraints_fe);\r
+  }\r
+\r
+  {\r
+    FESystem<3> fe(FE_Nedelec<3>(order),1);\r
+    test_boundary_values (fe, constraints_fes);\r
+  }\r
+\r
+  if(constraints_fes.n_constraints() == constraints_fe.n_constraints())\r
+  {\r
+    const IndexSet& lines = constraints_fes.get_local_lines ();\r
+\r
+    for(unsigned i = 0; i < lines.n_elements(); ++i)\r
+    {\r
+      if(!constraints_fe.is_constrained(lines.nth_index_in_set(i)))\r
+      {\r
+        deallog << "Failed" << std::endl;\r
+        return;\r
+      }\r
+\r
+      const std::vector<std::pair<types::global_dof_index,double> >& c1\r
+              = *constraints_fes.get_constraint_entries(lines.nth_index_in_set(i));\r
+      const std::vector<std::pair<types::global_dof_index,double> >& c2\r
+              = *constraints_fe.get_constraint_entries(lines.nth_index_in_set(i));\r
+\r
+      for(size_t j = 0; j < c1.size(); ++j)\r
+        if((c1[j].first != c2[j].first) || (fabs(c1[j].second - c2[j].second) > 1e-14))\r
+        {\r
+          deallog << "Failed" << std::endl;\r
+          return;\r
+        }\r
+    }\r
+  }\r
+  else\r
+    deallog << "Failed" << std::endl;\r
+\r
+  deallog << "OK" << std::endl;\r
+}\r
+\r
+int main ()\r
+{\r
+  deallog << std::setprecision (2);\r
+  deallog.attach (logfile);\r
+  deallog.depth_console (0);\r
+  deallog.threshold_double (1e-12);\r
+\r
+  test<2>(0);\r
+  test<2>(1);\r
+  test<2>(2);\r
+\r
+  test<3>(0);\r
+  test<3>(1);\r
+  //test<3>(2);\r
+}\r

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.