]> https://gitweb.dealii.org/ - dealii.git/commitdiff
make sure multivector_inner_product() updates Property of LAPACKFullMatrix accordingly
authorDenis Davydov <davydden@gmail.com>
Mon, 2 Oct 2017 18:00:26 +0000 (20:00 +0200)
committerDenis Davydov <davydden@gmail.com>
Tue, 3 Oct 2017 06:26:36 +0000 (08:26 +0200)
include/deal.II/lac/la_parallel_block_vector.templates.h
tests/mpi/parallel_block_vector_06.cc [new file with mode: 0644]
tests/mpi/parallel_block_vector_06.with_p4est=true.with_lapack=true.mpirun=3.output [new file with mode: 0644]

index 64c1afca1b44cedfac02d332a346ba0c51c8fab1..8902adcd7a20f98ba2a83bc7fca87f270f459ff6 100644 (file)
 #include <deal.II/lac/la_parallel_block_vector.h>
 #include <deal.II/lac/petsc_parallel_block_vector.h>
 #include <deal.II/lac/trilinos_parallel_block_vector.h>
+#include <deal.II/lac/lapack_support.h>
 
 
 DEAL_II_NAMESPACE_OPEN
 
+//Forward type declaration to have special treatment of LAPACKFullMatrix<number>
+// in multivector_inner_product()
+template <typename Number> class LAPACKFullMatrix;
 
 namespace LinearAlgebra
 {
@@ -781,6 +785,27 @@ namespace LinearAlgebra
 
 
 
+    namespace
+    {
+      template <typename FullMatrixType>
+      inline
+      void set_symmetric(FullMatrixType &, const bool)
+      {
+      }
+
+      template <typename NumberType>
+      inline
+      void set_symmetric(LAPACKFullMatrix<NumberType> &matrix, const bool symmetric)
+      {
+        if (symmetric)
+          matrix.set_property (LAPACKSupport::symmetric);
+        else
+          matrix.set_property (LAPACKSupport::general);
+      }
+    }
+
+
+
     template <typename Number>
     template <typename FullMatrixType>
     void
@@ -806,6 +831,7 @@ namespace LinearAlgebra
       // reset the matrix
       matrix = Number();
 
+      set_symmetric(matrix, symmetric);
       if (symmetric)
         {
           Assert (m == n,
diff --git a/tests/mpi/parallel_block_vector_06.cc b/tests/mpi/parallel_block_vector_06.cc
new file mode 100644 (file)
index 0000000..b8d77da
--- /dev/null
@@ -0,0 +1,204 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 BlockVector<Number>::scalar_product().
+// same as parallel_block_vector_05, but uses LAPACKFUllMatrix and tests
+// that internally we set it to be symmetric and so, for example, one can
+// immediately run Cholesky factorization
+
+#include "../tests.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/matrix_free/operators.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/numerics/vector_tools.h>
+#include <deal.II/lac/la_parallel_block_vector.h>
+#include <deal.II/lac/lapack_full_matrix.h>
+
+#include <iostream>
+
+
+
+
+template <int dim, int fe_degree>
+void test (const unsigned int n_blocks = 5)
+{
+  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, Functions::ZeroFunction<dim>(),
+                                            constraints);
+  constraints.close();
+
+  std::shared_ptr<MatrixFree<dim,number> > mf_data(new MatrixFree<dim,number> ());
+  {
+    const QGauss<1> quad (fe_degree+2);
+    typename MatrixFree<dim,number>::AdditionalData data;
+    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, fe_degree+2, 1, LinearAlgebra::distributed::Vector<number> > mf;
+  mf.initialize(mf_data);
+  mf.compute_diagonal();
+
+  LinearAlgebra::distributed::BlockVector<number> left(n_blocks), right(n_blocks);
+  for (unsigned int b = 0; b < n_blocks; ++b)
+    {
+      mf_data->initialize_dof_vector (left.block(b));
+      mf_data->initialize_dof_vector (right.block(b));
+
+      for (unsigned int i=0; i<right.block(b).local_size(); ++i)
+        {
+          const unsigned int glob_index =
+            owned_set.nth_index_in_set (i);
+          if (constraints.is_constrained(glob_index))
+            continue;
+          right.block(b).local_element(i) = (double)Testing::rand()/RAND_MAX;
+        }
+
+      mf.vmult (left.block(b), right.block(b));
+    }
+
+  {
+    // general lapack matrices:
+    LAPACKFullMatrix<number> product(n_blocks,n_blocks), diff(n_blocks,n_blocks);
+
+    deallog << "Symmetric:" << std::endl;
+
+    // inside LAPACK matrix should have its Property set as symmetric
+    // and so we can call compute_cholesky_factorization() below.
+    left.multivector_inner_product(product, right, true);
+
+    for (unsigned int i = 0; i < n_blocks; ++i)
+      for (unsigned int j = 0; j < n_blocks; ++j)
+        diff(i,j) = left.block(i) * right.block(j) - product(i,j);
+
+    const double diff_norm = diff.frobenius_norm();
+
+    deallog << "Norm of difference: " << diff_norm << std::endl;
+    deallog << "Cholesky: " << std::flush;
+    product.compute_cholesky_factorization();
+    deallog << "Ok" << std::endl;
+  }
+
+  // Non-symmetric case
+  {
+    const unsigned int n_blocks_2 = n_blocks+3;
+    LinearAlgebra::distributed::BlockVector<number> left2(n_blocks_2);
+    for (unsigned int b = 0; b < left2.n_blocks(); ++b)
+      {
+        mf_data->initialize_dof_vector (left2.block(b));
+        for (unsigned int i=0; i<left2.block(b).local_size(); ++i)
+          {
+            const unsigned int glob_index =
+              owned_set.nth_index_in_set (i);
+            if (constraints.is_constrained(glob_index))
+              continue;
+            left2.block(b).local_element(i) = (double)Testing::rand()/RAND_MAX;
+          }
+      }
+
+    LAPACKFullMatrix<number> product2(n_blocks_2,n_blocks), diff2(n_blocks_2,n_blocks);
+
+    deallog << "Non-Symmetric:" << std::endl;
+
+    left2.multivector_inner_product(product2, right, false);
+
+    for (unsigned int i = 0; i < n_blocks_2; ++i)
+      for (unsigned int j = 0; j < n_blocks; ++j)
+        diff2(i,j) = left2.block(i) * right.block(j) - product2(i,j);
+
+    const double diff_norm2 = diff2.frobenius_norm();
+
+    deallog << "Norm of difference: " << diff_norm2 << 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)
+    {
+      initlog();
+      deallog << std::setprecision(4);
+
+      test<2,1>();
+    }
+  else
+    {
+      test<2,1>();
+    }
+}
diff --git a/tests/mpi/parallel_block_vector_06.with_p4est=true.with_lapack=true.mpirun=3.output b/tests/mpi/parallel_block_vector_06.with_p4est=true.with_lapack=true.mpirun=3.output
new file mode 100644 (file)
index 0000000..3a27e7c
--- /dev/null
@@ -0,0 +1,6 @@
+
+DEAL:0::Symmetric:
+DEAL:0::Norm of difference: 0.000
+DEAL:0::Cholesky: Ok
+DEAL:0::Non-Symmetric:
+DEAL:0::Norm of difference: 0.000

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.