]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use `uint64` for `n_nonzero_elements`. 13599/head
authorMarc Fehling <mafehling.git@gmail.com>
Tue, 5 Apr 2022 02:26:49 +0000 (20:26 -0600)
committerMarc Fehling <mafehling.git@gmail.com>
Thu, 7 Apr 2022 02:32:25 +0000 (20:32 -0600)
done: PETSc MatrixBase, Trilinos SparseMatrix, PETSc&Trilinos BlockSparseMatrix
todo: dealii SparseMatrix, all SparsityPattern objects

Should be backwards compatible with implicit type conversion.

12 files changed:
doc/news/changes/incompatibilities/20220404Fehling [new file with mode: 0644]
include/deal.II/base/types.h
include/deal.II/lac/petsc_block_sparse_matrix.h
include/deal.II/lac/petsc_matrix_base.h
include/deal.II/lac/trilinos_block_sparse_matrix.h
include/deal.II/lac/trilinos_index_access.h
include/deal.II/lac/trilinos_sparse_matrix.h
include/deal.II/lac/trilinos_sparsity_pattern.h
source/lac/petsc_matrix_base.cc
source/lac/petsc_parallel_block_sparse_matrix.cc
source/lac/trilinos_block_sparse_matrix.cc
source/lac/trilinos_sparsity_pattern.cc

diff --git a/doc/news/changes/incompatibilities/20220404Fehling b/doc/news/changes/incompatibilities/20220404Fehling
new file mode 100644 (file)
index 0000000..3dcfb29
--- /dev/null
@@ -0,0 +1,5 @@
+Changed: All functions `n_nonzero_elements()`, which are members of the
+SparseMatrix and SparsityPattern family of classes, now always return
+`std::uint64_t` instead of types::global_dof_index.
+<br>
+(Marc Fehling, 2022/04/04)
index ac262ad545b759e6cbce1a0654e3f0303398036f..57d5a1a42d9efa2abc869a4a1459c180eba3b682 100644 (file)
@@ -167,11 +167,16 @@ namespace TrilinosWrappers
 {
   namespace types
   {
+    /**
+     * Declare type of 64 bit integer used in the Epetra package of Trilinos.
+     */
+    using int64_type = long long int;
+
 #ifdef DEAL_II_WITH_64BIT_INDICES
     /**
      * Declare type of integer used in the Epetra package of Trilinos.
      */
-    using int_type = long long int;
+    using int_type = int64_type;
 #else
     /**
      * Declare type of integer used in the Epetra package of Trilinos.
index 60ac2c3e93f5925f2424ca34a09b66a48372e63c..87b2499f42a0e7132699c6f3018c1fbc1f9ac8a2 100644 (file)
@@ -255,7 +255,7 @@ namespace PETScWrappers
        * returns the number of entries in the sparsity pattern; if any of the
        * entries should happen to be zero, it is counted anyway.
        */
-      size_type
+      std::uint64_t
       n_nonzero_elements() const;
 
       /**
index f0e94b0531ae9b271be642bd4c6997e8455c27ab..ac6c5473742a0dfcf8cd38b26826bf75522e7a33 100644 (file)
@@ -653,7 +653,7 @@ namespace PETScWrappers
      * returns the number of entries in the sparsity pattern; if any of the
      * entries should happen to be zero, it is counted anyway.
      */
-    size_type
+    std::uint64_t
     n_nonzero_elements() const;
 
     /**
index 4a165f5f6936b00486335b2a37375ac4ef2e8796..440fe662eea02c5292456882cda1c1c1cadec190 100644 (file)
@@ -218,7 +218,7 @@ namespace TrilinosWrappers
      * Return the total number of nonzero elements of this matrix (summed
      * over all MPI processes).
      */
-    size_type
+    std::uint64_t
     n_nonzero_elements() const;
 
     /**
index 73396a86bc0cefe583f07570beab94c346c51014..c6467298b5f03153595e2c1ec8216c72dfc4721a 100644 (file)
@@ -32,18 +32,14 @@ DEAL_II_NAMESPACE_OPEN
 namespace TrilinosWrappers
 {
   /**
-   * A helper function that queries the size of an Epetra_BlockMap object
-   * and calls either the 32 or 64 bit function to get the number of global
-   * elements in the map.
+   * A helper function that queries the size of an Epetra_BlockMap object. We
+   * always call the 64 bit function to get the number of global elements in the
+   * map.
    */
-  inline TrilinosWrappers::types::int_type
+  inline TrilinosWrappers::types::int64_type
   n_global_elements(const Epetra_BlockMap &map)
   {
-#  ifdef DEAL_II_WITH_64BIT_INDICES
     return map.NumGlobalElements64();
-#  else
-    return map.NumGlobalElements();
-#  endif
   }
 
   /**
@@ -133,17 +129,13 @@ namespace TrilinosWrappers
   }
 
   /**
-   * A helper function that finds the number of global entries by calling
-   * either the 32 or 64 bit function.
+   * A helper function that finds the number of global entries. We always call
+   * the 64 bit function.
    */
-  inline TrilinosWrappers::types::int_type
+  inline TrilinosWrappers::types::int64_type
   n_global_entries(const Epetra_CrsGraph &graph)
   {
-#  ifdef DEAL_II_WITH_64BIT_INDICES
     return graph.NumGlobalEntries64();
-#  else
-    return graph.NumGlobalEntries();
-#  endif
   }
 
   /**
index 8eec00f32ddd6db6a2d6183a8a7a69e75380db9c..88d6f536a6e479c32630d10b8d0f1d2650751fa4 100644 (file)
@@ -947,7 +947,7 @@ namespace TrilinosWrappers
      * Return the total number of nonzero elements of this matrix (summed
      * over all MPI processes).
      */
-    size_type
+    std::uint64_t
     n_nonzero_elements() const;
 
     /**
@@ -3015,14 +3015,14 @@ namespace TrilinosWrappers
 
 
 
-  inline SparseMatrix::size_type
+  inline std::uint64_t
   SparseMatrix::n_nonzero_elements() const
   {
-#      ifndef DEAL_II_WITH_64BIT_INDICES
-    return matrix->NumGlobalNonzeros();
-#      else
-    return matrix->NumGlobalNonzeros64();
-#      endif
+    // Trilinos uses 64bit functions internally for attribute access, which
+    // return `long long`. They also offer 32bit variants that return `int`,
+    // however those call the 64bit version and convert the values to 32bit.
+    // There is no necessity in using the 32bit versions at all.
+    return static_cast<std::uint64_t>(matrix->NumGlobalNonzeros64());
   }
 
 
index eebc438d45bb30ff4deba5d96d976a477882e76e..5598833b23290504f70d4ff9e38d49cf80c14f8a 100644 (file)
@@ -712,7 +712,7 @@ namespace TrilinosWrappers
     /**
      * Return the number of nonzero elements of this sparsity pattern.
      */
-    size_type
+    std::uint64_t
     n_nonzero_elements() const;
 
     /**
index 303b0ac0c26cfd59347f63a71f470c43e4b90cce..aa05073f28423b29893b9afd1c6279b0f016c080 100644 (file)
@@ -286,14 +286,16 @@ namespace PETScWrappers
 
 
 
-  MatrixBase::size_type
+  std::uint64_t
   MatrixBase::n_nonzero_elements() const
   {
     MatInfo              mat_info;
     const PetscErrorCode ierr = MatGetInfo(matrix, MAT_GLOBAL_SUM, &mat_info);
     AssertThrow(ierr == 0, ExcPETScError(ierr));
 
-    return static_cast<size_type>(mat_info.nz_used);
+    // MatInfo logs quantities as PetscLogDouble. So we need to cast it to match
+    // out interface.
+    return static_cast<std::uint64_t>(mat_info.nz_used);
   }
 
 
index 27ca8cb3bec6efe63ff1c6cb2840f4cc9cf082c1..9a35d250a4a1434a05c319792ea4b3409fa02c6c 100644 (file)
@@ -135,10 +135,10 @@ namespace PETScWrappers
       return index_sets;
     }
 
-    BlockSparseMatrix::size_type
+    std::uint64_t
     BlockSparseMatrix::n_nonzero_elements() const
     {
-      size_type n_nonzero = 0;
+      std::uint64_t n_nonzero = 0;
       for (size_type rows = 0; rows < this->n_block_rows(); ++rows)
         for (size_type cols = 0; cols < this->n_block_cols(); ++cols)
           n_nonzero += this->block(rows, cols).n_nonzero_elements();
index bd6b591eb9bc87156d88c26abfc0fbcc95562d3d..586b9bb8c686cb7edea6c4573a1767e6af3fff49 100644 (file)
@@ -232,10 +232,10 @@ namespace TrilinosWrappers
 
 
 
-  BlockSparseMatrix::size_type
+  std::uint64_t
   BlockSparseMatrix::n_nonzero_elements() const
   {
-    size_type n_nonzero = 0;
+    std::uint64_t n_nonzero = 0;
     for (size_type rows = 0; rows < this->n_block_rows(); ++rows)
       for (size_type cols = 0; cols < this->n_block_cols(); ++cols)
         n_nonzero += this->block(rows, cols).n_nonzero_elements();
index 19710f75adb37a51c6866e6cfff38a72caaa879f..7109b5a8e551d12b13a61efbb55f5b7605cb1633 100644 (file)
@@ -932,12 +932,12 @@ namespace TrilinosWrappers
 
 
 
-  SparsityPattern::size_type
+  std::uint64_t
   SparsityPattern::n_nonzero_elements() const
   {
-    TrilinosWrappers::types::int_type nnz = n_global_entries(*graph);
+    TrilinosWrappers::types::int64_type nnz = n_global_entries(*graph);
 
-    return static_cast<size_type>(nnz);
+    return static_cast<std::uint64_t>(nnz);
   }
 
 

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.