source: XIOS/dev/dev_trunk_omp/src/dht_auto_indexing.cpp @ 1670

Last change on this file since 1670 was 1646, checked in by yushan, 5 years ago

branch merged with trunk @1645. arch file (ep&mpi) added for ADA

File size: 3.4 KB
RevLine 
[892]1/*!
2   \file dht_auto_indexing.cpp
3   \author Ha NGUYEN
4   \since 6 Jul 2016
5   \date 6 Jul 2016
6
7   \brief Auto assign global index across processes.
8 */
9#include "dht_auto_indexing.hpp"
[1646]10#ifdef _usingEP
[1601]11using namespace ep_lib;
[1646]12#endif
[892]13
14namespace xios
15{
16
[924]17  CDHTAutoIndexing::~CDHTAutoIndexing()
18  {
19  }
[892]20
[924]21  /*!
22    \param [in]
23    \param [in]
24  */
[892]25
[924]26  CDHTAutoIndexing::CDHTAutoIndexing(const CArray<size_t,1>& hashValue,
27                                     const MPI_Comm& clientIntraComm)
28    : CClientClientDHTTemplate<size_t>(clientIntraComm)
[892]29  {
[924]30
31    nbIndexOnProc_ = hashValue.size();
32    size_t nbIndexAccum;
33    MPI_Scan(&nbIndexOnProc_, &nbIndexAccum, 1, MPI_UNSIGNED_LONG, MPI_SUM, clientIntraComm);
34
35    // Broadcasting the total number of indexes
36    int rank, size;
37    MPI_Comm_rank(clientIntraComm, &rank);
38    MPI_Comm_size(clientIntraComm, &size);
39    if (rank == (size-1)) nbIndexesGlobal_ = nbIndexAccum;
40    MPI_Bcast(&nbIndexesGlobal_, 1, MPI_UNSIGNED_LONG, size-1, clientIntraComm);
41
42    CArray<size_t,1>::const_iterator itbIdx = hashValue.begin(), itIdx,
43                                     iteIdx = hashValue.end();
44
45    size_t idx = 0;
46    beginIndexOnProc_ = nbIndexAccum - nbIndexOnProc_;
47    globalIndex_.resize(nbIndexOnProc_);
48    for (itIdx = itbIdx; itIdx != iteIdx; ++itIdx)
49    {
50      globalIndex_[idx] = beginIndexOnProc_ + idx;
51      ++idx ;
52    }
[892]53  }
54
[924]55  /*!
56   * \fn void CDHTAutoIndexing::CDHTAutoIndexing(Index2VectorInfoTypeMap& hashInitMap, const MPI_Comm& clientIntraComm)
57   * Assigns a global index to unique input indexes.
58   * The returned map has unique indexes as a key and global indexes as a mapped value.
59   * \param [in] hashInitMap map<size_t, vector<size_t>> is a map of unique indexes.
60   * \param [in] clientIntraComm
61  */
62  CDHTAutoIndexing::CDHTAutoIndexing(Index2VectorInfoTypeMap& hashInitMap,
63                                     const MPI_Comm& clientIntraComm)
64    : CClientClientDHTTemplate<size_t>(clientIntraComm)
65  {
[892]66
[924]67    nbIndexOnProc_ = hashInitMap.size();
68    size_t nbIndexAccum;
69    MPI_Scan(&nbIndexOnProc_, &nbIndexAccum, 1, MPI_UNSIGNED_LONG, MPI_SUM, clientIntraComm);
[892]70
[924]71    int rank, size;
72    MPI_Comm_rank(clientIntraComm, &rank);
73    MPI_Comm_size(clientIntraComm, &size);
74    if (rank == (size-1)) nbIndexesGlobal_ = nbIndexAccum;
75    MPI_Bcast(&nbIndexesGlobal_, 1, MPI_UNSIGNED_LONG, size-1, clientIntraComm);
76
77    Index2VectorInfoTypeMap::iterator itbIdx = hashInitMap.begin(), itIdx,
78                                      iteIdx = hashInitMap.end();
79    size_t idx = 0;
80    beginIndexOnProc_ = nbIndexAccum - nbIndexOnProc_;
81    globalIndex_.resize(nbIndexOnProc_);
82    for (itIdx = itbIdx; itIdx != iteIdx; ++itIdx)
83    {
[1158]84//      (itIdx->second)[0] = beginIndexOnProc_ + idx;
85      (itIdx->second)[1] = beginIndexOnProc_ + idx;
[924]86      globalIndex_[idx] = beginIndexOnProc_ + idx;
87      ++idx ;
88    }
89  }
90
91  /*!
92   * \fn size_t CDHTAutoIndexing::getNbIndexesGlobal() const
93   * Returns the total number of global indexes.
94  */
95  size_t CDHTAutoIndexing::getNbIndexesGlobal() const
[892]96  {
[924]97    return nbIndexesGlobal_;
[892]98  }
99
[924]100  /*!
101   * \fn size_t CDHTAutoIndexing::getIndexStart() const
102   * Returns the starting global index for a proc.
103  */
104  size_t CDHTAutoIndexing::getIndexStart() const
105  {
106    return beginIndexOnProc_;
107  }
108
109  /*!
110   * \fn size_t CDHTAutoIndexing::getIndexCount() const
111   * Returns the number of global indexes on a proc.
112  */
113  size_t CDHTAutoIndexing::getIndexCount() const
114  {
115    return nbIndexOnProc_;
116  }
117
[892]118}
Note: See TracBrowser for help on using the repository browser.