source: XIOS/trunk/src/transformation/grid_transformation.cpp @ 895

Last change on this file since 895 was 895, checked in by mhnguyen, 8 years ago

Adding a new transformation: Reduce a domain to an axis

Test
+) On Curie
+) Tests pass and are correct

  • Property svn:executable set to *
File size: 27.8 KB
Line 
1/*!
2   \file grid_transformation.cpp
3   \author Ha NGUYEN
4   \since 14 May 2015
5   \date 02 Jul 2015
6
7   \brief Interface for all transformations.
8 */
9#include "grid_transformation.hpp"
10#include "algo_types.hpp"
11#include "context.hpp"
12#include "context_client.hpp"
13#include "distribution_client.hpp"
14#include "mpi_tag.hpp"
15#include "grid.hpp"
16#include <boost/unordered_map.hpp>
17
18namespace xios {
19CGridTransformation::CGridTransformation(CGrid* destination, CGrid* source)
20 : CGridTransformationSelector(destination, source),
21  tmpGridDestination_(destination), originalGridSource_(source),
22  tempGridSrcs_(), tempGridDests_(),
23  dynamicalTransformation_(false), timeStamp_()
24{
25}
26
27CGridTransformation::~CGridTransformation()
28{
29}
30
31/*!
32  Select algorithm of a scalar correspoding to its transformation type and its position in each element
33  \param [in] elementPositionInGrid position of element in grid. e.g: a grid has 1 domain and 1 axis, then position of domain is 0 and position of axis is 1
34  \param [in] transType transformation type, for now we have
35  \param [in] transformationOrder position of the transformation in an element (an element can have several transformation)
36*/
37void CGridTransformation::selectScalarAlgo(int elementPositionInGrid, ETranformationType transType, int transformationOrder)
38{
39  int scalarSrcIndex = -1, axisSrcIndex = -1, domainSrcIndex = -1;
40  std::vector<CScalar*> scaListDestP = gridDestination_->getScalars();
41  std::vector<CScalar*> scaListSrcP  = gridSource_->getScalars();
42  std::vector<CAxis*> axisListSrcP   = gridSource_->getAxis();
43  std::vector<CDomain*> domainListSrcP = gridSource_->getDomains();
44
45  int scalarDstIndex =  elementPositionInGridDst2ScalarPosition_[elementPositionInGrid];
46  CScalar::TransMapTypes trans = scaListDestP[scalarDstIndex]->getAllTransformations();
47  CScalar::TransMapTypes::const_iterator it = trans.begin();
48
49  for (int i = 0; i < transformationOrder; ++i, ++it) {}  // Find the correct transformation
50
51  CReduceAxisToScalar* reduceAxis = 0;
52  CGenericAlgorithmTransformation* algo = 0;
53  switch (transType)
54  {
55    case TRANS_REDUCE_AXIS_TO_SCALAR:
56      reduceAxis = dynamic_cast<CReduceAxisToScalar*> (it->second);
57      axisSrcIndex = elementPositionInGridSrc2AxisPosition_[elementPositionInGrid];
58      algo = new CScalarAlgorithmReduceScalar(scaListDestP[scalarDstIndex], axisListSrcP[axisSrcIndex], reduceAxis);
59      break;
60    default:
61      break;
62  }
63  algoTransformation_.push_back(algo);
64}
65
66/*!
67  Select algorithm of an axis correspoding to its transformation type and its position in each element
68  \param [in] elementPositionInGrid position of element in grid. e.g: a grid has 1 domain and 1 axis, then position of domain is 0 and position of axis is 1
69  \param [in] transType transformation type, for now we have zoom_axis, inverse_axis, interpolate_axis
70  \param [in] transformationOrder position of the transformation in an element (an element can have several transformation)
71*/
72void CGridTransformation::selectAxisAlgo(int elementPositionInGrid, ETranformationType transType, int transformationOrder)
73{
74  int axisSrcIndex = -1, domainSrcIndex = -1;
75  std::vector<CAxis*> axisListDestP = gridDestination_->getAxis();
76  std::vector<CAxis*> axisListSrcP = gridSource_->getAxis();
77  std::vector<CDomain*> domainListSrcP = gridSource_->getDomains();
78
79  int axisDstIndex =  elementPositionInGridDst2AxisPosition_[elementPositionInGrid];
80  CAxis::TransMapTypes trans = axisListDestP[axisDstIndex]->getAllTransformations();
81  CAxis::TransMapTypes::const_iterator it = trans.begin();
82
83  for (int i = 0; i < transformationOrder; ++i, ++it) {}  // Find the correct transformation
84
85  CZoomAxis* zoomAxis = 0;
86  CInterpolateAxis* interpAxis = 0;
87  CReduceDomainToAxis* reduceDomain = 0;
88  CExtractDomainToAxis* extractDomain = 0;
89  CGenericAlgorithmTransformation* algo = 0;
90  switch (transType)
91  {
92    case TRANS_INTERPOLATE_AXIS:
93      interpAxis = dynamic_cast<CInterpolateAxis*> (it->second);
94      axisSrcIndex = elementPositionInGridSrc2AxisPosition_[elementPositionInGrid];
95      algo = new CAxisAlgorithmInterpolate(axisListDestP[axisDstIndex], axisListSrcP[axisSrcIndex], interpAxis);
96      break;
97    case TRANS_ZOOM_AXIS:
98      zoomAxis = dynamic_cast<CZoomAxis*> (it->second);
99      axisSrcIndex = elementPositionInGridSrc2AxisPosition_[elementPositionInGrid];
100      algo = new CAxisAlgorithmZoom(axisListDestP[axisDstIndex], axisListSrcP[axisSrcIndex], zoomAxis);
101      break;
102    case TRANS_INVERSE_AXIS:
103      axisSrcIndex = elementPositionInGridSrc2AxisPosition_[elementPositionInGrid];
104      algo = new CAxisAlgorithmInverse(axisListDestP[axisDstIndex], axisListSrcP[axisSrcIndex]);
105      break;
106    case TRANS_REDUCE_DOMAIN_TO_AXIS:
107      reduceDomain = dynamic_cast<CReduceDomainToAxis*> (it->second);
108      domainSrcIndex = elementPositionInGridSrc2DomainPosition_[elementPositionInGrid];
109      algo = new CAxisAlgorithmReduceDomain(axisListDestP[axisDstIndex], domainListSrcP[domainSrcIndex], reduceDomain);
110      break;
111    case TRANS_EXTRACT_DOMAIN_TO_AXIS:
112      extractDomain = dynamic_cast<CExtractDomainToAxis*> (it->second);
113      domainSrcIndex = elementPositionInGridSrc2DomainPosition_[elementPositionInGrid];
114      algo = new CAxisAlgorithmExtractDomain(axisListDestP[axisDstIndex], domainListSrcP[domainSrcIndex], extractDomain);
115      break;
116    default:
117      break;
118  }
119  algoTransformation_.push_back(algo);
120}
121
122/*!
123  Select algorithm of a domain correspoding to its transformation type and its position in each element
124  \param [in] elementPositionInGrid position of element in grid. e.g: a grid has 1 domain and 1 axis, then position of domain is 0 and position of axis is 1
125  \param [in] transType transformation type, for now we have zoom_domain, interpolate_domain
126  \param [in] transformationOrder position of the transformation in an element (an element can have several transformation)
127*/
128void CGridTransformation::selectDomainAlgo(int elementPositionInGrid, ETranformationType transType, int transformationOrder)
129{
130  std::vector<CDomain*> domainListDestP = gridDestination_->getDomains();
131  std::vector<CDomain*> domainListSrcP = gridSource_->getDomains();
132
133  int domainIndex =  elementPositionInGridDst2DomainPosition_[elementPositionInGrid];
134  CDomain::TransMapTypes trans = domainListDestP[domainIndex]->getAllTransformations();
135  CDomain::TransMapTypes::const_iterator it = trans.begin();
136
137  for (int i = 0; i < transformationOrder; ++i, ++it) {}  // Find the correct transformation
138
139  CZoomDomain* zoomDomain = 0;
140  CInterpolateDomain* interpFileDomain = 0;
141  CGenericAlgorithmTransformation* algo = 0;
142  switch (transType)
143  {
144    case TRANS_INTERPOLATE_DOMAIN:
145      interpFileDomain = dynamic_cast<CInterpolateDomain*> (it->second);
146      algo = new CDomainAlgorithmInterpolate(domainListDestP[domainIndex], domainListSrcP[domainIndex],interpFileDomain);
147      break;
148    case TRANS_ZOOM_DOMAIN:
149      zoomDomain = dynamic_cast<CZoomDomain*> (it->second);
150      algo = new CDomainAlgorithmZoom(domainListDestP[domainIndex], domainListSrcP[domainIndex], zoomDomain);
151      break;
152    default:
153      break;
154  }
155  algoTransformation_.push_back(algo);
156}
157
158/*!
159  If there are more than one transformation, a new temporary grid will be created and it will play the role of grid destination.
160This new created one keeps a pointer to the real transformed element of grid destination and generate new copies of other elements from grid source.
161  \param [in] elementPositionInGrid position of element in grid
162  \param [in] transType transformation type
163*/
164void CGridTransformation::setUpGridDestination(int elementPositionInGrid, ETranformationType transType, int nbTransformation)
165{
166  if (isSpecialTransformation(transType)) return;
167
168  if (!tempGridDests_.empty() && (getNbAlgo() == tempGridDests_.size()))
169  {
170    tempGridDests_.resize(0);
171  }
172
173  std::vector<CScalar*> scalarListDestP = gridDestination_->getScalars();
174  std::vector<CScalar*> scalarListSrcP = gridSource_->getScalars(), scalarDst;
175
176  std::vector<CAxis*> axisListDestP = gridDestination_->getAxis();
177  std::vector<CAxis*> axisListSrcP = gridSource_->getAxis(), axisDst;
178
179  std::vector<CDomain*> domListDestP = gridDestination_->getDomains();
180  std::vector<CDomain*> domListSrcP = gridSource_->getDomains(), domainDst;
181
182  CArray<int,1> axisDomainOrderSrc = gridSource_->axis_domain_order;
183  CArray<int,1> axisDomainOrderDst = gridDestination_->axis_domain_order;
184
185  int scalarIndex = -1, axisIndex = -1, domainIndex = -1;
186  switch (transType)
187  {
188    case TRANS_INTERPOLATE_DOMAIN:
189    case TRANS_ZOOM_DOMAIN:
190      domainIndex = elementPositionInGridDst2DomainPosition_[elementPositionInGrid];
191      break;
192
193    case TRANS_INTERPOLATE_AXIS:
194    case TRANS_ZOOM_AXIS:
195    case TRANS_INVERSE_AXIS:
196    case TRANS_REDUCE_DOMAIN_TO_AXIS:
197    case TRANS_EXTRACT_DOMAIN_TO_AXIS:
198      axisIndex =  elementPositionInGridDst2AxisPosition_[elementPositionInGrid];
199      break;
200
201    case TRANS_REDUCE_AXIS_TO_SCALAR:
202      scalarIndex = elementPositionInGridDst2ScalarPosition_[elementPositionInGrid];
203      break;
204    default:
205      break;
206  }
207
208  for (int idx = 0; idx < axisDomainOrderDst.numElements(); ++idx)
209  {
210    int dimElementDst = axisDomainOrderDst(idx);
211    if (2 == dimElementDst)
212    {
213      if (elementPositionInGrid == idx)
214        domainDst.push_back(domListDestP[domainIndex]);
215      else
216        domainDst.push_back(domListSrcP[elementPositionInGridSrc2DomainPosition_[elementPositionInGrid]]);
217    }
218    else if (1 == dimElementDst)
219    {
220      if (elementPositionInGrid == idx)
221        axisDst.push_back(axisListDestP[axisIndex]);
222      else
223        axisDst.push_back(axisListSrcP[elementPositionInGridSrc2AxisPosition_[elementPositionInGrid]]);
224    }
225    else
226    {
227      if (elementPositionInGrid == idx)
228        scalarDst.push_back(scalarListDestP[scalarIndex]);
229      else
230        scalarDst.push_back(scalarListSrcP[elementPositionInGridSrc2ScalarPosition_[elementPositionInGrid]]);
231    }
232  }
233
234  tmpGridDestination_ = CGrid::createGrid(domainDst, axisDst, scalarDst, gridDestination_->axis_domain_order);
235  tmpGridDestination_->computeGridGlobalDimension(domainDst, axisDst, scalarDst, gridDestination_->axis_domain_order);
236  tempGridDests_.push_back(tmpGridDestination_);
237}
238
239/*!
240  Assign the current grid destination to the grid source in the new transformation.
241The current grid destination plays the role of grid source in next transformation (if any).
242Only element on which the transformation is performed is modified
243  \param [in] elementPositionInGrid position of element in grid
244  \param [in] transType transformation type
245*/
246void CGridTransformation::setUpGridSource(int elementPositionInGrid, ETranformationType transType, int nbTransformation)
247{
248  if (!tempGridSrcs_.empty() && (getNbAlgo()-1) == tempGridSrcs_.size())
249  {
250    tempGridSrcs_.resize(0);
251  }
252
253  std::vector<CScalar*> scalarListDestP = tmpGridDestination_->getScalars();
254  std::vector<CScalar*> scalarListSrcP = gridSource_->getScalars(), scalarSrc;
255
256  std::vector<CAxis*> axisListDestP = tmpGridDestination_->getAxis();
257  std::vector<CAxis*> axisListSrcP = gridSource_->getAxis(), axisSrc;
258
259  std::vector<CDomain*> domListDestP = tmpGridDestination_->getDomains();
260  std::vector<CDomain*> domListSrcP = gridSource_->getDomains(), domainSrc;
261
262  CArray<int,1> axisDomainOrderDst = gridDestination_->axis_domain_order;
263
264  int axisIndex = -1, domainIndex = -1, scalarIndex = -1;
265  int axisListIndex = 0, domainListIndex = 0, scalarListIndex = 0;
266  switch (transType)
267  {
268    case TRANS_INTERPOLATE_DOMAIN:
269    case TRANS_ZOOM_DOMAIN:
270      domainIndex = elementPositionInGridDst2DomainPosition_[elementPositionInGrid];
271      break;
272
273    case TRANS_INTERPOLATE_AXIS:
274    case TRANS_ZOOM_AXIS:
275    case TRANS_INVERSE_AXIS:
276    case TRANS_REDUCE_DOMAIN_TO_AXIS:
277    case TRANS_EXTRACT_DOMAIN_TO_AXIS:
278      axisIndex =  elementPositionInGridDst2AxisPosition_[elementPositionInGrid];
279      break;
280
281    case TRANS_REDUCE_AXIS_TO_SCALAR:
282      scalarIndex = elementPositionInGridDst2ScalarPosition_[elementPositionInGrid];
283      break;
284    default:
285      break;
286  }
287
288  for (int idx = 0; idx < axisDomainOrderDst.numElements(); ++idx)
289  {
290    int dimElementDst = axisDomainOrderDst(idx);
291    if (2 == dimElementDst)
292    {
293      if (elementPositionInGrid == idx)
294        domainSrc.push_back(domListDestP[domainIndex]);
295      else
296      {
297        CDomain* domain = CDomain::createDomain();
298        domain->domain_ref.setValue(domListDestP[domainListIndex]->getId());
299        domain->solveRefInheritance(true);
300        domain->checkAttributesOnClient();
301        domainSrc.push_back(domain);
302      }
303      ++domainListIndex;
304    }
305    else if (1 == dimElementDst)
306    {
307      if (elementPositionInGrid == idx)
308        axisSrc.push_back(axisListDestP[axisIndex]);
309      else
310      {
311        CAxis* axis = CAxis::createAxis();
312        axis->axis_ref.setValue(axisListDestP[axisListIndex]->getId());
313        axis->solveRefInheritance(true);
314        axis->checkAttributesOnClient();
315        axisSrc.push_back(axis);
316      }
317      ++axisListIndex;
318    }
319    else
320    {
321      if (elementPositionInGrid == idx)
322        scalarSrc.push_back(scalarListDestP[scalarIndex]);
323      else
324      {
325        CScalar* scalar = CScalar::createScalar();
326        scalar->scalar_ref.setValue(scalarListDestP[scalarListIndex]->getId());
327        scalar->solveRefInheritance(true);
328        scalar->checkAttributesOnClient();
329        scalarSrc.push_back(scalar);
330      }
331      ++scalarListIndex;
332    }
333  }
334
335  gridSource_ = CGrid::createGrid(domainSrc, axisSrc, scalarSrc, tmpGridDestination_->axis_domain_order);
336  gridSource_->computeGridGlobalDimension(domainSrc, axisSrc, scalarSrc, tmpGridDestination_->axis_domain_order);
337
338  tempGridSrcs_.push_back(gridSource_);
339}
340
341/*!
342  Perform all transformations
343  For each transformation, there are some things to do:
344  -) Chose the correct algorithm by transformation type and position of element
345  -) Calculate the mapping of global index between the current grid source and grid destination
346  -) Calculate the mapping of global index between current grid DESTINATION and grid SOURCE
347  -) Make current grid destination become grid source in the next transformation
348*/
349void CGridTransformation::computeAll(const std::vector<CArray<double,1>* >& dataAuxInputs, Time timeStamp)
350{
351  if (nbNormalAlgos_ < 1) return;
352  if (!auxInputs_.empty() && !dynamicalTransformation_) { dynamicalTransformation_ = true; return; }
353  if (dynamicalTransformation_)
354  {
355    if (timeStamp_.insert(timeStamp).second)   //Reset map
356    {
357      std::list<SendingIndexGridSourceMap>().swap(localIndexToSendFromGridSource_);
358      std::list<RecvIndexGridDestinationMap>().swap(localIndexToReceiveOnGridDest_);
359      std::list<size_t>().swap(nbLocalIndexOnGridDest_);
360      std::list<std::vector<bool> >().swap(localMaskOnGridDest_);
361    }
362    else
363      return;
364  }
365
366  CContext* context = CContext::getCurrent();
367  CContextClient* client = context->client;
368
369  ListAlgoType::const_iterator itb = listAlgos_.begin(),
370                               ite = listAlgos_.end(), it;
371
372  CGenericAlgorithmTransformation* algo = 0;
373  int nbAgloTransformation = 0; // Only count for executed transformation. Generate domain is a special one, not executed in the list
374  for (it = itb; it != ite; ++it)
375  {
376    int elementPositionInGrid = it->first;
377    ETranformationType transType = (it->second).first;
378    int transformationOrder = (it->second).second;
379    SourceDestinationIndexMap globaIndexWeightFromSrcToDst;
380
381    // Create a temporary grid destination which contains transformed element of grid destination and
382    // non-transformed elements fo grid source
383    setUpGridDestination(elementPositionInGrid, transType, nbAgloTransformation);
384
385    // First of all, select an algorithm
386    if (!dynamicalTransformation_ || (algoTransformation_.size() < listAlgos_.size()))
387    {
388      selectAlgo(elementPositionInGrid, transType, transformationOrder, algoTypes_[std::distance(itb, it)]);
389      algo = algoTransformation_.back();
390    }
391    else
392      algo = algoTransformation_[std::distance(itb, it)];
393
394    if (0 != algo) // Only registered transformation can be executed
395    {
396      algo->computeIndexSourceMapping(dataAuxInputs);
397
398      // ComputeTransformation of global index of each element
399      int elementPosition = it->first;
400      algo->computeGlobalSourceIndex(elementPosition,
401                                     gridSource_,
402                                     tmpGridDestination_,
403                                     globaIndexWeightFromSrcToDst);
404
405      // Compute transformation of global indexes among grids
406      computeTransformationMapping(globaIndexWeightFromSrcToDst);
407
408      if (1 < nbNormalAlgos_)
409      {
410        // Now grid destination becomes grid source in a new transformation
411        if (nbAgloTransformation != (nbNormalAlgos_-1)) setUpGridSource(elementPositionInGrid, transType, nbAgloTransformation);
412      }
413      ++nbAgloTransformation;
414    }
415  }
416}
417
418/*!
419  Compute exchange index between grid source and grid destination
420  \param [in] globalIndexWeightFromDestToSource global index mapping between grid destination and grid source
421*/
422void CGridTransformation::computeTransformationMapping(const SourceDestinationIndexMap& globaIndexWeightFromSrcToDst)
423{
424  CContext* context = CContext::getCurrent();
425  CContextClient* client = context->client;
426  int nbClient = client->clientSize;
427  int clientRank = client->clientRank;
428
429  // Recalculate the distribution of grid destination
430  CDistributionClient distributionClientDest(client->clientRank, tmpGridDestination_);
431  CDistributionClient::GlobalLocalDataMap& globalLocalIndexGridDestSendToServer = distributionClientDest.getGlobalLocalDataSendToServer();
432
433  // Update number of local index on each transformation
434  size_t nbLocalIndex = globalLocalIndexGridDestSendToServer.size();
435  nbLocalIndexOnGridDest_.push_back(nbLocalIndex);
436  localMaskOnGridDest_.push_back(std::vector<bool>());
437  std::vector<bool>& tmpMask = localMaskOnGridDest_.back();
438  tmpMask.resize(nbLocalIndex,false);
439
440  // Find out number of index sent from grid source and number of index received on grid destination
441  SourceDestinationIndexMap::const_iterator itbIndex = globaIndexWeightFromSrcToDst.begin(),
442                                            iteIndex = globaIndexWeightFromSrcToDst.end(), itIndex;
443  typedef boost::unordered_map<size_t, std::vector<std::pair<size_t,double> > > SendIndexMap;
444  std::map<int,int> sendRankSizeMap,recvRankSizeMap;
445  int connectedClient = globaIndexWeightFromSrcToDst.size();
446  int* recvCount=new int[nbClient];
447  int* displ=new int[nbClient];
448  int* sendRankBuff=new int[connectedClient];
449  int* sendSizeBuff=new int[connectedClient];
450  int n = 0;
451  for (itIndex = itbIndex; itIndex != iteIndex; ++itIndex, ++n)
452  {
453    sendRankBuff[n] = itIndex->first;
454    const SendIndexMap& sendIndexMap = itIndex->second;
455    SendIndexMap::const_iterator itbSend = sendIndexMap.begin(), iteSend = sendIndexMap.end(), itSend;
456    int sendSize = 0;
457    for (itSend = itbSend; itSend != iteSend; ++itSend)
458    {
459      sendSize += itSend->second.size();
460    }
461    sendSizeBuff[n] = sendSize;
462    sendRankSizeMap[itIndex->first] = sendSize;
463  }
464  MPI_Allgather(&connectedClient,1,MPI_INT,recvCount,1,MPI_INT,client->intraComm);
465
466  displ[0]=0 ;
467  for(int n=1;n<nbClient;n++) displ[n]=displ[n-1]+recvCount[n-1];
468  int recvSize=displ[nbClient-1]+recvCount[nbClient-1];
469  int* recvRankBuff=new int[recvSize];
470  int* recvSizeBuff=new int[recvSize];
471  MPI_Allgatherv(sendRankBuff,connectedClient,MPI_INT,recvRankBuff,recvCount,displ,MPI_INT,client->intraComm);
472  MPI_Allgatherv(sendSizeBuff,connectedClient,MPI_INT,recvSizeBuff,recvCount,displ,MPI_INT,client->intraComm);
473  for (int i = 0; i < nbClient; ++i)
474  {
475    int currentPos = displ[i];
476    for (int j = 0; j < recvCount[i]; ++j)
477      if (recvRankBuff[currentPos+j] == clientRank)
478      {
479        recvRankSizeMap[i] = recvSizeBuff[currentPos+j];
480      }
481  }
482
483  // Sending global index of grid source to corresponding process as well as the corresponding mask
484  std::vector<MPI_Request> requests;
485  std::vector<MPI_Status> status;
486  boost::unordered_map<int, unsigned char* > recvMaskDst;
487  boost::unordered_map<int, unsigned long* > recvGlobalIndexSrc;
488  for (std::map<int,int>::const_iterator itRecv = recvRankSizeMap.begin(); itRecv != recvRankSizeMap.end(); ++itRecv)
489  {
490    int recvRank = itRecv->first;
491    int recvSize = itRecv->second;
492    recvMaskDst[recvRank] = new unsigned char [recvSize];
493    recvGlobalIndexSrc[recvRank] = new unsigned long [recvSize];
494
495    requests.push_back(MPI_Request());
496    MPI_Irecv(recvGlobalIndexSrc[recvRank], recvSize, MPI_UNSIGNED_LONG, recvRank, 46, client->intraComm, &requests.back());
497    requests.push_back(MPI_Request());
498    MPI_Irecv(recvMaskDst[recvRank], recvSize, MPI_UNSIGNED_CHAR, recvRank, 47, client->intraComm, &requests.back());
499  }
500
501  boost::unordered_map<int, CArray<size_t,1> > globalIndexDst;
502  boost::unordered_map<int, CArray<double,1> > weightDst;
503  boost::unordered_map<int, unsigned char* > sendMaskDst;
504  boost::unordered_map<int, unsigned long* > sendGlobalIndexSrc;
505  for (itIndex = itbIndex; itIndex != iteIndex; ++itIndex)
506  {
507    int sendRank = itIndex->first;
508    int sendSize = sendRankSizeMap[sendRank];
509    const SendIndexMap& sendIndexMap = itIndex->second;
510    SendIndexMap::const_iterator itbSend = sendIndexMap.begin(), iteSend = sendIndexMap.end(), itSend;
511    globalIndexDst[sendRank].resize(sendSize);
512    weightDst[sendRank].resize(sendSize);
513    sendMaskDst[sendRank] = new unsigned char [sendSize];
514    sendGlobalIndexSrc[sendRank] = new unsigned long [sendSize];
515    int countIndex = 0;
516    for (itSend = itbSend; itSend != iteSend; ++itSend)
517    {
518      const std::vector<std::pair<size_t,double> >& dstWeight = itSend->second;
519      for (int idx = 0; idx < dstWeight.size(); ++idx)
520      {
521        globalIndexDst[sendRank](countIndex) = dstWeight[idx].first;
522        weightDst[sendRank](countIndex) = dstWeight[idx].second;
523        if (0 < globalLocalIndexGridDestSendToServer.count(dstWeight[idx].first))
524          sendMaskDst[sendRank][countIndex] = 1;
525        else
526          sendMaskDst[sendRank][countIndex] = 0;
527        sendGlobalIndexSrc[sendRank][countIndex] = itSend->first;
528        ++countIndex;
529      }
530    }
531
532    // Send global index source and mask
533    requests.push_back(MPI_Request());
534    MPI_Isend(sendGlobalIndexSrc[sendRank], sendSize, MPI_UNSIGNED_LONG, sendRank, 46, client->intraComm, &requests.back());
535    requests.push_back(MPI_Request());
536    MPI_Isend(sendMaskDst[sendRank], sendSize, MPI_UNSIGNED_CHAR, sendRank, 47, client->intraComm, &requests.back());
537  }
538
539  status.resize(requests.size());
540  MPI_Waitall(requests.size(), &requests[0], &status[0]);
541
542  // Okie, now use the mask to identify which index source we need to send, then also signal the destination which masked index we will return
543  std::vector<MPI_Request>().swap(requests);
544  std::vector<MPI_Status>().swap(status);
545  // Okie, on destination side, we will wait for information of masked index of source
546  for (std::map<int,int>::const_iterator itSend = sendRankSizeMap.begin(); itSend != sendRankSizeMap.end(); ++itSend)
547  {
548    int recvRank = itSend->first;
549    int recvSize = itSend->second;
550
551    requests.push_back(MPI_Request());
552    MPI_Irecv(sendMaskDst[recvRank], recvSize, MPI_UNSIGNED_CHAR, recvRank, 48, client->intraComm, &requests.back());
553  }
554
555  // Ok, now we fill in local index of grid source (we even count for masked index)
556  CDistributionClient distributionClientSrc(client->clientRank, gridSource_);
557  CDistributionClient::GlobalLocalDataMap& globalLocalIndexGridSrcSendToServer = distributionClientSrc.getGlobalLocalDataSendToServer();
558  localIndexToSendFromGridSource_.push_back(SendingIndexGridSourceMap());
559  SendingIndexGridSourceMap& tmpSend = localIndexToSendFromGridSource_.back();
560  for (std::map<int,int>::const_iterator itRecv = recvRankSizeMap.begin(); itRecv != recvRankSizeMap.end(); ++itRecv)
561  {
562    int recvRank = itRecv->first;
563    int recvSize = itRecv->second;
564    unsigned char* recvMask = recvMaskDst[recvRank];
565    unsigned long* recvIndexSrc = recvGlobalIndexSrc[recvRank];
566    int realSendSize = 0;
567    for (int idx = 0; idx < recvSize; ++idx)
568    {
569      if (0 != (*(recvMask+idx))) // OKie, now we have a demand from non-masked index destination
570        if (0 < globalLocalIndexGridSrcSendToServer.count(*(recvIndexSrc+idx))) // check whether index source is masked
571         ++realSendSize;
572        else // inform the destination that this index is masked
573         *(recvMask+idx) = 0;
574    }
575
576    tmpSend[recvRank].resize(realSendSize);
577    realSendSize = 0;
578    for (int idx = 0; idx < recvSize; ++idx)
579    {
580      if (0 != (*(recvMask+idx))) // OKie, now we have a demand from non-masked index destination
581      {
582        tmpSend[recvRank](realSendSize) = globalLocalIndexGridSrcSendToServer[*(recvIndexSrc+idx)];
583         ++realSendSize;
584      }
585    }
586
587    // Okie, now inform the destination which source index are masked
588    requests.push_back(MPI_Request());
589    MPI_Isend(recvMaskDst[recvRank], recvSize, MPI_UNSIGNED_CHAR, recvRank, 48, client->intraComm, &requests.back());
590  }
591  status.resize(requests.size());
592  MPI_Waitall(requests.size(), &requests[0], &status[0]);
593
594  // Cool, now we can fill in local index of grid destination (counted for masked index)
595  localIndexToReceiveOnGridDest_.push_back(RecvIndexGridDestinationMap());
596  RecvIndexGridDestinationMap& recvTmp = localIndexToReceiveOnGridDest_.back();
597  for (std::map<int,int>::const_iterator itSend = sendRankSizeMap.begin(); itSend != sendRankSizeMap.end(); ++itSend)
598  {
599    int recvRank = itSend->first;
600    int recvSize = itSend->second;
601    unsigned char* recvMask = sendMaskDst[recvRank];
602
603    CArray<size_t,1>& recvIndexDst = globalIndexDst[recvRank];
604    CArray<double,1>& recvWeightDst = weightDst[recvRank];
605    int realRecvSize = 0;
606    for (int idx = 0; idx < recvSize; ++idx)
607    {
608      if (0 != *(recvMask+idx)) // OKie, now we have a non-masked index destination
609         ++realRecvSize;
610    }
611
612    int localIndexDst;
613    recvTmp[recvRank].resize(realRecvSize);
614    realRecvSize = 0;
615    for (int idx = 0; idx < recvSize; ++idx)
616    {
617      if (0 != *(recvMask+idx)) // OKie, now we have a demand from non-masked index destination
618      {
619        recvTmp[recvRank][realRecvSize].first = globalLocalIndexGridDestSendToServer[recvIndexDst(idx)];
620        recvTmp[recvRank][realRecvSize].second = recvWeightDst(idx);
621        tmpMask[globalLocalIndexGridDestSendToServer[recvIndexDst(idx)]] = true;
622         ++realRecvSize;
623      }
624    }
625  }
626
627  delete [] recvCount;
628  delete [] displ;
629  delete [] sendRankBuff;
630  delete [] recvRankBuff;
631  delete [] sendSizeBuff;
632  delete [] recvSizeBuff;
633
634  boost::unordered_map<int, unsigned char* >::const_iterator itChar;
635  for (itChar = sendMaskDst.begin(); itChar != sendMaskDst.end(); ++itChar)
636    delete [] itChar->second;
637  for (itChar = recvMaskDst.begin(); itChar != recvMaskDst.end(); ++itChar)
638    delete [] itChar->second;
639  boost::unordered_map<int, unsigned long* >::const_iterator itLong;
640  for (itLong = sendGlobalIndexSrc.begin(); itLong != sendGlobalIndexSrc.end(); ++itLong)
641    delete [] itLong->second;
642  for (itLong = recvGlobalIndexSrc.begin(); itLong != recvGlobalIndexSrc.end(); ++itLong)
643    delete [] itLong->second;
644
645}
646
647/*!
648  Local index of data which need sending from the grid source
649  \return local index of data
650*/
651const std::list<CGridTransformation::SendingIndexGridSourceMap>& CGridTransformation::getLocalIndexToSendFromGridSource() const
652{
653  return localIndexToSendFromGridSource_;
654}
655
656/*!
657  Local index of data which will be received on the grid destination
658  \return local index of data
659*/
660const std::list<CGridTransformation::RecvIndexGridDestinationMap>& CGridTransformation::getLocalIndexToReceiveOnGridDest() const
661{
662  return localIndexToReceiveOnGridDest_;
663}
664
665/*!
666 Number of index will be received on the grid destination
667  \return number of index of data
668*/
669const std::list<size_t>& CGridTransformation::getNbLocalIndexToReceiveOnGridDest() const
670{
671  return nbLocalIndexOnGridDest_;
672}
673
674/*!
675  Local mask of data which will be received on the grid destination
676  \return local mask of data
677*/
678const std::list<std::vector<bool> >& CGridTransformation::getLocalMaskIndexOnGridDest() const
679{
680  return localMaskOnGridDest_;
681}
682
683}
Note: See TracBrowser for help on using the repository browser.