parUtils.txx 32 KB

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