]> https://gitweb.dealii.org/ - dealii.git/commitdiff
MatrixFreeTools::compute_diagonal: enable for hp 11308/head
authorPeter Munch <peterrmuench@gmail.com>
Thu, 3 Dec 2020 17:48:32 +0000 (18:48 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Fri, 11 Dec 2020 11:34:22 +0000 (12:34 +0100)
include/deal.II/matrix_free/tools.h
tests/matrix_free/compute_diagonal_04.cc [new file with mode: 0644]
tests/matrix_free/compute_diagonal_04.with_mpi=true.with_p4est=true.mpirun=1.output [new file with mode: 0644]
tests/matrix_free/compute_diagonal_util.h

index 1cd9540ab9299dde7e996860a3cfc9d250b03c55..8bd19d735e9c669265ec767bd8baf41ee303ac95 100644 (file)
@@ -287,7 +287,7 @@ namespace MatrixFreeTools
                      n_components,
                      Number,
                      VectorizedArrayType>
-          phi(matrix_free, dof_no, quad_no, first_selected_component);
+          phi(matrix_free, range, dof_no, quad_no, first_selected_component);
 
         const unsigned int n_lanes = VectorizedArrayType::size();
 
diff --git a/tests/matrix_free/compute_diagonal_04.cc b/tests/matrix_free/compute_diagonal_04.cc
new file mode 100644 (file)
index 0000000..0f7b986
--- /dev/null
@@ -0,0 +1,92 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Same as compute_diagonal_01 but for hp.
+
+#include "compute_diagonal_util.h"
+
+template <int dim,
+          typename Number              = double,
+          typename VectorizedArrayType = VectorizedArray<Number>>
+void
+test()
+{
+  const unsigned int n_components = 1;
+
+  parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+
+  GridGenerator::hyper_cube(tria, -numbers::PI / 2, numbers::PI / 2);
+
+  tria.refine_global();
+  for (auto &cell : tria.active_cell_iterators())
+    if (cell->active() && cell->is_locally_owned() && cell->center()[0] < 0.0)
+      cell->set_refine_flag();
+  tria.execute_coarsening_and_refinement();
+
+
+  hp::FECollection<dim> fes{FE_Q<dim>(1), FE_Q<dim>(2)};
+
+  // setup dof-handlers
+  DoFHandler<dim> dof_handler(tria);
+
+  unsigned int counter = 0;
+
+  for (const auto &cell : dof_handler.active_cell_iterators())
+    cell->set_active_fe_index(++counter % 2);
+
+  dof_handler.distribute_dofs(fes);
+
+  AffineConstraints<Number> constraint;
+  DoFTools::make_hanging_node_constraints(dof_handler, constraint);
+
+  VectorTools::interpolate_boundary_values(
+    dof_handler, 0, Functions::ZeroFunction<dim>(n_components), constraint);
+
+  constraint.close();
+
+  typename MatrixFree<dim, Number, VectorizedArrayType>::AdditionalData
+    additional_data;
+  additional_data.mapping_update_flags = update_values | update_gradients;
+
+  MappingQ<dim> mapping(1);
+  QGauss<1>     quad(2);
+
+  MatrixFree<dim, Number, VectorizedArrayType> matrix_free;
+  matrix_free.reinit(mapping, dof_handler, constraint, quad, additional_data);
+
+
+  Test<dim, -1, 0, n_components, Number, VectorizedArrayType> test(
+    matrix_free,
+    constraint,
+    [](FEEvaluation<dim, -1, 0, n_components, Number, VectorizedArrayType>
+         &phi) {
+      phi.evaluate(false, true);
+      for (unsigned int q = 0; q < phi.n_q_points; ++q)
+        phi.submit_gradient(phi.get_gradient(q), q);
+      phi.integrate(false, true);
+    });
+
+  test.do_test();
+}
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  test<2>();
+}
diff --git a/tests/matrix_free/compute_diagonal_04.with_mpi=true.with_p4est=true.mpirun=1.output b/tests/matrix_free/compute_diagonal_04.with_mpi=true.with_p4est=true.mpirun=1.output
new file mode 100644 (file)
index 0000000..f04f613
--- /dev/null
@@ -0,0 +1,16 @@
+
+Process #0
+Local range: [0, 40), global size: 40
+Vector data:
+0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 4.741e+00 3.296e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 0.000e+00 0.000e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 4.741e+00 0.000e+00 
+DEAL:0::13.5700
+Process #0
+Local range: [0, 40), global size: 40
+Vector data:
+0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 4.741e+00 3.296e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 0.000e+00 0.000e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 4.741e+00 0.000e+00 
+DEAL:0::13.5700
+Process #0
+Local range: [0, 40), global size: 40
+Vector data:
+0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 4.741e+00 3.296e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 0.000e+00 0.000e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 3.259e+00 4.741e+00 3.111e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 4.741e+00 0.000e+00 
+DEAL:0::13.5700
index 2dd90473743c6d19f86ec4d1812551f0261e5106..6d2a3c5ddda641a034667e40ac22c24a7978a620 100644 (file)
@@ -208,7 +208,7 @@ public:
                  n_components,
                  Number,
                  VectorizedArrayType>
-      phi(data);
+      phi(data, pair);
     for (auto cell = pair.first; cell < pair.second; cell++)
       {
         phi.reinit(cell);

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.