]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a version of CUDAWrappers::MatrixFree::reinit() that takes a predicate for cell_loop 11780/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 19 Feb 2021 20:57:55 +0000 (20:57 +0000)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 19 Feb 2021 20:57:55 +0000 (20:57 +0000)
doc/news/changes/minor/20210219Turcksin [new file with mode: 0644]
include/deal.II/matrix_free/cuda_matrix_free.h
include/deal.II/matrix_free/cuda_matrix_free.templates.h
tests/cuda/matrix_free_matrix_vector_25.cu [new file with mode: 0644]
tests/cuda/matrix_free_matrix_vector_25.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output [new file with mode: 0644]
tests/cuda/matrix_free_matrix_vector_25.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=2.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20210219Turcksin b/doc/news/changes/minor/20210219Turcksin
new file mode 100644 (file)
index 0000000..4bb58a6
--- /dev/null
@@ -0,0 +1,5 @@
+Added a template version of CUDAWrappers::MatrixFree::reinit() that takes an
+IteratorFilters object. This allows to perform the operator evaluation on part
+of the domain.
+<br>
+(Bruno Turcksin, 2021/02/19)
index 9d3ff0c5983e34ccf7141fbf238d8a2779f939c1..d8c0b63cc326c6175412021d7af2e85ad90c1c27 100644 (file)
@@ -241,7 +241,21 @@ namespace CUDAWrappers
      * degrees of freedom, the DoFHandler and the mapping describe the
      * transformation from unit to real cell, and the finite element
      * underlying the DoFHandler together with the quadrature formula
-     * describe the local operations.
+     * describe the local operations. This function takes an IteratorFilters
+     * object (predicate) to loop over a subset of the active cells. When using
+     * MPI, the predicate should filter out non locally owned cells.
+     */
+    template <typename IteratorFiltersType>
+    void
+    reinit(const Mapping<dim> &             mapping,
+           const DoFHandler<dim> &          dof_handler,
+           const AffineConstraints<Number> &constraints,
+           const Quadrature<1> &            quad,
+           const IteratorFiltersType &      iterator_filter,
+           const AdditionalData &           additional_data = AdditionalData());
+
+    /**
+     * Same as above using Iterators::LocallyOwnedCell() as predicate.
      */
     void
     reinit(const Mapping<dim> &             mapping,
@@ -384,11 +398,13 @@ namespace CUDAWrappers
     /**
      * Initializes the data structures.
      */
+    template <typename IteratorFiltersType>
     void
     internal_reinit(const Mapping<dim> &             mapping,
                     const DoFHandler<dim> &          dof_handler,
                     const AffineConstraints<Number> &constraints,
                     const Quadrature<1> &            quad,
+                    const IteratorFiltersType &      iterator_filter,
                     std::shared_ptr<const MPI_Comm>  comm,
                     const AdditionalData             additional_data);
 
index af8904b1d1baf522798352c21586ea7593ae80d5..af6c89adc78c2b6fc2adde69ec6e4f2eee12d1f1 100644 (file)
@@ -624,12 +624,14 @@ namespace CUDAWrappers
 
 
   template <int dim, typename Number>
+  template <typename IteratorFiltersType>
   void
   MatrixFree<dim, Number>::reinit(const Mapping<dim> &             mapping,
                                   const DoFHandler<dim> &          dof_handler,
                                   const AffineConstraints<Number> &constraints,
                                   const Quadrature<1> &            quad,
-                                  const AdditionalData &additional_data)
+                                  const IteratorFiltersType &iterator_filter,
+                                  const AdditionalData &     additional_data)
   {
     const auto &triangulation = dof_handler.get_triangulation();
     if (const auto parallel_triangulation =
@@ -639,12 +641,37 @@ namespace CUDAWrappers
                       dof_handler,
                       constraints,
                       quad,
+                      iterator_filter,
                       std::make_shared<const MPI_Comm>(
                         parallel_triangulation->get_communicator()),
                       additional_data);
     else
-      internal_reinit(
-        mapping, dof_handler, constraints, quad, nullptr, additional_data);
+      internal_reinit(mapping,
+                      dof_handler,
+                      constraints,
+                      quad,
+                      iterator_filter,
+                      nullptr,
+                      additional_data);
+  }
+
+
+
+  template <int dim, typename Number>
+  void
+  MatrixFree<dim, Number>::reinit(const Mapping<dim> &             mapping,
+                                  const DoFHandler<dim> &          dof_handler,
+                                  const AffineConstraints<Number> &constraints,
+                                  const Quadrature<1> &            quad,
+                                  const AdditionalData &additional_data)
+  {
+    IteratorFilters::LocallyOwnedCell locally_owned_cell_filter;
+    reinit(mapping,
+           dof_handler,
+           constraints,
+           quad,
+           locally_owned_cell_filter,
+           additional_data);
   }
 
 
@@ -841,12 +868,14 @@ namespace CUDAWrappers
 
 
   template <int dim, typename Number>
+  template <typename IteratorFiltersType>
   void
   MatrixFree<dim, Number>::internal_reinit(
     const Mapping<dim> &             mapping,
     const DoFHandler<dim> &          dof_handler_,
     const AffineConstraints<Number> &constraints,
     const Quadrature<1> &            quad,
+    const IteratorFiltersType &      iterator_filter,
     std::shared_ptr<const MPI_Comm>  comm,
     const AdditionalData             additional_data)
   {
@@ -950,9 +979,8 @@ namespace CUDAWrappers
       this, mapping, fe, quad, shape_info, *dof_handler, update_flags);
 
     // Create a graph coloring
-    CellFilter begin(IteratorFilters::LocallyOwnedCell(),
-                     dof_handler->begin_active());
-    CellFilter end(IteratorFilters::LocallyOwnedCell(), dof_handler->end());
+    CellFilter begin(iterator_filter, dof_handler->begin_active());
+    CellFilter end(iterator_filter, dof_handler->end());
 
     if (begin != end)
       {
diff --git a/tests/cuda/matrix_free_matrix_vector_25.cu b/tests/cuda/matrix_free_matrix_vector_25.cu
new file mode 100644 (file)
index 0000000..475270f
--- /dev/null
@@ -0,0 +1,252 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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 matrix_free_matrix_vector_10 test but we only perform the operator
+// evaluation on the cells that have a material ID equal to zero.
+
+#include <deal.II/base/function.h>
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_values.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/manifold_lib.h>
+
+#include <deal.II/lac/affine_constraints.h>
+#include <deal.II/lac/trilinos_sparse_matrix.h>
+#include <deal.II/lac/trilinos_sparsity_pattern.h>
+
+#include <deal.II/numerics/vector_tools.h>
+
+#include <iostream>
+
+#include "../tests.h"
+
+#include "matrix_vector_mf.h"
+
+
+
+template <int dim, int fe_degree>
+void
+test()
+{
+  typedef double Number;
+
+  parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(1);
+  typename Triangulation<dim>::active_cell_iterator cell = tria.begin_active(),
+                                                    endc = tria.end();
+  for (; cell != endc; ++cell)
+    if (cell->is_locally_owned())
+      if (cell->center().norm() < 0.2)
+        cell->set_refine_flag();
+  tria.execute_coarsening_and_refinement();
+  if (dim < 3 && fe_degree < 2)
+    tria.refine_global(2);
+  else
+    tria.refine_global(1);
+  if (tria.begin(tria.n_levels() - 1)->is_locally_owned())
+    tria.begin(tria.n_levels() - 1)->set_refine_flag();
+  if (tria.last()->is_locally_owned())
+    tria.last()->set_refine_flag();
+  tria.execute_coarsening_and_refinement();
+  cell = tria.begin_active();
+  for (unsigned int i = 0; i < 10 - 3 * dim; ++i)
+    {
+      cell                 = tria.begin_active();
+      unsigned int counter = 0;
+      for (; cell != endc; ++cell, ++counter)
+        if (cell->is_locally_owned())
+          if (counter % (7 - i) == 0)
+            cell->set_refine_flag();
+      tria.execute_coarsening_and_refinement();
+    }
+  // Set material ID
+  for (cell = tria.begin_active(); cell < endc; ++cell)
+    {
+      if (cell->center()[0] < 0.5)
+        cell->set_material_id(0);
+      else
+        cell->set_material_id(1);
+    }
+
+  FE_Q<dim>       fe(fe_degree);
+  DoFHandler<dim> dof(tria);
+  dof.distribute_dofs(fe);
+
+  IndexSet owned_set = dof.locally_owned_dofs();
+  IndexSet relevant_set;
+  DoFTools::extract_locally_relevant_dofs(dof, relevant_set);
+
+  AffineConstraints<double> constraints(relevant_set);
+  DoFTools::make_hanging_node_constraints(dof, constraints);
+  VectorTools::interpolate_boundary_values(dof,
+                                           0,
+                                           Functions::ZeroFunction<dim>(),
+                                           constraints);
+  constraints.close();
+
+  deallog << "Testing " << dof.get_fe().get_name() << std::endl;
+
+  MappingQGeneric<dim>                  mapping(fe_degree);
+  CUDAWrappers::MatrixFree<dim, Number> mf_data;
+  const QGauss<1>                       quad(fe_degree + 1);
+  typename CUDAWrappers::MatrixFree<dim, Number>::AdditionalData
+    additional_data;
+  additional_data.mapping_update_flags = update_values | update_gradients |
+                                         update_JxW_values |
+                                         update_quadrature_points;
+  IteratorFilters::MaterialIdEqualTo filter(0, true);
+  mf_data.reinit(mapping, dof, constraints, quad, filter, additional_data);
+
+  const unsigned int coef_size =
+    tria.n_locally_owned_active_cells() * std::pow(fe_degree + 1, dim);
+  MatrixFreeTest<dim,
+                 fe_degree,
+                 Number,
+                 LinearAlgebra::distributed::Vector<Number, MemorySpace::CUDA>>
+                                                                mf(mf_data, coef_size);
+  LinearAlgebra::distributed::Vector<Number, MemorySpace::CUDA> in_dev;
+  LinearAlgebra::distributed::Vector<Number, MemorySpace::CUDA> out_dev;
+  mf_data.initialize_dof_vector(in_dev);
+  mf_data.initialize_dof_vector(out_dev);
+
+  LinearAlgebra::ReadWriteVector<Number> rw_in(owned_set);
+  for (unsigned int i = 0; i < in_dev.local_size(); ++i)
+    {
+      const unsigned int glob_index = owned_set.nth_index_in_set(i);
+      if (constraints.is_constrained(glob_index))
+        continue;
+      rw_in.local_element(i) = random_value<double>();
+    }
+  in_dev.import(rw_in, VectorOperation::insert);
+  mf.vmult(out_dev, in_dev);
+
+  LinearAlgebra::distributed::Vector<Number, MemorySpace::Host> out_host(
+    owned_set, MPI_COMM_WORLD);
+  LinearAlgebra::ReadWriteVector<Number> rw_out(owned_set);
+  rw_out.import(out_dev, VectorOperation::insert);
+  out_host.import(rw_out, VectorOperation::insert);
+
+  // assemble trilinos sparse matrix with
+  // (\nabla v, \nabla u) + (v, 10 * u) for
+  // reference
+  LinearAlgebra::distributed::Vector<Number, MemorySpace::Host> in_host(
+    owned_set, MPI_COMM_WORLD);
+  in_host.import(rw_in, VectorOperation::insert);
+  LinearAlgebra::distributed::Vector<Number, MemorySpace::Host> ref(
+    owned_set, MPI_COMM_WORLD);
+  TrilinosWrappers::SparseMatrix sparse_matrix;
+  {
+    TrilinosWrappers::SparsityPattern csp(owned_set, MPI_COMM_WORLD);
+    DoFTools::make_sparsity_pattern(dof,
+                                    csp,
+                                    constraints,
+                                    true,
+                                    Utilities::MPI::this_mpi_process(
+                                      MPI_COMM_WORLD));
+    csp.compress();
+    sparse_matrix.reinit(csp);
+  }
+  {
+    QGauss<dim> quadrature_formula(fe_degree + 1);
+
+    FEValues<dim> fe_values(dof.get_fe(),
+                            quadrature_formula,
+                            update_values | update_gradients |
+                              update_JxW_values);
+
+    const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
+    const unsigned int n_q_points    = quadrature_formula.size();
+
+    FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
+    std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
+
+    typename DoFHandler<dim>::active_cell_iterator cell = dof.begin_active(),
+                                                   endc = dof.end();
+    for (; cell != endc; ++cell)
+      if (cell->is_locally_owned() && cell->material_id() == 0)
+        {
+          cell_matrix = 0;
+          fe_values.reinit(cell);
+
+          for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
+            for (unsigned int i = 0; i < dofs_per_cell; ++i)
+              {
+                for (unsigned int j = 0; j < dofs_per_cell; ++j)
+                  cell_matrix(i, j) +=
+                    ((fe_values.shape_grad(i, q_point) *
+                        fe_values.shape_grad(j, q_point) +
+                      10. * fe_values.shape_value(i, q_point) *
+                        fe_values.shape_value(j, q_point)) *
+                     fe_values.JxW(q_point));
+              }
+
+          cell->get_dof_indices(local_dof_indices);
+          constraints.distribute_local_to_global(cell_matrix,
+                                                 local_dof_indices,
+                                                 sparse_matrix);
+        }
+  }
+  sparse_matrix.compress(VectorOperation::add);
+
+  sparse_matrix.vmult(ref, in_host);
+  out_host -= ref;
+
+  const double diff_norm = out_host.linfty_norm();
+
+  deallog << "Norm of difference: " << diff_norm << std::endl << std::endl;
+}
+
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+  unsigned int myid = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+  deallog.push(Utilities::int_to_string(myid));
+
+  init_cuda(true);
+
+  if (myid == 0)
+    {
+      initlog();
+      deallog << std::setprecision(4);
+
+      deallog.push("2d");
+      test<2, 1>();
+      deallog.pop();
+
+      deallog.push("3d");
+      test<3, 1>();
+      test<3, 2>();
+      deallog.pop();
+    }
+  else
+    {
+      test<2, 1>();
+      test<3, 1>();
+      test<3, 2>();
+    }
+}
diff --git a/tests/cuda/matrix_free_matrix_vector_25.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output b/tests/cuda/matrix_free_matrix_vector_25.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output
new file mode 100644 (file)
index 0000000..3b90478
--- /dev/null
@@ -0,0 +1,10 @@
+
+DEAL:0:2d::Testing FE_Q<2>(1)
+DEAL:0:2d::Norm of difference: 0
+DEAL:0:2d::
+DEAL:0:3d::Testing FE_Q<3>(1)
+DEAL:0:3d::Norm of difference: 0
+DEAL:0:3d::
+DEAL:0:3d::Testing FE_Q<3>(2)
+DEAL:0:3d::Norm of difference: 0
+DEAL:0:3d::
diff --git a/tests/cuda/matrix_free_matrix_vector_25.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=2.output b/tests/cuda/matrix_free_matrix_vector_25.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=2.output
new file mode 100644 (file)
index 0000000..3b90478
--- /dev/null
@@ -0,0 +1,10 @@
+
+DEAL:0:2d::Testing FE_Q<2>(1)
+DEAL:0:2d::Norm of difference: 0
+DEAL:0:2d::
+DEAL:0:3d::Testing FE_Q<3>(1)
+DEAL:0:3d::Norm of difference: 0
+DEAL:0:3d::
+DEAL:0:3d::Testing FE_Q<3>(2)
+DEAL:0:3d::Norm of difference: 0
+DEAL:0:3d::

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.