]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Test case for MappingManifold 9769/head
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Sun, 29 Mar 2020 18:08:12 +0000 (20:08 +0200)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Sun, 29 Mar 2020 18:13:00 +0000 (20:13 +0200)
tests/fe/cell_similarity_13.cc [new file with mode: 0644]
tests/fe/cell_similarity_13.output [new file with mode: 0644]

diff --git a/tests/fe/cell_similarity_13.cc b/tests/fe/cell_similarity_13.cc
new file mode 100644 (file)
index 0000000..6e47fb7
--- /dev/null
@@ -0,0 +1,107 @@
+// ---------------------------------------------------------------------
+//
+// 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 the CellSimilarity with MappingManifold on a somewhat arbitrarily
+// deformed shape - a cube with one refinement in flat coordinates where we
+// afterwards apply a spherical manifold to deform things.
+
+#include <deal.II/base/multithread_info.h>
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_values.h>
+#include <deal.II/fe/mapping_manifold.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/manifold_lib.h>
+
+#include "../tests.h"
+
+
+
+template <int dim>
+void
+test()
+{
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria, 0.2, 1.);
+  tria.refine_global(1);
+  SphericalManifold<dim> manifold;
+  tria.set_all_manifold_ids(0);
+  tria.set_manifold(0, manifold);
+
+  FE_Q<dim> fe(1);
+  deallog << "FE=" << fe.get_name() << std::endl;
+
+  MappingManifold<dim> mapping;
+  deallog << "Mapping=MappingManifold" << std::endl;
+
+  DoFHandler<dim> dof(tria);
+  dof.distribute_dofs(fe);
+
+  std::vector<Point<dim>> points(2);
+  for (unsigned int d = 0; d < dim; ++d)
+    points[0][d] = 0.1;
+  for (unsigned int d = 0; d < dim; ++d)
+    points[1][d] = 0.85;
+
+  const Quadrature<dim> quadrature(points);
+  FEValues<dim>         fe_values(mapping,
+                          fe,
+                          quadrature,
+                          update_gradients | update_jacobians);
+
+  for (const auto &cell : dof.active_cell_iterators())
+    {
+      fe_values.reinit(cell);
+
+      deallog << "Jacobians: ";
+      for (unsigned int q = 0; q < fe_values.n_quadrature_points; ++q)
+        {
+          deallog << "[ ";
+          for (unsigned int d = 0; d < dim; ++d)
+            for (unsigned int e = 0; e < dim; ++e)
+              deallog << fe_values.jacobian(q)[d][e] << " ";
+          deallog << " ] ";
+        }
+      deallog << std::endl;
+      deallog << "Derivatives of shape function: ";
+      for (unsigned int q = 0; q < fe_values.n_quadrature_points; ++q)
+        {
+          deallog << "[ ";
+          for (unsigned int d = 0; d < dim; ++d)
+            deallog << fe_values.shape_grad(fe.dofs_per_cell / 2, q)[d] << " ";
+          deallog << " ] ";
+        }
+      deallog << std::endl;
+    }
+}
+
+
+
+int
+main()
+{
+  initlog();
+  deallog << std::setprecision(8);
+  MultithreadInfo::set_thread_limit(1);
+
+  test<2>();
+  test<3>();
+}
diff --git a/tests/fe/cell_similarity_13.output b/tests/fe/cell_similarity_13.output
new file mode 100644 (file)
index 0000000..7331fb7
--- /dev/null
@@ -0,0 +1,29 @@
+
+DEAL::FE=FE_Q<2>(1)
+DEAL::Mapping=MappingManifold
+DEAL::Jacobians: [ 0.35298659 0.12255517 0.12255517 0.35298659  ] [ 0.42290226 -0.089001296 -0.089001296 0.42290226  ] 
+DEAL::Derivatives of shape function: [ -1.3286971 3.0109889  ] [ -2.0249615 -0.071468512  ] 
+DEAL::Jacobians: [ 0.38950228 0.087621384 0.043483822 0.36569083  ] [ 0.40834664 -0.061635082 -0.051345108 0.42090842  ] 
+DEAL::Derivatives of shape function: [ -0.54610111 2.5919440  ] [ -2.0749599 0.052528468  ] 
+DEAL::Jacobians: [ 0.36569083 0.043483822 0.087621384 0.38950228  ] [ 0.42090842 -0.051345108 -0.061635082 0.40834664  ] 
+DEAL::Derivatives of shape function: [ -0.84982882 2.4055156  ] [ -2.0025228 0.11553971  ] 
+DEAL::Jacobians: [ 0.37756894 0.061824859 0.061824859 0.37756894  ] [ 0.41458323 -0.049052802 -0.049052802 0.41458323  ] 
+DEAL::Derivatives of shape function: [ -0.67321583 2.4939061  ] [ -2.0359448 0.12091952  ] 
+DEAL::FE=FE_Q<3>(1)
+DEAL::Mapping=MappingManifold
+DEAL::Jacobians: [ 0.35158966 0.080459803 0.080459803 0.080459803 0.35158966 0.080459803 0.080459803 0.080459803 0.35158966  ] [ 0.42674598 -0.056388829 -0.056388829 -0.056388829 0.42674598 -0.056388829 -0.056388829 -0.056388829 0.42674598  ] 
+DEAL::Derivatives of shape function: [ -0.69673218 -0.69673218 2.6227104  ] [ -0.35033098 -0.35033098 -0.039858624  ] 
+DEAL::Jacobians: [ 0.38310458 0.075411822 0.075411822 0.036485993 0.37114052 0.028998614 0.036485993 0.028998614 0.37114052  ] [ 0.41316565 -0.048485032 -0.048485032 -0.040048398 0.41969227 -0.028552442 -0.040048398 -0.028552442 0.41969227  ] 
+DEAL::Derivatives of shape function: [ -0.42141902 -0.33613437 2.2943531  ] [ -0.34285258 -0.34404184 -0.0094032054  ] 
+DEAL::Jacobians: [ 0.37114052 0.036485993 0.028998614 0.075411822 0.38310458 0.075411822 0.028998614 0.036485993 0.37114052  ] [ 0.41969227 -0.040048398 -0.028552442 -0.048485032 0.41316565 -0.048485032 -0.028552442 -0.040048398 0.41969227  ] 
+DEAL::Derivatives of shape function: [ -0.33613437 -0.42141902 2.2943531  ] [ -0.34404184 -0.34285258 -0.0094032054  ] 
+DEAL::Jacobians: [ 0.37739297 0.055035159 0.045988840 0.055035159 0.37739297 0.045988840 0.021139738 0.021139738 0.38042405  ] [ 0.41524234 -0.041749005 -0.031973219 -0.041749005 0.41524234 -0.031973219 -0.024644034 -0.024644034 0.41488987  ] 
+DEAL::Derivatives of shape function: [ -0.31594994 -0.31594994 2.2055922  ] [ -0.34126383 -0.34126383 0.0016327071  ] 
+DEAL::Jacobians: [ 0.37114052 0.028998614 0.036485993 0.028998614 0.37114052 0.036485993 0.075411822 0.075411822 0.38310458  ] [ 0.41969227 -0.028552442 -0.040048398 -0.028552442 0.41969227 -0.040048398 -0.048485032 -0.048485032 0.41316565  ] 
+DEAL::Derivatives of shape function: [ -0.64660380 -0.64660380 2.2374673  ] [ -0.32707987 -0.32707987 -0.0089505256  ] 
+DEAL::Jacobians: [ 0.37739297 0.045988840 0.055035159 0.021139738 0.38042405 0.021139738 0.055035159 0.045988840 0.37739297  ] [ 0.41524234 -0.031973219 -0.041749005 -0.024644034 0.41488987 -0.024644034 -0.041749005 -0.031973219 0.41524234  ] 
+DEAL::Derivatives of shape function: [ -0.54179475 -0.44309641 2.2501338  ] [ -0.32661131 -0.33235552 0.0016224788  ] 
+DEAL::Jacobians: [ 0.38042405 0.021139738 0.021139738 0.045988840 0.37739297 0.055035159 0.045988840 0.055035159 0.37739297  ] [ 0.41488987 -0.024644034 -0.024644034 -0.031973219 0.41524234 -0.041749005 -0.031973219 -0.041749005 0.41524234  ] 
+DEAL::Derivatives of shape function: [ -0.44309641 -0.54179475 2.2501338  ] [ -0.33235552 -0.32661131 0.0016224788  ] 
+DEAL::Jacobians: [ 0.37890545 0.038260015 0.038260015 0.038260015 0.37890545 0.038260015 0.038260015 0.038260015 0.37890545  ] [ 0.41423674 -0.029655906 -0.029655906 -0.029655906 0.41423674 -0.029655906 -0.029655906 -0.029655906 0.41423674  ] 
+DEAL::Derivatives of shape function: [ -0.41957384 -0.41957384 2.2224695  ] [ -0.33099588 -0.33099588 0.0069236608  ] 

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.