]> https://gitweb.dealii.org/ - dealii.git/commitdiff
New test case
authorMartin Kronbichler <martin.kronbichler@uni-a.de>
Thu, 7 Sep 2023 18:46:09 +0000 (20:46 +0200)
committerMartin Kronbichler <martin.kronbichler@uni-a.de>
Thu, 7 Sep 2023 18:49:58 +0000 (20:49 +0200)
tests/fe/rtn_4.cc [new file with mode: 0644]
tests/fe/rtn_4.output [new file with mode: 0644]

diff --git a/tests/fe/rtn_4.cc b/tests/fe/rtn_4.cc
new file mode 100644 (file)
index 0000000..422f397
--- /dev/null
@@ -0,0 +1,62 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2023 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 that FEValues can be reinit'ed with a FE_RaviartThomasNodal on the
+// mesh of a 2d ball when no finite element entries are requested
+
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_raviart_thomas.h>
+#include <deal.II/fe/fe_values.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <string>
+#include <vector>
+
+#include "../tests.h"
+
+int
+main()
+{
+  initlog();
+
+  const int          dim = 2;
+  Triangulation<dim> tria;
+  GridGenerator::hyper_ball(tria, Point<dim>(), 1.);
+
+  FE_RaviartThomasNodal<dim> fe(0);
+  DoFHandler<dim>            dof(tria);
+  dof.distribute_dofs(fe);
+
+  FEValues<dim> fe_values(fe,
+                          QGauss<dim>(1),
+                          update_inverse_jacobians | update_quadrature_points);
+
+  for (const auto &cell : dof.active_cell_iterators())
+    {
+      fe_values.reinit(cell);
+      deallog << "Data on cell " << cell->id()
+              << ": p=" << fe_values.quadrature_point(0)
+              << " jac=" << Tensor<2, dim>(fe_values.inverse_jacobian(0))
+              << std::endl;
+    }
+  deallog << std::endl;
+
+  return 0;
+}
diff --git a/tests/fe/rtn_4.output b/tests/fe/rtn_4.output
new file mode 100644 (file)
index 0000000..5691eac
--- /dev/null
@@ -0,0 +1,7 @@
+
+DEAL::Data on cell 0_0:: p=0.00000 -0.500000 jac=1.00000 0.00000 0.00000 2.41421
+DEAL::Data on cell 1_0:: p=-0.500000 0.00000 jac=2.41421 0.00000 0.00000 1.00000
+DEAL::Data on cell 2_0:: p=0.00000 0.00000 jac=1.70711 0.00000 0.00000 1.70711
+DEAL::Data on cell 3_0:: p=0.500000 0.00000 jac=0.00000 1.00000 -2.41421 0.00000
+DEAL::Data on cell 4_0:: p=0.00000 0.500000 jac=0.00000 -2.41421 1.00000 0.00000
+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.