parUtils.txx 32 KB

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