source: XIOS/dev/dev_ym/XIOS_COUPLING/src/distribution/grid_remote_connector.cpp @ 1918

Last change on this file since 1918 was 1918, checked in by ymipsl, 4 years ago

Big update on on going work related to data distribution and transfer between clients and servers.

  • move all related file into distribution directorie
  • implement the concept of data "View"
  • implement the concept of "connector" which make the data transfer between 2 differents "Views"

YM

  • Property svn:eol-style set to native
  • Property svn:executable set to *
File size: 3.0 KB
Line 
1#include "grid_remote_connector.hpp"
2#include "client_client_dht_template.hpp"
3
4
5
6namespace xios
7{
8 
9  CGridRemoteConnector::CGridRemoteConnector(vector<CLocalView*>& srcView, vector<CDistributedView*>& dstView, MPI_Comm localComm) 
10                       : srcView_(srcView), dstView_(dstView), localComm_(localComm) 
11  {}
12
13  void CGridRemoteConnector::computeConnector(void)
14  {
15    computeGenericMethod() ;
16  }
17
18  void CGridRemoteConnector::computeGenericMethod(void)
19  {
20    // generic method, every element can be distributed
21    CClientClientDHTTemplate<int>::Index2VectorInfoTypeMap info ; // info map
22       
23// generate list of global index for dst view, and insert it into DHT map
24    int nDst = dstView_.size() ;
25    vector<size_t> dstGlobalIndex ;
26    vector<size_t> dstSliceSize(nDst) ;
27    dstSliceSize[nDst-1] = 1 ; 
28    for(int i=nDst-2; i>=0; i--)  dstSliceSize[i] = dstView_[i+1]->getGlobalSize()*dstSliceSize[i+1] ;
29    for(auto& ranks : dstView_[0]->getLocalSize())
30    {
31      dstGlobalIndex.clear() ;
32      int rank=ranks.first ;
33      size_t sliceIndex=0 ;
34      dstView_[nDst-1]->getGlobalIndex(rank, dstGlobalIndex, sliceIndex, dstSliceSize.data(), dstView_.data(), nDst-1) ;
35      for(auto globalIndex : dstGlobalIndex) info[globalIndex].push_back(rank) ; // insert into DHT
36    }
37   
38    CClientClientDHTTemplate<int> dataRanks(info, localComm_) ;
39   
40    // generate list of global index for src view
41    int nSrc = srcView_.size() ;
42    vector<size_t> srcSliceSize(nSrc) ;
43    srcSliceSize[nSrc-1] = 1 ;
44    for(int i=nSrc-2; i>=0; i--)  srcSliceSize[i] = srcView_[i+1]->getGlobalSize()*srcSliceSize[i+1] ;
45   
46    vector<size_t> srcGlobalIndex ;
47    size_t sliceIndex=0 ;
48    srcView_[nSrc-1]->getGlobalIndex(srcGlobalIndex, sliceIndex, srcSliceSize.data(), srcView_.data(), nSrc-1) ;
49   
50    CArray<size_t,1> srcGlobalIndexArray(srcGlobalIndex.data(), shape(srcGlobalIndex.size()),neverDeleteData) ;
51    dataRanks.computeIndexInfoMapping(srcGlobalIndexArray) ;
52    const auto& returnInfo = dataRanks.getInfoIndexMap() ;
53
54    vector<map<int, set<size_t>>> elements(nSrc) ; // internal representation of elements composing the grid
55
56    srcSliceSize[nSrc-1] = srcView_[nSrc-1]->getGlobalSize() ;
57    for(int i=nSrc-2 ; i>=0 ; i--) srcSliceSize[i] = srcView_[i]->getGlobalSize()*srcSliceSize[i+1] ;
58
59    for(auto& indRanks : returnInfo)
60    {
61      size_t gridIndexGlo=indRanks.first ;
62      auto& ranks = indRanks.second ;
63      for(int i=0;i<nSrc;i++)
64      {
65        auto& element = elements[i] ;
66        size_t localIndGlo = gridIndexGlo % srcSliceSize[i] ;
67        for(int rank : ranks) element[rank].insert(localIndGlo) ;
68      }
69    }
70
71    elements_.resize(nSrc) ;
72    for(int i=0 ; i<nSrc; i++)
73    {
74      auto& element=elements[i] ;
75      for(auto& rankInd : element)
76      {
77        int rank=rankInd.first ;
78        set<size_t>& indGlo = rankInd.second ;
79        CArray<size_t,1>& indGloArray = elements_[i][rank] ;
80        indGloArray.resize(indGlo.size()) ;
81        int j=0 ;
82        for (auto index : indGlo) { indGloArray(j) = index ; j++; }
83      }
84    }
85  }
86
87}
Note: See TracBrowser for help on using the repository browser.