--- /dev/null
+New: The method hp::FEFaceValues::reinit() can now also accept face iterators.
+<br>
+(Peter Munch, 2020/07/16)
#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>
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
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);
};
+ 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(
}
+
+ 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 -------------------------
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>>,
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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>();
+}
--- /dev/null
+
+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::