--- /dev/null
+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)
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();
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]
#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>
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
+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
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();
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]