]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add test case 11045/head
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Fri, 16 Oct 2020 10:24:19 +0000 (12:24 +0200)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Fri, 16 Oct 2020 15:23:43 +0000 (17:23 +0200)
tests/mappings/mapping_q_points_real_to_unit.cc [new file with mode: 0644]
tests/mappings/mapping_q_points_real_to_unit.output [new file with mode: 0644]

diff --git a/tests/mappings/mapping_q_points_real_to_unit.cc b/tests/mappings/mapping_q_points_real_to_unit.cc
new file mode 100644 (file)
index 0000000..cf747e3
--- /dev/null
@@ -0,0 +1,67 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 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.
+//
+// ---------------------------------------------------------------------
+
+
+// check MappingQGeneric::transform_real_to_unit_points on a set of
+// challenging points, especially with vectorization because we have nearby
+// points that succeed and others in the regime of negative Jacobian
+// determinants, respectively
+
+#include <deal.II/fe/mapping_q_generic.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+test(const unsigned int degree)
+{
+  Triangulation<dim> tria;
+  GridGenerator::hyper_shell(tria, Point<dim>(), 0.5, 1., 2 * dim);
+
+  std::vector<Point<dim>> real_points;
+  Point<dim>              p;
+  for (unsigned int d = 0; d < dim; ++d)
+    p[d] = std::sqrt(1. / dim) + 0.01 * d;
+  for (unsigned i = 1; i < 21; ++i)
+    real_points.push_back(
+      (i % 2 ? 0.05 * (static_cast<double>(i) - 10.) : 0.05 * i) * p);
+  std::vector<Point<dim>> unit_points(real_points.size());
+  MappingQGeneric<dim>    mapping(degree);
+  mapping.transform_points_real_to_unit_cell(tria.begin(),
+                                             real_points,
+                                             unit_points);
+
+  deallog << "Transform on cell with center: " << tria.begin()->center(true)
+          << " with mapping degree " << degree << std::endl;
+  for (unsigned int i = 0; i < real_points.size(); ++i)
+    deallog << "Transform " << real_points[i] << " gives " << unit_points[i]
+            << std::endl;
+  deallog << "OK" << std::endl;
+}
+
+
+
+int
+main()
+{
+  initlog();
+  deallog << std::setprecision(10);
+
+  test<2>(5);
+}
diff --git a/tests/mappings/mapping_q_points_real_to_unit.output b/tests/mappings/mapping_q_points_real_to_unit.output
new file mode 100644 (file)
index 0000000..2ff7f6f
--- /dev/null
@@ -0,0 +1,23 @@
+
+DEAL::Transform on cell with center: 0.5303300859 0.5303300859 with mapping degree 5
+DEAL::Transform -0.3181980515 -0.3226980515 gives inf 0.000000000
+DEAL::Transform 0.07071067812 0.07171067812 gives 0.5044698424 1.798577755
+DEAL::Transform -0.2474873734 -0.2509873734 gives inf 0.000000000
+DEAL::Transform 0.1414213562 0.1434213562 gives 0.5044698424 1.597155511
+DEAL::Transform -0.1767766953 -0.1792766953 gives inf 0.000000000
+DEAL::Transform 0.2121320344 0.2151320344 gives 0.5044698424 1.395733266
+DEAL::Transform -0.1060660172 -0.1075660172 gives inf 0.000000000
+DEAL::Transform 0.2828427125 0.2868427125 gives 0.5044698424 1.194311022
+DEAL::Transform -0.03535533906 -0.03585533906 gives inf 0.000000000
+DEAL::Transform 0.3535533906 0.3585533906 gives 0.5044698424 0.9928887775
+DEAL::Transform 0.03535533906 0.03585533906 gives 0.5044698424 1.899288878
+DEAL::Transform 0.4242640687 0.4302640687 gives 0.5044698424 0.7914665330
+DEAL::Transform 0.1060660172 0.1075660172 gives 0.5044698424 1.697866633
+DEAL::Transform 0.4949747468 0.5019747468 gives 0.5044698424 0.5900442885
+DEAL::Transform 0.1767766953 0.1792766953 gives 0.5044698424 1.496444389
+DEAL::Transform 0.5656854249 0.5736854249 gives 0.5044698424 0.3886220440
+DEAL::Transform 0.2474873734 0.2509873734 gives 0.5044698424 1.295022144
+DEAL::Transform 0.6363961031 0.6453961031 gives 0.5044698424 0.1871997995
+DEAL::Transform 0.3181980515 0.3226980515 gives 0.5044698424 1.093599900
+DEAL::Transform 0.7071067812 0.7171067812 gives 0.5044698424 -0.01422244502
+DEAL::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.