C++ mpi module for stochmagnet_main Package
MPI_WorldToWorldMessage.h
1 #ifndef MPI_WorldToWorldMessage_H
2 #define MPI_WorldToWorldMessage_H
3 
4 
5 //super class header
6 #include "MPI_Object.h"
7 
8 //environment header
9 #include "MPI_Environment.h"
10 
11 //type header
12 #include "MPI_Type.h"
13 
14 
63 
64  //attributes
65 private :
66 
67  //UML inheritage classes
69  typedef MPI_Object SuperClass;
70 
71 public:
72  // CONSTRUCTORS
76 
77  }
78 
79  // DESTRUCTORS
82  virtual ~MPI_WorldToWorldMessage(void) {
83 
84  }
85 
86 
87 public:
88 
89  //MEMORY
90 
104  virtual tMemSize getMemorySize() const override {
105  return sizeof(*this)+getContentsMemorySize();
106  }
107 
111  inline static CORE_UniquePointer<SelfClass> New() {
112  CORE_UniquePointer<SelfClass> p(new SelfClass(),
114  return p;
115  };
116 
117 
118  //sender data methods
119  //====================
120 
121  //primary type data
122  //------------------
123 
132  template<typename T>
133  inline static tBoolean Send(const MPI_Environment& env,const tMPICoreId& srcCore,
134  T* data,const tMPICount& nData) {
135  return (MPI_Bcast(data,nData,MPI_Type::GetPrimaryType<T>(),srcCore,env.getWorld())==MPI_SUCCESS);
136  }
137 
147  template<typename T>
148  inline static tBoolean BCast(const MPI_Environment& env,const tMPICoreId& srcCore,
149  T* data,const tMPICount& nData,const tMPIType& dataType) {
150  return (MPI_Bcast(data,nData,dataType,srcCore,env.getWorld())==MPI_SUCCESS);
151  }
160  template<typename T>
161  inline static tBoolean BCast(const MPI_Environment& env,
162  const tMPICoreId& srcCore,
163  std::valarray<T>& data) {
164  return (MPI_Bcast(((data.size()==0)?null:&data[0]),data.size(),MPI_Type::GetPrimaryType<T>(),srcCore,env.getWorld())==MPI_SUCCESS);
165  }
166 
167 
177  template<typename T>
178  inline static tBoolean IBCast(const MPI_Environment& env,const tMPICoreId& srcCore,
179  T* data,const tMPICount& nData,
180  tMPIRequest& request) {
181  return (MPI_Ibcast(data,nData,MPI_Type::GetPrimaryType<T>(),srcCore,env.getWorld(),&request)==MPI_SUCCESS);
182  }
183 
194  template<typename T>
195  inline static tBoolean IBCast(const MPI_Environment& env,const tMPICoreId& srcCore,
196  T* data,const tMPICount& nData,const tMPIType& dataType,
197  tMPIRequest& request) {
198  return (MPI_Ibcast(data,nData,dataType,srcCore,env.getWorld(),request)==MPI_SUCCESS);
199  }
209  template<typename T>
210  inline static tBoolean IBCast(const MPI_Environment& env,
211  const tMPICoreId& srcCore,
212  std::valarray<T>& data,
213  tMPIRequest& request) {
214  return (MPI_Ibcast(((data.size()==0)?null:&data[0]),data.size(),MPI_Type::GetPrimaryType<T>(),srcCore,env.getWorld(),request)==MPI_SUCCESS);
215  }
216 
217  //diffuse data methods
218  //====================
219 
220  //muliple primary data type
221  //-------------------------
222 
234  template<typename T>
235  inline static tBoolean Scatter(const MPI_Environment& env,const tMPICoreId& srcCore,
236  const std::valarray<T>& srcData,
237  std::valarray<T>& dstData) {
238 
239  //verify the lenth to send if enough big
240  return (MPI_Scatter(((srcData.size()==0)?null:&srcData[0]),dstData.size(),MPI_Type::GetPrimaryType<T>(),
241  ((dstData.size()==0)?null:&dstData[0]),dstData.size(),MPI_Type::GetPrimaryType<T>(),
242  srcCore,env.getWorld())==MPI_SUCCESS);
243  return true;
244  }
257  template<typename T>
258  inline static tBoolean Scatter(const MPI_Environment& env,const tMPICoreId& srcCore,
259  const T* srcData,const tMPICount& length,const tMPIType& srcDataType,
260  T* dstData,const tMPICount& dstDataSize,const tMPIType& dstDataType) {
261 
262  return (MPI_Scatter(srcData,length,srcDataType,
263  dstData,dstDataSize,dstDataType,
264  srcCore,env.getWorld())==MPI_SUCCESS);
265 
266  }
279  template<typename T>
280  inline static tBoolean Scatter(const MPI_Environment& env,const tMPICoreId& srcCore,
281  const T* srcData,
282  const std::valarray<tMPICount>& dataLengthPerCore,
283  const std::valarray<tMPICount>& dataLengthPerCoreOffset,
284  T* dstData,const tMPICount& dstDataSize) {
285 
286  return (MPI_Scatterv(srcData,
287  ((dataLengthPerCore.size()==0)?null:&dataLengthPerCore[0]),
288  ((dataLengthPerCoreOffset.size()==0)?null:&dataLengthPerCoreOffset[0]),
289  MPI_Type::GetPrimaryType<T>(),
290  dstData,
291  dstDataSize,
292  MPI_Type::GetPrimaryType<T>(),
293  srcCore,env.getWorld())==MPI_SUCCESS);
294 
295  }
296 
305  template<typename T>
306  inline static tBoolean Gather(const MPI_Environment& env,const tMPICoreId& dstCore,
307  const T& srcData,
308  std::valarray<T>& dstData) {
309  return (MPI_Gather(&srcData,1,MPI_Type::GetPrimaryType<T>(),
310  ((dstData.size()==0)?null:&dstData[0]),1,MPI_Type::GetPrimaryType<T>(),
311  dstCore,env.getWorld())==MPI_SUCCESS);
312 
313  }
314 
323  template<typename T>
324  inline static tBoolean Gather(const MPI_Environment& env,const tMPICoreId& dstCore,
325  const std::valarray<T>& srcData,
326  std::valarray<T>& dstData) {
327  if (dstData.size()<srcData.size()*env.getCoresNumber()) return false;
328  return (MPI_Gather(((srcData.size()==0)?null:&srcData[0]),(tMPIIndex)srcData.size(),MPI_Type::GetPrimaryType<T>(),
329  ((dstData.size()==0)?null:&dstData[0]),(tMPIIndex)srcData.size(),MPI_Type::GetPrimaryType<T>(),
330  dstCore,env.getWorld())==MPI_SUCCESS);
331 
332  }
345  template<typename T>
346  inline static tBoolean Gather(const MPI_Environment& env,const tMPICoreId& dstCore,
347  const T* srcData,const tMPICount& srcDataSize,const tMPIType& srcDataType,
348  T* dstData,const tMPICount& dstDataSize,const tMPIType& dstDataType) {
349  return (MPI_Gather(srcData,srcDataSize,srcDataType,
350  dstData,dstDataSize,dstDataType,
351  dstCore,env.getWorld())==MPI_SUCCESS);
352  }
353 
354 
365  template<typename T>
366  inline static tBoolean Gather(const MPI_Environment& env,const tMPICoreId& dstCore,
367  const std::valarray<T>& srcData,
368  std::valarray<T>& dstData,
369  const std::valarray<tMPICount>& dataLengthPerCore,
370  const std::valarray<tMPICount>& dataIndexPerCore) {
371 
372  return (MPI_Gatherv(((srcData.size()==0)?null:&srcData[0]),(tMPIIndex)srcData.size(),MPI_Type::GetPrimaryType<T>(),
373  ((dstData.size()==0)?null:&dstData[0]),
374  ((dataLengthPerCore.size()==0)?null:&dataLengthPerCore[0]),
375  ((dataIndexPerCore.size()==0)?null:&dataIndexPerCore[0]),
376  MPI_Type::GetPrimaryType<T>(),
377  dstCore,env.getWorld())==MPI_SUCCESS);
378 
379  }
392  template<typename T>
393  inline static tBoolean Gather(const MPI_Environment& env,const tMPICoreId& dstCore,
394  const T* srcData,const tMPICount& srcDataSize,
395  T* dstData,
396  const std::valarray<tMPICount>& dataLengthPerCore,
397  const std::valarray<tMPICount>& dataIndexPerCore) {
398  return (MPI_Gatherv(srcData,srcDataSize,MPI_Type::GetPrimaryType<T>(),
399  dstData,
400  ((dataLengthPerCore.size()==0)?null:&dataLengthPerCore[0]),
401  ((dataIndexPerCore.size()==0)?null:&dataIndexPerCore[0]),
402  MPI_Type::GetPrimaryType<T>(),
403  dstCore,env.getWorld())==MPI_SUCCESS);
404  }
405 
413  template<typename T>
414  inline static tBoolean AllGather(const MPI_Environment& env,
415  const T& srcData,
416  std::valarray<T>& dstData) {
417  //std::cout<<"begin MPI_WorldToWorldMessage::AllGather(env,srcData,dstData) coreId="<<env.getCoreId()<<"\n";
418  return (MPI_Allgather(&srcData,1,MPI_Type::GetPrimaryType<T>(),
419  ((dstData.size()==0)?null:&dstData[0]),1,MPI_Type::GetPrimaryType<T>(),
420  env.getWorld())==MPI_SUCCESS);
421  //std::cout<<"end MPI_WorldToWorldMessage::AllGather(env,srcData,dstData) coreId="<<env.getCoreId()<<"\n";
422  }
423 
431  template<typename T>
432  inline static tBoolean AllGather(const MPI_Environment& env,
433  const std::valarray<T>& srcData,
434  std::valarray<T>& dstData) {
435  if (dstData.size()!=srcData.size()*env.getCoresNumber()) return false;
436  return (MPI_Allgather(((srcData.size()==0)?null:&srcData[0]),(tMPIIndex)srcData.size(),MPI_Type::GetPrimaryType<T>(),
437  ((dstData.size()==0)?null:&dstData[0]),(tMPIIndex)srcData.size(),MPI_Type::GetPrimaryType<T>(),
438  env.getWorld())==MPI_SUCCESS);
439 
440  }
449  template<typename T>
450  inline static tBoolean AllGather(const MPI_Environment& env,
451  const T* srcData,const tIndex& nSrcData,
452  std::valarray<T>& dstData) {
453  //std::cout<<"begin MPI_WorldToWorldMessage::AllGather(env,srcData,nSrcData,dstData) coreId="<<env.getCoreId()<<"\n";
454  if (dstData.size()!=srcData.size()*env.getCoresNumber()) return false;
455  return (MPI_Allgather(srcData,nSrcData,MPI_Type::GetPrimaryType<T>(),
456  ((dstData.size()==0)?null:&dstData[0]),(tMPIIndex)nSrcData,MPI_Type::GetPrimaryType<T>(),
457  env.getWorld())==MPI_SUCCESS);
458  //std::cout<<"end MPI_WorldToWorldMessage::AllGather(env,srcData,nSrcData,dstData) coreId="<<env.getCoreId()<<"\n";
459 
460  }
461 
467  inline static int BuildOffset(const std::valarray<tMPIIndex>& dataLength,
468  std::valarray<tMPIIndex>& dataOffset) {
469  if (dataOffset.size()<=dataLength.size()) dataOffset.resize(dataLength.size()+1);
470 
471  //total number of received elements
472  tMPIIndex n=0;
473  //iterator on offset
474  tMPIIndex * iDataOffset=&dataOffset[0];
475  for(const auto& dLength:dataLength) {
476  (*iDataOffset)=n;
477  n+=dLength;
478  iDataOffset++;
479  }
480  (*iDataOffset)=n;
481  return n;
482  }
483 
493  template<typename T>
494  inline static tBoolean AllGather(const MPI_Environment& env,
495  const std::valarray<T>& srcData,
496  T* dstData,
497  const std::valarray<tMPIIndex>& srcDataLength,
498  const std::valarray<tMPIIndex>& displacement) {
499 
500  return AllGather(env,((srcData.size()==0)?null:&srcData[0]),(tMPIIndex)srcData.size(),srcDataLength,displacement);
501  }
502 
512  template<typename T>
513  inline static tBoolean AllGather(const MPI_Environment& env,
514  const T* srcData,const tMPIIndex& nSrcData,
515  T* dstData,
516  const std::valarray<tMPIIndex>& srcDataLength,
517  const std::valarray<tMPIIndex>& displacement) {
518  //std::cout<<"begin MPI_WorldToWorldMessage::AllGatherV(env,srcData,nSrcData,length,offset) coreId="<<env.getCoreId()<<"\n";
519  return (MPI_Allgatherv(srcData,(tMPIIndex)nSrcData,
520  MPI_Type::GetPrimaryType<T>(),
521  dstData,
522  ((srcDataLength.size()==0)?null:&srcDataLength[0]),
523  ((displacement.size()==0)?null:&displacement[0]),
524  MPI_Type::GetPrimaryType<T>(),
525  env.getWorld())==MPI_SUCCESS);
526  //std::cout<<"end MPI_WorldToWorldMessage::AllGatherV(env,srcData,nSrcData,length,offset) coreId="<<env.getCoreId()<<"\n";
527  }
528 
538  template<typename T>
539  inline static tBoolean AllGather(const MPI_Environment& env,
540  const std::valarray<T>& srcData,
541  std::valarray<T>& dstData,
542  const std::valarray<tMPIIndex>& srcDataLength,
543  const std::valarray<tMPIIndex>& displacement) {
544  //std::cout<<"begin MPI_WorldToWorldMessage::AllGatherV(env,srcData,nSrcData,length,offset) coreId="<<env.getCoreId()<<"\n";
545 
546  return AllGather(env,
547  ((srcData.size()==0)?null:&srcData[0]),srcData.size(),
548  ((dstData.size()==0)?null:&dstData[0]),srcDataLength,displacement);
549  //std::cout<<"end MPI_WorldToWorldMessage::AllGatherV(env,srcData,nSrcData,length,offset) coreId="<<env.getCoreId()<<"\n";
550  }
551 
552 
553 
554 
574  template<typename T>
575  inline static tBoolean Transpose(const MPI_Environment& env,
576  const std::valarray<T>& srcData,const tMPICount& packSize,
577  std::valarray<T>& dstData) {
578  return (MPI_Alltoall(((srcData.size()==0)?null:&srcData[0]),packSize,MPI_Type::GetPrimaryType<T>(),
579  ((dstData.size()==0)?null:&dstData[0]),packSize,MPI_Type::GetPrimaryType<T>(),
580  env.getWorld())==MPI_SUCCESS);
581  }
582 
583 
584  //reduction methods
585  //=================
586  //primary data reduction
587  //-----------------------
588 
598  template<typename T>
599  inline static tBoolean Reduce(const MPI_Environment& env,const tMPICoreId& dstCore,
600  const T& srcData,
601  T& dstData,const tMPIOperation& op) {
602 
603  return (MPI_Reduce(&srcData,&dstData,1,MPI_Type::GetPrimaryType<T>(),op,
604  dstCore,env.getWorld())==MPI_SUCCESS);
605  }
614  template<typename T>
615  inline static tBoolean Reduce(const MPI_Environment& env,const tMPICoreId& dstCore,
616  T& data,const tMPIOperation& op) {
617  tMPIInteger ret;
618  if (env.getCoreId()==dstCore)
619  ret=MPI_Reduce(MPI_IN_PLACE,&data,1,MPI_Type::GetPrimaryType<T>(),op,
620  dstCore,env.getWorld());
621  else
622  ret=MPI_Reduce(&data,&data,1,MPI_Type::GetPrimaryType<T>(),op,
623  dstCore,env.getWorld());
624  return (ret==MPI_SUCCESS);
625  }
626 
627  //multiple reduction
628  //===================
638  template<typename T>
639  inline static tBoolean Reduce(const MPI_Environment& env,const tMPICoreId& dstCore,
640  const std::valarray<T>& srcData,
641  std::valarray<T>& dstData,const tMPIOperation& op) {
642  if (srcData.size()<dstData.size()) {
643  throw CORE_Exception("mpi/core","MPI_WorldToWorldMessage::Reduce()",
644  "incompatible size of array");
645  }
646  return (MPI_Reduce(((srcData.size()==0)?null:&srcData[0]),
647  ((dstData.size()==0)?null:&dstData[0]),(tMPIIndex)dstData.size(),MPI_Type::GetPrimaryType<T>(),op,
648  dstCore,env.getWorld())==MPI_SUCCESS);
649  }
658  template<typename T>
659  inline static tBoolean Reduce(const MPI_Environment& env,const tMPICoreId& dstCore,
660  std::valarray<T>& data,const tMPIOperation& op) {
661 
662  tMPIInteger ret;
663  if (env.getCoreId()==dstCore)
664  ret=MPI_Reduce(MPI_IN_PLACE,(data.size()==0)?null:&data[0],(tMPIIndex)data.size(),MPI_Type::GetPrimaryType<T>(),op,
665  dstCore,env.getWorld());
666  else
667  ret=MPI_Reduce((data.size()==0)?null:&data[0],(data.size()==0)?null:&data[0],(tMPIIndex)data.size(),MPI_Type::GetPrimaryType<T>(),op,
668  dstCore,env.getWorld());
669  return (ret=MPI_SUCCESS);
670  }
679  template<typename T>
680  inline static tBoolean Reduce(const MPI_Environment& env,const tMPICoreId& dstCore,
681  T* data,const tMPIIndex& nData,const tMPIOperation& op) {
682 
683  tMPIInteger ret;
684  if (env.getCoreId()==dstCore)
685  ret=MPI_Reduce(MPI_IN_PLACE,data,nData,MPI_Type::GetPrimaryType<T>(),op,
686  dstCore,env.getWorld());
687  else
688  ret=MPI_Reduce(data,data,nData,MPI_Type::GetPrimaryType<T>(),op,
689  dstCore,env.getWorld());
690  return (ret=MPI_SUCCESS);
691  }
692 
693 
694  //collective reduction
695  //====================
702  template<typename T>
703  inline static tBoolean AllReduce(const MPI_Environment& env,
704  T& data,const tMPIOperation& op) {
705  //std::cout<<"begin MPI_WorldToWorldMessage::AllReduce(env,data,op) coreId="<<env.getCoreId()<<"\n";
706  return (MPI_Allreduce(MPI_IN_PLACE,&data,1,MPI_Type::GetPrimaryType<T>(),op,
707  env.getWorld())==MPI_SUCCESS);
708  //std::cout<<"end MPI_WorldToWorldMessage::AllReduce(env,srcData,op) coreId="<<env.getCoreId()<<"\n";
709  }
718  template<typename T>
719  inline static tBoolean AllReduce(const MPI_Environment& env,
720  const T& srcData,T& dstData,const tMPIOperation& op) {
721  //std::cout<<"begin MPI_WorldToWorldMessage::AllReduce(env,srcData,dstData,op) coreId="<<env.getCoreId()<<"\n";
722  MPI_Allreduce(&srcData,&dstData,1,MPI_Type::GetPrimaryType<T>(),op,
723  env.getWorld());
724  //std::cout<<"end MPI_WorldToWorldMessage::AllReduce(env,srcData,dstData,op) coreId="<<env.getCoreId()<<"\n";
725  return true;
726  }
727 
728 
729 
730  //collective multiple operator
731  //============================
732 
740  template<typename T>
741  inline static tBoolean AllReduce(const MPI_Environment& env,
742  T* data,const tMPIIndex& nData,const tMPIOperation& op) {
743  //std::cout<<"begin MPI_WorldToWorldMessage::AllReduce(env,data,op) coreId="<<env.getCoreId()<<"\n";
744  return (MPI_Allreduce(MPI_IN_PLACE,data,nData,MPI_Type::GetPrimaryType<T>(),op,
745  env.getWorld())==MPI_SUCCESS);
746  //std::cout<<"end MPI_WorldToWorldMessage::AllReduce(env,srcData,op) coreId="<<env.getCoreId()<<"\n";
747  }
748 
757  template<typename T>
758  inline static tBoolean AllReduce(const MPI_Environment& env,
759  const std::valarray<T>& srcData,std::valarray<T>& dstData,const tMPIOperation& op) {
760  if (srcData.size()<dstData.size()) {
761  throw CORE_Exception("mpi/core","MPI_WorldToWorldMessage::AllReduce()",
762  "incompatible size of array");
763  }
764  return (MPI_Allreduce((srcData.size()==0)?null:&srcData[0],
765  (dstData.size()==0)?null:&dstData[0],(tMPIIndex)dstData.size(),MPI_Type::GetPrimaryType<T>(),op,
766  env.getWorld())==MPI_SUCCESS);
767 
768  }
776  template<typename T,size_t N>
777  inline static tBoolean AllReduce(const MPI_Environment& env,
778  const std::array<T,N>& srcData,std::array<T,N>& dstData,const tMPIOperation& op) {
779 
780 
781  return (MPI_Allreduce(srcData.data(),dstData.data(),(tMPIIndex)dstData.size(),MPI_Type::GetPrimaryType<T>(),op,
782  env.getWorld())==MPI_SUCCESS);
783 
784  }
785 
786 
793  template<typename T>
794  inline static tBoolean AllReduce(const MPI_Environment& env,
795  std::valarray<T>& data,const tMPIOperation& op) {
796 
797  return (MPI_Allreduce(MPI_IN_PLACE,
798  (data.size()==0)?null:&data[0],
799  (tMPIIndex)data.size(),MPI_Type::GetPrimaryType<T>(),op,
800  env.getWorld())==MPI_SUCCESS);
801 
802  }
810  template<typename T,size_t N>
811  inline static tBoolean AllReduce(const MPI_Environment& env,
812  std::array<T,N>& data,const tMPIOperation& op) {
813 
814  return (MPI_Allreduce(MPI_IN_PLACE,data.data(),(tMPIIndex)data.size(),MPI_Type::GetPrimaryType<T>(),op,
815  env.getWorld())==MPI_SUCCESS);
816  }
817 
818 
819 
820 
821 };
822 
823 
824 #endif
this class describes the exceptions raised for CORE package
Definition: CORE_Exception.h:17
class Free introduced for deleting a smart pointer
Definition: CORE_Object.h:113
virtual tMemSize getContentsMemorySize() const
return nthe memory size of the included associations
Definition: CORE_Object.h:278
This class is a Environment class to define MPI world.
Definition: MPI_Environment.h:36
const tMPICoreId & getCoreId() const
get the id of the current process of this environment
Definition: MPI_Environment.h:200
const tMPICoreId & getCoresNumber() const
get the number of cores of this environment of common environment
Definition: MPI_Environment.h:180
const tMPIComm & getWorld() const
get the world of the environment for reading
Definition: MPI_Environment.h:165
This class is a base class of E-MicromM core package.
Definition: MPI_Object.h:32
This class is a class to send / receive from all cores of a cores to all cores.
Definition: MPI_WorldToWorldMessage.h:62
static tBoolean Gather(const MPI_Environment &env, const tMPICoreId &dstCore, const T *srcData, const tMPICount &srcDataSize, T *dstData, const std::valarray< tMPICount > &dataLengthPerCore, const std::valarray< tMPICount > &dataIndexPerCore)
blocking gathering a data from all cores of enviroments dest core
Definition: MPI_WorldToWorldMessage.h:393
static tBoolean AllGather(const MPI_Environment &env, const T &srcData, std::valarray< T > &dstData)
blocking gathering a data from all cores of enviroments to all cores
Definition: MPI_WorldToWorldMessage.h:414
static tBoolean AllReduce(const MPI_Environment &env, T *data, const tMPIIndex &nData, const tMPIOperation &op)
compute a data from data of others cores and to data of core 0 and copy the values to all cores
Definition: MPI_WorldToWorldMessage.h:741
static tBoolean AllGather(const MPI_Environment &env, const std::valarray< T > &srcData, std::valarray< T > &dstData)
blocking gathering a data from all cores of enviroments to all cores
Definition: MPI_WorldToWorldMessage.h:432
static tBoolean AllReduce(const MPI_Environment &env, std::valarray< T > &data, const tMPIOperation &op)
compute a data from data of others cores and to data of core 0 and copy the values to all cores
Definition: MPI_WorldToWorldMessage.h:794
virtual tMemSize getMemorySize() const override
return the memory size of the class and the memory size of all its attributes/associations
Definition: MPI_WorldToWorldMessage.h:104
static tBoolean BCast(const MPI_Environment &env, const tMPICoreId &srcCore, T *data, const tMPICount &nData, const tMPIType &dataType)
blocking send a data with flag to all cores of enviroments except the core with srcCore id
Definition: MPI_WorldToWorldMessage.h:148
static tBoolean Gather(const MPI_Environment &env, const tMPICoreId &dstCore, const T *srcData, const tMPICount &srcDataSize, const tMPIType &srcDataType, T *dstData, const tMPICount &dstDataSize, const tMPIType &dstDataType)
blocking gathering a data from all cores of enviroments dest core
Definition: MPI_WorldToWorldMessage.h:346
static tBoolean Gather(const MPI_Environment &env, const tMPICoreId &dstCore, const T &srcData, std::valarray< T > &dstData)
blocking gathering a data from all cores of enviroments to the dest core
Definition: MPI_WorldToWorldMessage.h:306
MPI_WorldToWorldMessage()
create
Definition: MPI_WorldToWorldMessage.h:75
virtual ~MPI_WorldToWorldMessage(void)
destroy
Definition: MPI_WorldToWorldMessage.h:82
static tBoolean Scatter(const MPI_Environment &env, const tMPICoreId &srcCore, const std::valarray< T > &srcData, std::valarray< T > &dstData)
blocking scattering a data to all cores of environment with srcCore id
Definition: MPI_WorldToWorldMessage.h:235
static tBoolean AllGather(const MPI_Environment &env, const T *srcData, const tMPIIndex &nSrcData, T *dstData, const std::valarray< tMPIIndex > &srcDataLength, const std::valarray< tMPIIndex > &displacement)
blocking gathering a variable data lengthfrom all cores of enviroments to all cores
Definition: MPI_WorldToWorldMessage.h:513
static tBoolean IBCast(const MPI_Environment &env, const tMPICoreId &srcCore, T *data, const tMPICount &nData, tMPIRequest &request)
blocking send a data with flag to all cores of enviroments except the core with srcCore id
Definition: MPI_WorldToWorldMessage.h:178
static tBoolean BCast(const MPI_Environment &env, const tMPICoreId &srcCore, std::valarray< T > &data)
blocking send a data with flag to all cores of envirments except the core with srcCore id
Definition: MPI_WorldToWorldMessage.h:161
static tBoolean Reduce(const MPI_Environment &env, const tMPICoreId &dstCore, std::valarray< T > &data, const tMPIOperation &op)
compute a data from data of others cores to data of dstCore
Definition: MPI_WorldToWorldMessage.h:659
static tBoolean Reduce(const MPI_Environment &env, const tMPICoreId &dstCore, T *data, const tMPIIndex &nData, const tMPIOperation &op)
compute a data from data of others cores to data of dstCore
Definition: MPI_WorldToWorldMessage.h:680
static tBoolean IBCast(const MPI_Environment &env, const tMPICoreId &srcCore, std::valarray< T > &data, tMPIRequest &request)
blocking send a data with flag to all cores of envirments except the core with srcCore id
Definition: MPI_WorldToWorldMessage.h:210
static tBoolean Transpose(const MPI_Environment &env, const std::valarray< T > &srcData, const tMPICount &packSize, std::valarray< T > &dstData)
blocking transpose a data of diffrent size from all cores of enviroments to dst core with pack Size
Definition: MPI_WorldToWorldMessage.h:575
static tBoolean Reduce(const MPI_Environment &env, const tMPICoreId &dstCore, const std::valarray< T > &srcData, std::valarray< T > &dstData, const tMPIOperation &op)
compute a data from data of others cores to data of dstCore
Definition: MPI_WorldToWorldMessage.h:639
static tBoolean AllGather(const MPI_Environment &env, const T *srcData, const tIndex &nSrcData, std::valarray< T > &dstData)
blocking gathering a data from all cores of enviroments to all cores
Definition: MPI_WorldToWorldMessage.h:450
static tBoolean AllGather(const MPI_Environment &env, const std::valarray< T > &srcData, std::valarray< T > &dstData, const std::valarray< tMPIIndex > &srcDataLength, const std::valarray< tMPIIndex > &displacement)
blocking gathering a variable data lengthfrom all cores of enviroments to all cores
Definition: MPI_WorldToWorldMessage.h:539
static tBoolean AllReduce(const MPI_Environment &env, const T &srcData, T &dstData, const tMPIOperation &op)
compute a data from data of others cores and to data of core 0 and copy the values to all cores
Definition: MPI_WorldToWorldMessage.h:719
static tBoolean Scatter(const MPI_Environment &env, const tMPICoreId &srcCore, const T *srcData, const tMPICount &length, const tMPIType &srcDataType, T *dstData, const tMPICount &dstDataSize, const tMPIType &dstDataType)
blocking scattering a data to all cores of environment with srcCore id
Definition: MPI_WorldToWorldMessage.h:258
static tBoolean AllReduce(const MPI_Environment &env, const std::array< T, N > &srcData, std::array< T, N > &dstData, const tMPIOperation &op)
compute a data from data of others cores and to data of core 0 and copy the values to all cores
Definition: MPI_WorldToWorldMessage.h:777
static tBoolean Reduce(const MPI_Environment &env, const tMPICoreId &dstCore, T &data, const tMPIOperation &op)
compute a data from data of others cores to data of dstCore
Definition: MPI_WorldToWorldMessage.h:615
static tBoolean Gather(const MPI_Environment &env, const tMPICoreId &dstCore, const std::valarray< T > &srcData, std::valarray< T > &dstData, const std::valarray< tMPICount > &dataLengthPerCore, const std::valarray< tMPICount > &dataIndexPerCore)
blocking gathering a data of diffzrent size from all cores of environments to dst core
Definition: MPI_WorldToWorldMessage.h:366
static tBoolean Gather(const MPI_Environment &env, const tMPICoreId &dstCore, const std::valarray< T > &srcData, std::valarray< T > &dstData)
blocking gathering a data from all cores of enviroments to the dest core
Definition: MPI_WorldToWorldMessage.h:324
static tBoolean IBCast(const MPI_Environment &env, const tMPICoreId &srcCore, T *data, const tMPICount &nData, const tMPIType &dataType, tMPIRequest &request)
blocking send a data with flag to all cores of enviroments except the core with srcCore id
Definition: MPI_WorldToWorldMessage.h:195
static int BuildOffset(const std::valarray< tMPIIndex > &dataLength, std::valarray< tMPIIndex > &dataOffset)
compute the displacement from the data length per core
Definition: MPI_WorldToWorldMessage.h:467
static tBoolean AllReduce(const MPI_Environment &env, T &data, const tMPIOperation &op)
compute a data from data of others cores and to data of core 0 and copy the values to all cores
Definition: MPI_WorldToWorldMessage.h:703
static CORE_UniquePointer< SelfClass > New()
create an unique instance of the class This
Definition: MPI_WorldToWorldMessage.h:111
static tBoolean AllGather(const MPI_Environment &env, const std::valarray< T > &srcData, T *dstData, const std::valarray< tMPIIndex > &srcDataLength, const std::valarray< tMPIIndex > &displacement)
blocking gathering a variable data lengthfrom all cores of enviroments to all cores
Definition: MPI_WorldToWorldMessage.h:494
static tBoolean AllReduce(const MPI_Environment &env, const std::valarray< T > &srcData, std::valarray< T > &dstData, const tMPIOperation &op)
compute a data from data of others cores and to data of core 0 and copy the values to all cores
Definition: MPI_WorldToWorldMessage.h:758
static tBoolean Scatter(const MPI_Environment &env, const tMPICoreId &srcCore, const T *srcData, const std::valarray< tMPICount > &dataLengthPerCore, const std::valarray< tMPICount > &dataLengthPerCoreOffset, T *dstData, const tMPICount &dstDataSize)
blocking scattering a data to all cores of environment with srcCore id
Definition: MPI_WorldToWorldMessage.h:280
static tBoolean Reduce(const MPI_Environment &env, const tMPICoreId &dstCore, const T &srcData, T &dstData, const tMPIOperation &op)
compute a data from data of others cores to data of dstCore
Definition: MPI_WorldToWorldMessage.h:599
static tBoolean AllReduce(const MPI_Environment &env, std::array< T, N > &data, const tMPIOperation &op)
compute a data from data of others cores and to data of core 0 and copy the values to all cores
Definition: MPI_WorldToWorldMessage.h:811
static tBoolean Send(const MPI_Environment &env, const tMPICoreId &srcCore, T *data, const tMPICount &nData)
blocking send a data with flag to all cores of enviroments except the core with srcCore id
Definition: MPI_WorldToWorldMessage.h:133