]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add MatrixFreeOperators::MassOperator
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 27 Oct 2016 13:41:58 +0000 (15:41 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 27 Oct 2016 20:37:47 +0000 (22:37 +0200)
doc/news/changes.h
include/deal.II/matrix_free/operators.h
tests/matrix_free/mass_operator_01.cc [new file with mode: 0644]
tests/matrix_free/mass_operator_01.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output [new file with mode: 0644]
tests/matrix_free/mass_operator_01.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=4.output [new file with mode: 0644]

index f4fb4d43d733fbe5c2ac0ea801fceb20d1d280a2..4d4a07a74a897649966468d63e366fcf38056726 100644 (file)
@@ -391,6 +391,10 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+ <li> New: Add MatrixFreeOperators::MassOperator representing a mass matrix.
+ <br>
+ (Daniel Arndt, 2016/10/27)
+ <li>
 
  <li> Fixed: GridIn::read_vtk() accidentally only read material ids of
  input cells correctly if the file listed them as integers. If they were
index c6dfc50a7081c683eed09dabd61086d656a2b679..24e6e510b6c65723e753cad376b07a8544f0d4f8 100644 (file)
@@ -158,7 +158,7 @@ namespace MatrixFreeOperators
     virtual void compute_diagonal () = 0;
 
     /**
-     * Get read access to diagonal of this operator.
+     * Get read access to the inverse diagonal of this operator.
      */
     const LinearAlgebra::distributed::Vector<number> &get_matrix_diagonal_inverse() const;
 
@@ -288,6 +288,46 @@ namespace MatrixFreeOperators
 
 
 
+  /**
+   * This class implements the operation of the action of a mass matrix.
+   *
+   * @author Daniel Arndt, 2016
+   */
+  template <int dim, int fe_degree, int n_components = 1, typename number = double>
+  class MassOperator : public Base<dim, number>
+  {
+  public:
+
+    /**
+     * Constructor. Initializes the MatrixFree object.
+     */
+    MassOperator ();
+
+    /**
+     * For preconditioning, we store a lumped mass matrix at the diagonal entries.
+     */
+    virtual void compute_diagonal ();
+
+  private:
+    /**
+     * Applies the mass matrix operation on an input vector. It is
+     * assumed that the passed input and output vector are correctly initialized
+     * using initialize_dof_vector().
+     */
+    virtual void apply_add (LinearAlgebra::distributed::Vector<number>       &dst,
+                            const LinearAlgebra::distributed::Vector<number> &src) const;
+
+    /**
+     * For this operator, there is just a cell contribution.
+     */
+    void local_apply_cell (const MatrixFree<dim,number>                     &data,
+                           LinearAlgebra::distributed::Vector<number>       &dst,
+                           const LinearAlgebra::distributed::Vector<number> &src,
+                           const std::pair<unsigned int,unsigned int>  &cell_range) const;
+  };
+
+
+
   // ------------------------------------ inline functions ---------------------
 
   template <int dim, int fe_degree, int n_components, typename Number>
@@ -717,6 +757,90 @@ namespace MatrixFreeOperators
     dst*= omega;
   }
 
+
+
+  //-----------------------------MassOperator----------------------------------
+
+  template <int dim, int fe_degree, int n_components, typename number>
+  MassOperator<dim, fe_degree, n_components, number>::
+  MassOperator ()
+    :
+    Base<dim, number>()
+  {}
+
+
+
+  template <int dim, int fe_degree, int n_components, typename number>
+  void
+  MassOperator<dim, fe_degree, n_components, number>::
+  compute_diagonal()
+  {
+    Assert((Base<dim, number>::data != NULL), ExcNotInitialized());
+
+    LinearAlgebra::distributed::Vector<number> ones;
+    Base<dim, number>::initialize_dof_vector(Base<dim, number>::inverse_diagonal_entries);
+    Base<dim, number>::initialize_dof_vector(ones);
+    Base<dim, number>::inverse_diagonal_entries = 1.;
+    ones = 1.;
+    ones.update_ghost_values();
+    apply_add(Base<dim, number>::inverse_diagonal_entries, ones);
+
+    const IndexSet &locally_owned_elements
+      = Base<dim, number>::inverse_diagonal_entries.locally_owned_elements();
+    locally_owned_elements.print(std::cout);
+    const std::vector<unsigned int> constrained_dofs
+      = Base<dim, number>::data->get_constrained_dofs();
+    for (unsigned int i=0; i<constrained_dofs.size(); ++i)
+      if (locally_owned_elements.is_element(constrained_dofs[i]))
+        Base<dim, number>::inverse_diagonal_entries(constrained_dofs[i]) = 1.;
+
+    const unsigned int local_size = Base<dim, number>::inverse_diagonal_entries.local_size();
+    for (unsigned int i=0; i<local_size; ++i)
+      {
+        const double tmp = Base<dim, number>::inverse_diagonal_entries.local_element(i);
+        Base<dim, number>::inverse_diagonal_entries.local_element(i)
+          =1./Base<dim, number>::inverse_diagonal_entries.local_element(i);
+      }
+
+    Base<dim, number>::inverse_diagonal_entries.compress(VectorOperation::insert);
+    Base<dim, number>::inverse_diagonal_entries.update_ghost_values();
+  }
+
+
+
+  template <int dim, int fe_degree, int n_components, typename number>
+  void
+  MassOperator<dim, fe_degree, n_components, number>::
+  apply_add (LinearAlgebra::distributed::Vector<number>       &dst,
+             const LinearAlgebra::distributed::Vector<number> &src) const
+  {
+    Base<dim, number>::data->cell_loop (&MassOperator::local_apply_cell,
+                                        this, dst, src);
+  }
+
+
+
+  template <int dim, int fe_degree, int n_components, typename number>
+  void
+  MassOperator<dim, fe_degree, n_components, number>::
+  local_apply_cell (const MatrixFree<dim,number>                     &data,
+                    LinearAlgebra::distributed::Vector<number>       &dst,
+                    const LinearAlgebra::distributed::Vector<number> &src,
+                    const std::pair<unsigned int,unsigned int>  &cell_range) const
+  {
+    FEEvaluation<dim, fe_degree, fe_degree+1, n_components, number> phi(*Base<dim, number>::data);
+    for (unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
+      {
+        phi.reinit (cell);
+        phi.read_dof_values(src);
+        phi.evaluate (true,false,false);
+        for (unsigned int q=0; q<phi.n_q_points; ++q)
+          phi.submit_value (phi.get_value(q), q);
+        phi.integrate (true,false);
+        phi.distribute_local_to_global (dst);
+      }
+  }
+
 } // end of namespace MatrixFreeOperators
 
 
diff --git a/tests/matrix_free/mass_operator_01.cc b/tests/matrix_free/mass_operator_01.cc
new file mode 100644 (file)
index 0000000..575b0bc
--- /dev/null
@@ -0,0 +1,228 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// this tests the correctness of MPI-parallel matrix free matrix-vector
+// products using MatrixFreeOperators::MassOperator by comparing the result
+// with a Trilinos sparse matrix assembled in the usual way. The mesh is
+// distributed among processors (hypercube) and has both hanging nodes
+// (by randomly refining some cells, so the mesh is going to be different
+// when run with different numbers of processors) and Dirichlet boundary
+// conditions
+//
+// this test does not use multithreading within the MPI nodes.
+
+#include "../tests.h"
+
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/function.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/dofs/dof_tools.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/lac/constraint_matrix.h>
+#include <deal.II/lac/trilinos_sparse_matrix.h>
+#include <deal.II/lac/trilinos_sparsity_pattern.h>
+#include <deal.II/matrix_free/operators.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_values.h>
+#include <deal.II/numerics/vector_tools.h>
+
+#include <iostream>
+
+
+
+
+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();
+  cell = tria.begin_active ();
+  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();
+    }
+
+  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);
+
+  ConstraintMatrix constraints (relevant_set);
+  DoFTools::make_hanging_node_constraints(dof, constraints);
+  VectorTools::interpolate_boundary_values (dof, 0, ZeroFunction<dim>(),
+                                            constraints);
+  constraints.close();
+
+  deallog << "Testing " << dof.get_fe().get_name() << std::endl;
+  //std::cout << "Number of cells: " << tria.n_global_active_cells() << std::endl;
+  //std::cout << "Number of degrees of freedom: " << dof.n_dofs() << std::endl;
+  //std::cout << "Number of constraints: " << constraints.n_constraints() << std::endl;
+
+  MatrixFree<dim,number> mf_data;
+  {
+    const QGauss<1> quad (fe_degree+1);
+    typename MatrixFree<dim,number>::AdditionalData data;
+    data.mpi_communicator = MPI_COMM_WORLD;
+    data.tasks_parallel_scheme =
+      MatrixFree<dim,number>::AdditionalData::none;
+    data.tasks_block_size = 7;
+    mf_data.reinit (dof, constraints, quad, data);
+  }
+
+  MatrixFreeOperators::MassOperator<dim,fe_degree, 1, number> mf;
+  mf.initialize(mf_data);
+  mf.compute_diagonal();
+  LinearAlgebra::distributed::Vector<number> in, out, ref;
+  mf_data.initialize_dof_vector (in);
+  out.reinit (in);
+  ref.reinit (in);
+
+  for (unsigned int i=0; i<in.local_size(); ++i)
+    {
+      const unsigned int glob_index =
+        owned_set.nth_index_in_set (i);
+      if (constraints.is_constrained(glob_index))
+        continue;
+      in.local_element(i) = (double)Testing::rand()/RAND_MAX;
+    }
+
+  mf.vmult (out, in);
+
+
+  // assemble trilinos sparse matrix with
+  // (v, u) for reference
+  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_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_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);
+  out -= ref;
+  const double diff_norm = out.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));
+
+  if (myid == 0)
+    {
+      std::ofstream logfile("output");
+      deallog.attach(logfile);
+      deallog << std::setprecision(4);
+      deallog.threshold_double(1.e-10);
+
+      deallog.push("2d");
+      test<2,1>();
+      test<2,2>();
+      deallog.pop();
+
+      deallog.push("3d");
+      test<3,1>();
+      test<3,2>();
+      deallog.pop();
+    }
+  else
+    {
+      test<2,1>();
+      test<2,2>();
+      test<3,1>();
+      test<3,2>();
+    }
+}
+
diff --git a/tests/matrix_free/mass_operator_01.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output b/tests/matrix_free/mass_operator_01.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=1.output
new file mode 100644 (file)
index 0000000..25e13d5
--- /dev/null
@@ -0,0 +1,13 @@
+
+DEAL:0:2d::Testing FE_Q<2>(1)
+DEAL:0:2d::Norm of difference: 0
+DEAL:0:2d::
+DEAL:0:2d::Testing FE_Q<2>(2)
+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/matrix_free/mass_operator_01.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=4.output b/tests/matrix_free/mass_operator_01.with_trilinos=true.with_mpi=true.with_p4est=true.mpirun=4.output
new file mode 100644 (file)
index 0000000..25e13d5
--- /dev/null
@@ -0,0 +1,13 @@
+
+DEAL:0:2d::Testing FE_Q<2>(1)
+DEAL:0:2d::Norm of difference: 0
+DEAL:0:2d::
+DEAL:0:2d::Testing FE_Q<2>(2)
+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.