parUtils.txx 32 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915
  1. /**
  2. @file parUtils.txx
  3. @brief Definitions of the templated functions in the par module.
  4. @author Rahul S. Sampath, rahul.sampath@gmail.com
  5. @author Hari Sundar, hsundar@gmail.com
  6. @author Shravan Veerapaneni, shravan@seas.upenn.edu
  7. @author Santi Swaroop Adavani, santis@gmail.com
  8. */
  9. #include "dtypes.h"
  10. #include <cassert>
  11. #include <iostream>
  12. #include <algorithm>
  13. #include <cstring>
  14. #include "ompUtils.h"
  15. #include <mpi.h>
  16. namespace pvfmm{
  17. namespace par{
  18. template <typename T>
  19. int Mpi_Alltoallv_sparse(T* sendbuf, int* sendcnts, int* sdispls,
  20. T* recvbuf, int* recvcnts, int* rdispls, const MPI_Comm &comm) {
  21. #ifndef ALLTOALLV_FIX
  22. return Mpi_Alltoallv
  23. (sendbuf, sendcnts, sdispls,
  24. recvbuf, recvcnts, rdispls, comm);
  25. #else
  26. int npes, rank;
  27. MPI_Comm_size(comm, &npes);
  28. MPI_Comm_rank(comm, &rank);
  29. int commCnt = 0;
  30. #pragma omp parallel for reduction(+:commCnt)
  31. for(int i = 0; i < rank; i++) {
  32. if(sendcnts[i] > 0) {
  33. commCnt++;
  34. }
  35. if(recvcnts[i] > 0) {
  36. commCnt++;
  37. }
  38. }
  39. #pragma omp parallel for reduction(+:commCnt)
  40. for(int i = (rank+1); i < npes; i++) {
  41. if(sendcnts[i] > 0) {
  42. commCnt++;
  43. }
  44. if(recvcnts[i] > 0) {
  45. commCnt++;
  46. }
  47. }
  48. MPI_Request* requests = new MPI_Request[commCnt];
  49. assert(requests);
  50. MPI_Status* statuses = new MPI_Status[commCnt];
  51. assert(statuses);
  52. commCnt = 0;
  53. //First place all recv requests. Do not recv from self.
  54. for(int i = 0; i < rank; i++) {
  55. if(recvcnts[i] > 0) {
  56. MPI_Irecv( &(recvbuf[rdispls[i]]) , recvcnts[i], par::Mpi_datatype<T>::value(), i, 1,
  57. comm, &(requests[commCnt]) );
  58. commCnt++;
  59. }
  60. }
  61. for(int i = (rank + 1); i < npes; i++) {
  62. if(recvcnts[i] > 0) {
  63. MPI_Irecv( &(recvbuf[rdispls[i]]) , recvcnts[i], par::Mpi_datatype<T>::value(), i, 1,
  64. comm, &(requests[commCnt]) );
  65. commCnt++;
  66. }
  67. }
  68. //Next send the messages. Do not send to self.
  69. for(int i = 0; i < rank; i++) {
  70. if(sendcnts[i] > 0) {
  71. MPI_Issend( &(sendbuf[sdispls[i]]), sendcnts[i], par::Mpi_datatype<T>::value(), i, 1,
  72. comm, &(requests[commCnt]) );
  73. commCnt++;
  74. }
  75. }
  76. for(int i = (rank + 1); i < npes; i++) {
  77. if(sendcnts[i] > 0) {
  78. MPI_Issend( &(sendbuf[sdispls[i]]), sendcnts[i], par::Mpi_datatype<T>::value(), i, 1,
  79. comm, &(requests[commCnt]) );
  80. commCnt++;
  81. }
  82. }
  83. //Now copy local portion.
  84. #pragma omp parallel for
  85. for(int i = 0; i < sendcnts[rank]; i++) {
  86. recvbuf[rdispls[rank] + i] = sendbuf[sdispls[rank] + i];
  87. }
  88. MPI_Waitall(commCnt, requests, statuses);
  89. delete [] requests;
  90. delete [] statuses;
  91. return 0;
  92. #endif
  93. }
  94. template <typename T>
  95. int Mpi_Alltoallv_dense(T* sbuff_, int* s_cnt_, int* sdisp_,
  96. T* rbuff_, int* r_cnt_, int* rdisp_, const MPI_Comm& comm){
  97. #ifndef ALLTOALLV_FIX
  98. return Mpi_Alltoallv
  99. (sbuff_, s_cnt_, sdisp_,
  100. rbuff_, r_cnt_, rdisp_, c);
  101. #else
  102. int np, pid;
  103. MPI_Comm_size(comm,&np);
  104. MPI_Comm_rank(comm,&pid);
  105. int range[2]={0,np-1};
  106. int split_id, partner;
  107. std::vector<int> s_cnt(np);
  108. #pragma omp parallel for
  109. for(int i=0;i<np;i++){
  110. s_cnt[i]=s_cnt_[i]*sizeof(T)+2*sizeof(int);
  111. }
  112. std::vector<int> sdisp(np); sdisp[0]=0;
  113. omp_par::scan(&s_cnt[0],&sdisp[0],np);
  114. char* sbuff=new char[sdisp[np-1]+s_cnt[np-1]];
  115. #pragma omp parallel for
  116. for(int i=0;i<np;i++){
  117. ((int*)&sbuff[sdisp[i]])[0]=s_cnt[i];
  118. ((int*)&sbuff[sdisp[i]])[1]=pid;
  119. memcpy(&sbuff[sdisp[i]]+2*sizeof(int),&sbuff_[sdisp_[i]],s_cnt[i]-2*sizeof(int));
  120. }
  121. while(range[0]<range[1]){
  122. split_id=(range[0]+range[1])/2;
  123. int new_range[2]={(pid<=split_id?range[0]:split_id+1),
  124. (pid<=split_id?split_id:range[1] )};
  125. int cmp_range[2]={(pid> split_id?range[0]:split_id+1),
  126. (pid> split_id?split_id:range[1] )};
  127. int new_np=new_range[1]-new_range[0]+1;
  128. int cmp_np=cmp_range[1]-cmp_range[0]+1;
  129. partner=pid+cmp_range[0]-new_range[0];
  130. if(partner>range[1]) partner=range[1];
  131. assert(partner>=range[0]);
  132. bool extra_partner=( (range[1]-range[0])%2==0 &&
  133. range[1] ==pid );
  134. //Communication.
  135. {
  136. int* s_lengths=&s_cnt[cmp_range[0]-range[0]];
  137. std::vector<int> s_len_ext(cmp_np,0);
  138. std::vector<int> r_cnt (new_np,0);
  139. std::vector<int> r_cnt_ext(new_np,0);
  140. MPI_Status status;
  141. //Exchange send sizes.
  142. MPI_Sendrecv (&s_lengths[0],cmp_np,MPI_INT, partner,0, &r_cnt [0],new_np,MPI_INT, partner, 0,comm,&status);
  143. if(extra_partner) MPI_Sendrecv(&s_len_ext[0],cmp_np,MPI_INT,split_id,0, &r_cnt_ext[0],new_np,MPI_INT,split_id, 0,comm,&status);
  144. //Allocate receive buffer.
  145. std::vector<int> rdisp (new_np,0);
  146. std::vector<int> rdisp_ext(new_np,0);
  147. omp_par::scan(&r_cnt [0],&rdisp [0],new_np);
  148. omp_par::scan(&r_cnt_ext[0],&rdisp_ext[0],new_np);
  149. int rbuff_size =rdisp [new_np-1]+r_cnt [new_np-1];
  150. int rbuff_size_ext=rdisp_ext[new_np-1]+r_cnt_ext[new_np-1];
  151. char* rbuff = new char[rbuff_size ];
  152. char* rbuffext=(extra_partner? new char[rbuff_size_ext]: NULL);
  153. //Sendrecv data.
  154. {
  155. int * s_cnt_tmp=&s_cnt[cmp_range[0]-range[0]] ;
  156. int * sdisp_tmp=&sdisp[cmp_range[0]-range[0]];
  157. char* sbuff_tmp=&sbuff[sdisp_tmp[0]];
  158. int sbuff_size=sdisp_tmp[cmp_np-1]+s_cnt_tmp[cmp_np-1]-sdisp_tmp[0];
  159. MPI_Sendrecv (sbuff_tmp,sbuff_size,MPI_BYTE, partner,0, &rbuff [0],rbuff_size ,MPI_BYTE, partner, 0,comm,&status);
  160. if(extra_partner) MPI_Sendrecv( NULL, 0,MPI_BYTE,split_id,0, &rbuffext[0],rbuff_size_ext,MPI_BYTE,split_id, 0,comm,&status);
  161. }
  162. //Rearrange received data.
  163. {
  164. //assert(!extra_partner);
  165. int * s_cnt_old=&s_cnt[new_range[0]-range[0]];
  166. int * sdisp_old=&sdisp[new_range[0]-range[0]];
  167. std::vector<int> s_cnt_new(&s_cnt_old[0],&s_cnt_old[new_np]);
  168. std::vector<int> sdisp_new(new_np ,0 );
  169. #pragma omp parallel for
  170. for(int i=0;i<new_np;i++){
  171. s_cnt_new[i]+=r_cnt[i]+r_cnt_ext[i];
  172. }
  173. omp_par::scan(&s_cnt_new[0],&sdisp_new[0],new_np);
  174. //Copy data to sbuff_new.
  175. char* sbuff_new=new char[sdisp_new[new_np-1]+s_cnt_new[new_np-1]];
  176. #pragma omp parallel for
  177. for(int i=0;i<new_np;i++){
  178. memcpy(&sbuff_new[sdisp_new[i] ],&sbuff [sdisp_old[i]],s_cnt_old[i]);
  179. memcpy(&sbuff_new[sdisp_new[i]+s_cnt_old[i] ],&rbuff [rdisp [i]],r_cnt [i]);
  180. memcpy(&sbuff_new[sdisp_new[i]+s_cnt_old[i]+r_cnt[i]],&rbuffext[rdisp_ext[i]],r_cnt_ext[i]);
  181. }
  182. //Free memory.
  183. if(sbuff !=NULL) delete[] sbuff ;
  184. if(rbuff !=NULL) delete[] rbuff ;
  185. if(rbuffext!=NULL) delete[] rbuffext;
  186. //Substitute data for next iteration.
  187. s_cnt=s_cnt_new;
  188. sdisp=sdisp_new;
  189. sbuff=sbuff_new;
  190. }
  191. }
  192. range[0]=new_range[0];
  193. range[1]=new_range[1];
  194. }
  195. //Copy data to rbuff_.
  196. std::vector<char*> buff_ptr(np);
  197. char* tmp_ptr=sbuff;
  198. for(int i=0;i<np;i++){
  199. int& blk_size=((int*)tmp_ptr)[0];
  200. buff_ptr[i]=tmp_ptr;
  201. tmp_ptr+=blk_size;
  202. }
  203. #pragma omp parallel for
  204. for(int i=0;i<np;i++){
  205. int& blk_size=((int*)buff_ptr[i])[0];
  206. int& src_pid=((int*)buff_ptr[i])[1];
  207. assert(blk_size-2*sizeof(int)<=r_cnt_[src_pid]*sizeof(T));
  208. memcpy(&rbuff_[rdisp_[src_pid]],buff_ptr[i]+2*sizeof(int),blk_size-2*sizeof(int));
  209. }
  210. //Free memory.
  211. if(sbuff !=NULL) delete[] sbuff;
  212. return 0;
  213. #endif
  214. }
  215. template<typename T>
  216. unsigned int defaultWeight(const T *a){
  217. return 1;
  218. }
  219. template<typename T>
  220. int partitionW(Vector<T>& nodeList, unsigned int (*getWeight)(const T *), const MPI_Comm& comm){
  221. int npes, rank;
  222. MPI_Comm_size(comm, &npes);
  223. MPI_Comm_rank(comm, &rank);
  224. long long npesLong = npes;
  225. if(getWeight == NULL) {
  226. getWeight = par::defaultWeight<T>;
  227. }
  228. long long nlSize = nodeList.Dim();
  229. long long off1= 0, off2= 0, localWt= 0, totalWt = 0;
  230. // First construct arrays of wts.
  231. Vector<long long> wts(nlSize);
  232. #pragma omp parallel for reduction(+:localWt)
  233. for (long long i = 0; i < nlSize; i++){
  234. wts[i] = (*getWeight)( &(nodeList[i]) );
  235. localWt+=wts[i];
  236. }
  237. // compute the total weight of the problem ...
  238. MPI_Allreduce(&localWt, &totalWt, 1, par::Mpi_datatype<long long>::value(), MPI_SUM, comm);
  239. MPI_Scan(&localWt, &off2, 1, par::Mpi_datatype<long long>::value(), MPI_SUM, comm );
  240. off1=off2-localWt;
  241. // perform a local scan on the weights first ...
  242. Vector<long long> lscn(nlSize);
  243. if(nlSize) {
  244. lscn[0]=off1;
  245. omp_par::scan(&wts[0],&lscn[0],nlSize);
  246. }
  247. Vector<int> int_buff(npesLong*4);
  248. Vector<int> sendSz (npesLong,&int_buff[0]+npesLong*0,false);
  249. Vector<int> recvSz (npesLong,&int_buff[0]+npesLong*1,false);
  250. Vector<int> sendOff(npesLong,&int_buff[0]+npesLong*2,false);
  251. Vector<int> recvOff(npesLong,&int_buff[0]+npesLong*3,false);
  252. // compute the partition offsets and sizes so that All2Allv can be performed.
  253. // initialize ...
  254. #pragma omp parallel for
  255. for (size_t i = 0; i < npesLong; i++) {
  256. sendSz[i] = 0;
  257. }
  258. //The Heart of the algorithm....
  259. if(nlSize>0 && totalWt>0) {
  260. long long pid1=( off1 *npesLong)/totalWt;
  261. long long pid2=((off2+1)*npesLong)/totalWt+1;
  262. assert((totalWt*pid2)/npesLong>=off2);
  263. pid1=(pid1< 0? 0:pid1);
  264. pid2=(pid2>npesLong?npesLong:pid2);
  265. #pragma omp parallel for
  266. for(int i=pid1;i<pid2;i++){
  267. long long wt1=(totalWt*(i ))/npesLong;
  268. long long wt2=(totalWt*(i+1))/npesLong;
  269. long long start = std::lower_bound(&lscn[0], &lscn[0]+nlSize, wt1, std::less<long long>())-&lscn[0];
  270. long long end = std::lower_bound(&lscn[0], &lscn[0]+nlSize, wt2, std::less<long long>())-&lscn[0];
  271. if(i== 0) start=0 ;
  272. if(i==npesLong-1) end =nlSize;
  273. sendSz[i]=end-start;
  274. }
  275. }else sendSz[0]=nlSize;
  276. // communicate with other procs how many you shall be sending and get how
  277. // many to recieve from whom.
  278. MPI_Alltoall(&sendSz[0], 1, par::Mpi_datatype<int>::value(),
  279. &recvSz[0], 1, par::Mpi_datatype<int>::value(), comm);
  280. // compute offsets ...
  281. sendOff[0] = 0; omp_par::scan(&sendSz[0],&sendOff[0],npesLong);
  282. recvOff[0] = 0; omp_par::scan(&recvSz[0],&recvOff[0],npesLong);
  283. // new value of nlSize, ie the local nodes.
  284. long long nn = recvSz[npesLong-1] + recvOff[npes-1];
  285. // allocate memory for the new arrays ...
  286. Vector<T> newNodes;
  287. {
  288. if(nodeList.Capacity()>nn+std::max(nn,nlSize)){
  289. newNodes.ReInit(nn,&nodeList[0]+std::max(nn,nlSize),false);
  290. //}else if(buff!=NULL && buff->Dim()>nn*sizeof(T)){
  291. // newNodes.ReInit(nn,(T*)&(*buff)[0],false);
  292. }else newNodes.Resize(nn);
  293. }
  294. // perform All2All ...
  295. par::Mpi_Alltoallv_sparse<T>(&nodeList[0], &sendSz[0], &sendOff[0],
  296. &newNodes[0], &recvSz[0], &recvOff[0], comm);
  297. // reset the pointer ...
  298. nodeList=newNodes;
  299. return 0;
  300. }//end function
  301. template<typename T>
  302. int partitionW(std::vector<T>& nodeList, unsigned int (*getWeight)(const T *), const MPI_Comm& comm){
  303. Vector<T> nodeList_=nodeList;
  304. int ret = par::partitionW<T>(nodeList_, getWeight, comm);
  305. nodeList.assign(&nodeList_[0],&nodeList_[0]+nodeList_.Dim());
  306. return ret;
  307. }
  308. template<typename T>
  309. int HyperQuickSort(const Vector<T>& arr_, Vector<T>& SortedElem, const MPI_Comm& comm_){ // O( ((N/p)+log(p))*(log(N/p)+log(p)) )
  310. // Copy communicator.
  311. MPI_Comm comm=comm_;
  312. // Get comm size and rank.
  313. int npes, myrank;
  314. MPI_Comm_size(comm, &npes);
  315. MPI_Comm_rank(comm, &myrank);
  316. int omp_p=omp_get_max_threads();
  317. srand(myrank);
  318. // Local and global sizes. O(log p)
  319. long long totSize, nelem = arr_.Dim();
  320. //assert(nelem); // TODO: Check if this is needed.
  321. MPI_Allreduce(&nelem, &totSize, 1, par::Mpi_datatype<long long>::value(), MPI_SUM, comm);
  322. // Local sort.
  323. Vector<T> arr=arr_;
  324. omp_par::merge_sort(&arr[0], &arr[0]+nelem);
  325. // Allocate memory.
  326. Vector<T> nbuff;
  327. Vector<T> nbuff_ext;
  328. Vector<T> rbuff ;
  329. Vector<T> rbuff_ext;
  330. // Binary split and merge in each iteration.
  331. while(npes>1 && totSize>0){ // O(log p) iterations.
  332. //Determine splitters. O( log(N/p) + log(p) )
  333. T split_key;
  334. long long totSize_new;
  335. //while(true)
  336. {
  337. // Take random splitters. O( 1 ) -- Let p * splt_count = glb_splt_count = const = 100~1000
  338. int splt_count=(1000*nelem)/totSize;
  339. if(npes>1000) splt_count=(((float)rand()/(float)RAND_MAX)*totSize<(1000*nelem)?1:0);
  340. if(splt_count>nelem) splt_count=nelem;
  341. std::vector<T> splitters(splt_count);
  342. for(int i=0;i<splt_count;i++)
  343. splitters[i]=arr[rand()%nelem];
  344. // Gather all splitters. O( log(p) )
  345. int glb_splt_count;
  346. std::vector<int> glb_splt_cnts(npes);
  347. std::vector<int> glb_splt_disp(npes,0);
  348. MPI_Allgather(&splt_count , 1, par::Mpi_datatype<int>::value(),
  349. &glb_splt_cnts[0], 1, par::Mpi_datatype<int>::value(), comm);
  350. omp_par::scan(&glb_splt_cnts[0],&glb_splt_disp[0],npes);
  351. glb_splt_count=glb_splt_cnts[npes-1]+glb_splt_disp[npes-1];
  352. std::vector<T> glb_splitters(glb_splt_count);
  353. MPI_Allgatherv(& splitters[0], splt_count, par::Mpi_datatype<T>::value(),
  354. &glb_splitters[0], &glb_splt_cnts[0], &glb_splt_disp[0],
  355. par::Mpi_datatype<T>::value(), comm);
  356. // Determine split key. O( log(N/p) + log(p) )
  357. std::vector<long long> disp(glb_splt_count,0);
  358. if(nelem>0){
  359. #pragma omp parallel for
  360. for(int i=0;i<glb_splt_count;i++){
  361. disp[i]=std::lower_bound(&arr[0], &arr[0]+nelem, glb_splitters[i])-&arr[0];
  362. }
  363. }
  364. std::vector<long long> glb_disp(glb_splt_count,0);
  365. MPI_Allreduce(&disp[0], &glb_disp[0], glb_splt_count, par::Mpi_datatype<long long>::value(), MPI_SUM, comm);
  366. long long* split_disp=&glb_disp[0];
  367. for(int i=0;i<glb_splt_count;i++)
  368. if(labs(glb_disp[i]-totSize/2)<labs(*split_disp-totSize/2)) split_disp=&glb_disp[i];
  369. split_key=glb_splitters[split_disp-&glb_disp[0]];
  370. totSize_new=(myrank<=(npes-1)/2?*split_disp:totSize-*split_disp);
  371. //double err=(((double)*split_disp)/(totSize/2))-1.0;
  372. //if(fabs(err)<0.01 || npes<=16) break;
  373. //else if(!myrank) std::cout<<err<<'\n';
  374. }
  375. // Split problem into two. O( N/p )
  376. int split_id=(npes-1)/2;
  377. {
  378. int new_p0=(myrank<=split_id?0:split_id+1);
  379. int cmp_p0=(myrank> split_id?0:split_id+1);
  380. int partner = myrank+cmp_p0-new_p0;
  381. if(partner>=npes) partner=npes-1;
  382. assert(partner>=0);
  383. bool extra_partner=( npes%2==1 && npes-1==myrank );
  384. // Exchange send sizes.
  385. char *sbuff, *lbuff;
  386. int rsize=0, ssize=0, lsize=0;
  387. int ext_rsize=0, ext_ssize=0;
  388. size_t split_indx=(nelem>0?std::lower_bound(&arr[0], &arr[0]+nelem, split_key)-&arr[0]:0);
  389. ssize= (myrank> split_id? split_indx: nelem -split_indx)*sizeof(T);
  390. sbuff=(char*)(myrank> split_id? &arr[0] : &arr[0]+split_indx);
  391. lsize= (myrank<=split_id? split_indx: nelem -split_indx)*sizeof(T);
  392. lbuff=(char*)(myrank<=split_id? &arr[0] : &arr[0]+split_indx);
  393. MPI_Status status;
  394. MPI_Sendrecv (& ssize,1,MPI_INT, partner,0, & rsize,1,MPI_INT, partner, 0,comm,&status);
  395. if(extra_partner) MPI_Sendrecv(&ext_ssize,1,MPI_INT,split_id,0, &ext_rsize,1,MPI_INT,split_id, 0,comm,&status);
  396. // Exchange data.
  397. rbuff .Resize( rsize/sizeof(T));
  398. rbuff_ext.Resize(ext_rsize/sizeof(T));
  399. MPI_Sendrecv (sbuff,ssize,MPI_BYTE, partner,0, &rbuff [0], rsize,MPI_BYTE, partner, 0,comm,&status);
  400. if(extra_partner) MPI_Sendrecv( NULL, 0,MPI_BYTE,split_id,0, &rbuff_ext[0],ext_rsize,MPI_BYTE,split_id, 0,comm,&status);
  401. int nbuff_size=lsize+rsize+ext_rsize;
  402. nbuff.Resize(nbuff_size/sizeof(T));
  403. omp_par::merge<T*>((T*)lbuff, (T*)(lbuff+lsize), &rbuff[0], &rbuff[0]+(rsize/sizeof(T)), &nbuff[0], omp_p, std::less<T>());
  404. if(ext_rsize>0 && nbuff.Dim()>0){
  405. nbuff_ext.Resize(nbuff_size/sizeof(T));
  406. omp_par::merge<T*>(&nbuff[0], &nbuff[0]+((lsize+rsize)/sizeof(T)), &rbuff_ext[0], &rbuff_ext[0]+(ext_rsize/sizeof(T)), &nbuff_ext[0], omp_p, std::less<T>());
  407. nbuff.Swap(nbuff_ext);
  408. nbuff_ext.Resize(0);
  409. }
  410. // Copy new data.
  411. totSize=totSize_new;
  412. nelem = nbuff_size/sizeof(T);
  413. arr.Swap(nbuff);
  414. nbuff.Resize(0);
  415. }
  416. {// Split comm. O( log(p) ) ??
  417. MPI_Comm scomm;
  418. MPI_Comm_split(comm, myrank<=split_id, myrank, &scomm );
  419. comm=scomm;
  420. npes =(myrank<=split_id? split_id+1: npes -split_id-1);
  421. myrank=(myrank<=split_id? myrank : myrank-split_id-1);
  422. }
  423. }
  424. SortedElem.Resize(nelem);
  425. memcpy(&SortedElem[0], &arr[0], nelem*sizeof(T));
  426. par::partitionW<T>(SortedElem, NULL , comm_);
  427. return 0;
  428. }//end function
  429. template<typename T>
  430. int HyperQuickSort(const std::vector<T>& arr_, std::vector<T>& SortedElem_, const MPI_Comm& comm_){
  431. Vector<T> SortedElem;
  432. const Vector<T> arr(arr_.size(),(T*)&arr_[0],false);
  433. int ret = HyperQuickSort(arr, SortedElem, comm_);
  434. SortedElem_.assign(&SortedElem[0],&SortedElem[0]+SortedElem.Dim());
  435. return ret;
  436. }
  437. template<typename A, typename B>
  438. struct SortPair{
  439. int operator<(const SortPair<A,B>& p1) const{ return key<p1.key;}
  440. A key;
  441. B data;
  442. };
  443. template<typename T>
  444. int SortScatterIndex(const Vector<T>& key, Vector<size_t>& scatter_index, const MPI_Comm& comm, const T* split_key_){
  445. typedef SortPair<T,size_t> Pair_t;
  446. int npes, rank;
  447. MPI_Comm_size(comm, &npes);
  448. MPI_Comm_rank(comm, &rank);
  449. long long npesLong = npes;
  450. //Vector<char> buff;
  451. //if(buff_!=NULL && buff_->Dim()>0){
  452. // buff.ReInit(buff_->Dim(),&(*buff_)[0],false);
  453. //}
  454. Vector<Pair_t> parray;
  455. { // Allocate memory
  456. //size_t parray_size=key.Dim()*sizeof(Pair_t);
  457. //if(buff.Dim()>parray_size){
  458. // parray.ReInit(key.Dim(),(Pair_t*)&buff[0],false);
  459. // buff.ReInit(buff.Dim()-parray_size,&buff[0]+parray_size,false);
  460. //}else
  461. parray.Resize(key.Dim());
  462. }
  463. { // Build global index.
  464. long long glb_dsp=0;
  465. long long loc_size=key.Dim();
  466. MPI_Scan(&loc_size, &glb_dsp, 1, par::Mpi_datatype<long long>::value(), MPI_SUM, comm);
  467. glb_dsp-=loc_size;
  468. #pragma omp parallel for
  469. for(size_t i=0;i<loc_size;i++){
  470. parray[i].key=key[i];
  471. parray[i].data=glb_dsp+i;
  472. }
  473. }
  474. Vector<Pair_t> psorted;
  475. { // Allocate memory
  476. //if(buff.Dim()>0){
  477. // psorted.ReInit(buff.Dim()/sizeof(Pair_t), (Pair_t*)&buff[0], false);
  478. //}
  479. }
  480. HyperQuickSort(parray, psorted, comm);
  481. if(split_key_!=NULL){ // Partition data
  482. Vector<T> split_key(npesLong);
  483. MPI_Allgather((void*)split_key_ , 1, par::Mpi_datatype<T>::value(),
  484. &split_key[0], 1, par::Mpi_datatype<T>::value(), comm);
  485. Vector<int> int_buff(npesLong*4);
  486. Vector<int> sendSz (npesLong,&int_buff[0]+npesLong*0,false);
  487. Vector<int> recvSz (npesLong,&int_buff[0]+npesLong*1,false);
  488. Vector<int> sendOff(npesLong,&int_buff[0]+npesLong*2,false);
  489. Vector<int> recvOff(npesLong,&int_buff[0]+npesLong*3,false);
  490. long long nlSize = psorted.Dim();
  491. // compute the partition offsets and sizes so that All2Allv can be performed.
  492. // initialize ...
  493. #pragma omp parallel for
  494. for (size_t i = 0; i < npesLong; i++) {
  495. sendSz[i] = 0;
  496. }
  497. //The Heart of the algorithm....
  498. if(nlSize>0) {
  499. // Determine processor range.
  500. long long pid1=std::lower_bound(&split_key[0], &split_key[0]+npesLong, psorted[ 0].key)-&split_key[0]-1;
  501. long long pid2=std::upper_bound(&split_key[0], &split_key[0]+npesLong, psorted[nlSize-1].key)-&split_key[0]+0;
  502. pid1=(pid1< 0? 0:pid1);
  503. pid2=(pid2>npesLong?npesLong:pid2);
  504. #pragma omp parallel for
  505. for(int i=pid1;i<pid2;i++){
  506. Pair_t p1; p1.key=split_key[ i];
  507. Pair_t p2; p2.key=split_key[i+1<npesLong?i+1:i];
  508. long long start = std::lower_bound(&psorted[0], &psorted[0]+nlSize, p1, std::less<Pair_t>())-&psorted[0];
  509. long long end = std::lower_bound(&psorted[0], &psorted[0]+nlSize, p2, std::less<Pair_t>())-&psorted[0];
  510. if(i== 0) start=0 ;
  511. if(i==npesLong-1) end =nlSize;
  512. sendSz[i]=end-start;
  513. }
  514. }
  515. // communicate with other procs how many you shall be sending and get how
  516. // many to recieve from whom.
  517. MPI_Alltoall(&sendSz[0], 1, par::Mpi_datatype<int>::value(),
  518. &recvSz[0], 1, par::Mpi_datatype<int>::value(), comm);
  519. // compute offsets ...
  520. sendOff[0] = 0; omp_par::scan(&sendSz[0],&sendOff[0],npesLong);
  521. recvOff[0] = 0; omp_par::scan(&recvSz[0],&recvOff[0],npesLong);
  522. // new value of nlSize, ie the local nodes.
  523. long long nn = recvSz[npesLong-1] + recvOff[npesLong-1];
  524. // allocate memory for the new arrays ...
  525. Vector<Pair_t> newNodes;
  526. {
  527. if(psorted.Capacity()>nn+std::max(nn,nlSize)){
  528. newNodes.ReInit(nn,&psorted[0]+std::max(nn,nlSize),false);
  529. }else newNodes.Resize(nn);
  530. }
  531. // perform All2All ...
  532. par::Mpi_Alltoallv_sparse<Pair_t>(&psorted[0], &sendSz[0], &sendOff[0],
  533. &newNodes[0], &recvSz[0], &recvOff[0], comm);
  534. // reset the pointer ...
  535. psorted=newNodes;
  536. }
  537. scatter_index.Resize(psorted.Dim());
  538. #pragma omp parallel for
  539. for(size_t i=0;i<psorted.Dim();i++){
  540. scatter_index[i]=psorted[i].data;
  541. }
  542. return 0;
  543. }
  544. template<typename T>
  545. int ScatterForward(Vector<T>& data_, const Vector<size_t>& scatter_index, const MPI_Comm& comm){
  546. typedef SortPair<size_t,size_t> Pair_t;
  547. int npes, rank;
  548. MPI_Comm_size(comm, &npes);
  549. MPI_Comm_rank(comm, &rank);
  550. long long npesLong = npes;
  551. size_t data_dim=0;
  552. long long send_size=0;
  553. long long recv_size=0;
  554. {
  555. recv_size=scatter_index.Dim();
  556. long long glb_size[2]={0,0};
  557. long long loc_size[2]={data_.Dim()*sizeof(T), recv_size};
  558. MPI_Allreduce(&loc_size, &glb_size, 2, par::Mpi_datatype<long long>::value(), MPI_SUM, comm);
  559. if(glb_size[0]==0 || glb_size[1]==0) return 0; //Nothing to be done.
  560. data_dim=glb_size[0]/glb_size[1];
  561. assert(glb_size[0]==data_dim*glb_size[1]);
  562. send_size=(data_.Dim()*sizeof(T))/data_dim;
  563. }
  564. Vector<char> recv_buff(recv_size*data_dim);
  565. Vector<char> send_buff(send_size*data_dim);
  566. // Global scan of data size.
  567. Vector<long long> glb_scan(npesLong);
  568. {
  569. long long glb_rank=0;
  570. MPI_Scan(&send_size, &glb_rank, 1, par::Mpi_datatype<long long>::value(), MPI_SUM, comm);
  571. glb_rank-=send_size;
  572. MPI_Allgather(&glb_rank , 1, par::Mpi_datatype<long long>::value(),
  573. &glb_scan[0], 1, par::Mpi_datatype<long long>::value(), comm);
  574. }
  575. // Sort scatter_index.
  576. Vector<Pair_t> psorted(recv_size);
  577. {
  578. #pragma omp parallel for
  579. for(size_t i=0;i<recv_size;i++){
  580. psorted[i].key=scatter_index[i];
  581. psorted[i].data=i;
  582. }
  583. omp_par::merge_sort(&psorted[0], &psorted[0]+recv_size);
  584. }
  585. // Exchange send, recv indices.
  586. Vector<size_t> recv_indx(recv_size);
  587. Vector<size_t> send_indx(send_size);
  588. Vector<int> sendSz(npesLong);
  589. Vector<int> sendOff(npesLong);
  590. Vector<int> recvSz(npesLong);
  591. Vector<int> recvOff(npesLong);
  592. {
  593. #pragma omp parallel for
  594. for(size_t i=0;i<recv_size;i++){
  595. recv_indx[i]=psorted[i].key;
  596. }
  597. #pragma omp parallel for
  598. for(size_t i=0;i<npesLong;i++){
  599. size_t start= std::lower_bound(&recv_indx[0], &recv_indx[0]+recv_size, glb_scan[ i])-&recv_indx[0];
  600. size_t end =(i+1<npesLong?std::lower_bound(&recv_indx[0], &recv_indx[0]+recv_size, glb_scan[i+1])-&recv_indx[0]:recv_size);
  601. recvSz[i]=end-start;
  602. recvOff[i]=start;
  603. }
  604. MPI_Alltoall(&recvSz[0], 1, par::Mpi_datatype<int>::value(),
  605. &sendSz[0], 1, par::Mpi_datatype<int>::value(), comm);
  606. sendOff[0] = 0; omp_par::scan(&sendSz[0],&sendOff[0],npesLong);
  607. assert(sendOff[npesLong-1]+sendSz[npesLong-1]==send_size);
  608. par::Mpi_Alltoallv_dense<size_t>(&recv_indx[0], &recvSz[0], &recvOff[0],
  609. &send_indx[0], &sendSz[0], &sendOff[0], comm);
  610. #pragma omp parallel for
  611. for(size_t i=0;i<send_size;i++){
  612. assert(send_indx[i]>=glb_scan[rank]);
  613. send_indx[i]-=glb_scan[rank];
  614. assert(send_indx[i]<send_size);
  615. }
  616. }
  617. // Prepare send buffer
  618. {
  619. char* data=(char*)&data_[0];
  620. #pragma omp parallel for
  621. for(size_t i=0;i<send_size;i++){
  622. size_t src_indx=send_indx[i]*data_dim;
  623. size_t trg_indx=i*data_dim;
  624. for(size_t j=0;j<data_dim;j++)
  625. send_buff[trg_indx+j]=data[src_indx+j];
  626. }
  627. }
  628. // All2Allv
  629. {
  630. #pragma omp parallel for
  631. for(size_t i=0;i<npesLong;i++){
  632. sendSz [i]*=data_dim;
  633. sendOff[i]*=data_dim;
  634. recvSz [i]*=data_dim;
  635. recvOff[i]*=data_dim;
  636. }
  637. par::Mpi_Alltoallv_dense<char>(&send_buff[0], &sendSz[0], &sendOff[0],
  638. &recv_buff[0], &recvSz[0], &recvOff[0], comm);
  639. }
  640. // Build output data.
  641. {
  642. data_.Resize(recv_size*data_dim/sizeof(T));
  643. char* data=(char*)&data_[0];
  644. #pragma omp parallel for
  645. for(size_t i=0;i<recv_size;i++){
  646. size_t src_indx=i*data_dim;
  647. size_t trg_indx=psorted[i].data*data_dim;
  648. for(size_t j=0;j<data_dim;j++)
  649. data[trg_indx+j]=recv_buff[src_indx+j];
  650. }
  651. }
  652. return 0;
  653. }
  654. template<typename T>
  655. int ScatterReverse(Vector<T>& data_, const Vector<size_t>& scatter_index, const MPI_Comm& comm, size_t loc_size){
  656. typedef SortPair<size_t,size_t> Pair_t;
  657. int npes, rank;
  658. MPI_Comm_size(comm, &npes);
  659. MPI_Comm_rank(comm, &rank);
  660. long long npesLong = npes;
  661. size_t data_dim=0;
  662. long long send_size=0;
  663. long long recv_size=0;
  664. {
  665. send_size=scatter_index.Dim();
  666. recv_size=loc_size;
  667. long long glb_size[3]={0,0};
  668. long long loc_size[3]={data_.Dim()*sizeof(T), send_size, recv_size};
  669. MPI_Allreduce(&loc_size, &glb_size, 3, par::Mpi_datatype<long long>::value(), MPI_SUM, comm);
  670. if(glb_size[0]==0 || glb_size[1]==0) return 0; //Nothing to be done.
  671. data_dim=glb_size[0]/glb_size[1];
  672. assert(glb_size[0]==data_dim*glb_size[1]);
  673. if(glb_size[1]!=glb_size[2]){
  674. recv_size=(((rank+1)*glb_size[1])/npesLong)-
  675. (( rank *glb_size[1])/npesLong);
  676. }
  677. }
  678. Vector<char> recv_buff(recv_size*data_dim);
  679. Vector<char> send_buff(send_size*data_dim);
  680. // Global data size.
  681. Vector<long long> glb_scan(npesLong);
  682. {
  683. long long glb_rank=0;
  684. MPI_Scan(&recv_size, &glb_rank, 1, par::Mpi_datatype<long long>::value(), MPI_SUM, comm);
  685. glb_rank-=recv_size;
  686. MPI_Allgather(&glb_rank , 1, par::Mpi_datatype<long long>::value(),
  687. &glb_scan[0], 1, par::Mpi_datatype<long long>::value(), comm);
  688. }
  689. // Sort scatter_index.
  690. Vector<Pair_t> psorted(send_size);
  691. {
  692. #pragma omp parallel for
  693. for(size_t i=0;i<send_size;i++){
  694. psorted[i].key=scatter_index[i];
  695. psorted[i].data=i;
  696. }
  697. omp_par::merge_sort(&psorted[0], &psorted[0]+send_size);
  698. }
  699. // Exchange send, recv indices.
  700. Vector<size_t> recv_indx(recv_size);
  701. Vector<size_t> send_indx(send_size);
  702. Vector<int> sendSz(npesLong);
  703. Vector<int> sendOff(npesLong);
  704. Vector<int> recvSz(npesLong);
  705. Vector<int> recvOff(npesLong);
  706. {
  707. #pragma omp parallel for
  708. for(size_t i=0;i<send_size;i++){
  709. send_indx[i]=psorted[i].key;
  710. }
  711. #pragma omp parallel for
  712. for(size_t i=0;i<npesLong;i++){
  713. size_t start= std::lower_bound(&send_indx[0], &send_indx[0]+send_size, glb_scan[ i])-&send_indx[0];
  714. size_t end =(i+1<npesLong?std::lower_bound(&send_indx[0], &send_indx[0]+send_size, glb_scan[i+1])-&send_indx[0]:send_size);
  715. sendSz[i]=end-start;
  716. sendOff[i]=start;
  717. }
  718. MPI_Alltoall(&sendSz[0], 1, par::Mpi_datatype<int>::value(),
  719. &recvSz[0], 1, par::Mpi_datatype<int>::value(), comm);
  720. recvOff[0] = 0; omp_par::scan(&recvSz[0],&recvOff[0],npesLong);
  721. assert(recvOff[npesLong-1]+recvSz[npesLong-1]==recv_size);
  722. par::Mpi_Alltoallv_dense<size_t>(&send_indx[0], &sendSz[0], &sendOff[0],
  723. &recv_indx[0], &recvSz[0], &recvOff[0], comm);
  724. #pragma omp parallel for
  725. for(size_t i=0;i<recv_size;i++){
  726. assert(recv_indx[i]>=glb_scan[rank]);
  727. recv_indx[i]-=glb_scan[rank];
  728. assert(recv_indx[i]<recv_size);
  729. }
  730. }
  731. // Prepare send buffer
  732. {
  733. char* data=(char*)&data_[0];
  734. #pragma omp parallel for
  735. for(size_t i=0;i<send_size;i++){
  736. size_t src_indx=psorted[i].data*data_dim;
  737. size_t trg_indx=i*data_dim;
  738. for(size_t j=0;j<data_dim;j++)
  739. send_buff[trg_indx+j]=data[src_indx+j];
  740. }
  741. }
  742. // All2Allv
  743. {
  744. #pragma omp parallel for
  745. for(size_t i=0;i<npesLong;i++){
  746. sendSz [i]*=data_dim;
  747. sendOff[i]*=data_dim;
  748. recvSz [i]*=data_dim;
  749. recvOff[i]*=data_dim;
  750. }
  751. par::Mpi_Alltoallv_dense<char>(&send_buff[0], &sendSz[0], &sendOff[0],
  752. &recv_buff[0], &recvSz[0], &recvOff[0], comm);
  753. }
  754. // Build output data.
  755. {
  756. data_.Resize(recv_size*data_dim/sizeof(T));
  757. char* data=(char*)&data_[0];
  758. #pragma omp parallel for
  759. for(size_t i=0;i<recv_size;i++){
  760. size_t src_indx=i*data_dim;
  761. size_t trg_indx=recv_indx[i]*data_dim;
  762. for(size_t j=0;j<data_dim;j++)
  763. data[trg_indx+j]=recv_buff[src_indx+j];
  764. }
  765. }
  766. return 0;
  767. }
  768. }//end namespace
  769. }//end namespace