Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Distribution.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 // Build maps for 1D or 2D matrix distribution
11 // Assumes square matrix
12 // Karen Devine, SNL
13 //
14 
15 #ifndef __TPETRA_DISTRIBUTION_HPP
16 #define __TPETRA_DISTRIBUTION_HPP
17 
18 #include <cstdio>
19 #include <cstdlib>
20 #include <iostream>
21 #include <fstream>
22 #include <set>
23 #ifndef __cplusplus
24 #define __cplusplus
25 #endif
26 
27 #include "Teuchos_Comm.hpp"
28 
29 namespace Tpetra {
30 
31 enum DistributionType {
32  TwoDRandom, // 2D randomly permuted distribution
33  TwoDLinear, // 2D linear distribution
34  TwoDVec, // 2D distribution based on vector assignment in file
35  OneDRandom, // 1D randomly permuted distribution
36  OneDLinear, // 1D linear distribution
37  OneDVec, // 1D distribution based on vector assignment in file
38  LowerTriangularBlock, // Seher Acer's lower-triangular block distrib
39  // for triangle counting
40  MMFile // Use values in matrix-market file as part assignment
41 };
42 
44 template <typename gno_t, typename scalar_t>
45 class Distribution {
46  public:
47  Distribution(size_t nrows_,
48  const Teuchos::RCP<const Teuchos::Comm<int> > &comm_,
49  const Teuchos::ParameterList &params)
50  : comm(comm_)
51  , me(comm_->getRank())
52  , np(comm_->getSize())
53  , nrows(nrows_) {}
54 
55  virtual ~Distribution(){};
56 
57  // Return the DistributionType for this distribution.
58  virtual enum DistributionType DistType() = 0;
59 
60  // Return whether this rank owns nonzero (i,j)
61  virtual bool Mine(gno_t i, gno_t j) = 0;
62  virtual bool Mine(gno_t i, gno_t j, int p) = 0;
63 
64  // Return whether this rank owns vector entry i
65  virtual bool VecMine(gno_t i) = 0;
66 
67  // Map of nonzeros needed for redistribution, handy for other things
68  using NZindex_t = std::pair<gno_t, gno_t>;
69  struct compareNzIndex { // sort nonzeros by row, then column
70  bool operator()(const NZindex_t &lhs, const NZindex_t &rhs) const {
71  if (lhs.first < rhs.first) return true;
72  if ((lhs.first == rhs.first) && (lhs.second < rhs.second)) return true;
73  return false;
74  }
75  };
76 
77  using LocalNZmap_t = std::map<NZindex_t, scalar_t, compareNzIndex>;
78 
79  // Redistribute nonzeros according to the needs of the Distribution
80  // Needed only when the final distribution cannot be determined until
81  // all nonzeros are known (e.g., final distribution depends on the number
82  // of nonzeros in a row).
83  // If the final distribution can be determined before all nonzeros (e.g.,
84  // Trilinos' traditional row map), the redistribution routine is a no-op.
85  // Thus, the default implementation is a no-op.
86  virtual void Redistribute(LocalNZmap_t &localNZ){};
87 
88  protected:
89  const Teuchos::RCP<const Teuchos::Comm<int> > comm;
90  int me; // my rank
91  int np; // number of ranks
92  size_t nrows; // global number of rows in the input matrix
93 
94  int HashToProc(gno_t i) {
95  // TODO PUT A GOOD HASH FUNCTION HERE!!!
96  // TODO FOR INITIAL TESTING, USE A CYCLIC HASH. LAME!
97  return (i % np);
98  }
99 };
100 
101 } // namespace Tpetra
102 
103 #include "Tpetra_Distribution2D.hpp"
104 #include "Tpetra_Distribution1D.hpp"
105 #include "Tpetra_DistributionMM.hpp"
106 #include "Tpetra_DistributionLowerTriangularBlock.hpp"
107 
108 #endif