]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add tests.
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 5 May 2022 20:32:06 +0000 (14:32 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 5 May 2022 20:32:06 +0000 (14:32 -0600)
tests/grid/reference_cell_type_02.cc [new file with mode: 0644]
tests/grid/reference_cell_type_02.output [new file with mode: 0644]
tests/grid/reference_cell_type_03.cc [new file with mode: 0644]
tests/grid/reference_cell_type_03.output [new file with mode: 0644]

diff --git a/tests/grid/reference_cell_type_02.cc b/tests/grid/reference_cell_type_02.cc
new file mode 100644 (file)
index 0000000..ebfd0e8
--- /dev/null
@@ -0,0 +1,118 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 - 2021 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Test ReferenceCell::volume()
+
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/fe/fe_pyramid_p.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_simplex_p.h>
+#include <deal.II/fe/fe_simplex_p_bubbles.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_values.h>
+#include <deal.II/fe/fe_wedge_p.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test(const ReferenceCell &reference_cell)
+{
+  Triangulation<dim> triangulation;
+  GridGenerator::reference_cell(triangulation, reference_cell);
+
+  std::unique_ptr<Quadrature<dim>> q;
+  if (reference_cell == ReferenceCells::Line)
+    q = std::make_unique<QGauss<dim>>(2);
+  else if (reference_cell == ReferenceCells::Triangle)
+    q = std::make_unique<QGaussSimplex<dim>>(2);
+  else if (reference_cell == ReferenceCells::Quadrilateral)
+    q = std::make_unique<QGauss<dim>>(2);
+  else if (reference_cell == ReferenceCells::Tetrahedron)
+    q = std::make_unique<QGaussSimplex<dim>>(2);
+  else if (reference_cell == ReferenceCells::Wedge)
+    q = std::make_unique<QGaussWedge<dim>>(2);
+  else if (reference_cell == ReferenceCells::Pyramid)
+    q = std::make_unique<QGaussPyramid<dim>>(2);
+  else if (reference_cell == ReferenceCells::Hexahedron)
+    q = std::make_unique<QGauss<dim>>(2);
+
+  std::unique_ptr<FiniteElement<dim>> fe;
+  if (reference_cell == ReferenceCells::Line)
+    fe = std::make_unique<FE_Q<dim>>(1);
+  else if (reference_cell == ReferenceCells::Triangle)
+    fe = std::make_unique<FE_SimplexP<dim>>(1);
+  else if (reference_cell == ReferenceCells::Quadrilateral)
+    fe = std::make_unique<FE_Q<dim>>(1);
+  else if (reference_cell == ReferenceCells::Tetrahedron)
+    fe = std::make_unique<FE_SimplexP<dim>>(1);
+  else if (reference_cell == ReferenceCells::Wedge)
+    fe = std::make_unique<FE_WedgeP<dim>>(1);
+  else if (reference_cell == ReferenceCells::Pyramid)
+    fe = std::make_unique<FE_PyramidP<dim>>(1);
+  else if (reference_cell == ReferenceCells::Hexahedron)
+    fe = std::make_unique<FE_Q<dim>>(1);
+
+  // Set up the objects to compute an integral on the reference cell
+  FEValues<dim> fe_values(*fe, *q, update_JxW_values);
+  fe_values.reinit(triangulation.begin_active());
+
+  double volume = 0;
+  for (unsigned int i = 0; i < q->size(); ++i)
+    volume += fe_values.JxW(i);
+
+  deallog << "ReferenceCell: " << reference_cell.to_string() << std::endl;
+  deallog << "  computed volume = " << volume << std::endl;
+  deallog << "  self-reported volume = " << reference_cell.volume()
+          << std::endl;
+
+  Assert(std::fabs(volume - reference_cell.volume()) < 1e-12,
+         ExcInternalError());
+}
+
+int
+main()
+{
+  initlog();
+
+  {
+    deallog.push("1D");
+    test<1>(ReferenceCells::Line);
+    deallog.pop();
+  }
+
+  {
+    deallog.push("2D");
+    test<2>(ReferenceCells::Quadrilateral);
+    test<2>(ReferenceCells::Triangle);
+    deallog.pop();
+  }
+
+  {
+    deallog.push("3D");
+    test<3>(ReferenceCells::Tetrahedron);
+    test<3>(ReferenceCells::Pyramid);
+    test<3>(ReferenceCells::Wedge);
+    test<3>(ReferenceCells::Hexahedron);
+    deallog.pop();
+  }
+}
diff --git a/tests/grid/reference_cell_type_02.output b/tests/grid/reference_cell_type_02.output
new file mode 100644 (file)
index 0000000..6bda28a
--- /dev/null
@@ -0,0 +1,22 @@
+
+DEAL:1D::ReferenceCell: Line
+DEAL:1D::  computed volume = 1.00000
+DEAL:1D::  self-reported volume = 1.00000
+DEAL:2D::ReferenceCell: Quad
+DEAL:2D::  computed volume = 1.00000
+DEAL:2D::  self-reported volume = 1.00000
+DEAL:2D::ReferenceCell: Tri
+DEAL:2D::  computed volume = 0.500000
+DEAL:2D::  self-reported volume = 0.500000
+DEAL:3D::ReferenceCell: Tet
+DEAL:3D::  computed volume = 0.166667
+DEAL:3D::  self-reported volume = 0.166667
+DEAL:3D::ReferenceCell: Pyramid
+DEAL:3D::  computed volume = 1.33333
+DEAL:3D::  self-reported volume = 1.33333
+DEAL:3D::ReferenceCell: Wedge
+DEAL:3D::  computed volume = 0.500000
+DEAL:3D::  self-reported volume = 0.500000
+DEAL:3D::ReferenceCell: Hex
+DEAL:3D::  computed volume = 1.00000
+DEAL:3D::  self-reported volume = 1.00000
diff --git a/tests/grid/reference_cell_type_03.cc b/tests/grid/reference_cell_type_03.cc
new file mode 100644 (file)
index 0000000..aa821e0
--- /dev/null
@@ -0,0 +1,125 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 - 2021 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Test ReferenceCell::barycenter()
+
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/fe/fe_pyramid_p.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_simplex_p.h>
+#include <deal.II/fe/fe_simplex_p_bubbles.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_values.h>
+#include <deal.II/fe/fe_wedge_p.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test(const ReferenceCell &reference_cell)
+{
+  Triangulation<dim> triangulation;
+  GridGenerator::reference_cell(triangulation, reference_cell);
+
+  std::unique_ptr<Quadrature<dim>> q;
+  if (reference_cell == ReferenceCells::Line)
+    q = std::make_unique<QGauss<dim>>(2);
+  else if (reference_cell == ReferenceCells::Triangle)
+    q = std::make_unique<QGaussSimplex<dim>>(2);
+  else if (reference_cell == ReferenceCells::Quadrilateral)
+    q = std::make_unique<QGauss<dim>>(2);
+  else if (reference_cell == ReferenceCells::Tetrahedron)
+    q = std::make_unique<QGaussSimplex<dim>>(2);
+  else if (reference_cell == ReferenceCells::Wedge)
+    q = std::make_unique<QGaussWedge<dim>>(2);
+  else if (reference_cell == ReferenceCells::Pyramid)
+    q = std::make_unique<QGaussPyramid<dim>>(2);
+  else if (reference_cell == ReferenceCells::Hexahedron)
+    q = std::make_unique<QGauss<dim>>(2);
+
+  std::unique_ptr<FiniteElement<dim>> fe;
+  if (reference_cell == ReferenceCells::Line)
+    fe = std::make_unique<FE_Q<dim>>(1);
+  else if (reference_cell == ReferenceCells::Triangle)
+    fe = std::make_unique<FE_SimplexP<dim>>(1);
+  else if (reference_cell == ReferenceCells::Quadrilateral)
+    fe = std::make_unique<FE_Q<dim>>(1);
+  else if (reference_cell == ReferenceCells::Tetrahedron)
+    fe = std::make_unique<FE_SimplexP<dim>>(1);
+  else if (reference_cell == ReferenceCells::Wedge)
+    fe = std::make_unique<FE_WedgeP<dim>>(1);
+  else if (reference_cell == ReferenceCells::Pyramid)
+    fe = std::make_unique<FE_PyramidP<dim>>(1);
+  else if (reference_cell == ReferenceCells::Hexahedron)
+    fe = std::make_unique<FE_Q<dim>>(1);
+
+  // Set up the objects to compute an integral on the reference cell
+  FEValues<dim> fe_values(*fe,
+                          *q,
+                          update_JxW_values | update_quadrature_points);
+  fe_values.reinit(triangulation.begin_active());
+
+  double     volume = 0;
+  Point<dim> barycenter;
+  for (unsigned int i = 0; i < q->size(); ++i)
+    {
+      volume += fe_values.JxW(i);
+      barycenter += fe_values.quadrature_point(i) * fe_values.JxW(i);
+    }
+  barycenter /= volume;
+
+  deallog << "ReferenceCell: " << reference_cell.to_string() << std::endl;
+  deallog << "  computed barycenter = " << barycenter << std::endl;
+  deallog << "  self-reported barycenter = " << reference_cell.barycenter<dim>()
+          << std::endl;
+
+  Assert((barycenter - reference_cell.barycenter<dim>()).norm() <= 1e-12,
+         ExcInternalError());
+}
+
+int
+main()
+{
+  initlog();
+
+  {
+    deallog.push("1D");
+    test<1>(ReferenceCells::Line);
+    deallog.pop();
+  }
+
+  {
+    deallog.push("2D");
+    test<2>(ReferenceCells::Quadrilateral);
+    test<2>(ReferenceCells::Triangle);
+    deallog.pop();
+  }
+
+  {
+    deallog.push("3D");
+    test<3>(ReferenceCells::Tetrahedron);
+    test<3>(ReferenceCells::Pyramid);
+    test<3>(ReferenceCells::Wedge);
+    test<3>(ReferenceCells::Hexahedron);
+    deallog.pop();
+  }
+}
diff --git a/tests/grid/reference_cell_type_03.output b/tests/grid/reference_cell_type_03.output
new file mode 100644 (file)
index 0000000..5929d52
--- /dev/null
@@ -0,0 +1,22 @@
+
+DEAL:1D::ReferenceCell: Line
+DEAL:1D::  computed barycenter = 0.500000
+DEAL:1D::  self-reported barycenter = 0.500000
+DEAL:2D::ReferenceCell: Quad
+DEAL:2D::  computed barycenter = 0.500000 0.500000
+DEAL:2D::  self-reported barycenter = 0.500000 0.500000
+DEAL:2D::ReferenceCell: Tri
+DEAL:2D::  computed barycenter = 0.333333 0.333333
+DEAL:2D::  self-reported barycenter = 0.333333 0.333333
+DEAL:3D::ReferenceCell: Tet
+DEAL:3D::  computed barycenter = 0.250000 0.250000 0.250000
+DEAL:3D::  self-reported barycenter = 0.250000 0.250000 0.250000
+DEAL:3D::ReferenceCell: Pyramid
+DEAL:3D::  computed barycenter = 2.08167e-17 6.24500e-17 0.250000
+DEAL:3D::  self-reported barycenter = 0.00000 0.00000 0.250000
+DEAL:3D::ReferenceCell: Wedge
+DEAL:3D::  computed barycenter = 0.333333 0.333333 0.500000
+DEAL:3D::  self-reported barycenter = 0.333333 0.333333 0.500000
+DEAL:3D::ReferenceCell: Hex
+DEAL:3D::  computed barycenter = 0.500000 0.500000 0.500000
+DEAL:3D::  self-reported barycenter = 0.500000 0.500000 0.500000

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.