From: Peter Munch Date: Fri, 8 May 2020 13:22:37 +0000 (+0200) Subject: Disable hold_all_faces_to_owned_cells for FE_Q X-Git-Tag: v9.2.0-rc1~66^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=351678d2235bae5a5ca02daa6790f3d1e6db793d;p=dealii.git Disable hold_all_faces_to_owned_cells for FE_Q --- diff --git a/include/deal.II/matrix_free/matrix_free.templates.h b/include/deal.II/matrix_free/matrix_free.templates.h index 6201082ce3..1b05a4712c 100644 --- a/include/deal.II/matrix_free/matrix_free.templates.h +++ b/include/deal.II/matrix_free/matrix_free.templates.h @@ -746,6 +746,8 @@ MatrixFree::initialize_indices( const unsigned int n_fe = dof_handlers.n_dof_handlers; const unsigned int n_active_cells = cell_level_index.size(); + std::vector is_fe_dg(n_fe, false); + AssertDimension(n_active_cells, cell_level_index.size()); AssertDimension(n_fe, locally_owned_set.size()); AssertDimension(n_fe, constraint.size()); @@ -769,8 +771,11 @@ MatrixFree::initialize_indices( fes.push_back(&fe[f]); if (fe.size() > 1) - dof_info[no].cell_active_fe_index.resize( - n_active_cells, numbers::invalid_unsigned_int); + { + dof_info[no].cell_active_fe_index.resize( + n_active_cells, numbers::invalid_unsigned_int); + is_fe_dg[no] = fe[0].dofs_per_vertex == 0; + } Assert(additional_data.cell_vectorization_category.empty(), ExcNotImplemented()); @@ -782,6 +787,7 @@ MatrixFree::initialize_indices( if (cell_categorization_enabled == true) dof_info[no].cell_active_fe_index.resize( n_active_cells, numbers::invalid_unsigned_int); + is_fe_dg[no] = dofh->get_fe().dofs_per_vertex == 0; } lexicographic[no].resize(fes.size()); @@ -1428,8 +1434,10 @@ MatrixFree::initialize_indices( } // compute tighter index sets for various sets of face integrals - for (internal::MatrixFreeFunctions::DoFInfo &di : dof_info) + for (unsigned int no = 0; no < n_fe; ++no) { + internal::MatrixFreeFunctions::DoFInfo &di = dof_info[no]; + const Utilities::MPI::Partitioner &part = *di.vector_partitioner; // partitioner 0: no face integrals, simply use the indices present @@ -1788,7 +1796,7 @@ MatrixFree::initialize_indices( loop_over_faces); - if (additional_data.hold_all_faces_to_owned_cells) + if (additional_data.hold_all_faces_to_owned_cells && is_fe_dg[no]) { ghost_indices.clear(); // partitioner 3: values on all faces diff --git a/tests/matrix_free/partitioner_02.cc b/tests/matrix_free/partitioner_02.cc new file mode 100644 index 0000000000..8330ce6865 --- /dev/null +++ b/tests/matrix_free/partitioner_02.cc @@ -0,0 +1,130 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 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. +// +// --------------------------------------------------------------------- + + + +// tests matrix-free partitioners for update_ghost_values and compress(add) + + +#include +#include + +#include + +#include + +#include +#include +#include + +#include + +#include +#include + +#include + +#include "../tests.h" + + +template > +void +test(const std::vector *> &finite_elements) +{ + Triangulation tria; + GridGenerator::hyper_cube(tria); + tria.refine_global(2); + + std::vector>> dof_handlers( + finite_elements.size()); + std::vector *> dof_handlers_(dof_handlers.size()); + + for (unsigned int i = 0; i < dof_handlers.size(); ++i) + { + dof_handlers[i] = std::make_shared>(tria); + dof_handlers[i]->distribute_dofs(*finite_elements[i]); + + dof_handlers_[i] = dof_handlers[i].get(); + } + + MappingQ mapping(1); + + std::vector> quads{QGauss<1>(finite_elements[0]->degree + 1)}; + + AffineConstraints constraint; + constraint.close(); + + std::vector *> constraints( + dof_handlers.size()); + for (unsigned int i = 0; i < dof_handlers.size(); ++i) + constraints[i] = &constraint; + + + typename MatrixFree::AdditionalData + additional_data; + additional_data.tasks_parallel_scheme = + MatrixFree::AdditionalData::none; + additional_data.mapping_update_flags = + update_gradients | update_JxW_values | update_quadrature_points; + additional_data.mapping_update_flags_inner_faces = + update_gradients | update_JxW_values | update_quadrature_points; + additional_data.mapping_update_flags_boundary_faces = + update_gradients | update_JxW_values | update_quadrature_points; + additional_data.mapping_update_flags_faces_by_cells = + update_gradients | update_JxW_values | update_quadrature_points; + additional_data.hold_all_faces_to_owned_cells = true; + + MatrixFree matrix_free; + + matrix_free.reinit( + mapping, dof_handlers_, constraints, quads, additional_data); + + deallog << "OK!" << std::endl; +} + + + +int +main() +{ + initlog(); + + { + const unsigned int dim = 2; + deallog.push("2d-dgq"); + FE_DGQ fe_dgq(1); + test({&fe_dgq}); + deallog.pop(); + } + + { + const unsigned int dim = 2; + deallog.push("2d-q"); + FE_Q fe_q(1); + test({&fe_q}); + deallog.pop(); + } + + { + const unsigned int dim = 2; + deallog.push("2d-dgq-q"); + FE_Q fe_q(1); + FE_DGQ fe_dgq(1); + test({&fe_q, &fe_dgq}); + deallog.pop(); + } +} diff --git a/tests/matrix_free/partitioner_02.output b/tests/matrix_free/partitioner_02.output new file mode 100644 index 0000000000..919ab8832b --- /dev/null +++ b/tests/matrix_free/partitioner_02.output @@ -0,0 +1,4 @@ + +DEAL:2d-dgq::OK! +DEAL:2d-q::OK! +DEAL:2d-dgq-q::OK!