This is convenient when paired with the new cell->face_iterators() function.
--- /dev/null
+New: Added two new overloads of FEValues::reinit() that take face and subface
+iterators instead of face and subface numbers.
+<br>
+(David Wells, 2019/10/27)
level_dof_access>> &cell,
const unsigned int face_no);
+ /**
+ * Reinitialize the gradients, Jacobi determinants, etc for face @p face
+ * and cell @p cell.
+ *
+ * @note @p face must be one of @p cell's face iterators.
+ */
+ template <template <int, int> class DoFHandlerType, bool level_dof_access>
+ void
+ reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
+ level_dof_access>> & cell,
+ const typename Triangulation<dim, spacedim>::face_iterator &face);
+
/**
* Reinitialize the gradients, Jacobi determinants, etc for the given face
* on a given cell of type "iterator into a Triangulation object", and the
reinit(const typename Triangulation<dim, spacedim>::cell_iterator &cell,
const unsigned int face_no);
+ /*
+ * Reinitialize the gradients, Jacobi determinants, etc for the given face
+ * on a given cell of type "iterator into a Triangulation object", and the
+ * given finite element. Since iterators into a triangulation alone only
+ * convey information about the geometry of a cell, but not about degrees of
+ * freedom possibly associated with this cell, you will not be able to call
+ * some functions of this class if they need information about degrees of
+ * freedom. These functions are, above all, the
+ * <tt>get_function_value/gradients/hessians/third_derivatives</tt>
+ * functions. If you want to call these functions, you have to call the @p
+ * reinit variants that take iterators into DoFHandler or other DoF handler
+ * type objects.
+ *
+ * @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);
+
/**
* Return a reference to this very object.
*
const unsigned int face_no,
const unsigned int subface_no);
+ /**
+ * Alternative reinitialization function that takes, as arguments, iterators
+ * to the face and subface instead of their numbers.
+ */
+ template <template <int, int> class DoFHandlerType, bool level_dof_access>
+ void
+ reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
+ level_dof_access>> & cell,
+ const typename Triangulation<dim, spacedim>::face_iterator &face,
+ const typename Triangulation<dim, spacedim>::face_iterator &subface);
+
/**
* Reinitialize the gradients, Jacobi determinants, etc for the given
* subface on a given cell of type "iterator into a Triangulation object", and
const unsigned int face_no,
const unsigned int subface_no);
+ /**
+ * Reinitialize the gradients, Jacobi determinants, etc for the given
+ * subface on a given cell of type "iterator into a Triangulation object", and
+ * the given finite element. Since iterators into a triangulation alone only
+ * convey information about the geometry of a cell, but not about degrees of
+ * freedom possibly associated with this cell, you will not be able to call
+ * some functions of this class if they need information about degrees of
+ * freedom. These functions are, above all, the
+ * <tt>get_function_value/gradients/hessians/third_derivatives</tt>
+ * functions. If you want to call these functions, you have to call the @p
+ * reinit variants that take iterators into DoFHandler or other DoF handler
+ * type objects.
+ *
+ * This does the same thing as the previous function but takes iterators
+ * instead of numbers as arguments.
+ *
+ * @note @p face and @p subface must correspond to a face (and a subface of
+ * that face) of @p cell.
+ */
+ void
+ reinit(const typename Triangulation<dim, spacedim>::cell_iterator &cell,
+ const typename Triangulation<dim, spacedim>::face_iterator &face,
+ const typename Triangulation<dim, spacedim>::face_iterator &subface);
+
/**
* Return a reference to this very object.
*
+template <int dim, int spacedim>
+template <template <int, int> class DoFHandlerType, bool lda>
+void
+FEFaceValues<dim, spacedim>::reinit(
+ const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>> &cell,
+ const typename Triangulation<dim, spacedim>::face_iterator & face)
+{
+ const auto face_n = cell->face_iterator_to_index(face);
+ reinit(cell, face_n);
+}
+
+
+
template <int dim, int spacedim>
void
FEFaceValues<dim, spacedim>::reinit(
+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 auto face_n = cell->face_iterator_to_index(face);
+ reinit(cell, face_n);
+}
+
+
+
template <int dim, int spacedim>
void
FEFaceValues<dim, spacedim>::do_reinit(const unsigned int face_no)
+template <int dim, int spacedim>
+template <template <int, int> class DoFHandlerType, bool lda>
+void
+FESubfaceValues<dim, spacedim>::reinit(
+ const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>> &cell,
+ const typename Triangulation<dim, spacedim>::face_iterator & face,
+ const typename Triangulation<dim, spacedim>::face_iterator &subface)
+{
+ reinit(cell,
+ cell->face_iterator_to_index(face),
+ face->child_iterator_to_index(subface));
+}
+
+
+
template <int dim, int spacedim>
void
FESubfaceValues<dim, spacedim>::reinit(
+template <int dim, int spacedim>
+void
+FESubfaceValues<dim, spacedim>::reinit(
+ const typename Triangulation<dim, spacedim>::cell_iterator &cell,
+ const typename Triangulation<dim, spacedim>::face_iterator &face,
+ const typename Triangulation<dim, spacedim>::face_iterator &subface)
+{
+ reinit(cell,
+ cell->face_iterator_to_index(face),
+ face->child_iterator_to_index(subface));
+}
+
+
+
template <int dim, int spacedim>
void
FESubfaceValues<dim, spacedim>::do_reinit(const unsigned int face_no,
DoFCellAccessor<dof_handler<deal_II_dimension, deal_II_space_dimension>,
lda>> &,
unsigned int);
+
+ template void
+ FEFaceValues<deal_II_dimension, deal_II_space_dimension>::reinit(
+ const TriaIterator<
+ DoFCellAccessor<dof_handler<deal_II_dimension, deal_II_space_dimension>,
+ lda>> &,
+ const TriaIterator<TriaAccessor<deal_II_dimension - 1,
+ deal_II_dimension,
+ deal_II_space_dimension>> &);
+
template void
FESubfaceValues<deal_II_dimension, deal_II_space_dimension>::reinit(
const TriaIterator<
DoFCellAccessor<dof_handler<deal_II_dimension, deal_II_space_dimension>,
lda>> &,
- unsigned int,
- unsigned int);
+ const unsigned int,
+ const unsigned int);
+
+ template void
+ FESubfaceValues<deal_II_dimension, deal_II_space_dimension>::reinit(
+ const TriaIterator<
+ DoFCellAccessor<dof_handler<deal_II_dimension, deal_II_space_dimension>,
+ lda>> &,
+ const TriaIterator<TriaAccessor<deal_II_dimension - 1,
+ deal_II_dimension,
+ deal_II_space_dimension>> &,
+ const TriaIterator<TriaAccessor<deal_II_dimension - 1,
+ deal_II_dimension,
+ deal_II_space_dimension>> &);
#endif
}
q,
UpdateFlags(uflags | update_quadrature_points));
+ // Test the iterator-based initialization routines too:
+ FEFaceValues<dim> fe2(mapping,
+ finel,
+ q,
+ UpdateFlags(uflags | update_quadrature_points));
+ FESubfaceValues<dim> sub2(mapping,
+ finel,
+ q,
+ UpdateFlags(uflags | update_quadrature_points));
+
+
sprintf(fname, "Face%dd-%s", dim, name);
deallog.push(fname);
if (!c->face(f)->has_children())
{
fe.reinit(c, f);
+ fe2.reinit(c, c->face(f));
unsigned int k = 0;
for (unsigned int my = 0; my <= ((dim > 2) ? div : 0); ++my)
c)
{
if (uflags & update_values)
- AssertThrow((fe.shape_value(i, k) ==
- fe.shape_value_component(i, k, c)),
- ExcInternalError());
+ {
+ AssertThrow(
+ (fe.shape_value(i, k) ==
+ fe.shape_value_component(i, k, c)),
+ ExcInternalError());
+ AssertThrow((fe.shape_value(i, k) ==
+ fe2.shape_value(i, k)),
+ ExcInternalError());
+ }
if (uflags & update_gradients)
- AssertThrow((fe.shape_grad(i, k) ==
- fe.shape_grad_component(i, k, c)),
- ExcInternalError());
+ {
+ AssertThrow(
+ (fe.shape_grad(i, k) ==
+ fe.shape_grad_component(i, k, c)),
+ ExcInternalError());
+ AssertThrow((fe.shape_grad(i, k) ==
+ fe2.shape_grad(i, k)),
+ ExcInternalError());
+ }
if (uflags & update_hessians)
- AssertThrow(
- (fe.shape_hessian(i, k) ==
- fe.shape_hessian_component(i, k, c)),
- ExcInternalError());
+ {
+ AssertThrow(
+ (fe.shape_hessian(i, k) ==
+ fe.shape_hessian_component(i, k, c)),
+ ExcInternalError());
+ AssertThrow((fe.shape_hessian(i, k) ==
+ fe2.shape_hessian(i, k)),
+ ExcInternalError());
+ }
if (uflags & update_3rd_derivatives)
- AssertThrow(
- (fe.shape_3rd_derivative(i, k) ==
- fe.shape_3rd_derivative_component(i, k, c)),
- ExcInternalError());
+ {
+ AssertThrow(
+ (fe.shape_3rd_derivative(i, k) ==
+ fe.shape_3rd_derivative_component(i,
+ k,
+ c)),
+ ExcInternalError());
+ AssertThrow((fe.shape_3rd_derivative(i, k) ==
+ fe2.shape_3rd_derivative(i, k)),
+ ExcInternalError());
+ }
}
else
{
++s)
{
sub.reinit(c, f, s);
+ sub2.reinit(c, c->face(f), c->face(f)->child(s));
unsigned int k = 0;
for (unsigned int my = 0; my <= ((dim > 2) ? div : 0); ++my)
const double v1 = sub.shape_value(i, k),
v2 =
sub.shape_value_component(
- i, k, c);
+ i, k, c),
+ v3 = sub2.shape_value(i, k);
Assert(v1 == v2, ExcInternalError());
+ Assert(v1 == v3, ExcInternalError());
}
if (uflags & update_gradients)
{
const Tensor<1, dim>
g1 = sub.shape_grad(i, k),
- g2 = sub.shape_grad_component(i, k, c);
+ g2 = sub.shape_grad_component(i, k, c),
+ g3 = sub2.shape_grad(i, k);
Assert(g1 == g2, ExcInternalError());
+ Assert(g1 == g3, ExcInternalError());
}
if (uflags & update_hessians)
{
const Tensor<2, dim>
s1 = sub.shape_hessian(i, k),
s2 =
- sub.shape_hessian_component(i, k, c);
+ sub.shape_hessian_component(i, k, c),
+ s3 = sub2.shape_hessian(i, k);
Assert(s1 == s2, ExcInternalError());
+ Assert(s1 == s3, ExcInternalError());
}
if (uflags & update_3rd_derivatives)
{
t2 =
sub.shape_3rd_derivative_component(i,
k,
- c);
+ c),
+ t3 = sub2.shape_3rd_derivative(i, k);
+
Assert(t1 == t2, ExcInternalError());
+ Assert(t1 == t3, ExcInternalError());
}
}
else