comm.txx 39 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147
  1. #include <type_traits>
  2. #include SCTL_INCLUDE(ompUtils.hpp)
  3. #include SCTL_INCLUDE(vector.hpp)
  4. namespace SCTL_NAMESPACE {
  5. inline Comm::Comm() {
  6. #ifdef SCTL_HAVE_MPI
  7. Init(MPI_COMM_SELF);
  8. #endif
  9. }
  10. inline Comm::Comm(const Comm& c) {
  11. #ifdef SCTL_HAVE_MPI
  12. Init(c.mpi_comm_);
  13. #endif
  14. }
  15. inline Comm Comm::Self() {
  16. #ifdef SCTL_HAVE_MPI
  17. Comm comm_self(MPI_COMM_SELF);
  18. return comm_self;
  19. #else
  20. Comm comm_self;
  21. return comm_self;
  22. #endif
  23. }
  24. inline Comm Comm::World() {
  25. #ifdef SCTL_HAVE_MPI
  26. Comm comm_world(MPI_COMM_WORLD);
  27. return comm_world;
  28. #else
  29. Comm comm_self;
  30. return comm_self;
  31. #endif
  32. }
  33. inline Comm& Comm::operator=(const Comm& c) {
  34. #ifdef SCTL_HAVE_MPI
  35. Init(c.mpi_comm_);
  36. #endif
  37. }
  38. inline Comm::~Comm() {
  39. #ifdef SCTL_HAVE_MPI
  40. while (!req.empty()) {
  41. delete (Vector<MPI_Request>*)req.top();
  42. req.pop();
  43. }
  44. MPI_Comm_free(&mpi_comm_);
  45. #endif
  46. }
  47. inline Comm Comm::Split(Integer clr) const {
  48. #ifdef SCTL_HAVE_MPI
  49. MPI_Comm new_comm;
  50. MPI_Comm_split(mpi_comm_, clr, mpi_rank_, &new_comm);
  51. Comm c(new_comm);
  52. MPI_Comm_free(&new_comm);
  53. return c;
  54. #else
  55. Comm c;
  56. return c;
  57. #endif
  58. }
  59. inline Integer Comm::Rank() const {
  60. #ifdef SCTL_HAVE_MPI
  61. return mpi_rank_;
  62. #else
  63. return 0;
  64. #endif
  65. }
  66. inline Integer Comm::Size() const {
  67. #ifdef SCTL_HAVE_MPI
  68. return mpi_size_;
  69. #else
  70. return 1;
  71. #endif
  72. }
  73. inline void Comm::Barrier() const {
  74. #ifdef SCTL_HAVE_MPI
  75. MPI_Barrier(mpi_comm_);
  76. #endif
  77. }
  78. template <class SType> void* Comm::Isend(ConstIterator<SType> sbuf, Long scount, Integer dest, Integer tag) const {
  79. static_assert(std::is_trivially_copyable<SType>::value, "Data is not trivially copyable!");
  80. #ifdef SCTL_HAVE_MPI
  81. if (!scount) return nullptr;
  82. Vector<MPI_Request>& request = *NewReq();
  83. request.ReInit(1);
  84. sbuf[0];
  85. sbuf[scount - 1];
  86. #ifndef NDEBUG
  87. MPI_Issend(&sbuf[0], scount, CommDatatype<SType>::value(), dest, tag, mpi_comm_, &request[0]);
  88. #else
  89. MPI_Isend(&sbuf[0], scount, CommDatatype<SType>::value(), dest, tag, mpi_comm_, &request[0]);
  90. #endif
  91. return &request;
  92. #else
  93. auto it = recv_req.find(tag);
  94. if (it == recv_req.end())
  95. send_req.insert(std::pair<Integer, ConstIterator<char>>(tag, (ConstIterator<char>)sbuf));
  96. else
  97. memcopy(it->second, (ConstIterator<char>)sbuf, scount * sizeof(SType));
  98. return nullptr;
  99. #endif
  100. }
  101. template <class RType> void* Comm::Irecv(Iterator<RType> rbuf, Long rcount, Integer source, Integer tag) const {
  102. static_assert(std::is_trivially_copyable<RType>::value, "Data is not trivially copyable!");
  103. #ifdef SCTL_HAVE_MPI
  104. if (!rcount) return nullptr;
  105. Vector<MPI_Request>& request = *NewReq();
  106. request.ReInit(1);
  107. rbuf[0];
  108. rbuf[rcount - 1];
  109. MPI_Irecv(&rbuf[0], rcount, CommDatatype<RType>::value(), source, tag, mpi_comm_, &request[0]);
  110. return &request;
  111. #else
  112. auto it = send_req.find(tag);
  113. if (it == send_req.end())
  114. recv_req.insert(std::pair<Integer, Iterator<char>>(tag, (Iterator<char>)rbuf));
  115. else
  116. memcopy((Iterator<char>)rbuf, it->second, rcount * sizeof(RType));
  117. return nullptr;
  118. #endif
  119. }
  120. inline void Comm::Wait(void* req_ptr) const {
  121. #ifdef SCTL_HAVE_MPI
  122. if (req_ptr == nullptr) return;
  123. Vector<MPI_Request>& request = *(Vector<MPI_Request>*)req_ptr;
  124. // std::vector<MPI_Status> status(request.Dim());
  125. if (request.Dim()) MPI_Waitall(request.Dim(), &request[0], MPI_STATUSES_IGNORE); //&status[0]);
  126. DelReq(&request);
  127. #endif
  128. }
  129. template <class SType, class RType> void Comm::Allgather(ConstIterator<SType> sbuf, Long scount, Iterator<RType> rbuf, Long rcount) const {
  130. static_assert(std::is_trivially_copyable<SType>::value, "Data is not trivially copyable!");
  131. static_assert(std::is_trivially_copyable<RType>::value, "Data is not trivially copyable!");
  132. #ifdef SCTL_HAVE_MPI
  133. if (scount) {
  134. sbuf[0];
  135. sbuf[scount - 1];
  136. }
  137. if (rcount) {
  138. rbuf[0];
  139. rbuf[rcount * Size() - 1];
  140. }
  141. MPI_Allgather((scount ? &sbuf[0] : nullptr), scount, CommDatatype<SType>::value(), (rcount ? &rbuf[0] : nullptr), rcount, CommDatatype<RType>::value(), mpi_comm_);
  142. #else
  143. memcopy((Iterator<char>)rbuf, (ConstIterator<char>)sbuf, scount * sizeof(SType));
  144. #endif
  145. }
  146. template <class SType, class RType> void Comm::Allgatherv(ConstIterator<SType> sbuf, Long scount, Iterator<RType> rbuf, ConstIterator<Long> rcounts, ConstIterator<Long> rdispls) const {
  147. static_assert(std::is_trivially_copyable<SType>::value, "Data is not trivially copyable!");
  148. static_assert(std::is_trivially_copyable<RType>::value, "Data is not trivially copyable!");
  149. #ifdef SCTL_HAVE_MPI
  150. Vector<int> rcounts_(mpi_size_), rdispls_(mpi_size_);
  151. Long rcount_sum = 0;
  152. #pragma omp parallel for schedule(static) reduction(+ : rcount_sum)
  153. for (Integer i = 0; i < mpi_size_; i++) {
  154. rcounts_[i] = rcounts[i];
  155. rdispls_[i] = rdispls[i];
  156. rcount_sum += rcounts[i];
  157. }
  158. if (scount) {
  159. sbuf[0];
  160. sbuf[scount - 1];
  161. }
  162. if (rcount_sum) {
  163. rbuf[0];
  164. rbuf[rcount_sum - 1];
  165. }
  166. MPI_Allgatherv((scount ? &sbuf[0] : nullptr), scount, CommDatatype<SType>::value(), (rcount_sum ? &rbuf[0] : nullptr), &rcounts_.begin()[0], &rdispls_.begin()[0], CommDatatype<RType>::value(), mpi_comm_);
  167. #else
  168. memcopy((Iterator<char>)(rbuf + rdispls[0]), (ConstIterator<char>)sbuf, scount * sizeof(SType));
  169. #endif
  170. }
  171. template <class SType, class RType> void Comm::Alltoall(ConstIterator<SType> sbuf, Long scount, Iterator<RType> rbuf, Long rcount) const {
  172. static_assert(std::is_trivially_copyable<SType>::value, "Data is not trivially copyable!");
  173. static_assert(std::is_trivially_copyable<RType>::value, "Data is not trivially copyable!");
  174. #ifdef SCTL_HAVE_MPI
  175. if (scount) {
  176. sbuf[0];
  177. sbuf[scount * Size() - 1];
  178. }
  179. if (rcount) {
  180. rbuf[0];
  181. rbuf[rcount * Size() - 1];
  182. }
  183. MPI_Alltoall((scount ? &sbuf[0] : nullptr), scount, CommDatatype<SType>::value(), (rcount ? &rbuf[0] : nullptr), rcount, CommDatatype<RType>::value(), mpi_comm_);
  184. #else
  185. memcopy((Iterator<char>)rbuf, (ConstIterator<char>)sbuf, scount * sizeof(SType));
  186. #endif
  187. }
  188. template <class SType, class RType> void* Comm::Ialltoallv_sparse(ConstIterator<SType> sbuf, ConstIterator<Long> scounts, ConstIterator<Long> sdispls, Iterator<RType> rbuf, ConstIterator<Long> rcounts, ConstIterator<Long> rdispls, Integer tag) const {
  189. static_assert(std::is_trivially_copyable<SType>::value, "Data is not trivially copyable!");
  190. static_assert(std::is_trivially_copyable<RType>::value, "Data is not trivially copyable!");
  191. #ifdef SCTL_HAVE_MPI
  192. Integer request_count = 0;
  193. for (Integer i = 0; i < mpi_size_; i++) {
  194. if (rcounts[i]) request_count++;
  195. if (scounts[i]) request_count++;
  196. }
  197. if (!request_count) return nullptr;
  198. Vector<MPI_Request>& request = *NewReq();
  199. request.ReInit(request_count);
  200. Integer request_iter = 0;
  201. for (Integer i = 0; i < mpi_size_; i++) {
  202. if (rcounts[i]) {
  203. rbuf[rdispls[i]];
  204. rbuf[rdispls[i] + rcounts[i] - 1];
  205. MPI_Irecv(&rbuf[rdispls[i]], rcounts[i], CommDatatype<RType>::value(), i, tag, mpi_comm_, &request[request_iter]);
  206. request_iter++;
  207. }
  208. }
  209. for (Integer i = 0; i < mpi_size_; i++) {
  210. if (scounts[i]) {
  211. sbuf[sdispls[i]];
  212. sbuf[sdispls[i] + scounts[i] - 1];
  213. MPI_Isend(&sbuf[sdispls[i]], scounts[i], CommDatatype<SType>::value(), i, tag, mpi_comm_, &request[request_iter]);
  214. request_iter++;
  215. }
  216. }
  217. return &request;
  218. #else
  219. memcopy((Iterator<char>)(rbuf + rdispls[0]), (ConstIterator<char>)(sbuf + sdispls[0]), scounts[0] * sizeof(SType));
  220. return nullptr;
  221. #endif
  222. }
  223. template <class Type> void Comm::Alltoallv(ConstIterator<Type> sbuf, ConstIterator<Long> scounts, ConstIterator<Long> sdispls, Iterator<Type> rbuf, ConstIterator<Long> rcounts, ConstIterator<Long> rdispls) const {
  224. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  225. #ifdef SCTL_HAVE_MPI
  226. { // Use Alltoallv_sparse of average connectivity<64
  227. Long connectivity = 0, glb_connectivity = 0;
  228. #pragma omp parallel for schedule(static) reduction(+ : connectivity)
  229. for (Integer i = 0; i < mpi_size_; i++) {
  230. if (rcounts[i]) connectivity++;
  231. }
  232. Allreduce(Ptr2ConstItr<Long>(&connectivity, 1), Ptr2Itr<Long>(&glb_connectivity, 1), 1, CommOp::SUM);
  233. if (glb_connectivity < 64 * Size()) {
  234. void* mpi_req = Ialltoallv_sparse(sbuf, scounts, sdispls, rbuf, rcounts, rdispls, 0);
  235. Wait(mpi_req);
  236. return;
  237. }
  238. }
  239. { // Use vendor MPI_Alltoallv
  240. //#ifndef ALLTOALLV_FIX
  241. Vector<int> scnt, sdsp, rcnt, rdsp;
  242. scnt.ReInit(mpi_size_);
  243. sdsp.ReInit(mpi_size_);
  244. rcnt.ReInit(mpi_size_);
  245. rdsp.ReInit(mpi_size_);
  246. Long stotal = 0, rtotal = 0;
  247. #pragma omp parallel for schedule(static) reduction(+ : stotal, rtotal)
  248. for (Integer i = 0; i < mpi_size_; i++) {
  249. scnt[i] = scounts[i];
  250. sdsp[i] = sdispls[i];
  251. rcnt[i] = rcounts[i];
  252. rdsp[i] = rdispls[i];
  253. stotal += scounts[i];
  254. rtotal += rcounts[i];
  255. }
  256. MPI_Alltoallv((stotal ? &sbuf[0] : nullptr), &scnt[0], &sdsp[0], CommDatatype<Type>::value(), (rtotal ? &rbuf[0] : nullptr), &rcnt[0], &rdsp[0], CommDatatype<Type>::value(), mpi_comm_);
  257. return;
  258. //#endif
  259. }
  260. // TODO: implement hypercube scheme
  261. #else
  262. memcopy((Iterator<char>)(rbuf + rdispls[0]), (ConstIterator<char>)(sbuf + sdispls[0]), scounts[0] * sizeof(Type));
  263. #endif
  264. }
  265. template <class Type> void Comm::Allreduce(ConstIterator<Type> sbuf, Iterator<Type> rbuf, Long count, CommOp op) const {
  266. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  267. #ifdef SCTL_HAVE_MPI
  268. if (!count) return;
  269. MPI_Op mpi_op;
  270. switch (op) {
  271. case CommOp::SUM:
  272. mpi_op = CommDatatype<Type>::sum();
  273. break;
  274. case CommOp::MIN:
  275. mpi_op = CommDatatype<Type>::min();
  276. break;
  277. case CommOp::MAX:
  278. mpi_op = CommDatatype<Type>::max();
  279. break;
  280. default:
  281. mpi_op = MPI_OP_NULL;
  282. break;
  283. }
  284. sbuf[0];
  285. sbuf[count - 1];
  286. rbuf[0];
  287. rbuf[count - 1];
  288. MPI_Allreduce(&sbuf[0], &rbuf[0], count, CommDatatype<Type>::value(), mpi_op, mpi_comm_);
  289. #else
  290. memcopy((Iterator<char>)rbuf, (ConstIterator<char>)sbuf, count * sizeof(Type));
  291. #endif
  292. }
  293. template <class Type> void Comm::Scan(ConstIterator<Type> sbuf, Iterator<Type> rbuf, int count, CommOp op) const {
  294. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  295. #ifdef SCTL_HAVE_MPI
  296. if (!count) return;
  297. MPI_Op mpi_op;
  298. switch (op) {
  299. case CommOp::SUM:
  300. mpi_op = CommDatatype<Type>::sum();
  301. break;
  302. case CommOp::MIN:
  303. mpi_op = CommDatatype<Type>::min();
  304. break;
  305. case CommOp::MAX:
  306. mpi_op = CommDatatype<Type>::max();
  307. break;
  308. default:
  309. mpi_op = MPI_OP_NULL;
  310. break;
  311. }
  312. sbuf[0];
  313. sbuf[count - 1];
  314. rbuf[0];
  315. rbuf[count - 1];
  316. MPI_Scan(&sbuf[0], &rbuf[0], count, CommDatatype<Type>::value(), mpi_op, mpi_comm_);
  317. #else
  318. memcopy((Iterator<char>)rbuf, (ConstIterator<char>)sbuf, count * sizeof(Type));
  319. #endif
  320. }
  321. template <class Type> void Comm::PartitionW(Vector<Type>& nodeList, const Vector<Long>* wts_) const {
  322. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  323. Integer npes = Size();
  324. if (npes == 1) return;
  325. Long nlSize = nodeList.Dim();
  326. Vector<Long> wts;
  327. Long localWt = 0;
  328. if (wts_ == nullptr) { // Construct arrays of wts.
  329. wts.ReInit(nlSize);
  330. #pragma omp parallel for schedule(static)
  331. for (Long i = 0; i < nlSize; i++) {
  332. wts[i] = 1;
  333. }
  334. localWt = nlSize;
  335. } else {
  336. wts.ReInit(nlSize, (Iterator<Long>)wts_->begin(), false);
  337. #pragma omp parallel for reduction(+ : localWt)
  338. for (Long i = 0; i < nlSize; i++) {
  339. localWt += wts[i];
  340. }
  341. }
  342. Long off1 = 0, off2 = 0, totalWt = 0;
  343. { // compute the total weight of the problem ...
  344. Allreduce<Long>(Ptr2ConstItr<Long>(&localWt, 1), Ptr2Itr<Long>(&totalWt, 1), 1, CommOp::SUM);
  345. Scan<Long>(Ptr2ConstItr<Long>(&localWt, 1), Ptr2Itr<Long>(&off2, 1), 1, CommOp::SUM);
  346. off1 = off2 - localWt;
  347. }
  348. Vector<Long> lscn;
  349. if (nlSize) { // perform a local scan on the weights first ...
  350. lscn.ReInit(nlSize);
  351. lscn[0] = off1;
  352. omp_par::scan(wts.begin(), lscn.begin(), nlSize);
  353. }
  354. Vector<Long> sendSz, recvSz, sendOff, recvOff;
  355. sendSz.ReInit(npes);
  356. recvSz.ReInit(npes);
  357. sendOff.ReInit(npes);
  358. recvOff.ReInit(npes);
  359. sendSz.SetZero();
  360. if (nlSize > 0 && totalWt > 0) { // Compute sendSz
  361. Long pid1 = (off1 * npes) / totalWt;
  362. Long pid2 = ((off2 + 1) * npes) / totalWt + 1;
  363. assert((totalWt * pid2) / npes >= off2);
  364. pid1 = (pid1 < 0 ? 0 : pid1);
  365. pid2 = (pid2 > npes ? npes : pid2);
  366. #pragma omp parallel for schedule(static)
  367. for (Integer i = pid1; i < pid2; i++) {
  368. Long wt1 = (totalWt * (i)) / npes;
  369. Long wt2 = (totalWt * (i + 1)) / npes;
  370. Long start = std::lower_bound(lscn.begin(), lscn.begin() + nlSize, wt1, std::less<Long>()) - lscn.begin();
  371. Long end = std::lower_bound(lscn.begin(), lscn.begin() + nlSize, wt2, std::less<Long>()) - lscn.begin();
  372. if (i == 0) start = 0;
  373. if (i == npes - 1) end = nlSize;
  374. sendSz[i] = end - start;
  375. }
  376. } else {
  377. sendSz[0] = nlSize;
  378. }
  379. // Exchange sendSz, recvSz
  380. Alltoall<Long>(sendSz.begin(), 1, recvSz.begin(), 1);
  381. { // Compute sendOff, recvOff
  382. sendOff[0] = 0;
  383. omp_par::scan(sendSz.begin(), sendOff.begin(), npes);
  384. recvOff[0] = 0;
  385. omp_par::scan(recvSz.begin(), recvOff.begin(), npes);
  386. assert(sendOff[npes - 1] + sendSz[npes - 1] == nlSize);
  387. }
  388. // perform All2All ...
  389. Vector<Type> newNodes;
  390. newNodes.ReInit(recvSz[npes - 1] + recvOff[npes - 1]);
  391. void* mpi_req = Ialltoallv_sparse<Type>(nodeList.begin(), sendSz.begin(), sendOff.begin(), newNodes.begin(), recvSz.begin(), recvOff.begin());
  392. Wait(mpi_req);
  393. // reset the pointer ...
  394. nodeList.Swap(newNodes);
  395. }
  396. template <class Type> void Comm::PartitionN(Vector<Type>& v, Long N) const {
  397. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  398. Integer rank = Rank();
  399. Integer np = Size();
  400. if (np == 1) return;
  401. Vector<Long> v_cnt(np), v_dsp(np + 1);
  402. Vector<Long> N_cnt(np), N_dsp(np + 1);
  403. { // Set v_cnt, v_dsp
  404. v_dsp[0] = 0;
  405. Long cnt = v.Dim();
  406. Allgather(Ptr2ConstItr<Long>(&cnt, 1), 1, v_cnt.begin(), 1);
  407. omp_par::scan(v_cnt.begin(), v_dsp.begin(), np);
  408. v_dsp[np] = v_cnt[np - 1] + v_dsp[np - 1];
  409. }
  410. { // Set N_cnt, N_dsp
  411. N_dsp[0] = 0;
  412. Long cnt = N;
  413. Allgather(Ptr2ConstItr<Long>(&cnt, 1), 1, N_cnt.begin(), 1);
  414. omp_par::scan(N_cnt.begin(), N_dsp.begin(), np);
  415. N_dsp[np] = N_cnt[np - 1] + N_dsp[np - 1];
  416. }
  417. { // Adjust for dof
  418. Long dof = (N_dsp[np] ? v_dsp[np] / N_dsp[np] : 0);
  419. assert(dof * N_dsp[np] == v_dsp[np]);
  420. if (dof == 0) return;
  421. if (dof != 1) {
  422. #pragma omp parallel for schedule(static)
  423. for (Integer i = 0; i < np; i++) N_cnt[i] *= dof;
  424. #pragma omp parallel for schedule(static)
  425. for (Integer i = 0; i <= np; i++) N_dsp[i] *= dof;
  426. }
  427. }
  428. Vector<Type> v_(N_cnt[rank]);
  429. { // Set v_
  430. Vector<Long> scnt(np), sdsp(np);
  431. Vector<Long> rcnt(np), rdsp(np);
  432. #pragma omp parallel for schedule(static)
  433. for (Integer i = 0; i < np; i++) {
  434. { // Set scnt
  435. Long n0 = N_dsp[i + 0];
  436. Long n1 = N_dsp[i + 1];
  437. if (n0 < v_dsp[rank + 0]) n0 = v_dsp[rank + 0];
  438. if (n1 < v_dsp[rank + 0]) n1 = v_dsp[rank + 0];
  439. if (n0 > v_dsp[rank + 1]) n0 = v_dsp[rank + 1];
  440. if (n1 > v_dsp[rank + 1]) n1 = v_dsp[rank + 1];
  441. scnt[i] = n1 - n0;
  442. }
  443. { // Set rcnt
  444. Long n0 = v_dsp[i + 0];
  445. Long n1 = v_dsp[i + 1];
  446. if (n0 < N_dsp[rank + 0]) n0 = N_dsp[rank + 0];
  447. if (n1 < N_dsp[rank + 0]) n1 = N_dsp[rank + 0];
  448. if (n0 > N_dsp[rank + 1]) n0 = N_dsp[rank + 1];
  449. if (n1 > N_dsp[rank + 1]) n1 = N_dsp[rank + 1];
  450. rcnt[i] = n1 - n0;
  451. }
  452. }
  453. sdsp[0] = 0;
  454. omp_par::scan(scnt.begin(), sdsp.begin(), np);
  455. rdsp[0] = 0;
  456. omp_par::scan(rcnt.begin(), rdsp.begin(), np);
  457. void* mpi_request = Ialltoallv_sparse(v.begin(), scnt.begin(), sdsp.begin(), v_.begin(), rcnt.begin(), rdsp.begin());
  458. Wait(mpi_request);
  459. }
  460. v.Swap(v_);
  461. }
  462. template <class Type> void Comm::PartitionS(Vector<Type>& nodeList, const Type& splitter) const {
  463. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  464. Integer npes = Size();
  465. if (npes == 1) return;
  466. Vector<Type> mins(npes);
  467. Allgather(Ptr2ConstItr<Type>(&splitter, 1), 1, mins.begin(), 1);
  468. Vector<Long> scnt(npes), sdsp(npes);
  469. Vector<Long> rcnt(npes), rdsp(npes);
  470. { // Compute scnt, sdsp
  471. #pragma omp parallel for schedule(static)
  472. for (Integer i = 0; i < npes; i++) {
  473. sdsp[i] = std::lower_bound(nodeList.begin(), nodeList.begin() + nodeList.Dim(), mins[i]) - nodeList.begin();
  474. }
  475. #pragma omp parallel for schedule(static)
  476. for (Integer i = 0; i < npes - 1; i++) {
  477. scnt[i] = sdsp[i + 1] - sdsp[i];
  478. }
  479. scnt[npes - 1] = nodeList.Dim() - sdsp[npes - 1];
  480. }
  481. { // Compute rcnt, rdsp
  482. rdsp[0] = 0;
  483. Alltoall(scnt.begin(), 1, rcnt.begin(), 1);
  484. omp_par::scan(rcnt.begin(), rdsp.begin(), npes);
  485. }
  486. { // Redistribute nodeList
  487. Vector<Type> nodeList_(rdsp[npes - 1] + rcnt[npes - 1]);
  488. void* mpi_request = Ialltoallv_sparse(nodeList.begin(), scnt.begin(), sdsp.begin(), nodeList_.begin(), rcnt.begin(), rdsp.begin());
  489. Wait(mpi_request);
  490. nodeList.Swap(nodeList_);
  491. }
  492. }
  493. template <class Type> void Comm::SortScatterIndex(const Vector<Type>& key, Vector<Long>& scatter_index, const Type* split_key_) const {
  494. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  495. typedef SortPair<Type, Long> Pair_t;
  496. Integer npes = Size(), rank = Rank();
  497. Vector<Pair_t> parray(key.Dim());
  498. { // Build global index.
  499. Long glb_dsp = 0;
  500. Long loc_size = key.Dim();
  501. Scan(Ptr2ConstItr<Long>(&loc_size, 1), Ptr2Itr<Long>(&glb_dsp, 1), 1, CommOp::SUM);
  502. glb_dsp -= loc_size;
  503. #pragma omp parallel for schedule(static)
  504. for (Long i = 0; i < loc_size; i++) {
  505. parray[i].key = key[i];
  506. parray[i].data = glb_dsp + i;
  507. }
  508. }
  509. Vector<Pair_t> psorted;
  510. HyperQuickSort(parray, psorted);
  511. if (npes > 1 && split_key_ != nullptr) { // Partition data
  512. Vector<Type> split_key(npes);
  513. Allgather(Ptr2ConstItr<Type>(split_key_, 1), 1, split_key.begin(), 1);
  514. Vector<Long> sendSz(npes);
  515. Vector<Long> recvSz(npes);
  516. Vector<Long> sendOff(npes);
  517. Vector<Long> recvOff(npes);
  518. Long nlSize = psorted.Dim();
  519. sendSz.SetZero();
  520. if (nlSize > 0) { // Compute sendSz
  521. // Determine processor range.
  522. Long pid1 = std::lower_bound(split_key.begin(), split_key.begin() + npes, psorted[0].key) - split_key.begin() - 1;
  523. Long pid2 = std::upper_bound(split_key.begin(), split_key.begin() + npes, psorted[nlSize - 1].key) - split_key.begin() + 1;
  524. pid1 = (pid1 < 0 ? 0 : pid1);
  525. pid2 = (pid2 > npes ? npes : pid2);
  526. #pragma omp parallel for schedule(static)
  527. for (Integer i = pid1; i < pid2; i++) {
  528. Pair_t p1;
  529. p1.key = split_key[i];
  530. Pair_t p2;
  531. p2.key = split_key[i + 1 < npes ? i + 1 : i];
  532. Long start = std::lower_bound(psorted.begin(), psorted.begin() + nlSize, p1, std::less<Pair_t>()) - psorted.begin();
  533. Long end = std::lower_bound(psorted.begin(), psorted.begin() + nlSize, p2, std::less<Pair_t>()) - psorted.begin();
  534. if (i == 0) start = 0;
  535. if (i == npes - 1) end = nlSize;
  536. sendSz[i] = end - start;
  537. }
  538. }
  539. // Exchange sendSz, recvSz
  540. Alltoall<Long>(sendSz.begin(), 1, recvSz.begin(), 1);
  541. // compute offsets ...
  542. { // Compute sendOff, recvOff
  543. sendOff[0] = 0;
  544. omp_par::scan(sendSz.begin(), sendOff.begin(), npes);
  545. recvOff[0] = 0;
  546. omp_par::scan(recvSz.begin(), recvOff.begin(), npes);
  547. assert(sendOff[npes - 1] + sendSz[npes - 1] == nlSize);
  548. }
  549. // perform All2All ...
  550. Vector<Pair_t> newNodes(recvSz[npes - 1] + recvOff[npes - 1]);
  551. void* mpi_req = Ialltoallv_sparse<Pair_t>(psorted.begin(), sendSz.begin(), sendOff.begin(), newNodes.begin(), recvSz.begin(), recvOff.begin());
  552. Wait(mpi_req);
  553. // reset the pointer ...
  554. psorted.Swap(newNodes);
  555. }
  556. scatter_index.ReInit(psorted.Dim());
  557. #pragma omp parallel for schedule(static)
  558. for (Long i = 0; i < psorted.Dim(); i++) {
  559. scatter_index[i] = psorted[i].data;
  560. }
  561. }
  562. template <class Type> void Comm::ScatterForward(Vector<Type>& data_, const Vector<Long>& scatter_index) const {
  563. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  564. typedef SortPair<Long, Long> Pair_t;
  565. Integer npes = Size(), rank = Rank();
  566. Long data_dim = 0;
  567. Long send_size = 0;
  568. Long recv_size = 0;
  569. { // Set data_dim, send_size, recv_size
  570. recv_size = scatter_index.Dim();
  571. StaticArray<Long, 2> glb_size;
  572. StaticArray<Long, 2> loc_size;
  573. loc_size[0] = data_.Dim();
  574. loc_size[1] = recv_size;
  575. Allreduce(loc_size, glb_size, 2, CommOp::SUM);
  576. if (glb_size[0] == 0 || glb_size[1] == 0) return; // Nothing to be done.
  577. data_dim = glb_size[0] / glb_size[1];
  578. SCTL_ASSERT(glb_size[0] == data_dim * glb_size[1]);
  579. send_size = data_.Dim() / data_dim;
  580. }
  581. if (npes == 1) { // Scatter directly
  582. Vector<Type> data;
  583. data.ReInit(recv_size * data_dim);
  584. #pragma omp parallel for schedule(static)
  585. for (Long i = 0; i < recv_size; i++) {
  586. Long src_indx = scatter_index[i] * data_dim;
  587. Long trg_indx = i * data_dim;
  588. for (Long j = 0; j < data_dim; j++) data[trg_indx + j] = data_[src_indx + j];
  589. }
  590. data_.Swap(data);
  591. return;
  592. }
  593. Vector<Long> glb_scan;
  594. { // Global scan of data size.
  595. glb_scan.ReInit(npes);
  596. Long glb_rank = 0;
  597. Scan(Ptr2ConstItr<Long>(&send_size, 1), Ptr2Itr<Long>(&glb_rank, 1), 1, CommOp::SUM);
  598. glb_rank -= send_size;
  599. Allgather(Ptr2ConstItr<Long>(&glb_rank, 1), 1, glb_scan.begin(), 1);
  600. }
  601. Vector<Pair_t> psorted;
  602. { // Sort scatter_index.
  603. psorted.ReInit(recv_size);
  604. #pragma omp parallel for schedule(static)
  605. for (Long i = 0; i < recv_size; i++) {
  606. psorted[i].key = scatter_index[i];
  607. psorted[i].data = i;
  608. }
  609. omp_par::merge_sort(psorted.begin(), psorted.begin() + recv_size);
  610. }
  611. Vector<Long> recv_indx(recv_size);
  612. Vector<Long> send_indx(send_size);
  613. Vector<Long> sendSz(npes);
  614. Vector<Long> sendOff(npes);
  615. Vector<Long> recvSz(npes);
  616. Vector<Long> recvOff(npes);
  617. { // Exchange send, recv indices.
  618. #pragma omp parallel for schedule(static)
  619. for (Long i = 0; i < recv_size; i++) {
  620. recv_indx[i] = psorted[i].key;
  621. }
  622. #pragma omp parallel for schedule(static)
  623. for (Integer i = 0; i < npes; i++) {
  624. Long start = std::lower_bound(recv_indx.begin(), recv_indx.begin() + recv_size, glb_scan[i]) - recv_indx.begin();
  625. Long end = (i + 1 < npes ? std::lower_bound(recv_indx.begin(), recv_indx.begin() + recv_size, glb_scan[i + 1]) - recv_indx.begin() : recv_size);
  626. recvSz[i] = end - start;
  627. recvOff[i] = start;
  628. }
  629. Alltoall(recvSz.begin(), 1, sendSz.begin(), 1);
  630. sendOff[0] = 0;
  631. omp_par::scan(sendSz.begin(), sendOff.begin(), npes);
  632. assert(sendOff[npes - 1] + sendSz[npes - 1] == send_size);
  633. Alltoallv(recv_indx.begin(), recvSz.begin(), recvOff.begin(), send_indx.begin(), sendSz.begin(), sendOff.begin());
  634. #pragma omp parallel for schedule(static)
  635. for (Long i = 0; i < send_size; i++) {
  636. assert(send_indx[i] >= glb_scan[rank]);
  637. send_indx[i] -= glb_scan[rank];
  638. assert(send_indx[i] < send_size);
  639. }
  640. }
  641. Vector<Type> send_buff;
  642. { // Prepare send buffer
  643. send_buff.ReInit(send_size * data_dim);
  644. ConstIterator<Type> data = data_.begin();
  645. #pragma omp parallel for schedule(static)
  646. for (Long i = 0; i < send_size; i++) {
  647. Long src_indx = send_indx[i] * data_dim;
  648. Long trg_indx = i * data_dim;
  649. for (Long j = 0; j < data_dim; j++) send_buff[trg_indx + j] = data[src_indx + j];
  650. }
  651. }
  652. Vector<Type> recv_buff;
  653. { // All2Allv
  654. recv_buff.ReInit(recv_size * data_dim);
  655. #pragma omp parallel for schedule(static)
  656. for (Integer i = 0; i < npes; i++) {
  657. sendSz[i] *= data_dim;
  658. sendOff[i] *= data_dim;
  659. recvSz[i] *= data_dim;
  660. recvOff[i] *= data_dim;
  661. }
  662. Alltoallv(send_buff.begin(), sendSz.begin(), sendOff.begin(), recv_buff.begin(), recvSz.begin(), recvOff.begin());
  663. }
  664. { // Build output data.
  665. data_.ReInit(recv_size * data_dim);
  666. Iterator<Type> data = data_.begin();
  667. #pragma omp parallel for schedule(static)
  668. for (Long i = 0; i < recv_size; i++) {
  669. Long src_indx = i * data_dim;
  670. Long trg_indx = psorted[i].data * data_dim;
  671. for (Long j = 0; j < data_dim; j++) data[trg_indx + j] = recv_buff[src_indx + j];
  672. }
  673. }
  674. }
  675. template <class Type> void Comm::ScatterReverse(Vector<Type>& data_, const Vector<Long>& scatter_index_, Long loc_size_) const {
  676. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  677. typedef SortPair<Long, Long> Pair_t;
  678. Integer npes = Size(), rank = Rank();
  679. Long data_dim = 0;
  680. Long send_size = 0;
  681. Long recv_size = 0;
  682. { // Set data_dim, send_size, recv_size
  683. recv_size = loc_size_;
  684. StaticArray<Long, 3> glb_size;
  685. StaticArray<Long, 3> loc_size;
  686. loc_size[0] = data_.Dim();
  687. loc_size[1] = scatter_index_.Dim();
  688. loc_size[2] = recv_size;
  689. Allreduce(loc_size, glb_size, 3, CommOp::SUM);
  690. if (glb_size[0] == 0 || glb_size[1] == 0) return; // Nothing to be done.
  691. SCTL_ASSERT(glb_size[0] % glb_size[1] == 0);
  692. data_dim = glb_size[0] / glb_size[1];
  693. SCTL_ASSERT(loc_size[0] % data_dim == 0);
  694. send_size = loc_size[0] / data_dim;
  695. if (glb_size[0] != glb_size[2] * data_dim) {
  696. recv_size = (((rank + 1) * (glb_size[0] / data_dim)) / npes) - ((rank * (glb_size[0] / data_dim)) / npes);
  697. }
  698. }
  699. if (npes == 1) { // Scatter directly
  700. Vector<Type> data;
  701. data.ReInit(recv_size * data_dim);
  702. #pragma omp parallel for schedule(static)
  703. for (Long i = 0; i < recv_size; i++) {
  704. Long src_indx = i * data_dim;
  705. Long trg_indx = scatter_index_[i] * data_dim;
  706. for (Long j = 0; j < data_dim; j++) data[trg_indx + j] = data_[src_indx + j];
  707. }
  708. data_.Swap(data);
  709. return;
  710. }
  711. Vector<Long> scatter_index;
  712. {
  713. StaticArray<Long, 2> glb_rank;
  714. StaticArray<Long, 3> glb_size;
  715. StaticArray<Long, 2> loc_size;
  716. loc_size[0] = data_.Dim() / data_dim;
  717. loc_size[1] = scatter_index_.Dim();
  718. Scan(loc_size, glb_rank, 2, CommOp::SUM);
  719. Allreduce(loc_size, glb_size, 2, CommOp::SUM);
  720. SCTL_ASSERT(glb_size[0] == glb_size[1]);
  721. glb_rank[0] -= loc_size[0];
  722. glb_rank[1] -= loc_size[1];
  723. Vector<Long> glb_scan0(npes + 1);
  724. Vector<Long> glb_scan1(npes + 1);
  725. Allgather(glb_rank + 0, 1, glb_scan0.begin(), 1);
  726. Allgather(glb_rank + 1, 1, glb_scan1.begin(), 1);
  727. glb_scan0[npes] = glb_size[0];
  728. glb_scan1[npes] = glb_size[1];
  729. if (loc_size[0] != loc_size[1] || glb_rank[0] != glb_rank[1]) { // Repartition scatter_index
  730. scatter_index.ReInit(loc_size[0]);
  731. Vector<Long> send_dsp(npes + 1);
  732. Vector<Long> recv_dsp(npes + 1);
  733. #pragma omp parallel for schedule(static)
  734. for (Integer i = 0; i <= npes; i++) {
  735. send_dsp[i] = std::min(std::max(glb_scan0[i], glb_rank[1]), glb_rank[1] + loc_size[1]) - glb_rank[1];
  736. recv_dsp[i] = std::min(std::max(glb_scan1[i], glb_rank[0]), glb_rank[0] + loc_size[0]) - glb_rank[0];
  737. }
  738. // Long commCnt=0;
  739. Vector<Long> send_cnt(npes + 0);
  740. Vector<Long> recv_cnt(npes + 0);
  741. #pragma omp parallel for schedule(static) // reduction(+:commCnt)
  742. for (Integer i = 0; i < npes; i++) {
  743. send_cnt[i] = send_dsp[i + 1] - send_dsp[i];
  744. recv_cnt[i] = recv_dsp[i + 1] - recv_dsp[i];
  745. // if(send_cnt[i] && i!=rank) commCnt++;
  746. // if(recv_cnt[i] && i!=rank) commCnt++;
  747. }
  748. void* mpi_req = Ialltoallv_sparse<Long>(scatter_index_.begin(), send_cnt.begin(), send_dsp.begin(), scatter_index.begin(), recv_cnt.begin(), recv_dsp.begin(), 0);
  749. Wait(mpi_req);
  750. } else {
  751. scatter_index.ReInit(scatter_index_.Dim(), (Iterator<Long>)scatter_index_.begin(), false);
  752. }
  753. }
  754. Vector<Long> glb_scan(npes);
  755. { // Global data size.
  756. Long glb_rank = 0;
  757. Scan(Ptr2ConstItr<Long>(&recv_size, 1), Ptr2Itr<Long>(&glb_rank, 1), 1, CommOp::SUM);
  758. glb_rank -= recv_size;
  759. Allgather(Ptr2ConstItr<Long>(&glb_rank, 1), 1, glb_scan.begin(), 1);
  760. }
  761. Vector<Pair_t> psorted(send_size);
  762. { // Sort scatter_index.
  763. #pragma omp parallel for schedule(static)
  764. for (Long i = 0; i < send_size; i++) {
  765. psorted[i].key = scatter_index[i];
  766. psorted[i].data = i;
  767. }
  768. omp_par::merge_sort(psorted.begin(), psorted.begin() + send_size);
  769. }
  770. Vector<Long> recv_indx(recv_size);
  771. Vector<Long> send_indx(send_size);
  772. Vector<Long> sendSz(npes);
  773. Vector<Long> sendOff(npes);
  774. Vector<Long> recvSz(npes);
  775. Vector<Long> recvOff(npes);
  776. { // Exchange send, recv indices.
  777. #pragma omp parallel for schedule(static)
  778. for (Long i = 0; i < send_size; i++) {
  779. send_indx[i] = psorted[i].key;
  780. }
  781. #pragma omp parallel for schedule(static)
  782. for (Integer i = 0; i < npes; i++) {
  783. Long start = std::lower_bound(send_indx.begin(), send_indx.begin() + send_size, glb_scan[i]) - send_indx.begin();
  784. Long end = (i + 1 < npes ? std::lower_bound(send_indx.begin(), send_indx.begin() + send_size, glb_scan[i + 1]) - send_indx.begin() : send_size);
  785. sendSz[i] = end - start;
  786. sendOff[i] = start;
  787. }
  788. Alltoall(sendSz.begin(), 1, recvSz.begin(), 1);
  789. recvOff[0] = 0;
  790. omp_par::scan(recvSz.begin(), recvOff.begin(), npes);
  791. assert(recvOff[npes - 1] + recvSz[npes - 1] == recv_size);
  792. Alltoallv(send_indx.begin(), sendSz.begin(), sendOff.begin(), recv_indx.begin(), recvSz.begin(), recvOff.begin());
  793. #pragma omp parallel for schedule(static)
  794. for (Long i = 0; i < recv_size; i++) {
  795. assert(recv_indx[i] >= glb_scan[rank]);
  796. recv_indx[i] -= glb_scan[rank];
  797. assert(recv_indx[i] < recv_size);
  798. }
  799. }
  800. Vector<Type> send_buff;
  801. { // Prepare send buffer
  802. send_buff.ReInit(send_size * data_dim);
  803. ConstIterator<Type> data = data_.begin();
  804. #pragma omp parallel for schedule(static)
  805. for (Long i = 0; i < send_size; i++) {
  806. Long src_indx = psorted[i].data * data_dim;
  807. Long trg_indx = i * data_dim;
  808. for (Long j = 0; j < data_dim; j++) send_buff[trg_indx + j] = data[src_indx + j];
  809. }
  810. }
  811. Vector<Type> recv_buff;
  812. { // All2Allv
  813. recv_buff.ReInit(recv_size * data_dim);
  814. #pragma omp parallel for schedule(static)
  815. for (Integer i = 0; i < npes; i++) {
  816. sendSz[i] *= data_dim;
  817. sendOff[i] *= data_dim;
  818. recvSz[i] *= data_dim;
  819. recvOff[i] *= data_dim;
  820. }
  821. Alltoallv(send_buff.begin(), sendSz.begin(), sendOff.begin(), recv_buff.begin(), recvSz.begin(), recvOff.begin());
  822. }
  823. { // Build output data.
  824. data_.ReInit(recv_size * data_dim);
  825. Iterator<Type> data = data_.begin();
  826. #pragma omp parallel for schedule(static)
  827. for (Long i = 0; i < recv_size; i++) {
  828. Long src_indx = i * data_dim;
  829. Long trg_indx = recv_indx[i] * data_dim;
  830. for (Long j = 0; j < data_dim; j++) data[trg_indx + j] = recv_buff[src_indx + j];
  831. }
  832. }
  833. }
  834. #ifdef SCTL_HAVE_MPI
  835. inline Vector<MPI_Request>* Comm::NewReq() const {
  836. if (req.empty()) req.push(new Vector<MPI_Request>);
  837. Vector<MPI_Request>& request = *(Vector<MPI_Request>*)req.top();
  838. req.pop();
  839. return &request;
  840. }
  841. inline void Comm::DelReq(Vector<MPI_Request>* req_ptr) const {
  842. if (req_ptr) req.push(req_ptr);
  843. }
  844. #define HS_MPIDATATYPE(CTYPE, MPITYPE) \
  845. template <> class Comm::CommDatatype<CTYPE> { \
  846. public: \
  847. static MPI_Datatype value() { return MPITYPE; } \
  848. static MPI_Op sum() { return MPI_SUM; } \
  849. static MPI_Op min() { return MPI_MIN; } \
  850. static MPI_Op max() { return MPI_MAX; } \
  851. }
  852. HS_MPIDATATYPE(short, MPI_SHORT);
  853. HS_MPIDATATYPE(int, MPI_INT);
  854. HS_MPIDATATYPE(long, MPI_LONG);
  855. HS_MPIDATATYPE(unsigned short, MPI_UNSIGNED_SHORT);
  856. HS_MPIDATATYPE(unsigned int, MPI_UNSIGNED);
  857. HS_MPIDATATYPE(unsigned long, MPI_UNSIGNED_LONG);
  858. HS_MPIDATATYPE(float, MPI_FLOAT);
  859. HS_MPIDATATYPE(double, MPI_DOUBLE);
  860. HS_MPIDATATYPE(long double, MPI_LONG_DOUBLE);
  861. HS_MPIDATATYPE(long long, MPI_LONG_LONG_INT);
  862. HS_MPIDATATYPE(char, MPI_CHAR);
  863. HS_MPIDATATYPE(unsigned char, MPI_UNSIGNED_CHAR);
  864. #undef HS_MPIDATATYPE
  865. #endif
  866. template <class Type> void Comm::HyperQuickSort(const Vector<Type>& arr_, Vector<Type>& SortedElem) const { // O( ((N/p)+log(p))*(log(N/p)+log(p)) )
  867. static_assert(std::is_trivially_copyable<Type>::value, "Data is not trivially copyable!");
  868. #ifdef SCTL_HAVE_MPI
  869. Integer npes, myrank, omp_p;
  870. { // Get comm size and rank.
  871. npes = Size();
  872. myrank = Rank();
  873. omp_p = omp_get_max_threads();
  874. }
  875. srand(myrank);
  876. Long totSize, nelem = arr_.Dim();
  877. { // Local and global sizes. O(log p)
  878. assert(nelem); // TODO: Check if this is needed.
  879. Allreduce<Long>(Ptr2ConstItr<Long>(&nelem, 1), Ptr2Itr<Long>(&totSize, 1), 1, CommOp::SUM);
  880. }
  881. if (npes == 1) { // SortedElem <--- local_sort(arr_)
  882. SortedElem = arr_;
  883. omp_par::merge_sort(SortedElem.begin(), SortedElem.begin() + nelem);
  884. return;
  885. }
  886. Vector<Type> arr;
  887. { // arr <-- local_sort(arr_)
  888. arr = arr_;
  889. omp_par::merge_sort(arr.begin(), arr.begin() + nelem);
  890. }
  891. Vector<Type> nbuff, nbuff_ext, rbuff, rbuff_ext; // Allocate memory.
  892. MPI_Comm comm = mpi_comm_; // Copy comm
  893. bool free_comm = false; // Flag to free comm.
  894. // Binary split and merge in each iteration.
  895. while (npes > 1 && totSize > 0) { // O(log p) iterations.
  896. Type split_key;
  897. Long totSize_new;
  898. { // Determine split_key. O( log(N/p) + log(p) )
  899. Integer glb_splt_count;
  900. Vector<Type> glb_splitters;
  901. { // Take random splitters. glb_splt_count = const = 100~1000
  902. Integer splt_count;
  903. { // Set splt_coun. O( 1 ) -- Let p * splt_count = t
  904. splt_count = (100 * nelem) / totSize;
  905. if (npes > 100) splt_count = (drand48() * totSize) < (100 * nelem) ? 1 : 0;
  906. if (splt_count > nelem) splt_count = nelem;
  907. }
  908. Vector<Type> splitters(splt_count);
  909. for (Integer i = 0; i < splt_count; i++) {
  910. splitters[i] = arr[rand() % nelem];
  911. }
  912. Vector<Integer> glb_splt_cnts(npes), glb_splt_disp(npes);
  913. { // Set glb_splt_count, glb_splt_cnts, glb_splt_disp
  914. MPI_Allgather(&splt_count, 1, CommDatatype<Integer>::value(), &glb_splt_cnts[0], 1, CommDatatype<Integer>::value(), comm);
  915. glb_splt_disp[0] = 0;
  916. omp_par::scan(glb_splt_cnts.begin(), glb_splt_disp.begin(), npes);
  917. glb_splt_count = glb_splt_cnts[npes - 1] + glb_splt_disp[npes - 1];
  918. SCTL_ASSERT(glb_splt_count);
  919. }
  920. { // Gather all splitters. O( log(p) )
  921. glb_splitters.ReInit(glb_splt_count);
  922. Vector<int> glb_splt_cnts_(npes), glb_splt_disp_(npes);
  923. for (Integer i = 0; i < npes; i++) {
  924. glb_splt_cnts_[i] = glb_splt_cnts[i];
  925. glb_splt_disp_[i] = glb_splt_disp[i];
  926. }
  927. MPI_Allgatherv((splt_count ? &splitters[0] : nullptr), splt_count, CommDatatype<Type>::value(), &glb_splitters[0], &glb_splt_cnts_[0], &glb_splt_disp_[0], CommDatatype<Type>::value(), comm);
  928. }
  929. }
  930. // Determine split key. O( log(N/p) + log(p) )
  931. Vector<Long> lrank(glb_splt_count);
  932. { // Compute local rank
  933. #pragma omp parallel for schedule(static)
  934. for (Integer i = 0; i < glb_splt_count; i++) {
  935. lrank[i] = std::lower_bound(arr.begin(), arr.begin() + nelem, glb_splitters[i]) - arr.begin();
  936. }
  937. }
  938. Vector<Long> grank(glb_splt_count);
  939. { // Compute global rank
  940. MPI_Allreduce(&lrank[0], &grank[0], glb_splt_count, CommDatatype<Long>::value(), CommDatatype<Long>::sum(), comm);
  941. }
  942. { // Determine split_key, totSize_new
  943. ConstIterator<Long> split_disp = grank.begin();
  944. for (Integer i = 0; i < glb_splt_count; i++) {
  945. if (labs(grank[i] - totSize / 2) < labs(*split_disp - totSize / 2)) {
  946. split_disp = grank.begin() + i;
  947. }
  948. }
  949. split_key = glb_splitters[split_disp - grank.begin()];
  950. if (myrank <= (npes - 1) / 2)
  951. totSize_new = split_disp[0];
  952. else
  953. totSize_new = totSize - split_disp[0];
  954. // double err=(((double)*split_disp)/(totSize/2))-1.0;
  955. // if(fabs<double>(err)<0.01 || npes<=16) break;
  956. // else if(!myrank) std::cout<<err<<'\n';
  957. }
  958. }
  959. Integer split_id = (npes - 1) / 2;
  960. { // Split problem into two. O( N/p )
  961. Integer new_p0 = (myrank <= split_id ? 0 : split_id + 1);
  962. Integer cmp_p0 = (myrank > split_id ? 0 : split_id + 1);
  963. Integer partner;
  964. { // Set partner
  965. partner = myrank + cmp_p0 - new_p0;
  966. if (partner >= npes) partner = npes - 1;
  967. assert(partner >= 0);
  968. }
  969. bool extra_partner = (npes % 2 == 1 && npes - 1 == myrank);
  970. Long ssize = 0, lsize = 0;
  971. ConstIterator<Type> sbuff, lbuff;
  972. { // Set ssize, lsize, sbuff, lbuff
  973. Long split_indx = std::lower_bound(arr.begin(), arr.begin() + nelem, split_key) - arr.begin();
  974. ssize = (myrank > split_id ? split_indx : nelem - split_indx);
  975. sbuff = (myrank > split_id ? arr.begin() : arr.begin() + split_indx);
  976. lsize = (myrank <= split_id ? split_indx : nelem - split_indx);
  977. lbuff = (myrank <= split_id ? arr.begin() : arr.begin() + split_indx);
  978. }
  979. Long rsize = 0, ext_rsize = 0;
  980. { // Get rsize, ext_rsize
  981. Long ext_ssize = 0;
  982. MPI_Status status;
  983. MPI_Sendrecv(&ssize, 1, CommDatatype<Long>::value(), partner, 0, &rsize, 1, CommDatatype<Long>::value(), partner, 0, comm, &status);
  984. if (extra_partner) MPI_Sendrecv(&ext_ssize, 1, CommDatatype<Long>::value(), split_id, 0, &ext_rsize, 1, CommDatatype<Long>::value(), split_id, 0, comm, &status);
  985. }
  986. { // Exchange data.
  987. rbuff.ReInit(rsize);
  988. rbuff_ext.ReInit(ext_rsize);
  989. MPI_Status status;
  990. MPI_Sendrecv((ssize ? &sbuff[0] : nullptr), ssize, CommDatatype<Type>::value(), partner, 0, (rsize ? &rbuff[0] : nullptr), rsize, CommDatatype<Type>::value(), partner, 0, comm, &status);
  991. if (extra_partner) MPI_Sendrecv(nullptr, 0, CommDatatype<Type>::value(), split_id, 0, (ext_rsize ? &rbuff_ext[0] : nullptr), ext_rsize, CommDatatype<Type>::value(), split_id, 0, comm, &status);
  992. }
  993. Long nbuff_size = lsize + rsize + ext_rsize;
  994. { // nbuff <-- merge(lbuff, rbuff, rbuff_ext)
  995. nbuff.ReInit(lsize + rsize);
  996. omp_par::merge<ConstIterator<Type>>(lbuff, (lbuff + lsize), rbuff.begin(), rbuff.begin() + rsize, nbuff.begin(), omp_p, std::less<Type>());
  997. if (ext_rsize > 0 && nbuff.Dim() > 0) {
  998. nbuff_ext.ReInit(nbuff_size);
  999. omp_par::merge(nbuff.begin(), nbuff.begin() + (lsize + rsize), rbuff_ext.begin(), rbuff_ext.begin() + ext_rsize, nbuff_ext.begin(), omp_p, std::less<Type>());
  1000. nbuff.Swap(nbuff_ext);
  1001. nbuff_ext.ReInit(0);
  1002. }
  1003. }
  1004. // Copy new data.
  1005. totSize = totSize_new;
  1006. nelem = nbuff_size;
  1007. arr.Swap(nbuff);
  1008. nbuff.ReInit(0);
  1009. }
  1010. { // Split comm. O( log(p) ) ??
  1011. MPI_Comm scomm;
  1012. MPI_Comm_split(comm, myrank <= split_id, myrank, &scomm);
  1013. if (free_comm) MPI_Comm_free(&comm);
  1014. comm = scomm;
  1015. free_comm = true;
  1016. npes = (myrank <= split_id ? split_id + 1 : npes - split_id - 1);
  1017. myrank = (myrank <= split_id ? myrank : myrank - split_id - 1);
  1018. }
  1019. }
  1020. if (free_comm) MPI_Comm_free(&comm);
  1021. SortedElem = arr;
  1022. PartitionW<Type>(SortedElem);
  1023. #else
  1024. SortedElem = arr_;
  1025. std::sort(SortedElem.begin(), SortedElem.begin() + SortedElem.Dim());
  1026. #endif
  1027. }
  1028. } // end namespace