]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add dof_indices() to FEInterfaceValues 12628/head
authorJean-Paul Pelteret <jppelteret@gmail.com>
Thu, 5 Aug 2021 17:57:34 +0000 (19:57 +0200)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Fri, 6 Aug 2021 05:34:38 +0000 (07:34 +0200)
doc/news/changes/minor/20210805Jean-PaulPelteret [new file with mode: 0644]
examples/step-74/step-74.cc
include/deal.II/fe/fe_interface_values.h
tests/simplex/step-74.cc

diff --git a/doc/news/changes/minor/20210805Jean-PaulPelteret b/doc/news/changes/minor/20210805Jean-PaulPelteret
new file mode 100644 (file)
index 0000000..d98b095
--- /dev/null
@@ -0,0 +1,5 @@
+New: The function FEInterfaceValues::dof_indices() has been added to make
+writing range-based for loops over all local interface degrees of freedom much
+simpler.
+<br>
+(Jean-Paul Pelteret, 2021/08/05)
index 8a0ee148afabd4cbc9159fe8466cc3f2383aab81..01623b5128b624e9895b5d8770d96ca5047cc941 100644 (file)
@@ -488,9 +488,6 @@ namespace Step74
       const FEInterfaceValues<dim> &fe_iv =
         scratch_data.reinit(cell, f, sf, ncell, nf, nsf);
 
-      const auto &       q_points   = fe_iv.get_quadrature_points();
-      const unsigned int n_q_points = q_points.size();
-
       copy_data.face_data.emplace_back();
       CopyDataFace &     copy_data_face = copy_data.face_data.back();
       const unsigned int n_dofs_face    = fe_iv.n_current_interface_dofs();
@@ -504,10 +501,10 @@ namespace Step74
       const double extent2 = ncell->measure() / ncell->face(nf)->measure();
       const double penalty = get_penalty_factor(degree, extent1, extent2);
 
-      for (unsigned int point = 0; point < n_q_points; ++point)
+      for (const unsigned int point : fe_iv.quadrature_point_indices())
         {
-          for (unsigned int i = 0; i < n_dofs_face; ++i)
-            for (unsigned int j = 0; j < n_dofs_face; ++j)
+          for (const unsigned int i : fe_iv.dof_indices())
+            for (const unsigned int j : fe_iv.dof_indices())
               copy_data_face.cell_matrix(i, j) +=
                 (-diffusion_coefficient *                     // - nu
                    fe_iv.jump_in_shape_values(i, point) *     // [v_h]
index f285cce413553c5a7d35164cf615bb9565316c08..93a74a1db1f11497ddb8a562cbf2df6b11c59317 100644 (file)
@@ -18,6 +18,8 @@
 
 #include <deal.II/base/config.h>
 
+#include <deal.II/base/std_cxx20/iota_view.h>
+
 #include <deal.II/fe/fe_values.h>
 #include <deal.II/fe/mapping_q1.h>
 
@@ -721,6 +723,36 @@ public:
   unsigned
   n_current_interface_dofs() const;
 
+  /**
+   * Return an object that can be thought of as an array containing all
+   * indices from zero (inclusive) to `n_current_interface_dofs()` (exclusive).
+   * This allows one to write code using range-based `for` loops of the
+   * following kind:
+   * @code
+   *   FEInterfaceValues<dim> fe_iv (...);
+   *   FullMatrix<double>     cell_matrix (...);
+   *
+   *   for (auto &cell : dof_handler.active_cell_iterators())
+   *     {
+   *       cell_matrix = 0;
+   *       for (const auto &face : cell->face_iterators())
+   *         {
+   *           fe_iv.values.reinit(cell, face, ...);
+   *           for (const auto q : fe_iv.quadrature_point_indices())
+   *             for (const auto i : fe_iv.dof_indices())
+   *               for (const auto j : fe_iv.dof_indices())
+   *                 cell_matrix(i,j) += ...; // Do something for DoF indices
+   *                                          // (i,j) at quadrature point q
+   *         }
+   *     }
+   * @endcode
+   * Here, we are looping over all degrees of freedom on all cell interfaces,
+   * with `i` and `j` taking on all valid indices for interface degrees of
+   * freedom, as defined by the finite element passed to `fe_iv`.
+   */
+  std_cxx20::ranges::iota_view<unsigned int, unsigned int>
+  dof_indices() const;
+
   /**
    * Return the set of joint DoF indices. This includes indices from both cells.
    * If reinit was called with an active cell iterator, the indices are based
@@ -1341,6 +1373,15 @@ FEInterfaceValues<dim, spacedim>::n_current_interface_dofs() const
 
 
 
+template <int dim, int spacedim>
+inline std_cxx20::ranges::iota_view<unsigned int, unsigned int>
+FEInterfaceValues<dim, spacedim>::dof_indices() const
+{
+  return {0U, n_current_interface_dofs()};
+}
+
+
+
 template <int dim, int spacedim>
 bool
 FEInterfaceValues<dim, spacedim>::at_boundary() const
index 13d9bde58952e19d7812e098a4f7199e543aed20..815c18aa387eb95c5f447c7b8a22c571e3cedac9 100644 (file)
@@ -458,9 +458,6 @@ namespace Step74
       const FEInterfaceValues<dim> &fe_iv =
         scratch_data.reinit(cell, f, sf, ncell, nf, nsf);
 
-      const auto &       q_points   = fe_iv.get_quadrature_points();
-      const unsigned int n_q_points = q_points.size();
-
       copy_data.face_data.emplace_back();
       CopyDataFace &     copy_data_face = copy_data.face_data.back();
       const unsigned int n_dofs_face    = fe_iv.n_current_interface_dofs();
@@ -474,10 +471,10 @@ namespace Step74
       const double extent2 = ncell->measure() / ncell->face(nf)->measure();
       const double penalty = compute_penalty(degree, extent1, extent2);
 
-      for (unsigned int point = 0; point < n_q_points; ++point)
+      for (const unsigned int point : fe_iv.quadrature_point_indices())
         {
-          for (unsigned int i = 0; i < n_dofs_face; ++i)
-            for (unsigned int j = 0; j < n_dofs_face; ++j)
+          for (const unsigned int i : fe_iv.dof_indices())
+            for (const unsigned int j : fe_iv.dof_indices())
               copy_data_face.cell_matrix(i, j) +=
                 (-diffusion_coefficient *                 // - nu
                    fe_iv.jump_in_shape_values(i, point) * // [v_h]

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.