source: XIOS/trunk/src/filter/grid_transformation.cpp @ 623

Last change on this file since 623 was 623, checked in by mhnguyen, 9 years ago

Implementing transformation algorithm: zoom axis (local commit)

+) Implement zoom axis: zoomed points are points not masked
+) Correct some minor bugs

Test
+) Ok with normal cases: zoom in the last of transformation list
+) There is still a bug in case of zoom then inverse

File size: 22.7 KB
Line 
1/*!
2   \file grid_transformation.cpp
3   \author Ha NGUYEN
4   \since 14 May 2015
5   \date 09 June 2015
6
7   \brief Interface for all transformations.
8 */
9#include "grid_transformation.hpp"
10#include "axis_algorithm_inverse.hpp"
11#include "axis_algorithm_zoom.hpp"
12#include "context.hpp"
13#include "context_client.hpp"
14#include "transformation_mapping.hpp"
15
16#include "axis_algorithm_transformation.hpp"
17
18namespace xios {
19CGridTransformation::CGridTransformation(CGrid* destination, CGrid* source)
20: gridSource_(source), gridDestination_(destination), originalGridSource_(source),
21  globalIndexOfCurrentGridSource_(0), globalIndexOfOriginalGridSource_(0)
22{
23  //Verify the compatibity between two grids
24  int numElement = gridDestination_->axis_domain_order.numElements();
25  if (numElement != gridSource_->axis_domain_order.numElements())
26    ERROR("CGridTransformation::CGridTransformation(CGrid* destination, CGrid* source)",
27       << "Two grids have different number of elements"
28       << "Number of elements of grid source " <<gridSource_->getId() << " is " << gridSource_->axis_domain_order.numElements()  << std::endl
29       << "Number of elements of grid destination " <<gridDestination_->getId() << " is " << numElement);
30
31  for (int i = 0; i < numElement; ++i)
32  {
33    if (gridDestination_->axis_domain_order(i) != gridSource_->axis_domain_order(i))
34      ERROR("CGridTransformation::CGridTransformation(CGrid* destination, CGrid* source)",
35         << "Transformed grid and its grid source have incompatible elements"
36         << "Grid source " <<gridSource_->getId() << std::endl
37         << "Grid destination " <<gridDestination_->getId());
38  }
39
40  std::vector<CAxis*> axisSrcTmp = gridSource_->getAxis(), axisSrc;
41  std::vector<CDomain*> domainSrcTmp = gridSource_->getDomains(), domainSrc;
42  for (int idx = 0; idx < axisSrcTmp.size(); ++idx)
43  {
44    CAxis* axis = CAxis::createAxis();
45    axis->setAttributes(axisSrcTmp[idx]);
46    axisSrc.push_back(axis);
47  }
48
49  for (int idx = 0; idx < domainSrcTmp.size(); ++idx)
50  {
51    CDomain* domain = CDomain::createDomain();
52    domain->setAttributes(domainSrcTmp[idx]);
53    domainSrc.push_back(domain);
54  }
55
56  gridSource_ = CGrid::createGrid(domainSrc, axisSrc, gridDestination_->axis_domain_order);
57  gridSourceDimensionSize_ = gridSource_->getGlobalDimension();
58  gridDestinationDimensionSize_ = gridDestination_->getGlobalDimension();
59
60  initializeMappingOfOriginalGridSource();
61  initializeAlgorithms();
62}
63
64void CGridTransformation::initializeMappingOfOriginalGridSource()
65{
66  CContext* context = CContext::getCurrent();
67  CContextClient* client=context->client;
68
69  CDistributionClient distribution(client->clientRank, originalGridSource_);
70  const CArray<size_t,1>& globalIndexGridDestSendToServer = distribution.getGlobalDataIndexSendToServer();
71
72  globalIndexOfCurrentGridSource_   = new CArray<size_t,1>(globalIndexGridDestSendToServer.numElements());
73  globalIndexOfOriginalGridSource_  = new CArray<size_t,1>(globalIndexGridDestSendToServer.numElements());
74  *globalIndexOfCurrentGridSource_  = globalIndexGridDestSendToServer;
75  *globalIndexOfOriginalGridSource_ = globalIndexGridDestSendToServer;
76}
77
78CGridTransformation::~CGridTransformation()
79{
80  std::list<CGenericAlgorithmTransformation*>::const_iterator itb = algoTransformation_.begin(), it,
81                                                              ite = algoTransformation_.end();
82  for (it = itb; it != ite; ++it) delete (*it);
83
84  std::map<int, std::vector<CArray<int,1>* > >::const_iterator itMapRecv, iteMapRecv;
85  itMapRecv = localIndexToReceiveOnGridDest_.begin();
86  iteMapRecv = localIndexToReceiveOnGridDest_.end();
87  for (; itMapRecv != iteMapRecv; ++itMapRecv)
88  {
89    int numVec = (itMapRecv->second).size();
90    for (int idx = 0; idx < numVec; ++idx) delete (itMapRecv->second)[idx];
91  }
92
93  std::map<int, CArray<int,1>* >::const_iterator itMap, iteMap;
94  itMap = localIndexToSendFromGridSource_.begin();
95  iteMap = localIndexToSendFromGridSource_.end();
96  for (; itMap != iteMap; ++itMap) delete (itMap->second);
97
98  if (0 != globalIndexOfCurrentGridSource_) delete globalIndexOfCurrentGridSource_;
99  if (0 != globalIndexOfOriginalGridSource_) delete globalIndexOfOriginalGridSource_;
100}
101
102void CGridTransformation::initializeAlgorithms()
103{
104  initializeAxisAlgorithms();
105  initializeDomainAlgorithms();
106}
107
108/*!
109  Initialize the algorithms corresponding to transformation info contained in each axis.
110If an axis has transformations, these transformations will be represented in form of vector of CTransformation pointers
111In general, each axis can have several transformations performed on itself. However, should they be done seperately or combinely (of course in order)?
112For now, one approach is to do these combinely but maybe this needs changing.
113*/
114void CGridTransformation::initializeAxisAlgorithms()
115{
116  std::vector<int> axisPositionInGrid;
117  std::vector<CAxis*> axisListDestP = gridDestination_->getAxis();
118  std::vector<CAxis*> axisListSrcP = gridSource_->getAxis();
119  if (!axisListDestP.empty())
120  {
121    int idx = 0;
122    for (int i = 0; i < gridDestination_->axis_domain_order.numElements(); ++i)
123    {
124      if (false == (gridDestination_->axis_domain_order)(i))
125      {
126        axisPositionInGrid.push_back(idx);
127        ++idx;
128      }
129      else idx += 2;
130    }
131
132    for (int i = 0; i < axisListDestP.size(); ++i)
133    {
134      elementPosition2AxisPositionInGrid_[axisPositionInGrid[i]] = i;
135      if (axisListDestP[i]->hasTransformation())
136      {
137        CAxis::TransMapTypes trans = axisListDestP[i]->getAllTransformations();
138        CAxis::TransMapTypes::const_iterator itb = trans.begin(), it,
139                                             ite = trans.end();
140        int transformationOrder = 0;
141        for (it = itb; it != ite; ++it)
142        {
143          listAlgos_.push_back(std::make_pair(axisPositionInGrid[i], std::make_pair(it->first, transformationOrder)));
144          ++transformationOrder;
145        }
146      }
147    }
148  }
149}
150
151void CGridTransformation::initializeDomainAlgorithms()
152{
153
154}
155
156void CGridTransformation::selectAlgo(int elementPositionInGrid, ETranformationType transType, int transformationOrder)
157{
158   selectAxisAlgo(elementPositionInGrid, transType, transformationOrder);
159}
160
161void CGridTransformation::selectAxisAlgo(int elementPositionInGrid, ETranformationType transType, int transformationOrder)
162{
163  std::vector<CAxis*> axisListDestP = gridDestination_->getAxis();
164  std::vector<CAxis*> axisListSrcP = gridSource_->getAxis();
165
166  int axisIndex =  elementPosition2AxisPositionInGrid_[elementPositionInGrid];
167  CAxis::TransMapTypes trans = axisListDestP[axisIndex]->getAllTransformations();
168  CAxis::TransMapTypes::const_iterator it = trans.begin();
169
170  for (int i = 0; i < transformationOrder; ++i, ++it) {}  // Find the correct transformation
171
172  CZoomAxis* zoomAxis = 0;
173  CGenericAlgorithmTransformation* algo = 0;
174  switch (transType)
175  {
176    case TRANS_ZOOM_AXIS:
177      zoomAxis = dynamic_cast<CZoomAxis*> (it->second);
178      algo = new CAxisAlgorithmZoom(axisListDestP[axisIndex], axisListSrcP[axisIndex], zoomAxis);
179      break;
180    case TRANS_INVERSE_AXIS:
181      algo = new CAxisAlgorithmInverse(axisListDestP[axisIndex], axisListSrcP[axisIndex]);
182      break;
183    default:
184      break;
185  }
186  algoTransformation_.push_back(algo);
187
188}
189
190void CGridTransformation::selectDomainAlgo(int elementPositionInGrid, ETranformationType transType, int transformationOrder)
191{
192
193}
194
195void CGridTransformation::setUpGrid(int elementPositionInGrid, ETranformationType transType)
196{
197  std::vector<CAxis*> axisListDestP = gridDestination_->getAxis();
198  std::vector<CAxis*> axisListSrcP = gridSource_->getAxis();
199
200  int axisIndex;
201  switch (transType)
202  {
203    case TRANS_ZOOM_AXIS:
204    case TRANS_INVERSE_AXIS:
205      axisIndex =  elementPosition2AxisPositionInGrid_[elementPositionInGrid];
206      axisListSrcP[axisIndex]->duplicateAttributes(axisListDestP[axisIndex]);
207      break;
208    default:
209      break;
210  }
211}
212
213void CGridTransformation::computeAll()
214{
215  CContext* context = CContext::getCurrent();
216  CContextClient* client=context->client;
217
218  ListAlgoType::const_iterator itb = listAlgos_.begin(),
219                               ite = listAlgos_.end(), it;
220  CGenericAlgorithmTransformation* algo = 0;
221  for (it = itb; it != ite; ++it)
222  {
223    std::map<size_t, std::set<size_t> > globaIndexMapFromDestToSource;
224    int elementPositionInGrid = it->first;
225    ETranformationType transType = (it->second).first;
226    int transformationOrder = (it->second).second;
227
228    // First of all, select an algorithm
229    selectAlgo(elementPositionInGrid, transType, transformationOrder);
230    algo = algoTransformation_.back();
231
232    // Recalculate the distribution of grid destination
233    CDistributionClient distributionClientDest(client->clientRank, gridDestination_);
234    const CArray<size_t,1>& globalIndexGridDestSendToServer = distributionClientDest.getGlobalDataIndexSendToServer();
235
236   std::cout << "global index grid  dest send to server " << globalIndexGridDestSendToServer << std::endl;
237    // ComputeTransformation of global index of each element
238    std::vector<int> gridDestinationDimensionSize = gridDestination_->getGlobalDimension();
239    int elementPosition = it->first;
240    algo->computeGlobalSourceIndex(elementPosition,
241                                   gridDestinationDimensionSize,
242                                   globalIndexGridDestSendToServer,
243                                   globaIndexMapFromDestToSource);
244
245    // Compute transformation of global indexes among grids
246    computeTransformationFromOriginalGridSource(globaIndexMapFromDestToSource);
247
248    // Now grid destination becomes grid source in a new transformation
249    setUpGrid(elementPositionInGrid, transType);
250  }
251
252 std::cout << "global index destination 0 final " << *globalIndexOfCurrentGridSource_ << std::endl;
253 std::cout << "global index destination 1 final " << *globalIndexOfOriginalGridSource_ << std::endl;
254  updateFinalGridDestination();
255  computeFinalTransformationMapping();
256}
257
258
259/*!
260  After applying the algorithms, there are some informations on grid destination needing change, for now, there are:
261   +) mask
262*/
263void CGridTransformation::updateFinalGridDestination()
264{
265  CContext* context = CContext::getCurrent();
266  CContextClient* client=context->client;
267
268  //First of all, retrieve info of local mask of grid destination
269  CDistributionClient distributionClientDest(client->clientRank, gridDestination_);
270  const CArray<int, 1>& localMaskIndexOnClientDest = distributionClientDest.getLocalMaskIndexOnClient();
271  std::cout << "local mask " << localMaskIndexOnClientDest << std::endl;
272
273  const CArray<size_t,1>& globalIndexOnClientDest = distributionClientDest.getGlobalDataIndexSendToServer();
274  std::cout << "global index " << globalIndexOnClientDest <<  std::endl;
275  CArray<size_t, 1>::const_iterator itbArr, itArr, iteArr;
276  itbArr = globalIndexOnClientDest.begin();
277  iteArr = globalIndexOnClientDest.end();
278
279  // Then find out which index became invalid (become masked after being applied the algorithms, or demande some masked points from grid source)
280  int num = globalIndexOfOriginalGridSource_->numElements();
281  const size_t sfmax = NumTraits<unsigned long>::sfmax();
282  int maskIndexNum = 0;
283  for (int idx = 0; idx < num; ++idx)
284  {
285    if (sfmax == (*globalIndexOfOriginalGridSource_)(idx))
286    {
287      size_t maskedGlobalIndex = (*globalIndexOfCurrentGridSource_)(idx);
288      itArr = std::find(itbArr, iteArr, maskedGlobalIndex);
289      if (iteArr != itArr) ++maskIndexNum;
290    }
291  }
292
293  CArray<int,1>* maskIndexToModify = new CArray<int,1>(maskIndexNum);
294  maskIndexNum = 0;
295  for (int idx = 0; idx < num; ++idx)
296  {
297    if (sfmax == (*globalIndexOfOriginalGridSource_)(idx))
298    {
299      size_t maskedGlobalIndex = (*globalIndexOfCurrentGridSource_)(idx);
300      itArr = std::find(itbArr, iteArr, maskedGlobalIndex);
301      if (iteArr != itArr)
302      {
303        int localIdx = std::distance(itbArr, itArr);
304        (*maskIndexToModify)(maskIndexNum) = (localMaskIndexOnClientDest)(localIdx);
305        ++maskIndexNum;
306      }
307    }
308  }
309
310  std::cout << "index to modify " << *maskIndexToModify << std::endl;
311  gridDestination_->modifyMask(*maskIndexToModify);
312
313  delete maskIndexToModify;
314}
315
316/*!
317  A transformation from a grid source to grid destination often passes several intermediate grids, which play a role of
318temporary grid source and/or grid destination. This function makes sure that global index of original grid source are mapped correctly to
319the final grid destination
320*/
321void CGridTransformation::computeTransformationFromOriginalGridSource(const std::map<size_t, std::set<size_t> >& globaIndexMapFromDestToSource)
322{
323  CContext* context = CContext::getCurrent();
324  CContextClient* client=context->client;
325
326  CTransformationMapping transformationMap(gridDestination_, gridSource_);
327
328    // Then compute transformation mapping among clients
329  transformationMap.computeTransformationMapping(globaIndexMapFromDestToSource);
330
331  const std::map<int,std::vector<std::vector<size_t> > >& globalIndexToReceive = transformationMap.getGlobalIndexReceivedOnGridDestMapping();
332  const std::map<int,std::vector<size_t> >& globalIndexToSend = transformationMap.getGlobalIndexSendToGridDestMapping();
333
334 // Sending global index of original grid source
335  std::map<int,std::vector<size_t> >::const_iterator itbSend = globalIndexToSend.begin(), itSend,
336                                                     iteSend = globalIndexToSend.end();
337  CArray<size_t,1>::const_iterator itbArr = globalIndexOfCurrentGridSource_->begin(), itArr,
338                                   iteArr = globalIndexOfCurrentGridSource_->end();
339 int sendBuffSize = 0;
340 for (itSend = itbSend; itSend != iteSend; ++itSend) sendBuffSize += (itSend->second).size();
341
342 std::cout << "global index destination 0 before" << *globalIndexOfCurrentGridSource_ << std::endl;
343 std::cout << "global index destination 1 before" << *globalIndexOfOriginalGridSource_ << std::endl;
344
345 typedef unsigned long Scalar;
346 unsigned long* sendBuff, *currentSendBuff;
347 if (0 != sendBuffSize) sendBuff = new unsigned long [sendBuffSize];
348 for (StdSize idx = 0; idx < sendBuffSize; ++idx) sendBuff[idx] = NumTraits<Scalar>::sfmax();
349
350 int currentBuffPosition = 0;
351 for (itSend = itbSend; itSend != iteSend; ++itSend)
352 {
353   int destRank = itSend->first;
354   const std::vector<size_t>& globalIndexOfCurrentGridSourceToSend = itSend->second;
355   int countSize = globalIndexOfCurrentGridSourceToSend.size();
356   for (int idx = 0; idx < (countSize); ++idx)
357   {
358     itArr = std::find(itbArr, iteArr, globalIndexOfCurrentGridSourceToSend[idx]);
359     if (iteArr != itArr)
360     {
361       int index = std::distance(itbArr, itArr);
362       sendBuff[idx+currentBuffPosition] = (*globalIndexOfOriginalGridSource_)(index);
363     }
364   }
365   currentSendBuff = sendBuff + currentBuffPosition;
366   MPI_Send(currentSendBuff, countSize, MPI_UNSIGNED_LONG, destRank, 14, client->intraComm);
367   currentBuffPosition += countSize;
368 }
369
370 // Receiving global index of grid source sending from current grid source
371 std::map<int,std::vector<std::vector<size_t> > >::const_iterator itbRecv = globalIndexToReceive.begin(), itRecv,
372                                                                  iteRecv = globalIndexToReceive.end();
373 int recvBuffSize = 0;
374 for (itRecv = itbRecv; itRecv != iteRecv; ++itRecv) recvBuffSize += (itRecv->second).size();
375
376 unsigned long* recvBuff, *currentRecvBuff;
377 if (0 != recvBuffSize) recvBuff = new unsigned long [recvBuffSize];
378 for (StdSize idx = 0; idx < recvBuffSize; ++idx) recvBuff[idx] = NumTraits<Scalar>::sfmax();
379
380 int currentRecvBuffPosition = 0;
381 for (itRecv = itbRecv; itRecv != iteRecv; ++itRecv)
382 {
383   MPI_Status status;
384   int srcRank = itRecv->first;
385   int countSize = (itRecv->second).size();
386   currentRecvBuff = recvBuff + currentRecvBuffPosition;
387   MPI_Recv(currentRecvBuff, countSize, MPI_UNSIGNED_LONG, srcRank, 14, client->intraComm, &status);
388   currentRecvBuffPosition += countSize;
389 }
390
391 int nbCurrentGridSource = 0;
392 for (itRecv = itbRecv; itRecv != iteRecv; ++itRecv)
393 {
394   int ssize = (itRecv->second).size();
395   for (int idx = 0; idx < ssize; ++idx)
396   {
397     nbCurrentGridSource += (itRecv->second)[idx].size();
398   }
399 }
400
401 if (globalIndexOfCurrentGridSource_->numElements()  != nbCurrentGridSource)
402 {
403   if ((0 != nbCurrentGridSource) && (0 != globalIndexOfCurrentGridSource_))
404   {
405     delete globalIndexOfCurrentGridSource_;
406     globalIndexOfCurrentGridSource_ = new CArray<size_t,1>(nbCurrentGridSource);
407   }
408 }
409
410 if (globalIndexOfOriginalGridSource_->numElements() != nbCurrentGridSource)
411 {
412   if ((0 != nbCurrentGridSource) && (0 != globalIndexOfOriginalGridSource_))
413   {
414     delete globalIndexOfOriginalGridSource_;
415     globalIndexOfOriginalGridSource_ = new CArray<size_t,1>(nbCurrentGridSource);
416   }
417 }
418
419 int k = 0;
420 currentRecvBuff = recvBuff;
421 for (itRecv = itbRecv; itRecv != iteRecv; ++itRecv)
422 {
423   int countSize = (itRecv->second).size();
424   for (int idx = 0; idx < countSize; ++idx, ++currentRecvBuff)
425   {
426     int ssize = (itRecv->second)[idx].size();
427     for (int i = 0; i < ssize; ++i)
428     {
429       (*globalIndexOfCurrentGridSource_)(k) = (itRecv->second)[idx][i];
430       (*globalIndexOfOriginalGridSource_)(k) = *currentRecvBuff;
431       ++k;
432     }
433   }
434 }
435
436 std::cout << "global index destination 0 after " << *globalIndexOfCurrentGridSource_ << std::endl;
437 std::cout << "global index destination 1 after " << *globalIndexOfOriginalGridSource_ << std::endl;
438 if (0 != sendBuffSize) delete [] sendBuff;
439 if (0 != recvBuffSize) delete [] recvBuff;
440}
441
442/*!
443  Compute transformation mapping between grid source and grid destination
444  The transformation between grid source and grid destination is represented in form of mapping between global index
445of two grids. Then local index mapping between data on each grid will be found out thanks to these global indexes
446*/
447void CGridTransformation::computeFinalTransformationMapping()
448{
449  CContext* context = CContext::getCurrent();
450  CContextClient* client=context->client;
451
452  CTransformationMapping transformationMap(gridDestination_, originalGridSource_);
453
454  std::map<size_t, std::set<size_t> > globaIndexMapFromDestToSource;
455
456  int nb = globalIndexOfCurrentGridSource_->numElements();
457  const size_t sfmax = NumTraits<unsigned long>::sfmax();
458  for (int idx = 0; idx < nb; ++idx)
459  {
460    if (sfmax != (*globalIndexOfOriginalGridSource_)(idx))
461      globaIndexMapFromDestToSource[(*globalIndexOfCurrentGridSource_)(idx)].insert((*globalIndexOfOriginalGridSource_)(idx));
462  }
463
464  // Then compute transformation mapping among clients
465  transformationMap.computeTransformationMapping(globaIndexMapFromDestToSource);
466
467  const std::map<int,std::vector<std::vector<size_t> > >& globalIndexToReceive = transformationMap.getGlobalIndexReceivedOnGridDestMapping();
468  const std::map<int,std::vector<size_t> >& globalIndexToSend = transformationMap.getGlobalIndexSendToGridDestMapping();
469
470  CDistributionClient distributionClientDest(client->clientRank, gridDestination_);
471  CDistributionClient distributionClientSrc(client->clientRank, originalGridSource_);
472
473//  const CArray<int, 1>& localIndexOnClientDest = distributionClientDest.getLocalDataIndexOnClient(); //gridDestination_->getDistributionClient()->getLocalDataIndexOnClient();
474  const CArray<int, 1>& localIndexOnClientDest = distributionClientDest.getLocalDataIndexSendToServer();
475  const CArray<size_t,1>& globalIndexOnClientDest = distributionClientDest.getGlobalDataIndexSendToServer(); //gridDestination_->getDistributionClient()->getGlobalDataIndexSendToServer();
476
477 std::cout << "dest: local index " << localIndexOnClientDest << std::endl;
478 std::cout << "dest: global index " << globalIndexOnClientDest << std::endl;
479  const CArray<int, 1>& localIndexOnClientSrc = distributionClientSrc.getLocalDataIndexOnClient(); //gridSource_->getDistributionClient()->getLocalDataIndexOnClient();
480  const CArray<size_t,1>& globalIndexOnClientSrc = distributionClientSrc.getGlobalDataIndexSendToServer(); //gridSource_->getDistributionClient()->getGlobalDataIndexSendToServer();
481 std::cout << "src: local index " << localIndexOnClientSrc << std::endl;
482 std::cout << "src: global index " << globalIndexOnClientSrc << std::endl;
483  std::vector<size_t>::const_iterator itbVec, itVec, iteVec;
484  CArray<size_t, 1>::const_iterator itbArr, itArr, iteArr;
485
486  std::map<int,std::vector<std::vector<size_t> > >::const_iterator itbMapRecv, itMapRecv, iteMapRecv;
487
488  // Find out local index on grid destination (received)
489  itbMapRecv = globalIndexToReceive.begin();
490  iteMapRecv = globalIndexToReceive.end();
491  itbArr = globalIndexOnClientDest.begin();
492  iteArr = globalIndexOnClientDest.end();
493  for (itMapRecv = itbMapRecv; itMapRecv != iteMapRecv; ++itMapRecv)
494  {
495    int sourceRank = itMapRecv->first;
496    int numGlobalIndex = (itMapRecv->second).size();
497    for (int i = 0; i < numGlobalIndex; ++i)
498    {
499      int vecSize = ((itMapRecv->second)[i]).size();
500      CArray<int,1>* ptr = new CArray<int,1>(vecSize);
501      localIndexToReceiveOnGridDest_[sourceRank].push_back(ptr);
502      for (int idx = 0; idx < vecSize; ++idx)
503      {
504        itArr = std::find(itbArr, iteArr, (itMapRecv->second)[i][idx]);
505        if (iteArr != itArr)
506        {
507          int localIdx = std::distance(itbArr, itArr);
508          (*localIndexToReceiveOnGridDest_[sourceRank][i])(idx) = localIndexOnClientDest(localIdx); // Local index of un-extracted data (only domain)
509//          (*localIndexToReceiveOnGridDest_[sourceRank][i])(idx) = (localIdx); // Local index of extracted data
510        }
511      }
512    }
513//    std::cout << "local index to receive from source Rank = " << sourceRank << (*localIndexToReceiveOnGridDest_[sourceRank][i]) << std::endl;
514  }
515
516  std::map<int,std::vector<size_t> >::const_iterator itbMap, itMap, iteMap;
517  // Find out local index on grid source (to send)
518  itbMap = globalIndexToSend.begin();
519  iteMap = globalIndexToSend.end();
520  itbArr = globalIndexOnClientSrc.begin();
521  iteArr = globalIndexOnClientSrc.end();
522  for (itMap = itbMap; itMap != iteMap; ++itMap)
523  {
524    CArray<int,1>* ptr = new CArray<int,1>((itMap->second).size());
525    localIndexToSendFromGridSource_[itMap->first] = ptr;
526    int destRank = itMap->first;
527    int vecSize = (itMap->second).size();
528    for (int idx = 0; idx < vecSize; ++idx)
529    {
530      itArr = std::find(itbArr, iteArr, (itMap->second)[idx]);
531      if (iteArr != itArr)
532      {
533        int localIdx = std::distance(itbArr, itArr);
534//        (*localIndexToSendFromGridSource_[destRank])(idx) = localIndexOnClientSrc(localIdx);
535        (*localIndexToSendFromGridSource_[destRank])(idx) = (localIdx);
536      }
537    }
538    std::cout << "local index to send to dest Rank = " << destRank << (*localIndexToSendFromGridSource_[destRank]) << std::endl;
539  }
540}
541
542/*!
543  Local index of data which need sending from the grid source
544  \return local index of data
545*/
546std::map<int, CArray<int,1>* > CGridTransformation::getLocalIndexToSendFromGridSource()
547{
548  return localIndexToSendFromGridSource_;
549}
550
551/*!
552  Local index of data which will be received on the grid destination
553  \return local index of data
554*/
555std::map<int, std::vector<CArray<int,1>* > > CGridTransformation::getLocalIndexToReceiveOnGridDest()
556{
557  return localIndexToReceiveOnGridDest_;
558}
559
560}
Note: See TracBrowser for help on using the repository browser.