]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Test inverse Jacobians.
authorkronbichler <kronbichler@0785d39b-7218-0410-832d-ea1e28bc413d>
Sun, 11 Mar 2012 13:11:13 +0000 (13:11 +0000)
committerkronbichler <kronbichler@0785d39b-7218-0410-832d-ea1e28bc413d>
Sun, 11 Mar 2012 13:11:13 +0000 (13:11 +0000)
git-svn-id: https://svn.dealii.org/trunk@25249 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
tests/fe/jacobians.cc [new file with mode: 0644]
tests/fe/jacobians/cmp/generic [new file with mode: 0644]

index c16b1b6e5f9745ab1dada629ba84c5737aac72e4..090b99eea9b3946181b514abf2f34f22652744d7 100644 (file)
@@ -152,6 +152,12 @@ enabled due to a missing include file in file
 <h3>Specific improvements</h3>
 
 <ol>
+<li> Fixed:
+The method FEValues<dim>::inverse_jacobian() previously returned the transpose of the
+inverse Jacobians instead of just the inverse Jacobian as documented. This is now fixed.
+<br>
+(Sebastian Pauletti, Katharina Kormann, Martin Kronbichler, 2012/03/11)
+
 <li> Extended:
 SolutionTransfer is now also able to transfer solutions between hp::DoFHandler where
 the finite element index changes during refinement.
diff --git a/tests/fe/jacobians.cc b/tests/fe/jacobians.cc
new file mode 100644 (file)
index 0000000..83eb479
--- /dev/null
@@ -0,0 +1,91 @@
+//----------------------------  jacobians.cc  ---------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2012 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  jacobians.cc  ---------------------------
+
+// Show the Jacobians and inverse Jacobians on hyperball with one quadrature
+// point
+
+#include "../tests.h"
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/fe/fe_nothing.h>
+#include <deal.II/fe/fe_values.h>
+#include <deal.II/fe/mapping_q.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+
+template<int dim>
+void test()
+{
+  Triangulation<dim> tria;
+  GridGenerator::hyper_ball (tria);
+  static const HyperBallBoundary<dim> boundary;
+  tria.set_boundary (0, boundary);
+
+  MappingQ<dim> mapping(3);
+  FE_Nothing<dim> dummy;
+                               // choose a point that is not right in the
+                               // middle of the cell so that the Jacobian
+                               // contains many nonzero entries
+  Point<dim> quad_p;
+  for (int d=0; d<dim; ++d)
+    quad_p(d) = 0.42 + 0.11 * d;
+  Quadrature<dim> quad(quad_p);
+  FEValues<dim> fe_val (mapping, dummy, quad,
+                       update_jacobians | update_inverse_jacobians);
+  deallog << dim << "d Jacobians:" << std::endl;
+  typename Triangulation<dim>::active_cell_iterator
+    cell = tria.begin_active(), endc = tria.end();
+  for ( ; cell != endc; ++cell)
+    {
+      fe_val.reinit (cell);
+      
+      for (unsigned int d=0; d<dim; ++d)
+       for (unsigned int e=0; e<dim; ++e)
+         deallog << fe_val.jacobian(0)[d][e] << " ";
+      deallog << std::endl;
+    }
+  deallog << std::endl;
+
+  deallog << dim << "d inverse Jacobians:" << std::endl;
+  cell = tria.begin_active();
+  endc = tria.end();
+  for ( ; cell != endc; ++cell)
+    {
+      fe_val.reinit (cell);
+      
+      for (unsigned int d=0; d<dim; ++d)
+       for (unsigned int e=0; e<dim; ++e)
+         deallog << fe_val.inverse_jacobian(0)[d][e] << " ";
+      deallog << std::endl;
+    }
+  deallog << std::endl;
+}
+
+
+int
+main()
+{
+  std::ofstream logfile ("jacobians/output");
+  deallog << std::setprecision(4) << std::fixed;
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  test<2>();
+  test<3>();
+
+  return 0;
+}
+
+
+
diff --git a/tests/fe/jacobians/cmp/generic b/tests/fe/jacobians/cmp/generic
new file mode 100644 (file)
index 0000000..56cc9df
--- /dev/null
@@ -0,0 +1,33 @@
+
+DEAL::2d Jacobians:
+DEAL::1.0000 0 0 0.4142 
+DEAL::0.4142 0 0 1.0000 
+DEAL::0.5858 0 0 0.5858 
+DEAL::0 -0.4142 1.0000 0 
+DEAL::0 1.0000 -0.4142 0 
+DEAL::
+DEAL::2d inverse Jacobians:
+DEAL::1.0000 0 0 2.4142 
+DEAL::2.4142 0 0 1.0000 
+DEAL::1.7071 0 0 1.7071 
+DEAL::0 1.0000 -2.4142 0 
+DEAL::0 -2.4142 1.0000 0 
+DEAL::
+DEAL::3d Jacobians:
+DEAL::0.4226 0 0 0 0.4226 0 0 0 0.4226 
+DEAL::0.7887 0 0 0 0.7887 0 0 0 0.3660 
+DEAL::0 -0.3660 0 0.7887 0 0 0 0 0.7887 
+DEAL::0.7887 0 0 0 0 0.7887 0 -0.3660 0 
+DEAL::0.3660 0 0 0 0.7887 0 0 0 0.7887 
+DEAL::0.7887 0 0 0 0.3660 0 0 0 0.7887 
+DEAL::0 0.7887 0 -0.3660 0 0 0 0 0.7887 
+DEAL::
+DEAL::3d inverse Jacobians:
+DEAL::2.3660 0 0 0 2.3660 0 0 0 2.3660 
+DEAL::1.2679 0 0 0 1.2679 0 0 0 2.7321 
+DEAL::0 1.2679 0 -2.7321 0 0 0 0 1.2679 
+DEAL::1.2679 0 0 0 0 -2.7321 0 1.2679 0 
+DEAL::2.7321 0 0 0 1.2679 0 0 0 1.2679 
+DEAL::1.2679 0 0 0 2.7321 0 0 0 1.2679 
+DEAL::0 -2.7321 0 1.2679 0 0 0 0 1.2679 
+DEAL::

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.