XIOS  1.0
Xml I/O Server
 Tout Classes Espaces de nommage Fichiers Fonctions Variables Définitions de type Énumérations Valeurs énumérées Amis Macros
dht_auto_indexing.cpp
Aller à la documentation de ce fichier.
1 
9 #include "dht_auto_indexing.hpp"
10 
11 namespace xios
12 {
13 
15  {
16  }
17 
24  const MPI_Comm& clientIntraComm)
25  : CClientClientDHTTemplate<size_t>(clientIntraComm)
26  {
27 
28  nbIndexOnProc_ = hashValue.size();
29  size_t nbIndexAccum;
30  MPI_Scan(&nbIndexOnProc_, &nbIndexAccum, 1, MPI_UNSIGNED_LONG, MPI_SUM, clientIntraComm);
31 
32  // Broadcasting the total number of indexes
33  int rank, size;
34  MPI_Comm_rank(clientIntraComm, &rank);
35  MPI_Comm_size(clientIntraComm, &size);
36  if (rank == (size-1)) nbIndexesGlobal_ = nbIndexAccum;
37  MPI_Bcast(&nbIndexesGlobal_, 1, MPI_UNSIGNED_LONG, size-1, clientIntraComm);
38 
39  CArray<size_t,1>::const_iterator itbIdx = hashValue.begin(), itIdx,
40  iteIdx = hashValue.end();
41 
42  size_t idx = 0;
43  beginIndexOnProc_ = nbIndexAccum - nbIndexOnProc_;
44  globalIndex_.resize(nbIndexOnProc_);
45  for (itIdx = itbIdx; itIdx != iteIdx; ++itIdx)
46  {
47  globalIndex_[idx] = beginIndexOnProc_ + idx;
48  ++idx ;
49  }
50  }
51 
60  const MPI_Comm& clientIntraComm)
61  : CClientClientDHTTemplate<size_t>(clientIntraComm)
62  {
63 
64  nbIndexOnProc_ = hashInitMap.size();
65  size_t nbIndexAccum;
66  MPI_Scan(&nbIndexOnProc_, &nbIndexAccum, 1, MPI_UNSIGNED_LONG, MPI_SUM, clientIntraComm);
67 
68  int rank, size;
69  MPI_Comm_rank(clientIntraComm, &rank);
70  MPI_Comm_size(clientIntraComm, &size);
71  if (rank == (size-1)) nbIndexesGlobal_ = nbIndexAccum;
72  MPI_Bcast(&nbIndexesGlobal_, 1, MPI_UNSIGNED_LONG, size-1, clientIntraComm);
73 
74  Index2VectorInfoTypeMap::iterator itbIdx = hashInitMap.begin(), itIdx,
75  iteIdx = hashInitMap.end();
76  size_t idx = 0;
77  beginIndexOnProc_ = nbIndexAccum - nbIndexOnProc_;
78  globalIndex_.resize(nbIndexOnProc_);
79  for (itIdx = itbIdx; itIdx != iteIdx; ++itIdx)
80  {
81 // (itIdx->second)[0] = beginIndexOnProc_ + idx;
82  (itIdx->second)[1] = beginIndexOnProc_ + idx;
83  globalIndex_[idx] = beginIndexOnProc_ + idx;
84  ++idx ;
85  }
86  }
87 
93  {
94  return nbIndexesGlobal_;
95  }
96 
102  {
103  return beginIndexOnProc_;
104  }
105 
111  {
112  return nbIndexOnProc_;
113  }
114 
115 }
virtual size_t size(void) const
Definition: array_new.hpp:548
Auto assign global index across processes.
std::vector< size_t > globalIndex_
This class provides the similar features like.
#define xios(arg)
size_t getNbIndexesGlobal() const
Returns the total number of global indexes.
std::unordered_map< size_t, std::vector< InfoType > > Index2VectorInfoTypeMap
size_t getIndexStart() const
Returns the starting global index for a proc.
CDHTAutoIndexing(const CArray< size_t, 1 > &hashValue, const MPI_Comm &clientIntraComm)
size_t getIndexCount() const
Returns the number of global indexes on a proc.