]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add tests for twisted grid generator input. 2591/head
authorDavid Wells <wellsd2@rpi.edu>
Sat, 14 May 2016 19:34:06 +0000 (15:34 -0400)
committerDavid Wells <wellsd2@rpi.edu>
Sat, 14 May 2016 19:39:40 +0000 (15:39 -0400)
These tests verify that the right exceptions are raised when input data
would result in cells with negative measure.

tests/grid/twisted_parallelepiped_01.cc [new file with mode: 0644]
tests/grid/twisted_parallelepiped_01.debug.output [new file with mode: 0644]
tests/grid/twisted_parallelepiped_02.cc [new file with mode: 0644]
tests/grid/twisted_parallelepiped_02.release.output [new file with mode: 0644]

diff --git a/tests/grid/twisted_parallelepiped_01.cc b/tests/grid/twisted_parallelepiped_01.cc
new file mode 100644 (file)
index 0000000..a6ca228
--- /dev/null
@@ -0,0 +1,87 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/base/exceptions.h>
+#include <deal.II/base/point.h>
+#include <deal.II/base/std_cxx11/array.h>
+#include <deal.II/base/tensor.h>
+
+#include <fstream>
+
+/*
+ * Verify that the correct exception (ExcInvalidInputOrientation) is thrown
+ * when this function is called for invalid input.
+ */
+
+template<int dim>
+void check_parallelepiped (std::ofstream &logfile)
+{
+  // Data structure defining dim coordinates that make up a
+  // parallelepiped.
+  std_cxx11::array<Tensor<1, dim>, dim> edges;
+
+  switch (dim)
+    {
+    case 1:
+      edges[0][0] = -0.5;
+      break;
+
+    case 2:
+      edges[0][1] = 0.5;
+      edges[1][0] = 0.5;
+      break;
+
+    case 3:
+      edges[0][0] = 1.0;
+      edges[1][1] = 1.0;
+      edges[2][2] = -1.0;
+      break;
+
+    default:
+      Assert (false, ExcInternalError ());
+    }
+
+  Point<dim> origin;
+
+  Triangulation<dim> triangulation;
+  std::vector<unsigned int> subdivisions;
+
+  try
+    {
+      GridGenerator::subdivided_parallelepiped<dim>(triangulation, origin, edges,
+                                                    subdivisions, false);
+    }
+  catch (ExceptionBase &exc)
+    {
+      logfile << exc.get_exc_name() << std::endl;
+    }
+}
+
+int main ()
+{
+  deal_II_exceptions::disable_abort_on_exception();
+  std::ofstream logfile ("output");
+
+  check_parallelepiped<1>(logfile);
+  check_parallelepiped<2>(logfile);
+  check_parallelepiped<3>(logfile);
+  logfile << "OK" << std::endl;
+}
diff --git a/tests/grid/twisted_parallelepiped_01.debug.output b/tests/grid/twisted_parallelepiped_01.debug.output
new file mode 100644 (file)
index 0000000..3a8ca65
--- /dev/null
@@ -0,0 +1,4 @@
+ExcInvalidInputOrientation ("The triangulation you are trying to create will consist of cells" " with negative measures. This is usually the result of input data" " that does not define a right-handed coordinate system. The usual" " fix for this is to ensure that in 1D the given point is to the" " right of the origin (or the given edge tensor is positive), in 2D" " that the two edges (and their cross product) obey the right-hand" " rule (which may usually be done by switching the order of the" " points or edge tensors), or in 3D that the edges form a" " right-handed coordinate system (which may also be accomplished by" " switching the order of the first two points or edge tensors).")
+ExcInvalidInputOrientation ("The triangulation you are trying to create will consist of cells" " with negative measures. This is usually the result of input data" " that does not define a right-handed coordinate system. The usual" " fix for this is to ensure that in 1D the given point is to the" " right of the origin (or the given edge tensor is positive), in 2D" " that the two edges (and their cross product) obey the right-hand" " rule (which may usually be done by switching the order of the" " points or edge tensors), or in 3D that the edges form a" " right-handed coordinate system (which may also be accomplished by" " switching the order of the first two points or edge tensors).")
+ExcInvalidInputOrientation ("The triangulation you are trying to create will consist of cells" " with negative measures. This is usually the result of input data" " that does not define a right-handed coordinate system. The usual" " fix for this is to ensure that in 1D the given point is to the" " right of the origin (or the given edge tensor is positive), in 2D" " that the two edges (and their cross product) obey the right-hand" " rule (which may usually be done by switching the order of the" " points or edge tensors), or in 3D that the edges form a" " right-handed coordinate system (which may also be accomplished by" " switching the order of the first two points or edge tensors).")
+OK
diff --git a/tests/grid/twisted_parallelepiped_02.cc b/tests/grid/twisted_parallelepiped_02.cc
new file mode 100644 (file)
index 0000000..c8a37f0
--- /dev/null
@@ -0,0 +1,87 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/base/exceptions.h>
+#include <deal.II/base/point.h>
+#include <deal.II/base/std_cxx11/array.h>
+#include <deal.II/base/tensor.h>
+
+#include <fstream>
+
+/*
+ * Verify that the correct exception (ExcGridHasInvalidCell) is thrown when
+ * this function is called for invalid input.
+ */
+
+template<int dim>
+void check_parallelepiped (std::ofstream &logfile)
+{
+  // Data structure defining dim coordinates that make up a
+  // parallelepiped.
+  std_cxx11::array<Tensor<1, dim>, dim> edges;
+
+  switch (dim)
+    {
+    case 1:
+      edges[0][0] = -0.5;
+      break;
+
+    case 2:
+      edges[0][1] = 0.5;
+      edges[1][0] = 0.5;
+      break;
+
+    case 3:
+      edges[0][0] = 1.0;
+      edges[1][1] = 1.0;
+      edges[2][2] = -1.0;
+      break;
+
+    default:
+      Assert (false, ExcInternalError ());
+    }
+
+  Point<dim> origin;
+
+  Triangulation<dim> triangulation;
+  std::vector<unsigned int> subdivisions;
+
+  try
+    {
+      GridGenerator::subdivided_parallelepiped<dim>(triangulation, origin, edges,
+                                                    subdivisions, false);
+    }
+  catch (ExceptionBase &exc)
+    {
+      logfile << exc.get_exc_name() << std::endl;
+    }
+}
+
+int main ()
+{
+  deal_II_exceptions::disable_abort_on_exception();
+  std::ofstream logfile ("output");
+
+  check_parallelepiped<1>(logfile);
+  check_parallelepiped<2>(logfile);
+  check_parallelepiped<3>(logfile);
+  logfile << "OK" << std::endl;
+}
diff --git a/tests/grid/twisted_parallelepiped_02.release.output b/tests/grid/twisted_parallelepiped_02.release.output
new file mode 100644 (file)
index 0000000..b4453d8
--- /dev/null
@@ -0,0 +1,4 @@
+ExcGridHasInvalidCell(cell_no)
+ExcGridHasInvalidCell(cell_no)
+ExcGridHasInvalidCell(cell_no)
+OK

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.