]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add new hp::FEFaceValues::reinit() functions 10705/head
authorPeter Munch <peterrmuench@gmail.com>
Wed, 15 Jul 2020 10:47:45 +0000 (12:47 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Thu, 16 Jul 2020 06:46:21 +0000 (08:46 +0200)
doc/news/changes/minor/20200716Munch [new file with mode: 0644]
include/deal.II/hp/fe_values.h
source/hp/fe_values.cc
source/hp/fe_values.inst.in
tests/hp/fe_face_values_reinit_01.cc [new file with mode: 0644]
tests/hp/fe_face_values_reinit_01.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20200716Munch b/doc/news/changes/minor/20200716Munch
new file mode 100644 (file)
index 0000000..4e5125b
--- /dev/null
@@ -0,0 +1,3 @@
+New: The method hp::FEFaceValues::reinit() can now also accept face iterators.
+<br>
+(Peter Munch, 2020/07/16)
index 2debb37cb5c4b84740cc1bcf3b430ba0dcd0ed85..9b2df92b936fc09086c253966ef87202f4c4c7bd 100644 (file)
@@ -23,6 +23,9 @@
 #include <deal.II/fe/fe.h>
 #include <deal.II/fe/fe_values.h>
 
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_iterator.h>
+
 #include <deal.II/hp/fe_collection.h>
 #include <deal.II/hp/mapping_collection.h>
 #include <deal.II/hp/q_collection.h>
@@ -476,6 +479,19 @@ namespace hp
            const unsigned int mapping_index = numbers::invalid_unsigned_int,
            const unsigned int fe_index      = numbers::invalid_unsigned_int);
 
+    /**
+     * Reinitialize the object for the given cell and face.
+     *
+     * @note @p face must be one of @p cell's face iterators.
+     */
+    template <bool lda>
+    void
+    reinit(const TriaIterator<DoFCellAccessor<dim, spacedim, lda>> &   cell,
+           const typename Triangulation<dim, spacedim>::face_iterator &face,
+           const unsigned int q_index       = numbers::invalid_unsigned_int,
+           const unsigned int mapping_index = numbers::invalid_unsigned_int,
+           const unsigned int fe_index      = numbers::invalid_unsigned_int);
+
     /**
      * Like the previous function, but for non-DoFHandler iterators. The reason
      * this function exists is so that one can use this class for
@@ -494,6 +510,18 @@ namespace hp
            const unsigned int q_index       = numbers::invalid_unsigned_int,
            const unsigned int mapping_index = numbers::invalid_unsigned_int,
            const unsigned int fe_index      = numbers::invalid_unsigned_int);
+
+    /**
+     * Reinitialize the object for the given cell and face.
+     *
+     * @note @p face must be one of @p cell's face iterators.
+     */
+    void
+    reinit(const typename Triangulation<dim, spacedim>::cell_iterator &cell,
+           const typename Triangulation<dim, spacedim>::face_iterator &face,
+           const unsigned int q_index       = numbers::invalid_unsigned_int,
+           const unsigned int mapping_index = numbers::invalid_unsigned_int,
+           const unsigned int fe_index      = numbers::invalid_unsigned_int);
   };
 
 
index aa7c60f98393756672b956308792ffb22049ba39..6c1db5ab089d758184f205af4ac09f275a421267 100644 (file)
@@ -381,6 +381,22 @@ namespace hp
 
 
 
+  template <int dim, int spacedim>
+  template <bool lda>
+  void
+  FEFaceValues<dim, spacedim>::reinit(
+    const TriaIterator<DoFCellAccessor<dim, spacedim, lda>> &   cell,
+    const typename Triangulation<dim, spacedim>::face_iterator &face,
+    const unsigned int                                          q_index,
+    const unsigned int                                          mapping_index,
+    const unsigned int                                          fe_index)
+  {
+    const auto face_n = cell->face_iterator_to_index(face);
+    reinit(cell, face_n, q_index, mapping_index, fe_index);
+  }
+
+
+
   template <int dim, int spacedim>
   void
   FEFaceValues<dim, spacedim>::reinit(
@@ -417,6 +433,21 @@ namespace hp
   }
 
 
+
+  template <int dim, int spacedim>
+  void
+  FEFaceValues<dim, spacedim>::reinit(
+    const typename Triangulation<dim, spacedim>::cell_iterator &cell,
+    const typename Triangulation<dim, spacedim>::face_iterator &face,
+    const unsigned int                                          q_index,
+    const unsigned int                                          mapping_index,
+    const unsigned int                                          fe_index)
+  {
+    const auto face_n = cell->face_iterator_to_index(face);
+    reinit(cell, face_n, q_index, mapping_index, fe_index);
+  }
+
+
   // -------------------------- FESubfaceValues -------------------------
 
 
index 592ab00cbd8db9524627e0261d3fc9a4e850372f..01acaeb46503f5bdce6144caff688d86c75a932c 100644 (file)
@@ -105,6 +105,15 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
         unsigned int,
         unsigned int);
       template void
+      FEFaceValues<deal_II_dimension, deal_II_space_dimension>::reinit(
+        const TriaIterator<
+          DoFCellAccessor<deal_II_dimension, deal_II_space_dimension, lda>> &,
+        const Triangulation<deal_II_dimension,
+                            deal_II_space_dimension>::face_iterator &,
+        unsigned int,
+        unsigned int,
+        unsigned int);
+      template void
       FESubfaceValues<deal_II_dimension, deal_II_space_dimension>::reinit(
         TriaIterator<
           DoFCellAccessor<deal_II_dimension, deal_II_space_dimension, lda>>,
diff --git a/tests/hp/fe_face_values_reinit_01.cc b/tests/hp/fe_face_values_reinit_01.cc
new file mode 100644 (file)
index 0000000..f725949
--- /dev/null
@@ -0,0 +1,131 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Test different versions of hp::FEFaceValues::reinit().
+
+
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_nothing.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_refinement.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+
+#include <deal.II/hp/dof_handler.h>
+#include <deal.II/hp/fe_collection.h>
+#include <deal.II/hp/fe_values.h>
+
+#include "../tests.h"
+
+
+
+template <int dim, int spacedim = dim>
+void
+test()
+{
+  Triangulation<dim> triangulation;
+  GridGenerator::hyper_cube(triangulation);
+
+  hp::FECollection<dim>      fe_collection(FE_Q<dim>(1));
+  hp::QCollection<dim - 1>   quad_collection(QGauss<dim - 1>(2));
+  hp::MappingCollection<dim> mapping_collection(MappingQ<dim, spacedim>(1));
+
+  hp::DoFHandler<dim> dof_handler(triangulation);
+
+  dof_handler.distribute_dofs(fe_collection);
+
+
+
+  hp::FEFaceValues<dim, spacedim> hp_fe_face_values(mapping_collection,
+                                                    fe_collection,
+                                                    quad_collection,
+                                                    update_quadrature_points);
+
+  for (const auto &cell : triangulation.active_cell_iterators())
+    {
+      for (const auto face : cell->face_indices())
+        {
+          hp_fe_face_values.reinit(cell, face);
+
+          for (const auto &p : hp_fe_face_values.get_present_fe_values()
+                                 .get_quadrature_points())
+            deallog << p << " ";
+          deallog << std::endl;
+        }
+      deallog << std::endl;
+    }
+
+  for (const auto &cell : triangulation.active_cell_iterators())
+    {
+      for (const auto face : cell->face_iterators())
+        {
+          hp_fe_face_values.reinit(cell, face);
+
+          for (const auto &p : hp_fe_face_values.get_present_fe_values()
+                                 .get_quadrature_points())
+            deallog << p << " ";
+          deallog << std::endl;
+        }
+      deallog << std::endl;
+    }
+
+  for (const auto &cell : dof_handler.active_cell_iterators())
+    {
+      for (const auto face : cell->face_indices())
+        {
+          hp_fe_face_values.reinit(cell, face);
+
+          for (const auto &p : hp_fe_face_values.get_present_fe_values()
+                                 .get_quadrature_points())
+            deallog << p << " ";
+          deallog << std::endl;
+        }
+      deallog << std::endl;
+    }
+
+  for (const auto &cell : dof_handler.active_cell_iterators())
+    {
+      for (const auto face : cell->face_iterators())
+        {
+          hp_fe_face_values.reinit(cell, face);
+
+          for (const auto &p : hp_fe_face_values.get_present_fe_values()
+                                 .get_quadrature_points())
+            deallog << p << " ";
+          deallog << std::endl;
+        }
+      deallog << std::endl;
+    }
+}
+
+
+
+int
+main()
+{
+  initlog();
+  deallog.get_file_stream().precision(2);
+
+  test<2>();
+  test<3>();
+}
diff --git a/tests/hp/fe_face_values_reinit_01.output b/tests/hp/fe_face_values_reinit_01.output
new file mode 100644 (file)
index 0000000..d7ee8d4
--- /dev/null
@@ -0,0 +1,49 @@
+
+DEAL::0.00000 0.211325 0.00000 0.788675 
+DEAL::1.00000 0.211325 1.00000 0.788675 
+DEAL::0.211325 0.00000 0.788675 0.00000 
+DEAL::0.211325 1.00000 0.788675 1.00000 
+DEAL::
+DEAL::0.00000 0.211325 0.00000 0.788675 
+DEAL::1.00000 0.211325 1.00000 0.788675 
+DEAL::0.211325 0.00000 0.788675 0.00000 
+DEAL::0.211325 1.00000 0.788675 1.00000 
+DEAL::
+DEAL::0.00000 0.211325 0.00000 0.788675 
+DEAL::1.00000 0.211325 1.00000 0.788675 
+DEAL::0.211325 0.00000 0.788675 0.00000 
+DEAL::0.211325 1.00000 0.788675 1.00000 
+DEAL::
+DEAL::0.00000 0.211325 0.00000 0.788675 
+DEAL::1.00000 0.211325 1.00000 0.788675 
+DEAL::0.211325 0.00000 0.788675 0.00000 
+DEAL::0.211325 1.00000 0.788675 1.00000 
+DEAL::
+DEAL::0.00000 0.211325 0.211325 0.00000 0.788675 0.211325 0.00000 0.211325 0.788675 0.00000 0.788675 0.788675 
+DEAL::1.00000 0.211325 0.211325 1.00000 0.788675 0.211325 1.00000 0.211325 0.788675 1.00000 0.788675 0.788675 
+DEAL::0.211325 0.00000 0.211325 0.211325 0.00000 0.788675 0.788675 0.00000 0.211325 0.788675 0.00000 0.788675 
+DEAL::0.211325 1.00000 0.211325 0.211325 1.00000 0.788675 0.788675 1.00000 0.211325 0.788675 1.00000 0.788675 
+DEAL::0.211325 0.211325 0.00000 0.788675 0.211325 0.00000 0.211325 0.788675 0.00000 0.788675 0.788675 0.00000 
+DEAL::0.211325 0.211325 1.00000 0.788675 0.211325 1.00000 0.211325 0.788675 1.00000 0.788675 0.788675 1.00000 
+DEAL::
+DEAL::0.00000 0.211325 0.211325 0.00000 0.788675 0.211325 0.00000 0.211325 0.788675 0.00000 0.788675 0.788675 
+DEAL::1.00000 0.211325 0.211325 1.00000 0.788675 0.211325 1.00000 0.211325 0.788675 1.00000 0.788675 0.788675 
+DEAL::0.211325 0.00000 0.211325 0.211325 0.00000 0.788675 0.788675 0.00000 0.211325 0.788675 0.00000 0.788675 
+DEAL::0.211325 1.00000 0.211325 0.211325 1.00000 0.788675 0.788675 1.00000 0.211325 0.788675 1.00000 0.788675 
+DEAL::0.211325 0.211325 0.00000 0.788675 0.211325 0.00000 0.211325 0.788675 0.00000 0.788675 0.788675 0.00000 
+DEAL::0.211325 0.211325 1.00000 0.788675 0.211325 1.00000 0.211325 0.788675 1.00000 0.788675 0.788675 1.00000 
+DEAL::
+DEAL::0.00000 0.211325 0.211325 0.00000 0.788675 0.211325 0.00000 0.211325 0.788675 0.00000 0.788675 0.788675 
+DEAL::1.00000 0.211325 0.211325 1.00000 0.788675 0.211325 1.00000 0.211325 0.788675 1.00000 0.788675 0.788675 
+DEAL::0.211325 0.00000 0.211325 0.211325 0.00000 0.788675 0.788675 0.00000 0.211325 0.788675 0.00000 0.788675 
+DEAL::0.211325 1.00000 0.211325 0.211325 1.00000 0.788675 0.788675 1.00000 0.211325 0.788675 1.00000 0.788675 
+DEAL::0.211325 0.211325 0.00000 0.788675 0.211325 0.00000 0.211325 0.788675 0.00000 0.788675 0.788675 0.00000 
+DEAL::0.211325 0.211325 1.00000 0.788675 0.211325 1.00000 0.211325 0.788675 1.00000 0.788675 0.788675 1.00000 
+DEAL::
+DEAL::0.00000 0.211325 0.211325 0.00000 0.788675 0.211325 0.00000 0.211325 0.788675 0.00000 0.788675 0.788675 
+DEAL::1.00000 0.211325 0.211325 1.00000 0.788675 0.211325 1.00000 0.211325 0.788675 1.00000 0.788675 0.788675 
+DEAL::0.211325 0.00000 0.211325 0.211325 0.00000 0.788675 0.788675 0.00000 0.211325 0.788675 0.00000 0.788675 
+DEAL::0.211325 1.00000 0.211325 0.211325 1.00000 0.788675 0.788675 1.00000 0.211325 0.788675 1.00000 0.788675 
+DEAL::0.211325 0.211325 0.00000 0.788675 0.211325 0.00000 0.211325 0.788675 0.00000 0.788675 0.788675 0.00000 
+DEAL::0.211325 0.211325 1.00000 0.788675 0.211325 1.00000 0.211325 0.788675 1.00000 0.788675 0.788675 1.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.