From: Peter Munch Date: Wed, 15 Jul 2020 10:47:45 +0000 (+0200) Subject: Add new hp::FEFaceValues::reinit() functions X-Git-Tag: v9.3.0-rc1~1271^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=9a610454ced61ab4d5e762569860faae8815f355;p=dealii.git Add new hp::FEFaceValues::reinit() functions --- diff --git a/doc/news/changes/minor/20200716Munch b/doc/news/changes/minor/20200716Munch new file mode 100644 index 0000000000..4e5125b76d --- /dev/null +++ b/doc/news/changes/minor/20200716Munch @@ -0,0 +1,3 @@ +New: The method hp::FEFaceValues::reinit() can now also accept face iterators. +
+(Peter Munch, 2020/07/16) diff --git a/include/deal.II/hp/fe_values.h b/include/deal.II/hp/fe_values.h index 2debb37cb5..9b2df92b93 100644 --- a/include/deal.II/hp/fe_values.h +++ b/include/deal.II/hp/fe_values.h @@ -23,6 +23,9 @@ #include #include +#include +#include + #include #include #include @@ -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 + void + reinit(const TriaIterator> & cell, + const typename Triangulation::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::cell_iterator &cell, + const typename Triangulation::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); }; diff --git a/source/hp/fe_values.cc b/source/hp/fe_values.cc index aa7c60f983..6c1db5ab08 100644 --- a/source/hp/fe_values.cc +++ b/source/hp/fe_values.cc @@ -381,6 +381,22 @@ namespace hp + template + template + void + FEFaceValues::reinit( + const TriaIterator> & cell, + const typename Triangulation::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 void FEFaceValues::reinit( @@ -417,6 +433,21 @@ namespace hp } + + template + void + FEFaceValues::reinit( + const typename Triangulation::cell_iterator &cell, + const typename Triangulation::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 ------------------------- diff --git a/source/hp/fe_values.inst.in b/source/hp/fe_values.inst.in index 592ab00cbd..01acaeb465 100644 --- a/source/hp/fe_values.inst.in +++ b/source/hp/fe_values.inst.in @@ -105,6 +105,15 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS; unsigned int, unsigned int); template void + FEFaceValues::reinit( + const TriaIterator< + DoFCellAccessor> &, + const Triangulation::face_iterator &, + unsigned int, + unsigned int, + unsigned int); + template void FESubfaceValues::reinit( TriaIterator< DoFCellAccessor>, diff --git a/tests/hp/fe_face_values_reinit_01.cc b/tests/hp/fe_face_values_reinit_01.cc new file mode 100644 index 0000000000..f725949283 --- /dev/null +++ b/tests/hp/fe_face_values_reinit_01.cc @@ -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 +#include + +#include +#include +#include + +#include +#include +#include +#include +#include + +#include +#include +#include + +#include "../tests.h" + + + +template +void +test() +{ + Triangulation triangulation; + GridGenerator::hyper_cube(triangulation); + + hp::FECollection fe_collection(FE_Q(1)); + hp::QCollection quad_collection(QGauss(2)); + hp::MappingCollection mapping_collection(MappingQ(1)); + + hp::DoFHandler dof_handler(triangulation); + + dof_handler.distribute_dofs(fe_collection); + + + + hp::FEFaceValues 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 index 0000000000..d7ee8d4a7c --- /dev/null +++ b/tests/hp/fe_face_values_reinit_01.output @@ -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::