]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add more parameters to LA::distributed::BlockVector::mmult() 5937/head
authorDenis Davydov <davydden@gmail.com>
Wed, 21 Feb 2018 13:09:11 +0000 (14:09 +0100)
committerDenis Davydov <davydden@gmail.com>
Fri, 23 Feb 2018 17:56:13 +0000 (18:56 +0100)
change meaning of arguments in BlockVector::mmult()

doc/news/changes/minor/20180210DenisDavydov
include/deal.II/lac/la_parallel_block_vector.h
include/deal.II/lac/la_parallel_block_vector.templates.h
source/lac/la_parallel_block_vector.inst.in
tests/mpi/parallel_block_vector_08.cc
tests/mpi/parallel_block_vector_09.cc
tests/mpi/parallel_block_vector_12.cc [new file with mode: 0644]
tests/mpi/parallel_block_vector_12.with_p4est=true.mpirun=4.output [new file with mode: 0644]

index ddda9265bbf59f5d0b9e2460b93925f7f5366ba4..27f476c8eac2634a5a1b699e499ebf6be05efe9c 100644 (file)
@@ -1,4 +1,4 @@
-New: add parallel::distributed::BlockVector::mmult(const BlockVector &, const FullMatrixType &)
+New: add parallel::distributed::BlockVector::mmult(BlockVector &, const FullMatrixType &, const NumberType, const NumberType) const
 and parallel::distributed::BlockVector::multivector_inner_product_with_metric(const FullMatrixType &, const BlockVector &V, const bool) const
 to operate on multivectors with a metric tensor.
 <br>
index 31d644831ac8466091e82eca1fcf613c08cec9a2..9d0c60b3ff88fdd05e45f67a081db30f7d7c337e 100644 (file)
@@ -491,16 +491,18 @@ namespace LinearAlgebra
 
       /**
        * Set each block of this vector as follows:
-       * $U^i = \sum_{j} V_j A^{ji}$ where $U^i$
-       * and $V_j$ indicate the $i$th block (not element) of $U$ and the
-       * $j$th block of $V$, respectively.
+       * $V^i = s V^i + b \sum_{j} U_j A^{ji}$ where $V^i$
+       * and $U_j$ indicate the $i$th block (not element) of $V$ and the
+       * $j$th block of $U$, respectively.
        *
        * Obviously, this function can only be used if all blocks of both vectors
        * are of the same size.
        */
       template <typename FullMatrixType>
-      void mmult(const BlockVector<Number> &V,
-                 const FullMatrixType &matrix);
+      void mmult(BlockVector<Number> &V,
+                 const FullMatrixType &matrix,
+                 const Number s = Number(0.),
+                 const Number b = Number(1.)) const;
 
       /**
        * Add @p a to all components. Note that @p a is a scalar not a vector.
index 0090e56c5c840a17b5f98635e3e4af96fc644ece..df6ad24706866e190ed4c22ce8c1b91877bc123e 100644 (file)
@@ -835,9 +835,9 @@ namespace LinearAlgebra
         return;
 
       Assert (matrix.m() == m,
-              dealii::ExcDimensionMismatch(matrix.m(),m));
+              ExcDimensionMismatch(matrix.m(),m));
       Assert (matrix.n() == n,
-              dealii::ExcDimensionMismatch(matrix.n(),n));
+              ExcDimensionMismatch(matrix.n(),n));
 
       // reset the matrix
       matrix = typename FullMatrixType::value_type(0.0);
@@ -846,7 +846,7 @@ namespace LinearAlgebra
       if (symmetric)
         {
           Assert (m == n,
-                  dealii::ExcDimensionMismatch(m,n));
+                  ExcDimensionMismatch(m,n));
 
           for (unsigned int i = 0; i < m; i++)
             for (unsigned int j = i; j < n; j++)
@@ -888,14 +888,14 @@ namespace LinearAlgebra
         return res;
 
       Assert (matrix.m() == m,
-              dealii::ExcDimensionMismatch(matrix.m(),m));
+              ExcDimensionMismatch(matrix.m(),m));
       Assert (matrix.n() == n,
-              dealii::ExcDimensionMismatch(matrix.n(),n));
+              ExcDimensionMismatch(matrix.n(),n));
 
       if (symmetric)
         {
           Assert (m == n,
-                  dealii::ExcDimensionMismatch(m,n));
+                  ExcDimensionMismatch(m,n));
 
           for (unsigned int i = 0; i < m; i++)
             {
@@ -919,11 +919,13 @@ namespace LinearAlgebra
     template <typename Number>
     template <typename FullMatrixType>
     void
-    BlockVector<Number>::mmult(const BlockVector<Number> &V,
-                               const FullMatrixType &matrix)
+    BlockVector<Number>::mmult(BlockVector<Number> &V,
+                               const FullMatrixType &matrix,
+                               const Number s,
+                               const Number b) const
     {
-      const unsigned int n = this->n_blocks();
-      const unsigned int m = V.n_blocks();
+      const unsigned int m = this->n_blocks();
+      const unsigned int n = V.n_blocks();
 
       // in case one vector is empty and the second one is not, the
       // FullMatrix resized to (m,n) will have 0 both in m() and n()
@@ -933,14 +935,21 @@ namespace LinearAlgebra
         return;
 
       Assert (matrix.m() == m,
-              dealii::ExcDimensionMismatch(matrix.m(),m));
+              ExcDimensionMismatch(matrix.m(),m));
       Assert (matrix.n() == n,
-              dealii::ExcDimensionMismatch(matrix.n(),n));
+              ExcDimensionMismatch(matrix.n(),n));
 
-      (*this) = Number();
       for (unsigned int i = 0; i < n; i++)
-        for (unsigned int j = 0; j < m; j++)
-          this->block(i).add (matrix(j,i), V.block(j));
+        {
+          // below we make this work gracefully for identity-like matrices in
+          // which case the two loops over j won't do any work as A(j,i)==0
+          const unsigned int k = std::min(i,m-1);
+          V.block(i).sadd(s,matrix(k,i)*b, this->block(k));
+          for (unsigned int j = 0  ; j < k; j++)
+            V.block(i).add (matrix(j,i)*b, this->block(j));
+          for (unsigned int j = k+1; j < m; j++)
+            V.block(i).add (matrix(j,i)*b, this->block(j));
+        }
     }
 
 
index 657e150417dcdf42058892aef91a1915d0fa5041..b87cb3e6044723a27053bc5486c2dfc42e2a58da 100644 (file)
@@ -36,10 +36,10 @@ for (SCALAR : REAL_AND_COMPLEX_SCALARS)
     BlockVector<SCALAR>::multivector_inner_product_with_metric(const LAPACKFullMatrix<SCALAR> &, const BlockVector<SCALAR> &V, const bool) const;
     template
     void
-    BlockVector<SCALAR>::mmult(const BlockVector<SCALAR> &V,const FullMatrix<SCALAR> &);
+    BlockVector<SCALAR>::mmult(BlockVector<SCALAR> &V,const FullMatrix<SCALAR> &, const SCALAR, const SCALAR) const;
     template
     void
-    BlockVector<SCALAR>::mmult(const BlockVector<SCALAR> &V,const LAPACKFullMatrix<SCALAR> &);
+    BlockVector<SCALAR>::mmult(BlockVector<SCALAR> &V,const LAPACKFullMatrix<SCALAR> &, const SCALAR, const SCALAR) const;
     \}
     \}
 }
index bccb635727d4a3d93dc7e2ce2c6fa6cf23f38da3..8c2d4d3d69e5808f2c5a1012eadb144306130c8c 100644 (file)
@@ -129,7 +129,7 @@ void test (const unsigned int n_blocks = 5)
     for (unsigned int j = 0; j < n_blocks; ++j)
       metric(i,j) = 0.3 + (i*3 + j*7);
 
-  left.mmult(right,metric);
+  right.mmult(left,metric);
 
   for (unsigned int i = 0; i < n_blocks; ++i)
     for (unsigned int j = 0; j < n_blocks; ++j)
index a66bc3780e959a3f555fcd5fec58532395a7f37a..f4bdf5e6838b18e1ab7193eed27432f270517ddd 100644 (file)
@@ -139,7 +139,7 @@ void test (const unsigned int n_blocks = 5)
       }
 
   const double res = left.multivector_inner_product_with_metric(metric, right);
-  left2.mmult(left, metric);
+  left.mmult(left2, metric);
   const double res2 = left2 * right;
 
   const double diff_norm = std::abs(res - res2);
diff --git a/tests/mpi/parallel_block_vector_12.cc b/tests/mpi/parallel_block_vector_12.cc
new file mode 100644 (file)
index 0000000..4518331
--- /dev/null
@@ -0,0 +1,181 @@
+// ---------------------------------------------------------------------
+//
+// 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>::mmult(const BlockVector<Number> &V,const FullMatrixType &matrix, s,b)
+// also for vectors of different number of blocks.
+// Triangulation and Mass operator are the same as in matrix_free/mass_operator_01.cc
+
+#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/manifold_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 <iostream>
+
+
+
+
+template <int dim, int fe_degree>
+void test (const unsigned int n = 5, const unsigned int m = 3)
+{
+  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), right(m), left2(n);
+  for (unsigned int b = 0; b < n; ++b)
+    {
+      mf_data->initialize_dof_vector (left.block(b));
+      mf_data->initialize_dof_vector (left2.block(b));
+      left.block(b) = 0.;
+      left2.block(b) = 0.;
+      for (unsigned int i=0; i<left.block(b).local_size(); ++i)
+        {
+          const unsigned int glob_index =
+            owned_set.nth_index_in_set (i);
+          if (constraints.is_constrained(glob_index))
+            continue;
+          left.block(b).local_element(i) = random_value<double>();
+        }
+    }
+
+  for (unsigned int b = 0; b < m; ++b)
+    {
+      mf_data->initialize_dof_vector (right.block(b));
+      right.block(b) = 0.;
+      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) = random_value<double>();
+        }
+    }
+
+
+  FullMatrix<number> metric(m,n);
+  for (unsigned int i = 0; i < m; ++i)
+    for (unsigned int j = 0; j < n; ++j)
+      metric(i,j) = 0.3 + (i*3 + j*7);
+
+  right.mmult(left,metric);       // L = RM
+  right.mmult(left,metric,2.,3.); // L = 2L + 3 RM = 5RM
+
+  for (unsigned int i = 0; i < n; ++i)
+    for (unsigned int j = 0; j < m; ++j)
+      left2.block(i).add(5.*metric(j,i), right.block(j));
+
+  left2.add(-1., left);
+
+  const double diff_norm = left2.linfty_norm();
+  deallog << "Norm of difference: " << diff_norm << 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>(5,3);
+      test<2,1>(3,3);
+      test<2,1>(3,5);
+    }
+  else
+    {
+      test<2,1>(5,3);
+      test<2,1>(3,3);
+      test<2,1>(3,5);
+    }
+}
diff --git a/tests/mpi/parallel_block_vector_12.with_p4est=true.mpirun=4.output b/tests/mpi/parallel_block_vector_12.with_p4est=true.mpirun=4.output
new file mode 100644 (file)
index 0000000..48c41a9
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL:0::Norm of difference: 0.000
+DEAL:0::Norm of difference: 0.000
+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.