]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Disable hold_all_faces_to_owned_cells for FE_Q 10087/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 8 May 2020 13:22:37 +0000 (15:22 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Fri, 8 May 2020 20:11:20 +0000 (22:11 +0200)
include/deal.II/matrix_free/matrix_free.templates.h
tests/matrix_free/partitioner_02.cc [new file with mode: 0644]
tests/matrix_free/partitioner_02.output [new file with mode: 0644]

index 6201082ce371abfd33cf71ce1d1270a1d1888d13..1b05a4712cce53e8477a2c1cd5646d25ed2e9799 100644 (file)
@@ -746,6 +746,8 @@ MatrixFree<dim, Number, VectorizedArrayType>::initialize_indices(
   const unsigned int n_fe           = dof_handlers.n_dof_handlers;
   const unsigned int n_active_cells = cell_level_index.size();
 
+  std::vector<bool> 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<dim, Number, VectorizedArrayType>::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<dim, Number, VectorizedArrayType>::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<dim, Number, VectorizedArrayType>::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<dim, Number, VectorizedArrayType>::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 (file)
index 0000000..8330ce6
--- /dev/null
@@ -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 <deal.II/base/mpi.h>
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_dgq.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/matrix_free/fe_evaluation.h>
+#include <deal.II/matrix_free/matrix_free.h>
+
+#include <deal.II/multigrid/mg_constrained_dofs.h>
+
+#include "../tests.h"
+
+
+template <int dim,
+          typename Number              = double,
+          typename VectorizedArrayType = VectorizedArray<Number>>
+void
+test(const std::vector<FiniteElement<dim, dim> *> &finite_elements)
+{
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(2);
+
+  std::vector<std::shared_ptr<DoFHandler<dim, dim>>> dof_handlers(
+    finite_elements.size());
+  std::vector<const DoFHandler<dim, dim> *> dof_handlers_(dof_handlers.size());
+
+  for (unsigned int i = 0; i < dof_handlers.size(); ++i)
+    {
+      dof_handlers[i] = std::make_shared<DoFHandler<dim, dim>>(tria);
+      dof_handlers[i]->distribute_dofs(*finite_elements[i]);
+
+      dof_handlers_[i] = dof_handlers[i].get();
+    }
+
+  MappingQ<dim> mapping(1);
+
+  std::vector<Quadrature<1>> quads{QGauss<1>(finite_elements[0]->degree + 1)};
+
+  AffineConstraints<Number> constraint;
+  constraint.close();
+
+  std::vector<const AffineConstraints<Number> *> constraints(
+    dof_handlers.size());
+  for (unsigned int i = 0; i < dof_handlers.size(); ++i)
+    constraints[i] = &constraint;
+
+
+  typename MatrixFree<dim, Number, VectorizedArrayType>::AdditionalData
+    additional_data;
+  additional_data.tasks_parallel_scheme =
+    MatrixFree<dim, Number, VectorizedArrayType>::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<dim, Number, VectorizedArrayType> 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<dim> fe_dgq(1);
+    test<dim>({&fe_dgq});
+    deallog.pop();
+  }
+
+  {
+    const unsigned int dim = 2;
+    deallog.push("2d-q");
+    FE_Q<dim> fe_q(1);
+    test<dim>({&fe_q});
+    deallog.pop();
+  }
+
+  {
+    const unsigned int dim = 2;
+    deallog.push("2d-dgq-q");
+    FE_Q<dim>   fe_q(1);
+    FE_DGQ<dim> fe_dgq(1);
+    test<dim>({&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 (file)
index 0000000..919ab88
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL:2d-dgq::OK!
+DEAL:2d-q::OK!
+DEAL:2d-dgq-q::OK!

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.