]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added get_mpi_communicator() to Trilinos Matrix.
authorLuca Heltai <luca.heltai@sissa.it>
Fri, 24 Jul 2015 18:23:19 +0000 (20:23 +0200)
committerTimo Heister <timo.heister@gmail.com>
Sun, 26 Jul 2015 21:42:23 +0000 (17:42 -0400)
include/deal.II/base/index_set.h
include/deal.II/lac/trilinos_sparse_matrix.h
include/deal.II/lac/trilinos_vector.h
source/lac/trilinos_sparse_matrix.cc

index 3549ce1d566cbca0d48b04f7e692a9722aacdfe4..a80e6b5ef860c5665ba50b27f0ea633768f27a77 100644 (file)
@@ -551,9 +551,13 @@ IndexSet::IndexSet (const Epetra_Map &map)
   index_space_size (map.NumGlobalElements()),
   largest_range (numbers::invalid_unsigned_int)
 {
-  int *map_indices = map.MyGlobalElements();
-  for (int i=0; i<map.NumMyElements(); ++i)
-    add_index((size_type)map_indices[i]);
+  const size_type n_indices = map.NumMyElements();
+#ifndef DEAL_II_WITH_64BIT_INDICES
+  unsigned int *indices = (unsigned int *)map.MyGlobalElements();
+#else
+  size_type *indices = (size_type *)map.MyGlobalElements64();
+#endif
+  add_indices(indices, indices+n_indices);
   compress();
 }
 #endif
index c12d11debca0c6f23fd843456328e730bb63bd48..9d7a7a99503f684d073f60cb9966e900dc62ebd8 100644 (file)
@@ -662,7 +662,7 @@ namespace TrilinosWrappers
      * the compress() step).
      */
     SparseMatrix (const Epetra_Map  &parallel_partitioning,
-                  const size_type    n_max_entries_per_row = 0);
+                  const size_type    n_max_entries_per_row = 0) DEAL_II_DEPRECATED;
 
     /**
      * Same as before, but now set a value of nonzeros for each matrix row.
@@ -672,7 +672,7 @@ namespace TrilinosWrappers
      * respective SparseMatrix::reinit call considerably faster.
      */
     SparseMatrix (const Epetra_Map                &parallel_partitioning,
-                  const std::vector<unsigned int> &n_entries_per_row);
+                  const std::vector<unsigned int> &n_entries_per_row) DEAL_II_DEPRECATED;
 
     /**
      * This constructor is similar to the one above, but it now takes two
@@ -692,7 +692,7 @@ namespace TrilinosWrappers
      */
     SparseMatrix (const Epetra_Map &row_parallel_partitioning,
                   const Epetra_Map &col_parallel_partitioning,
-                  const size_type   n_max_entries_per_row = 0);
+                  const size_type   n_max_entries_per_row = 0) DEAL_II_DEPRECATED;
 
     /**
      * This constructor is similar to the one above, but it now takes two
@@ -710,7 +710,7 @@ namespace TrilinosWrappers
      */
     SparseMatrix (const Epetra_Map                &row_parallel_partitioning,
                   const Epetra_Map                &col_parallel_partitioning,
-                  const std::vector<unsigned int> &n_entries_per_row);
+                  const std::vector<unsigned int> &n_entries_per_row) DEAL_II_DEPRECATED;
 
     /**
      * This function is initializes the Trilinos Epetra matrix according to
@@ -739,7 +739,7 @@ namespace TrilinosWrappers
     template<typename SparsityType>
     void reinit (const Epetra_Map    &parallel_partitioning,
                  const SparsityType  &sparsity_pattern,
-                 const bool          exchange_data = false);
+                 const bool          exchange_data = false) DEAL_II_DEPRECATED;
 
     /**
      * This function is similar to the other initialization function above,
@@ -757,7 +757,7 @@ namespace TrilinosWrappers
     void reinit (const Epetra_Map    &row_parallel_partitioning,
                  const Epetra_Map    &col_parallel_partitioning,
                  const SparsityType  &sparsity_pattern,
-                 const bool          exchange_data = false);
+                 const bool          exchange_data = false) DEAL_II_DEPRECATED;
 
     /**
      * This function initializes the Trilinos matrix using the deal.II sparse
@@ -780,7 +780,7 @@ namespace TrilinosWrappers
                  const ::dealii::SparseMatrix<number> &dealii_sparse_matrix,
                  const double                          drop_tolerance=1e-13,
                  const bool                            copy_values=true,
-                 const ::dealii::SparsityPattern      *use_this_sparsity=0);
+                 const ::dealii::SparsityPattern      *use_this_sparsity=0) DEAL_II_DEPRECATED;
 
     /**
      * This function is similar to the other initialization function with
@@ -801,7 +801,7 @@ namespace TrilinosWrappers
                  const ::dealii::SparseMatrix<number>  &dealii_sparse_matrix,
                  const double                           drop_tolerance=1e-13,
                  const bool                             copy_values=true,
-                 const ::dealii::SparsityPattern      *use_this_sparsity=0);
+                 const ::dealii::SparsityPattern      *use_this_sparsity=0) DEAL_II_DEPRECATED;
 //@}
     /**
      * @name Constructors and initialization using an IndexSet description
@@ -1030,6 +1030,11 @@ namespace TrilinosWrappers
      */
     size_type memory_consumption () const;
 
+    /**
+     * Return the MPI communicator object in use with this matrix.
+     */
+    MPI_Comm get_mpi_communicator () const;
+
 //@}
     /**
      * @name Modifying entries
@@ -1664,7 +1669,7 @@ namespace TrilinosWrappers
      * sets the partitioning of the domain space of this matrix, i.e., the
      * partitioning of the vectors this matrix has to be multiplied with.
      */
-    const Epetra_Map &domain_partitioner () const;
+    const Epetra_Map &domain_partitioner ()  const DEAL_II_DEPRECATED;
 
     /**
      * Return a const reference to the underlying Trilinos Epetra_Map that
@@ -1672,14 +1677,14 @@ namespace TrilinosWrappers
      * partitioning of the vectors that are result from matrix-vector
      * products.
      */
-    const Epetra_Map &range_partitioner () const;
+    const Epetra_Map &range_partitioner () const DEAL_II_DEPRECATED;
 
     /**
      * Return a const reference to the underlying Trilinos Epetra_Map that
      * sets the partitioning of the matrix rows. Equal to the partitioning of
      * the range.
      */
-    const Epetra_Map &row_partitioner () const;
+    const Epetra_Map &row_partitioner () const DEAL_II_DEPRECATED;
 
     /**
      * Return a const reference to the underlying Trilinos Epetra_Map that
@@ -1687,8 +1692,29 @@ namespace TrilinosWrappers
      * equal to the partitioner Epetra_Map for the domain because of overlap
      * in the matrix.
      */
-    const Epetra_Map &col_partitioner () const;
+    const Epetra_Map &col_partitioner () const DEAL_II_DEPRECATED;
 //@}
+
+    /**
+     * @name Partitioners
+     */
+//@{
+
+    /**
+     * Return the partitioning of the domain space of this matrix, i.e., the
+     * partitioning of the vectors this matrix has to be multiplied with.
+     */
+    IndexSet locally_owned_domain_indices() const;
+
+    /**
+     * Return the partitioning of the range space of this matrix, i.e., the
+     * partitioning of the vectors that are result from matrix-vector
+     * products.
+     */
+    IndexSet locally_owned_range_indices() const;
+
+//@}
+
     /**
      * @name Iterators
      */
@@ -2617,6 +2643,24 @@ namespace TrilinosWrappers
 
 
 
+  inline
+  IndexSet
+  SparseMatrix::locally_owned_domain_indices () const
+  {
+    return IndexSet(matrix->DomainMap());
+  }
+
+
+
+  inline
+  IndexSet
+  SparseMatrix::locally_owned_range_indices () const
+  {
+    return IndexSet(matrix->RangeMap());
+  }
+
+
+
   inline
   const Epetra_Map &
   SparseMatrix::row_partitioner () const
index a34869895f16853635edfd38f8a8527055bad256..afaa0940247c1642d9accc9e165595f6f24299f5 100644 (file)
@@ -988,7 +988,7 @@ namespace internal
                                 TrilinosWrappers::MPI::Vector &v,
                                 bool fast)
       {
-        v.reinit(IndexSet(matrix.range_partitioner()), MPI_COMM_WORLD, fast);
+        v.reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator(), fast);
       }
 
       template <typename Matrix>
@@ -997,7 +997,7 @@ namespace internal
                                 TrilinosWrappers::MPI::Vector &v,
                                 bool fast)
       {
-        v.reinit(IndexSet(matrix.domain_partitioner()), MPI_COMM_WORLD, fast);
+        v.reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator(), fast);
       }
     };
 
@@ -1015,7 +1015,7 @@ namespace internal
                                 TrilinosWrappers::Vector &v,
                                 bool fast)
       {
-        v.reinit(IndexSet(matrix.range_partitioner()), MPI_COMM_WORLD, fast);
+        v.reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator(), fast);
       }
 
       template <typename Matrix>
@@ -1024,7 +1024,7 @@ namespace internal
                                 TrilinosWrappers::Vector &v,
                                 bool fast)
       {
-        v.reinit(IndexSet(matrix.domain_partitioner()), MPI_COMM_WORLD, fast);
+        v.reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator(), fast);
       }
     };
 
index 7687b1a9e1fe53401497df902e3ba7ac60d0eaf1..66bc09313b0b771a230beee0f1df491ae21823ae 100644 (file)
@@ -2310,6 +2310,22 @@ namespace TrilinosWrappers
     return ((sizeof(TrilinosScalar)+sizeof(TrilinosWrappers::types::int_type))*
             matrix->NumMyNonzeros() + sizeof(int)*local_size() + static_memory);
   }
+
+  MPI_Comm SparseMatrix::get_mpi_communicator () const
+  {
+
+#ifdef DEAL_II_WITH_MPI
+
+    const Epetra_MpiComm *mpi_comm
+      = dynamic_cast<const Epetra_MpiComm *>(&range_partitioner().Comm());
+    return mpi_comm->Comm();
+#else
+
+    return MPI_COMM_SELF;
+
+#endif
+
+  }
 }
 
 

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.