00001
00002
00003
00009
00010
00011
00012
00013
00014
00015 #include "adio.h"
00016 #include "adio_cb_config_list.h"
00017 #include "ad_bgl.h"
00018 #include "ad_bgl_pset.h"
00019 #include "ad_bgl_aggrs.h"
00020 #ifdef AGGREGATION_PROFILE
00021 #include "mpe.h"
00022 #endif
00023
00024 #ifdef USE_DBG_LOGGING
00025 #define AGG_DEBUG 1
00026 #endif
00027
00028
00029
00030 static int aggrsInPsetSize=0;
00031 static int *aggrsInPset=NULL;
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065 static void
00066 ADIOI_BGL_compute_agg_ranklist_serial ( ADIO_File fd,
00067 const ADIOI_BGL_ConfInfo_t *confInfo,
00068 ADIOI_BGL_ProcInfo_t *all_procInfo,
00069 int *aggrsInPset );
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079 int
00080 ADIOI_BGL_gen_agg_ranklist(ADIO_File fd, int n_aggrs_per_pset)
00081 {
00082 int r, s;
00083 ADIOI_BGL_ProcInfo_t *procInfo, *all_procInfo;
00084 ADIOI_BGL_ConfInfo_t *confInfo;
00085
00086 MPI_Comm_size( fd->comm, &s );
00087 MPI_Comm_rank( fd->comm, &r );
00088
00089
00090 confInfo = ADIOI_BGL_ConfInfo_new ();
00091 procInfo = ADIOI_BGL_ProcInfo_new ();
00092 ADIOI_BGL_persInfo_init( confInfo, procInfo, s, r, n_aggrs_per_pset );
00093
00094
00095
00096 all_procInfo = ADIOI_BGL_ProcInfo_new_n (s);
00097 if(s > aggrsInPsetSize)
00098 {
00099 if(aggrsInPset) ADIOI_Free(aggrsInPset);
00100 aggrsInPset = (int *) ADIOI_Malloc (s *sizeof(int));
00101 aggrsInPsetSize = s;
00102 }
00103
00104
00105 MPI_Gather( (void *)procInfo, sizeof(ADIOI_BGL_ProcInfo_t), MPI_BYTE,
00106 (void *)all_procInfo, sizeof(ADIOI_BGL_ProcInfo_t), MPI_BYTE,
00107 0,
00108 fd->comm );
00109
00110
00111 if (r == 0) {
00112 ADIOI_BGL_compute_agg_ranklist_serial (fd, confInfo, all_procInfo, aggrsInPset);
00113
00114 }
00115 ADIOI_BGL_ProcInfo_free (all_procInfo);
00116
00117
00118
00119 ADIOI_cb_bcast_rank_map(fd);
00120
00121
00122 MPI_Bcast( (void *)aggrsInPset,
00123 fd->hints->cb_nodes * sizeof(int), MPI_BYTE,
00124 0,
00125 fd->comm );
00126
00127 ADIOI_BGL_persInfo_free( confInfo, procInfo );
00128 return 0;
00129 }
00130
00131
00132
00133
00134 static int
00135 ADIOI_BGL_select_agg_in_pset (const ADIOI_BGL_ConfInfo_t *confInfo,
00136 ADIOI_BGL_ProcInfo_t *pset_procInfo,
00137 int nCN_in_pset,
00138 int *tmp_ranklist)
00139 {
00140
00141
00142 int i, j, k;
00143
00144
00145 int nAggrs = nCN_in_pset * confInfo->aggRatio;
00146 if (nAggrs < ADIOI_BGL_NAGG_PSET_MIN) nAggrs = ADIOI_BGL_NAGG_PSET_MIN;
00147
00148
00149 if (!confInfo->isVNM)
00150 {
00151 for (i=0; i<nAggrs; i++) tmp_ranklist[i] = pset_procInfo[i].rank;
00152 }
00153
00154
00155 else
00156 {
00157
00158 j = 0;
00159 for (k=0; k < confInfo->cpuidSize; k++){
00160 for (i=0; i< nCN_in_pset ; i++) {
00161 if (pset_procInfo[i].cpuid == k)
00162 tmp_ranklist[j++] = pset_procInfo[i].rank;
00163 if ( j >= nAggrs) break;
00164 }
00165 if ( j >= nAggrs) break;
00166 }
00167 }
00168
00169 return nAggrs;
00170 }
00171
00172
00173
00174
00175
00176
00177 static int
00178 ADIOI_BGL_compute_agg_ranklist_serial_do (const ADIOI_BGL_ConfInfo_t *confInfo,
00179 ADIOI_BGL_ProcInfo_t *all_procInfo,
00180 int *aggrsInPset,
00181 int *tmp_ranklist)
00182 {
00183 int i, j;
00184
00185
00186 int *psetNumList = (int *) ADIOI_Malloc ( confInfo->nProcs * sizeof(int) );
00187
00188
00189
00190 int n_psets=0;
00191 for (i=0; i<confInfo->nProcs; i++) {
00192
00193 ADIOI_BGL_ProcInfo_t *info_p = all_procInfo+i;
00194
00195 int exist = 0;
00196 for (j=n_psets-1; j>=0; j--)
00197 if (info_p->psetNum == psetNumList[j]) { exist=1; break; }
00198
00199 if (!exist) {
00200 psetNumList [n_psets] = info_p->psetNum;
00201 n_psets ++;
00202 }
00203 }
00204
00205
00206
00207
00208 ADIOI_BGL_ProcInfo_t *sorted_procInfo = ADIOI_BGL_ProcInfo_new_n ( n_psets * confInfo->virtualPsetSize );
00209 int *PsetIdx = (int *) ADIOI_Malloc ( n_psets * sizeof(int) );
00210 AD_BGL_assert ( (PsetIdx != NULL) );
00211
00212
00213 for (i=0; i<n_psets; i++) {
00214 PsetIdx[i] = i*confInfo->virtualPsetSize;
00215 }
00216
00217
00218 for (i=0; i<confInfo->nProcs; i++) {
00219 int pset_id = all_procInfo[i].psetNum;
00220
00221 for (j=n_psets-1; j>=0; j--) if (pset_id == psetNumList[j]) break;
00222 AD_BGL_assert ( (j >= 0) );
00223
00224 sorted_procInfo[ PsetIdx[j] ++ ] = all_procInfo[i];
00225 }
00226
00227 ADIOI_Free(psetNumList);
00228
00229
00230 int naggs = 0;
00231 for (i=0; i<n_psets; i++) {
00232
00233
00234 int nCN_in_pset = PsetIdx[i] - i*confInfo->virtualPsetSize;
00235
00236
00237 int local_naggs = ADIOI_BGL_select_agg_in_pset( confInfo,
00238 sorted_procInfo + i*confInfo->virtualPsetSize,
00239 nCN_in_pset,
00240 tmp_ranklist + naggs);
00241 aggrsInPset[i+1] = local_naggs;
00242
00243 naggs += local_naggs;
00244 }
00245 aggrsInPset[0] = n_psets;
00246
00247
00248 ADIOI_Free ( PsetIdx );
00249 ADIOI_BGL_ProcInfo_free ( sorted_procInfo );
00250 return naggs;
00251 }
00252
00253
00254
00255
00256 static void
00257 ADIOI_BGL_compute_agg_ranklist_serial ( ADIO_File fd,
00258 const ADIOI_BGL_ConfInfo_t *confInfo,
00259 ADIOI_BGL_ProcInfo_t *all_procInfo,
00260 int *aggrsInPset )
00261 {
00262 # if AGG_DEBUG
00263 int i;
00264 # endif
00265 int naggs;
00266 int *tmp_ranklist;
00267
00268
00269 tmp_ranklist = (int *) ADIOI_Malloc (confInfo->nProcs * sizeof(int));
00270
00271 # if AGG_DEBUG
00272 for (i=0; i<confInfo->nProcs; i++) {
00273 DBG_FPRINTF(stderr, "\tcpuid %1d, rank = %6d\n", all_procInfo[i].cpuid, all_procInfo[i].rank );
00274 }
00275 # endif
00276
00277 naggs =
00278 ADIOI_BGL_compute_agg_ranklist_serial_do (confInfo, all_procInfo, aggrsInPset, tmp_ranklist);
00279
00280 # define VERIFY 0
00281 # if VERIFY
00282 DBG_FPRINTF(stderr, "\tconfInfo = %3d,%3d,%3d,%3d,%3d,%3d,%.4f; naggs = %d\n",
00283 confInfo->PsetSize ,
00284 confInfo->numPsets ,
00285 confInfo->isVNM ,
00286 confInfo->virtualPsetSize ,
00287 confInfo->nProcs ,
00288 confInfo->nAggrs ,
00289 confInfo->aggRatio ,
00290 naggs );
00291 # endif
00292
00293 # if AGG_DEBUG
00294 for (i=0; i<naggs; i++) {
00295 DBG_FPRINTF(stderr, "\taggr %-4d = %6d\n", i, tmp_ranklist[i] );
00296 }
00297 # endif
00298
00299
00300 if(fd->hints->ranklist != NULL) ADIOI_Free (fd->hints->ranklist);
00301
00302 fd->hints->cb_nodes = naggs;
00303 fd->hints->ranklist = (int *) ADIOI_Malloc (naggs * sizeof(int));
00304 memcpy( fd->hints->ranklist, tmp_ranklist, naggs*sizeof(int) );
00305
00306
00307 ADIOI_Free( tmp_ranklist );
00308 return;
00309 }
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343 int ADIOI_BGL_Calc_aggregator(ADIO_File fd,
00344 ADIO_Offset off,
00345 ADIO_Offset min_off,
00346 ADIO_Offset *len,
00347 ADIO_Offset fd_size,
00348 ADIO_Offset *fd_start,
00349 ADIO_Offset *fd_end)
00350 {
00351 int rank_index, rank;
00352 ADIO_Offset avail_bytes;
00353
00354 AD_BGL_assert ( (off <= fd_end[fd->hints->cb_nodes-1] && off >= min_off && fd_start[0] >= min_off ) );
00355
00356
00357 int ub = fd->hints->cb_nodes;
00358 int lb = 0;
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374 rank_index = fd->hints->cb_nodes / 2;
00375 while ( off < fd_start[rank_index] || off > fd_end[rank_index] ) {
00376 if ( off > fd_end [rank_index] ) {
00377 lb = rank_index;
00378 rank_index = (rank_index + ub) / 2;
00379 }
00380 else
00381 if ( off < fd_start[rank_index] ) {
00382 ub = rank_index;
00383 rank_index = (rank_index + lb) / 2;
00384 }
00385 }
00386
00387
00388
00389 if (rank_index >= fd->hints->cb_nodes || rank_index < 0) {
00390 FPRINTF(stderr, "Error in ADIOI_Calc_aggregator(): rank_index(%d) >= fd->hints->cb_nodes (%d) fd_size=%lld off=%lld\n",
00391 rank_index,fd->hints->cb_nodes,fd_size,off);
00392 MPI_Abort(MPI_COMM_WORLD, 1);
00393 }
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404 avail_bytes = fd_end[rank_index] + 1 - off;
00405 if (avail_bytes < *len && avail_bytes > 0) {
00406
00407
00408 *len = avail_bytes;
00409 }
00410
00411
00412
00413 rank = fd->hints->ranklist[rank_index];
00414
00415 return rank;
00416 }
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432 void ADIOI_BGL_GPFS_Calc_file_domains(ADIO_Offset *st_offsets,
00433 ADIO_Offset *end_offsets,
00434 int nprocs,
00435 int nprocs_for_coll,
00436 ADIO_Offset *min_st_offset_ptr,
00437 ADIO_Offset **fd_start_ptr,
00438 ADIO_Offset **fd_end_ptr,
00439 ADIO_Offset *fd_size_ptr,
00440 void *fs_ptr)
00441 {
00442 ADIO_Offset min_st_offset, max_end_offset, *fd_start, *fd_end, *fd_size;
00443 int i, aggr;
00444
00445 #ifdef AGGREGATION_PROFILE
00446 MPE_Log_event (5004, 0, NULL);
00447 #endif
00448
00449 # if AGG_DEBUG
00450 static char myname[] = "ADIOI_BGL_GPFS_Calc_file_domains";
00451 DBG_FPRINTF(stderr, "%s(%d): %d aggregator(s)\n",
00452 myname,__LINE__,nprocs_for_coll);
00453 # endif
00454 __blksize_t blksize = 1048576;
00455 if(fs_ptr && ((ADIOI_BGL_fs*)fs_ptr)->blksize)
00456 blksize = ((ADIOI_BGL_fs*)fs_ptr)->blksize;
00457 # if AGG_DEBUG
00458 DBG_FPRINTF(stderr,"%s(%d): Blocksize=%ld\n",myname,__LINE__,blksize);
00459 # endif
00460
00461 min_st_offset = st_offsets [0];
00462 max_end_offset = end_offsets[0];
00463 for (i=1; i<nprocs; i++) {
00464 min_st_offset = ADIOI_MIN(min_st_offset, st_offsets[i]);
00465 max_end_offset = ADIOI_MAX(max_end_offset, end_offsets[i]);
00466 }
00467
00468
00469
00470
00471
00472
00473 ADIO_Offset gpfs_ub = (max_end_offset +blksize-1) / blksize * blksize - 1;
00474 ADIO_Offset gpfs_lb = min_st_offset / blksize * blksize;
00475 ADIO_Offset gpfs_ub_rdoff = (max_end_offset +blksize-1) / blksize * blksize - 1 - max_end_offset;
00476 ADIO_Offset gpfs_lb_rdoff = min_st_offset - min_st_offset / blksize * blksize;
00477 ADIO_Offset fd_gpfs_range = gpfs_ub - gpfs_lb + 1;
00478
00479 int naggs = nprocs_for_coll;
00480
00481
00482
00483
00484
00485
00486
00487
00488
00489
00490
00491
00492
00493
00494
00495
00496
00497
00498 fd_size = (ADIO_Offset *) ADIOI_Malloc(nprocs_for_coll * sizeof(ADIO_Offset));
00499 *fd_start_ptr = (ADIO_Offset *) ADIOI_Malloc(nprocs_for_coll * sizeof(ADIO_Offset));
00500 *fd_end_ptr = (ADIO_Offset *) ADIOI_Malloc(nprocs_for_coll * sizeof(ADIO_Offset));
00501 fd_start = *fd_start_ptr;
00502 fd_end = *fd_end_ptr;
00503
00504 ADIO_Offset n_gpfs_blk = fd_gpfs_range / blksize;
00505 ADIO_Offset nb_cn_small = n_gpfs_blk/naggs;
00506 ADIO_Offset naggs_large = n_gpfs_blk - naggs * (n_gpfs_blk/naggs);
00507 ADIO_Offset naggs_small = naggs - naggs_large;
00508
00509
00510
00511
00512
00513 for (i=0; i<naggs; i++)
00514 if (i < naggs_small) fd_size[i] = nb_cn_small * blksize;
00515 else fd_size[i] = (nb_cn_small+1) * blksize;
00516
00517
00518
00519
00520 # if AGG_DEBUG
00521 DBG_FPRINTF(stderr,"%s(%d): "
00522 "gpfs_ub %llu, "
00523 "gpfs_lb %llu, "
00524 "gpfs_ub_rdoff %llu, "
00525 "gpfs_lb_rdoff %llu, "
00526 "fd_gpfs_range %llu, "
00527 "n_gpfs_blk %llu, "
00528 "nb_cn_small %llu, "
00529 "naggs_large %llu, "
00530 "naggs_small %llu, "
00531 "\n",
00532 myname,__LINE__,
00533 gpfs_ub ,
00534 gpfs_lb ,
00535 gpfs_ub_rdoff,
00536 gpfs_lb_rdoff,
00537 fd_gpfs_range,
00538 n_gpfs_blk ,
00539 nb_cn_small ,
00540 naggs_large ,
00541 naggs_small
00542 );
00543 # endif
00544
00545 fd_size[0] -= gpfs_lb_rdoff;
00546 fd_size[naggs-1] -= gpfs_ub_rdoff;
00547
00548
00549 ADIO_Offset offset = min_st_offset;
00550 for (aggr=0; aggr<naggs; aggr++) {
00551 fd_start[aggr] = offset;
00552 fd_end [aggr] = offset + fd_size[aggr] - 1;
00553 offset += fd_size[aggr];
00554 }
00555
00556 *fd_size_ptr = fd_size[0];
00557 *min_st_offset_ptr = min_st_offset;
00558
00559 #ifdef AGGREGATION_PROFILE
00560 MPE_Log_event (5005, 0, NULL);
00561 #endif
00562 ADIOI_Free (fd_size);
00563 }
00564
00565
00566
00567
00568
00569 int ADIOI_BGL_Aggrs_index( ADIO_File fd, int myrank )
00570 {
00571 int i;
00572 for (i=0; i<fd->hints->cb_nodes; i++)
00573 if (fd->hints->ranklist[i] == myrank) return i;
00574 return -1;
00575 }
00576
00577
00578
00579
00580
00581
00582
00583
00584
00585 void ADIOI_BGL_Calc_my_req(ADIO_File fd, ADIO_Offset *offset_list, ADIO_Offset *len_list,
00586 int contig_access_count, ADIO_Offset
00587 min_st_offset, ADIO_Offset *fd_start,
00588 ADIO_Offset *fd_end, ADIO_Offset fd_size,
00589 int nprocs,
00590 int *count_my_req_procs_ptr,
00591 int **count_my_req_per_proc_ptr,
00592 ADIOI_Access **my_req_ptr,
00593 int **buf_idx_ptr)
00594
00595
00596 {
00597 int *count_my_req_per_proc, count_my_req_procs, *buf_idx;
00598 int i, l, proc;
00599 ADIO_Offset fd_len, rem_len, curr_idx, off;
00600 ADIOI_Access *my_req;
00601
00602 #ifdef AGGREGATION_PROFILE
00603 MPE_Log_event (5024, 0, NULL);
00604 #endif
00605
00606 *count_my_req_per_proc_ptr = (int *) ADIOI_Calloc(nprocs,sizeof(int));
00607 count_my_req_per_proc = *count_my_req_per_proc_ptr;
00608
00609
00610
00611
00612
00613 buf_idx = (int *) ADIOI_Malloc(nprocs*sizeof(int));
00614
00615
00616
00617
00618
00619
00620 for (i=0; i < nprocs; i++) buf_idx[i] = -1;
00621
00622
00623
00624
00625 for (i=0; i < contig_access_count; i++) {
00626
00627
00628 if (len_list[i] == 0)
00629 continue;
00630 off = offset_list[i];
00631 fd_len = len_list[i];
00632
00633
00634
00635
00636
00637 proc = ADIOI_BGL_Calc_aggregator(fd, off, min_st_offset, &fd_len, fd_size,
00638 fd_start, fd_end);
00639 count_my_req_per_proc[proc]++;
00640
00641
00642
00643
00644
00645 rem_len = len_list[i] - fd_len;
00646
00647 while (rem_len > 0) {
00648 off += fd_len;
00649 fd_len = rem_len;
00650 proc = ADIOI_BGL_Calc_aggregator(fd, off, min_st_offset, &fd_len,
00651 fd_size, fd_start, fd_end);
00652
00653 count_my_req_per_proc[proc]++;
00654 rem_len -= fd_len;
00655 }
00656 }
00657
00658
00659
00660 *my_req_ptr = (ADIOI_Access *)
00661 ADIOI_Malloc(nprocs*sizeof(ADIOI_Access));
00662 my_req = *my_req_ptr;
00663
00664 count_my_req_procs = 0;
00665 for (i=0; i < nprocs; i++) {
00666 if (count_my_req_per_proc[i]) {
00667 my_req[i].offsets = (ADIO_Offset *)
00668 ADIOI_Malloc(count_my_req_per_proc[i] * sizeof(ADIO_Offset));
00669 my_req[i].lens = (int *)
00670 ADIOI_Malloc(count_my_req_per_proc[i] * sizeof(int));
00671 count_my_req_procs++;
00672 }
00673 my_req[i].count = 0;
00674
00675 }
00676
00677
00678 curr_idx = 0;
00679 for (i=0; i<contig_access_count; i++) {
00680
00681
00682 if (len_list[i] == 0)
00683 continue;
00684 off = offset_list[i];
00685 fd_len = len_list[i];
00686 proc = ADIOI_BGL_Calc_aggregator(fd, off, min_st_offset, &fd_len, fd_size,
00687 fd_start, fd_end);
00688
00689
00690 if (buf_idx[proc] == -1)
00691 {
00692 ADIOI_Assert(curr_idx == (int) curr_idx);
00693 buf_idx[proc] = (int) curr_idx;
00694 }
00695
00696 l = my_req[proc].count;
00697 curr_idx += fd_len;
00698
00699 rem_len = len_list[i] - fd_len;
00700
00701
00702
00703
00704
00705
00706 my_req[proc].offsets[l] = off;
00707 ADIOI_Assert(fd_len == (int) fd_len);
00708 my_req[proc].lens[l] = (int) fd_len;
00709 my_req[proc].count++;
00710
00711 while (rem_len > 0) {
00712 off += fd_len;
00713 fd_len = rem_len;
00714 proc = ADIOI_BGL_Calc_aggregator(fd, off, min_st_offset, &fd_len,
00715 fd_size, fd_start, fd_end);
00716
00717 if (buf_idx[proc] == -1)
00718 {
00719 ADIOI_Assert(curr_idx == (int) curr_idx);
00720 buf_idx[proc] = (int) curr_idx;
00721 }
00722
00723 l = my_req[proc].count;
00724 curr_idx += fd_len;
00725 rem_len -= fd_len;
00726
00727 my_req[proc].offsets[l] = off;
00728 ADIOI_Assert(fd_len == (int) fd_len);
00729 my_req[proc].lens[l] = (int) fd_len;
00730 my_req[proc].count++;
00731 }
00732 }
00733
00734 #ifdef AGG_DEBUG
00735 for (i=0; i<nprocs; i++) {
00736 if (count_my_req_per_proc[i] > 0) {
00737 DBG_FPRINTF(stderr, "data needed from %d (count = %d):\n", i,
00738 my_req[i].count);
00739 for (l=0; l < my_req[i].count; l++) {
00740 DBG_FPRINTF(stderr, " off[%d] = %lld, len[%d] = %d\n", l,
00741 my_req[i].offsets[l], l, my_req[i].lens[l]);
00742 }
00743 }
00744 DBG_FPRINTF(stderr, "buf_idx[%d] = 0x%x\n", i, buf_idx[i]);
00745 }
00746 #endif
00747
00748 *count_my_req_procs_ptr = count_my_req_procs;
00749 *buf_idx_ptr = buf_idx;
00750 #ifdef AGGREGATION_PROFILE
00751 MPE_Log_event (5025, 0, NULL);
00752 #endif
00753 }
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766
00767
00768
00769
00770
00771
00772 void ADIOI_BGL_Calc_others_req(ADIO_File fd, int count_my_req_procs,
00773 int *count_my_req_per_proc,
00774 ADIOI_Access *my_req,
00775 int nprocs, int myrank,
00776 int *count_others_req_procs_ptr,
00777 ADIOI_Access **others_req_ptr)
00778 {
00779
00780
00781
00782
00783
00784
00785
00786
00787 int *count_others_req_per_proc, count_others_req_procs;
00788 int i;
00789 ADIOI_Access *others_req;
00790
00791
00792 int *scounts, *sdispls, *rcounts, *rdispls;
00793
00794
00795
00796
00797
00798
00799 void *sendBufForOffsets=(void*)0xFFFFFFFF,
00800 *sendBufForLens =(void*)0xFFFFFFFF,
00801 *recvBufForOffsets=(void*)0xFFFFFFFF,
00802 *recvBufForLens =(void*)0xFFFFFFFF;
00803
00804
00805 #ifdef AGGREGATION_PROFILE
00806 MPE_Log_event (5026, 0, NULL);
00807 #endif
00808
00809
00810
00811
00812
00813 count_others_req_per_proc = (int *) ADIOI_Malloc(nprocs*sizeof(int));
00814
00815 MPI_Alltoall(count_my_req_per_proc, 1, MPI_INT,
00816 count_others_req_per_proc, 1, MPI_INT, fd->comm);
00817
00818
00819
00820
00821
00822
00823 *others_req_ptr = (ADIOI_Access *)
00824 ADIOI_Malloc(nprocs*sizeof(ADIOI_Access));
00825 others_req = *others_req_ptr;
00826
00827 scounts = ADIOI_Malloc(nprocs*sizeof(int));
00828 sdispls = ADIOI_Malloc(nprocs*sizeof(int));
00829 rcounts = ADIOI_Malloc(nprocs*sizeof(int));
00830 rdispls = ADIOI_Malloc(nprocs*sizeof(int));
00831
00832
00833
00834
00835
00836
00837 count_others_req_procs = 0;
00838 for (i=0; i<nprocs; i++) {
00839 if (count_others_req_per_proc[i]) {
00840 others_req[i].count = count_others_req_per_proc[i];
00841
00842 others_req[i].offsets = (ADIO_Offset *)
00843 ADIOI_Malloc(count_others_req_per_proc[i]*sizeof(ADIO_Offset));
00844 others_req[i].lens = (int *)
00845 ADIOI_Malloc(count_others_req_per_proc[i]*sizeof(int));
00846
00847 if ( (MPIR_Upint)others_req[i].offsets < (MPIR_Upint)recvBufForOffsets )
00848 recvBufForOffsets = others_req[i].offsets;
00849 if ( (MPIR_Upint)others_req[i].lens < (MPIR_Upint)recvBufForLens )
00850 recvBufForLens = others_req[i].lens;
00851
00852 others_req[i].mem_ptrs = (MPI_Aint *)
00853 ADIOI_Malloc(count_others_req_per_proc[i]*sizeof(MPI_Aint));
00854
00855 count_others_req_procs++;
00856 }
00857 else
00858 {
00859 others_req[i].count = 0;
00860 others_req[i].offsets = NULL;
00861 others_req[i].lens = NULL;
00862 }
00863 }
00864
00865 if ( recvBufForOffsets == (void*)0xFFFFFFFF) recvBufForOffsets = NULL;
00866 if ( recvBufForLens == (void*)0xFFFFFFFF) recvBufForLens = NULL;
00867
00868
00869
00870
00871
00872
00873
00874
00875 for (i=0; i<nprocs; i++)
00876 {
00877 if ( (my_req[i].count) &&
00878 ((MPIR_Upint)my_req[i].offsets <= (MPIR_Upint)sendBufForOffsets) )
00879 sendBufForOffsets = my_req[i].offsets;
00880
00881 if ( (my_req[i].count) &&
00882 ((MPIR_Upint)my_req[i].lens <= (MPIR_Upint)sendBufForLens) )
00883 sendBufForLens = my_req[i].lens;
00884 }
00885
00886
00887 if ( sendBufForOffsets == (void*)0xFFFFFFFF) sendBufForOffsets = NULL;
00888 if ( sendBufForLens == (void*)0xFFFFFFFF) sendBufForLens = NULL;
00889
00890
00891 for (i=0; i<nprocs; i++)
00892 {
00893
00894 scounts[i] = count_my_req_per_proc[i];
00895 if ( scounts[i] == 0 )
00896 sdispls[i] = 0;
00897 else
00898 sdispls[i] = (int)
00899 ( ( (MPIR_Upint)my_req[i].offsets -
00900 (MPIR_Upint)sendBufForOffsets ) /
00901 (MPIR_Upint)sizeof(ADIO_Offset) );
00902
00903
00904 rcounts[i] = count_others_req_per_proc[i];
00905 if ( rcounts[i] == 0 )
00906 rdispls[i] = 0;
00907 else
00908 rdispls[i] = (int)
00909 ( ( (MPIR_Upint)others_req[i].offsets -
00910 (MPIR_Upint)recvBufForOffsets ) /
00911 (MPIR_Upint)sizeof(ADIO_Offset) );
00912 }
00913
00914
00915 MPI_Alltoallv(sendBufForOffsets,
00916 scounts, sdispls, ADIO_OFFSET,
00917 recvBufForOffsets,
00918 rcounts, rdispls, ADIO_OFFSET,
00919 fd->comm);
00920
00921
00922
00923
00924
00925 for (i=0; i<nprocs; i++)
00926 {
00927
00928 scounts[i] = count_my_req_per_proc[i];
00929 if ( scounts[i] == 0 )
00930 sdispls[i] = 0;
00931 else
00932 sdispls[i] = (int)
00933 ( ( (MPIR_Upint)my_req[i].lens -
00934 (MPIR_Upint)sendBufForLens ) /
00935 (MPIR_Upint) sizeof(int) );
00936
00937
00938 rcounts[i] = count_others_req_per_proc[i];
00939 if ( rcounts[i] == 0 )
00940 rdispls[i] = 0;
00941 else
00942 rdispls[i] = (int)
00943 ( ( (MPIR_Upint)others_req[i].lens -
00944 (MPIR_Upint)recvBufForLens ) /
00945 (MPIR_Upint) sizeof(int) );
00946 }
00947
00948
00949 MPI_Alltoallv(sendBufForLens,
00950 scounts, sdispls, MPI_INT,
00951 recvBufForLens,
00952 rcounts, rdispls, MPI_INT,
00953 fd->comm);
00954
00955
00956 ADIOI_Free(count_others_req_per_proc);
00957 ADIOI_Free (scounts);
00958 ADIOI_Free (sdispls);
00959 ADIOI_Free (rcounts);
00960 ADIOI_Free (rdispls);
00961
00962 *count_others_req_procs_ptr = count_others_req_procs;
00963 #ifdef AGGREGATION_PROFILE
00964 MPE_Log_event (5027, 0, NULL);
00965 #endif
00966 }