]> https://gitweb.dealii.org/ - dealii.git/commitdiff
consistently use local_size().
authorDavid Wells <drwells@email.unc.edu>
Sat, 23 May 2020 01:52:17 +0000 (21:52 -0400)
committerDaniel Arndt <arndtd@ornl.gov>
Thu, 11 Feb 2021 16:48:49 +0000 (11:48 -0500)
This is the name given to the function by la_parallel_vector. Defining this
everywhere gives us a true check as to whether or not a vector is distributed -
i.e., if size() != local_size() then the vector is distributed and otherwise it
is not.

include/deal.II/lac/block_vector_base.h
include/deal.II/lac/la_parallel_vector.h
include/deal.II/lac/read_write_vector.h
include/deal.II/lac/trilinos_epetra_vector.h
include/deal.II/lac/trilinos_tpetra_vector.h
include/deal.II/lac/trilinos_tpetra_vector.templates.h
include/deal.II/lac/vector.h
source/lac/trilinos_epetra_vector.cc
tests/mpi/local_size.cc [new file with mode: 0644]
tests/mpi/local_size.with_trilinos_with_tpetra=on.with_petsc=on.mpirun=4.output [new file with mode: 0644]

index 870a6ae6355075515a8fa28a87ca87ab06f94519..fc5dca1fae70e495300ce5b7ac70b41e7587b7b5 100644 (file)
@@ -544,6 +544,14 @@ public:
   std::size_t
   size() const;
 
+  /**
+   * Return local dimension of the vector. This is the sum of the local
+   * dimensions (i.e., values stored on the current processor) of all
+   * components.
+   */
+  std::size_t
+  local_size() const;
+
   /**
    * Return an index set that describes which elements of this vector are
    * owned by the current processor. Note that this index set does not include
@@ -1430,6 +1438,18 @@ BlockVectorBase<VectorType>::size() const
 
 
 
+template <class VectorType>
+inline std::size_t
+BlockVectorBase<VectorType>::local_size() const
+{
+  std::size_t local_size = 0;
+  for (unsigned int b = 0; b < n_blocks(); ++b)
+    local_size += block(b).local_size();
+  return local_size;
+}
+
+
+
 template <class VectorType>
 inline IndexSet
 BlockVectorBase<VectorType>::locally_owned_elements() const
index 209de68f5c72611e04d4fca61d10de490bf5ce7b..05e2ad31fa701bdbdfea2f7ffb2a12386fe575df 100644 (file)
@@ -1514,7 +1514,7 @@ namespace LinearAlgebra
     inline typename Vector<Number, MemorySpace>::size_type
     Vector<Number, MemorySpace>::local_size() const
     {
-      return partitioner->local_size();
+      return locally_owned_size();
     }
 
 
index bb6f47e98d41696d5549ae4aadf51a3baa20583e..1df77dd799144e06a2b538bb0a7201c8bcf84f79 100644 (file)
@@ -399,10 +399,21 @@ namespace LinearAlgebra
      * This function returns the number of elements stored. It is smaller or
      * equal to the dimension of the vector space that is modeled by an object
      * of this kind. This dimension is return by size().
+     *
+     * @deprecated use local_size() instead.
      */
+    DEAL_II_DEPRECATED
     size_type
     n_elements() const;
 
+    /**
+     * Return the local size of the vector, i.e., the number of indices
+     * owned locally.
+     */
+    size_type
+    local_size() const;
+
+
     /**
      * Return the IndexSet that represents the indices of the elements stored.
      */
@@ -821,6 +832,15 @@ namespace LinearAlgebra
 
 
 
+  template <typename Number>
+  inline typename ReadWriteVector<Number>::size_type
+  ReadWriteVector<Number>::local_size() const
+  {
+    return stored_elements.n_elements();
+  }
+
+
+
   template <typename Number>
   inline const IndexSet &
   ReadWriteVector<Number>::get_stored_elements() const
index 34ab495165da96f6e9859b17836566d8641b8be5..421c7865c2641cee92c0fd5d5488c0fe2f0c1559 100644 (file)
@@ -284,6 +284,13 @@ namespace LinearAlgebra
       virtual size_type
       size() const override;
 
+      /**
+       * Return the local size of the vector, i.e., the number of indices
+       * owned locally.
+       */
+      size_type
+      local_size() const;
+
       /**
        * Return the MPI communicator object in use with this object.
        */
index 95fb600629af9015ff6bc43d5c2cd92dd847b0f0..729a6d480648e04c3775d85704d592a9008c5f91 100644 (file)
@@ -307,6 +307,13 @@ namespace LinearAlgebra
       virtual size_type
       size() const override;
 
+      /**
+       * Return the local size of the vector, i.e., the number of indices
+       * owned locally.
+       */
+      size_type
+      local_size() const;
+
       /**
        * Return the MPI communicator object in use with this object.
        */
index d75d263d5a59e63033afd1c55bedb701bf4bded2..55b04638a2351abbb4cb40735d42e541e983c213 100644 (file)
@@ -548,6 +548,15 @@ namespace LinearAlgebra
 
 
 
+    template <typename Number>
+    typename Vector<Number>::size_type
+    Vector<Number>::local_size() const
+    {
+      return vector->getLocalLength();
+    }
+
+
+
     template <typename Number>
     MPI_Comm
     Vector<Number>::get_mpi_communicator() const
index 221ff506a0a15f4a71cb97bd511ff31006d92844..de32e87ddc0eb8c7f2ccd390d0930a1052e40dc0 100644 (file)
@@ -953,6 +953,16 @@ public:
   size_type
   size() const;
 
+  /**
+   * Return local dimension of the vector. Since this vector does not support
+   * distributed data this is always the same value as size().
+   *
+   * @note This function exists for compatibility with
+   * LinearAlgebra::ReadWriteVector.
+   */
+  size_type
+  local_size() const;
+
   /**
    * Return whether the vector contains only elements with value zero. This
    * function is mainly for internal consistency checks and should seldom be
@@ -1088,6 +1098,16 @@ Vector<Number>::size() const
 }
 
 
+
+template <typename Number>
+inline typename Vector<Number>::size_type
+Vector<Number>::local_size() const
+{
+  return values.size();
+}
+
+
+
 template <typename Number>
 inline bool
 Vector<Number>::in_local_range(const size_type) const
index 54b45876421a54b76735494bb2636219d60de0b6..e45f395f88735ff699553d1949787990acf7706f 100644 (file)
@@ -540,6 +540,14 @@ namespace LinearAlgebra
 
 
 
+    Vector::size_type
+    Vector::local_size() const
+    {
+      return vector->MyLength();
+    }
+
+
+
     MPI_Comm
     Vector::get_mpi_communicator() const
     {
diff --git a/tests/mpi/local_size.cc b/tests/mpi/local_size.cc
new file mode 100644 (file)
index 0000000..1fcb51b
--- /dev/null
@@ -0,0 +1,167 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 - 2020 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.
+//
+// ---------------------------------------------------------------------
+
+// check Vector::local_size() for all supported vector types
+
+#include <deal.II/base/index_set.h>
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/lac/la_parallel_block_vector.h>
+#include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/la_vector.h>
+#include <deal.II/lac/petsc_block_vector.h>
+#include <deal.II/lac/petsc_vector.h>
+#include <deal.II/lac/trilinos_epetra_vector.h>
+#include <deal.II/lac/trilinos_parallel_block_vector.h>
+#include <deal.II/lac/trilinos_tpetra_vector.h>
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/vector.h>
+
+#include "../tests.h"
+
+
+template <typename VEC>
+void
+check_serial()
+{
+  const auto dofs_per_proc = 4;
+  VEC vec(dofs_per_proc);
+  deallog << "type: " << Utilities::type_to_string(vec) << std::endl;
+  deallog << "local size: " << vec.local_size() << std::endl;
+  deallog << "size: " << vec.size() << std::endl;
+}
+
+
+
+template <typename VEC>
+void
+check_unghosted_parallel()
+{
+  const auto n_procs = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+  const auto my_rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+
+  const auto dofs_per_proc = 4;
+  const auto n_dofs = dofs_per_proc*n_procs;
+
+  IndexSet local_indices(n_dofs);
+  const auto my_dofs_begin = dofs_per_proc*my_rank;
+  const auto my_dofs_end = dofs_per_proc*(my_rank + 1);
+  local_indices.add_range(my_dofs_begin, my_dofs_end);
+  local_indices.compress();
+
+  VEC vec(local_indices, MPI_COMM_WORLD);
+  deallog << "type: " << Utilities::type_to_string(vec) << std::endl;
+  deallog << "index set size: " << local_indices.n_elements() << std::endl;
+  deallog << "local size: " << vec.local_size() << std::endl;
+  deallog << "size: " << vec.size() << std::endl;
+}
+
+
+
+template <typename VEC>
+void
+check_ghosted_parallel()
+{
+  const auto n_procs = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+  const auto my_rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+
+  const auto dofs_per_proc = 4;
+  const auto n_dofs = dofs_per_proc*n_procs;
+
+  IndexSet local_indices(n_dofs);
+  IndexSet ghost_indices(n_dofs);
+  const auto my_dofs_begin = dofs_per_proc*my_rank;
+  const auto my_dofs_end = dofs_per_proc*(my_rank + 1);
+  local_indices.add_range(my_dofs_begin, my_dofs_end);
+  local_indices.compress();
+  if (my_rank == 0)
+    {
+      ghost_indices.add_index(my_dofs_end);
+    }
+  else
+    {
+      ghost_indices.add_index(my_dofs_begin - 1);
+      ghost_indices.add_index(my_dofs_end % n_dofs);
+    }
+
+  VEC vec(local_indices, ghost_indices, MPI_COMM_WORLD);
+  deallog << "type: " << Utilities::type_to_string(vec) << std::endl;
+  deallog << "index set size: " << local_indices.n_elements() << std::endl;
+  deallog << "local size: " << vec.local_size() << std::endl;
+  deallog << "size: " << vec.size() << std::endl;
+}
+
+
+
+template <typename VEC>
+void
+check_ghosted_parallel_block()
+{
+  const auto n_procs = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+  const auto my_rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+
+  const auto dofs_per_proc = 4;
+  const auto n_dofs = dofs_per_proc*n_procs;
+
+  IndexSet local_indices(n_dofs);
+  IndexSet ghost_indices(n_dofs);
+  const auto my_dofs_begin = dofs_per_proc*my_rank;
+  const auto my_dofs_end = dofs_per_proc*(my_rank + 1);
+  local_indices.add_range(my_dofs_begin, my_dofs_end);
+  local_indices.compress();
+  if (my_rank == 0)
+    {
+      ghost_indices.add_index(my_dofs_end);
+    }
+  else
+    {
+      ghost_indices.add_index(my_dofs_begin - 1);
+      ghost_indices.add_index(my_dofs_end % n_dofs);
+    }
+
+  std::vector<IndexSet> local_blocks {local_indices, local_indices};
+  // for variety do not ghost the second component
+  std::vector<IndexSet> ghost_blocks {ghost_indices, IndexSet()};
+
+  VEC vec(local_blocks, ghost_blocks, MPI_COMM_WORLD);
+  deallog << "type: " << Utilities::type_to_string(vec) << std::endl;
+  deallog << "local size: " << vec.local_size() << std::endl;
+  deallog << "size: " << vec.size() << std::endl;
+}
+
+
+
+int
+main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  MPILogInitAll                    log;
+
+  // non-block vectors:
+  check_serial<Vector<double>>();
+  check_serial<LinearAlgebra::Vector<double>>();
+
+  check_unghosted_parallel<LinearAlgebra::EpetraWrappers::Vector>();
+  check_unghosted_parallel<LinearAlgebra::TpetraWrappers::Vector<double>>();
+
+  check_ghosted_parallel<LinearAlgebra::distributed::Vector<double>>();
+  check_ghosted_parallel<PETScWrappers::MPI::Vector>();
+  check_ghosted_parallel<TrilinosWrappers::MPI::Vector>();
+
+  // block vectors:
+  check_ghosted_parallel_block<LinearAlgebra::distributed::BlockVector<double>>();
+  check_ghosted_parallel_block<PETScWrappers::MPI::BlockVector>();
+  check_ghosted_parallel_block<TrilinosWrappers::MPI::BlockVector>();
+}
diff --git a/tests/mpi/local_size.with_trilinos_with_tpetra=on.with_petsc=on.mpirun=4.output b/tests/mpi/local_size.with_trilinos_with_tpetra=on.with_petsc=on.mpirun=4.output
new file mode 100644 (file)
index 0000000..622b61f
--- /dev/null
@@ -0,0 +1,147 @@
+
+DEAL:0::type: dealii::Vector<double>
+DEAL:0::local size: 4
+DEAL:0::size: 4
+DEAL:0::type: dealii::LinearAlgebra::Vector<double>
+DEAL:0::local size: 4
+DEAL:0::size: 4
+DEAL:0::type: dealii::LinearAlgebra::EpetraWrappers::Vector
+DEAL:0::index set size: 4
+DEAL:0::local size: 4
+DEAL:0::size: 16
+DEAL:0::type: dealii::LinearAlgebra::TpetraWrappers::Vector<double>
+DEAL:0::index set size: 4
+DEAL:0::local size: 4
+DEAL:0::size: 16
+DEAL:0::type: dealii::LinearAlgebra::distributed::Vector<double, dealii::MemorySpace::Host>
+DEAL:0::index set size: 4
+DEAL:0::local size: 4
+DEAL:0::size: 16
+DEAL:0::type: dealii::PETScWrappers::MPI::Vector
+DEAL:0::index set size: 4
+DEAL:0::local size: 4
+DEAL:0::size: 16
+DEAL:0::type: dealii::TrilinosWrappers::MPI::Vector
+DEAL:0::index set size: 4
+DEAL:0::local size: 5
+DEAL:0::size: 16
+DEAL:0::type: dealii::LinearAlgebra::distributed::BlockVector<double>
+DEAL:0::local size: 8
+DEAL:0::size: 32
+DEAL:0::type: dealii::PETScWrappers::MPI::BlockVector
+DEAL:0::local size: 8
+DEAL:0::size: 32
+DEAL:0::type: dealii::TrilinosWrappers::MPI::BlockVector
+DEAL:0::local size: 9
+DEAL:0::size: 32
+
+DEAL:1::type: dealii::Vector<double>
+DEAL:1::local size: 4
+DEAL:1::size: 4
+DEAL:1::type: dealii::LinearAlgebra::Vector<double>
+DEAL:1::local size: 4
+DEAL:1::size: 4
+DEAL:1::type: dealii::LinearAlgebra::EpetraWrappers::Vector
+DEAL:1::index set size: 4
+DEAL:1::local size: 4
+DEAL:1::size: 16
+DEAL:1::type: dealii::LinearAlgebra::TpetraWrappers::Vector<double>
+DEAL:1::index set size: 4
+DEAL:1::local size: 4
+DEAL:1::size: 16
+DEAL:1::type: dealii::LinearAlgebra::distributed::Vector<double, dealii::MemorySpace::Host>
+DEAL:1::index set size: 4
+DEAL:1::local size: 4
+DEAL:1::size: 16
+DEAL:1::type: dealii::PETScWrappers::MPI::Vector
+DEAL:1::index set size: 4
+DEAL:1::local size: 4
+DEAL:1::size: 16
+DEAL:1::type: dealii::TrilinosWrappers::MPI::Vector
+DEAL:1::index set size: 4
+DEAL:1::local size: 6
+DEAL:1::size: 16
+DEAL:1::type: dealii::LinearAlgebra::distributed::BlockVector<double>
+DEAL:1::local size: 8
+DEAL:1::size: 32
+DEAL:1::type: dealii::PETScWrappers::MPI::BlockVector
+DEAL:1::local size: 8
+DEAL:1::size: 32
+DEAL:1::type: dealii::TrilinosWrappers::MPI::BlockVector
+DEAL:1::local size: 10
+DEAL:1::size: 32
+
+
+DEAL:2::type: dealii::Vector<double>
+DEAL:2::local size: 4
+DEAL:2::size: 4
+DEAL:2::type: dealii::LinearAlgebra::Vector<double>
+DEAL:2::local size: 4
+DEAL:2::size: 4
+DEAL:2::type: dealii::LinearAlgebra::EpetraWrappers::Vector
+DEAL:2::index set size: 4
+DEAL:2::local size: 4
+DEAL:2::size: 16
+DEAL:2::type: dealii::LinearAlgebra::TpetraWrappers::Vector<double>
+DEAL:2::index set size: 4
+DEAL:2::local size: 4
+DEAL:2::size: 16
+DEAL:2::type: dealii::LinearAlgebra::distributed::Vector<double, dealii::MemorySpace::Host>
+DEAL:2::index set size: 4
+DEAL:2::local size: 4
+DEAL:2::size: 16
+DEAL:2::type: dealii::PETScWrappers::MPI::Vector
+DEAL:2::index set size: 4
+DEAL:2::local size: 4
+DEAL:2::size: 16
+DEAL:2::type: dealii::TrilinosWrappers::MPI::Vector
+DEAL:2::index set size: 4
+DEAL:2::local size: 6
+DEAL:2::size: 16
+DEAL:2::type: dealii::LinearAlgebra::distributed::BlockVector<double>
+DEAL:2::local size: 8
+DEAL:2::size: 32
+DEAL:2::type: dealii::PETScWrappers::MPI::BlockVector
+DEAL:2::local size: 8
+DEAL:2::size: 32
+DEAL:2::type: dealii::TrilinosWrappers::MPI::BlockVector
+DEAL:2::local size: 10
+DEAL:2::size: 32
+
+
+DEAL:3::type: dealii::Vector<double>
+DEAL:3::local size: 4
+DEAL:3::size: 4
+DEAL:3::type: dealii::LinearAlgebra::Vector<double>
+DEAL:3::local size: 4
+DEAL:3::size: 4
+DEAL:3::type: dealii::LinearAlgebra::EpetraWrappers::Vector
+DEAL:3::index set size: 4
+DEAL:3::local size: 4
+DEAL:3::size: 16
+DEAL:3::type: dealii::LinearAlgebra::TpetraWrappers::Vector<double>
+DEAL:3::index set size: 4
+DEAL:3::local size: 4
+DEAL:3::size: 16
+DEAL:3::type: dealii::LinearAlgebra::distributed::Vector<double, dealii::MemorySpace::Host>
+DEAL:3::index set size: 4
+DEAL:3::local size: 4
+DEAL:3::size: 16
+DEAL:3::type: dealii::PETScWrappers::MPI::Vector
+DEAL:3::index set size: 4
+DEAL:3::local size: 4
+DEAL:3::size: 16
+DEAL:3::type: dealii::TrilinosWrappers::MPI::Vector
+DEAL:3::index set size: 4
+DEAL:3::local size: 6
+DEAL:3::size: 16
+DEAL:3::type: dealii::LinearAlgebra::distributed::BlockVector<double>
+DEAL:3::local size: 8
+DEAL:3::size: 32
+DEAL:3::type: dealii::PETScWrappers::MPI::BlockVector
+DEAL:3::local size: 8
+DEAL:3::size: 32
+DEAL:3::type: dealii::TrilinosWrappers::MPI::BlockVector
+DEAL:3::local size: 10
+DEAL:3::size: 32
+

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.