File: ad_iwrite_coll.c

package info (click to toggle)
openmpi 5.0.8-3
  • links: PTS, VCS
  • area: main
  • in suites:
  • size: 201,692 kB
  • sloc: ansic: 613,078; makefile: 42,353; sh: 11,194; javascript: 9,244; f90: 7,052; java: 6,404; perl: 5,179; python: 1,859; lex: 740; fortran: 61; cpp: 20; tcl: 12
file content (1502 lines) | stat: -rw-r--r-- 55,955 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
/*
 * Copyright (C) by Argonne National Laboratory
 *     See COPYRIGHT in top-level directory
 */

#include "adio.h"
#include "adio_extern.h"
#include "mpiu_greq.h"
#include "mpioimpl.h"

#ifdef AGGREGATION_PROFILE
#include "mpe.h"
#endif

#ifdef HAVE_MPI_GREQUEST_EXTENSIONS

/* ADIOI_GEN_IwriteStridedColl */
struct ADIOI_GEN_IwriteStridedColl_vars {
    /* requests */
    MPI_Request req_offset[2];  /* ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL */
    MPI_Request req_ind_io;     /* ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL_INDIO */
    MPI_Request req_err;        /* ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL_BCAST */

    /* parameters */
    ADIO_File fd;
    const void *buf;
    int count;
    MPI_Datatype datatype;
    int file_ptr_type;
    ADIO_Offset offset;

    /* stack variables */
    ADIOI_Access *my_req;
    /* array of nprocs access structures, one for each other process in
     * whose file domain this process's request lies */

    ADIOI_Access *others_req;
    /* array of nprocs access structures, one for each other process
     * whose request lies in this process's file domain. */

    int nprocs;
    int nprocs_for_coll;
    int myrank;
    int contig_access_count;
    int interleave_count;
    int buftype_is_contig;
    int *count_my_req_per_proc;
    int count_my_req_procs;
    int count_others_req_procs;
    ADIO_Offset start_offset;
    ADIO_Offset end_offset;
    ADIO_Offset orig_fp;
    ADIO_Offset fd_size;
    ADIO_Offset min_st_offset;
    ADIO_Offset *offset_list;
    ADIO_Offset *st_offsets;
    ADIO_Offset *fd_start;
    ADIO_Offset *fd_end;
    ADIO_Offset *end_offsets;
    MPI_Aint *buf_idx;
    ADIO_Offset *len_list;
    int old_error;
    int tmp_error;
    int error_code;
};

/* ADIOI_Iexch_and_write */
struct ADIOI_Iexch_and_write_vars {
    /* requests */
    MPI_Request req1;           /* ADIOI_IWC_STATE_IEXCH_AND_WRITE */
    MPI_Request req3;           /* ADIOI_IWC_STATE_IEXCH_AND_WRITE_L1_BODY */

    /* parameters */
    ADIO_File fd;
    void *buf;
    MPI_Datatype datatype;
    int nprocs;
    int myrank;
    ADIOI_Access *others_req;
    ADIO_Offset *offset_list;
    ADIO_Offset *len_list;
    int contig_access_count;
    ADIO_Offset min_st_offset;
    ADIO_Offset fd_size;
    ADIO_Offset *fd_start;
    ADIO_Offset *fd_end;
    MPI_Aint *buf_idx;

    /* stack variables */
    /* Not convinced end_loc-st_loc couldn't be > int, so make these offsets */
    ADIO_Offset size;
    int hole;
    int m;
    int ntimes;
    int max_ntimes;
    int buftype_is_contig;
    ADIO_Offset st_loc;
    ADIO_Offset end_loc;
    ADIO_Offset off;
    ADIO_Offset done;
    char *write_buf;
    int *curr_offlen_ptr;
    int *count;
    int *send_size;
    int *recv_size;
    int *partial_recv;
    int *sent_to_proc;
    int *start_pos;
    int *send_buf_idx;
    int *curr_to_proc;
    int *done_to_proc;
    ADIOI_Flatlist_node *flat_buf;
    MPI_Aint buftype_extent;
    int coll_bufsize;

    /* next function to be called */
    void (*next_fn) (ADIOI_NBC_Request *, int *);
};

/* ADIOI_W_Iexchange_data */
struct ADIOI_W_Iexchange_data_vars {
    /* requests */
    MPI_Request req1;           /* ADIOI_IWC_STATE_W_IEXCHANGE_DATA */
    MPI_Request req2;           /* ADIOI_IWC_STATE_W_IEXCHANGE_DATA_HOLE */
    MPI_Request *req3;          /* ADIOI_IWC_STATE_W_IEXCHANGE_DATA_SEND */

    /* parameters */
    ADIO_File fd;
    void *buf;
    char *write_buf;
    ADIOI_Flatlist_node *flat_buf;
    ADIO_Offset *offset_list;
    ADIO_Offset *len_list;
    int *send_size;
    int *recv_size;
    ADIO_Offset off;
    int size;
    int *count;
    int *start_pos;
    int *partial_recv;
    int *sent_to_proc;
    int nprocs;
    int myrank;
    int buftype_is_contig;
    int contig_access_count;
    ADIO_Offset min_st_offset;
    ADIO_Offset fd_size;
    ADIO_Offset *fd_start;
    ADIO_Offset *fd_end;
    ADIOI_Access *others_req;
    int *send_buf_idx;
    int *curr_to_proc;
    int *done_to_proc;
    int *hole;
    int iter;
    MPI_Aint buftype_extent;
    MPI_Aint *buf_idx;

    /* stack variables */
    int nprocs_recv;
    int nprocs_send;
    int err;
    char **send_buf;
    MPI_Request *requests;
    MPI_Request *send_req;
    MPI_Datatype *recv_types;
    int sum;
    ADIO_Offset *srt_off;

    /* next function to be called */
    void (*next_fn) (ADIOI_NBC_Request *, int *);
};


void ADIOI_Fill_send_buffer(ADIO_File fd, void *buf, ADIOI_Flatlist_node
                            * flat_buf, char **send_buf, ADIO_Offset
                            * offset_list, ADIO_Offset * len_list, int *send_size,
                            MPI_Request * requests, int *sent_to_proc,
                            int nprocs, int myrank,
                            int contig_access_count, ADIO_Offset
                            min_st_offset, ADIO_Offset fd_size,
                            ADIO_Offset * fd_start, ADIO_Offset * fd_end,
                            int *send_buf_idx, int *curr_to_proc,
                            int *done_to_proc, int iter, MPI_Aint buftype_extent);
void ADIOI_Heap_merge(ADIOI_Access * others_req, int *count,
                      ADIO_Offset * srt_off, int *srt_len, int *start_pos,
                      int nprocs, int nprocs_recv, int total_elements);


/* prototypes of functions used for nonblocking collective writes only. */
static void ADIOI_GEN_IwriteStridedColl_inter(ADIOI_NBC_Request *, int *);
static void ADIOI_GEN_IwriteStridedColl_indio(ADIOI_NBC_Request *, int *);
static void ADIOI_GEN_IwriteStridedColl_exch(ADIOI_NBC_Request *, int *);
static void ADIOI_GEN_IwriteStridedColl_bcast(ADIOI_NBC_Request *, int *);
static void ADIOI_GEN_IwriteStridedColl_free(ADIOI_NBC_Request *, int *);
static void ADIOI_GEN_IwriteStridedColl_fini(ADIOI_NBC_Request *, int *);

static void ADIOI_Iexch_and_write(ADIOI_NBC_Request *, int *);
static void ADIOI_Iexch_and_write_l1_begin(ADIOI_NBC_Request *, int *);
static void ADIOI_Iexch_and_write_l1_body(ADIOI_NBC_Request *, int *);
static void ADIOI_Iexch_and_write_l1_end(ADIOI_NBC_Request *, int *);
static void ADIOI_Iexch_and_write_reset(ADIOI_NBC_Request *, int *);
static void ADIOI_Iexch_and_write_l2_begin(ADIOI_NBC_Request *, int *);
static void ADIOI_Iexch_and_write_l2_end(ADIOI_NBC_Request *, int *);
static void ADIOI_Iexch_and_write_fini(ADIOI_NBC_Request *, int *);

static void ADIOI_W_Iexchange_data(ADIOI_NBC_Request *, int *);
static void ADIOI_W_Iexchange_data_hole(ADIOI_NBC_Request *, int *);
static void ADIOI_W_Iexchange_data_send(ADIOI_NBC_Request *, int *);
static void ADIOI_W_Iexchange_data_wait(ADIOI_NBC_Request *, int *);
static void ADIOI_W_Iexchange_data_fini(ADIOI_NBC_Request *, int *);

static MPIX_Grequest_class ADIOI_GEN_greq_class = 0;
static int ADIOI_GEN_iwc_query_fn(void *extra_state, MPI_Status * status);
static int ADIOI_GEN_iwc_free_fn(void *extra_state);
static int ADIOI_GEN_iwc_poll_fn(void *extra_state, MPI_Status * status);
static int ADIOI_GEN_iwc_wait_fn(int count, void **array_of_states,
                                 double timeout, MPI_Status * status);


/* Non-blocking version of ADIOI_GEN_WriteStridedColl() */
void ADIOI_GEN_IwriteStridedColl(ADIO_File fd, const void *buf, int count,
                                 MPI_Datatype datatype, int file_ptr_type,
                                 ADIO_Offset offset, MPI_Request * request, int *error_code)
{
    /* Uses a generalized version of the extended two-phase method described
     * in "An Extended Two-Phase Method for Accessing Sections of
     * Out-of-Core Arrays", Rajeev Thakur and Alok Choudhary,
     * Scientific Programming, (5)4:301--317, Winter 1996.
     * http://www.mcs.anl.gov/home/thakur/ext2ph.ps */

    ADIOI_NBC_Request *nbc_req = NULL;
    ADIOI_GEN_IwriteStridedColl_vars *vars = NULL;
    int nprocs, myrank;

#if 0
    /* FIXME: need an implementation of ADIOI_IOIstridedColl */
    if (fd->hints->cb_pfr != ADIOI_HINT_DISABLE) {
        /* Cast away const'ness as the below function is used for read
         * and write */
        ADIOI_IOIstridedColl(fd, (char *) buf, count, ADIOI_WRITE, datatype,
                             file_ptr_type, offset, request, error_code);
        return;
    }
#endif

    /* top-level struct keeping the status of function progress */
    nbc_req = (ADIOI_NBC_Request *) ADIOI_Calloc(1, sizeof(ADIOI_NBC_Request));
    nbc_req->rdwr = ADIOI_WRITE;

    /* create a generalized request */
    if (ADIOI_GEN_greq_class == 0) {
        MPIX_Grequest_class_create(ADIOI_GEN_iwc_query_fn,
                                   ADIOI_GEN_iwc_free_fn, MPIU_Greq_cancel_fn,
                                   ADIOI_GEN_iwc_poll_fn, ADIOI_GEN_iwc_wait_fn,
                                   &ADIOI_GEN_greq_class);
    }
    MPIX_Grequest_class_allocate(ADIOI_GEN_greq_class, nbc_req, request);
    memcpy(&nbc_req->req, request, sizeof(MPI_Request));

    /* create a struct for parameters and variables */
    vars =
        (ADIOI_GEN_IwriteStridedColl_vars *) ADIOI_Calloc(1,
                                                          sizeof(ADIOI_GEN_IwriteStridedColl_vars));
    nbc_req->data.wr.wsc_vars = vars;

    /* save the parameters */
    vars->fd = fd;
    vars->buf = buf;
    vars->count = count;
    vars->datatype = datatype;
    vars->file_ptr_type = file_ptr_type;
    vars->offset = offset;

    MPI_Comm_size(fd->comm, &nprocs);
    MPI_Comm_rank(fd->comm, &myrank);
    vars->nprocs = nprocs;
    vars->myrank = myrank;

    /* the number of processes that actually perform I/O, nprocs_for_coll,
     * is stored in the hints off the ADIO_File structure
     */
    vars->nprocs_for_coll = fd->hints->cb_nodes;
    vars->orig_fp = fd->fp_ind;

    /* only check for interleaving if cb_write isn't disabled */
    if (fd->hints->cb_write != ADIOI_HINT_DISABLE) {
        /* For this process's request, calculate the list of offsets and
         * lengths in the file and determine the start and end offsets. */

        /* Note: end_offset points to the last byte-offset that will be accessed.
         * e.g., if start_offset=0 and 100 bytes to be read, end_offset=99 */

        ADIOI_Calc_my_off_len(fd, count, datatype, file_ptr_type, offset,
                              &vars->offset_list, &vars->len_list,
                              &vars->start_offset, &vars->end_offset, &vars->contig_access_count);

        /* each process communicates its start and end offsets to other
         * processes. The result is an array each of start and end offsets
         * stored in order of process rank. */

        vars->st_offsets = (ADIO_Offset *) ADIOI_Malloc(nprocs * 2 * sizeof(ADIO_Offset));
        vars->end_offsets = vars->st_offsets + nprocs;

        *error_code = MPI_Iallgather(&vars->start_offset, 1, ADIO_OFFSET,
                                     vars->st_offsets, 1, ADIO_OFFSET,
                                     fd->comm, &vars->req_offset[0]);
        if (*error_code != MPI_SUCCESS)
            return;
        *error_code = MPI_Iallgather(&vars->end_offset, 1, ADIO_OFFSET,
                                     vars->end_offsets, 1, ADIO_OFFSET,
                                     fd->comm, &vars->req_offset[1]);

        nbc_req->data.wr.state = ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL;
        return;
    }

    ADIOI_GEN_IwriteStridedColl_indio(nbc_req, error_code);
}

static void ADIOI_GEN_IwriteStridedColl_inter(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_GEN_IwriteStridedColl_vars *vars = nbc_req->data.wr.wsc_vars;
    int nprocs = vars->nprocs;
    ADIO_Offset *st_offsets = vars->st_offsets;
    ADIO_Offset *end_offsets = vars->end_offsets;
    int i, interleave_count = 0;

    /* are the accesses of different processes interleaved? */
    for (i = 1; i < nprocs; i++)
        if ((st_offsets[i] < end_offsets[i - 1]) && (st_offsets[i] <= end_offsets[i]))
            interleave_count++;
    /* This is a rudimentary check for interleaving, but should suffice
     * for the moment. */

    vars->interleave_count = interleave_count;

    ADIOI_GEN_IwriteStridedColl_indio(nbc_req, error_code);
}

static void ADIOI_GEN_IwriteStridedColl_indio(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_GEN_IwriteStridedColl_vars *vars = nbc_req->data.wr.wsc_vars;
    ADIOI_Icalc_others_req_vars *cor_vars = NULL;
    ADIO_File fd = vars->fd;
    const void *buf;
    int count, file_ptr_type;
    MPI_Datatype datatype = vars->datatype;
    ADIO_Offset offset;
    int filetype_is_contig;
    ADIO_Offset off;
    int nprocs;

    ADIOI_Datatype_iscontig(datatype, &vars->buftype_is_contig);

    if (fd->hints->cb_write == ADIOI_HINT_DISABLE ||
        (!vars->interleave_count && (fd->hints->cb_write == ADIOI_HINT_AUTO))) {
        buf = vars->buf;
        count = vars->count;
        file_ptr_type = vars->file_ptr_type;
        offset = vars->offset;

        /* use independent accesses */
        if (fd->hints->cb_write != ADIOI_HINT_DISABLE) {
            ADIOI_Free(vars->offset_list);
            ADIOI_Free(vars->st_offsets);
        }

        fd->fp_ind = vars->orig_fp;
        ADIOI_Datatype_iscontig(fd->filetype, &filetype_is_contig);

#if defined(ROMIO_RUN_ON_LINUX) && !defined(HAVE_AIO_LITE_H)
        /* NOTE: This is currently a workaround to avoid weird errors, e.g.,
         * stack fault, occurred on Linux. When the host OS is Linux and
         * aio-lite is not used, a blocking ADIO function is used here.
         * See https://trac.mpich.org/projects/mpich/ticket/2201. */
        MPI_Status status;
        if (vars->buftype_is_contig && filetype_is_contig) {
            if (file_ptr_type == ADIO_EXPLICIT_OFFSET) {
                off = fd->disp + (ADIO_Offset) (fd->etype_size) * offset;
                ADIO_WriteContig(fd, buf, count, datatype,
                                 ADIO_EXPLICIT_OFFSET, off, &status, error_code);
            } else
                ADIO_WriteContig(fd, buf, count, datatype, ADIO_INDIVIDUAL, 0, &status, error_code);
        } else {
            ADIO_WriteStrided(fd, buf, count, datatype, file_ptr_type, offset, &status, error_code);
        }
        ADIOI_GEN_IwriteStridedColl_fini(nbc_req, error_code);
#else
        if (vars->buftype_is_contig && filetype_is_contig) {
            if (file_ptr_type == ADIO_EXPLICIT_OFFSET) {
                off = fd->disp + (ADIO_Offset) (fd->etype_size) * offset;
                ADIO_IwriteContig(fd, buf, count, datatype,
                                  ADIO_EXPLICIT_OFFSET, off, &vars->req_ind_io, error_code);
            } else
                ADIO_IwriteContig(fd, buf, count, datatype, ADIO_INDIVIDUAL,
                                  0, &vars->req_ind_io, error_code);
        } else {
            ADIO_IwriteStrided(fd, buf, count, datatype, file_ptr_type,
                               offset, &vars->req_ind_io, error_code);
        }

        nbc_req->data.wr.state = ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL_INDIO;
#endif
        return;
    }

    nprocs = vars->nprocs;

    /* Divide the I/O workload among "nprocs_for_coll" processes. This is
     * done by (logically) dividing the file into file domains (FDs); each
     * process may directly access only its own file domain. */

    ADIOI_Calc_file_domains(vars->st_offsets, vars->end_offsets, nprocs,
                            vars->nprocs_for_coll, &vars->min_st_offset,
                            &vars->fd_start, &vars->fd_end,
                            fd->hints->min_fdomain_size, &vars->fd_size, fd->hints->striping_unit);

    /* calculate what portions of the access requests of this process are
     * located in what file domains */

    ADIOI_Calc_my_req(fd, vars->offset_list, vars->len_list,
                      vars->contig_access_count, vars->min_st_offset,
                      vars->fd_start, vars->fd_end, vars->fd_size,
                      nprocs, &vars->count_my_req_procs,
                      &vars->count_my_req_per_proc, &vars->my_req, &vars->buf_idx);

    /* based on everyone's my_req, calculate what requests of other
     * processes lie in this process's file domain.
     * count_others_req_procs = number of processes whose requests lie in
     * this process's file domain (including this process itself)
     * count_others_req_per_proc[i] indicates how many separate contiguous
     * requests of proc. i lie in this process's file domain. */

    cor_vars = (ADIOI_Icalc_others_req_vars *) ADIOI_Calloc(1, sizeof(ADIOI_Icalc_others_req_vars));
    nbc_req->cor_vars = cor_vars;
    cor_vars->fd = vars->fd;
    cor_vars->count_my_req_procs = vars->count_my_req_procs;
    cor_vars->count_my_req_per_proc = vars->count_my_req_per_proc;
    cor_vars->my_req = vars->my_req;
    cor_vars->nprocs = vars->nprocs;
    cor_vars->myrank = vars->myrank;
    cor_vars->count_others_req_procs_ptr = &vars->count_others_req_procs;
    cor_vars->others_req_ptr = &vars->others_req;
    cor_vars->next_fn = ADIOI_GEN_IwriteStridedColl_exch;

    ADIOI_Icalc_others_req(nbc_req, error_code);
}

static void ADIOI_GEN_IwriteStridedColl_exch(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_GEN_IwriteStridedColl_vars *vars = nbc_req->data.wr.wsc_vars;
    ADIOI_Iexch_and_write_vars *eaw_vars = NULL;
    ADIOI_Access *my_req = vars->my_req;

    ADIOI_Free(vars->count_my_req_per_proc);
    ADIOI_Free(my_req[0].offsets);
    ADIOI_Free(my_req);

    /* exchange data and write in sizes of no more than coll_bufsize. */
    /* Cast away const'ness for the below function */
    eaw_vars = (ADIOI_Iexch_and_write_vars *) ADIOI_Calloc(1, sizeof(ADIOI_Iexch_and_write_vars));
    nbc_req->data.wr.eaw_vars = eaw_vars;
    eaw_vars->fd = vars->fd;
    eaw_vars->buf = (char *) vars->buf;
    eaw_vars->datatype = vars->datatype;
    eaw_vars->nprocs = vars->nprocs;
    eaw_vars->myrank = vars->myrank;
    eaw_vars->others_req = vars->others_req;
    eaw_vars->offset_list = vars->offset_list;
    eaw_vars->len_list = vars->len_list;
    eaw_vars->contig_access_count = vars->contig_access_count;
    eaw_vars->min_st_offset = vars->min_st_offset;
    eaw_vars->fd_size = vars->fd_size;
    eaw_vars->fd_start = vars->fd_start;
    eaw_vars->fd_end = vars->fd_end;
    eaw_vars->buf_idx = vars->buf_idx;
    eaw_vars->next_fn = ADIOI_GEN_IwriteStridedColl_bcast;

    ADIOI_Iexch_and_write(nbc_req, error_code);
}

static void ADIOI_GEN_IwriteStridedColl_bcast(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_GEN_IwriteStridedColl_vars *vars = nbc_req->data.wr.wsc_vars;
    ADIO_File fd = vars->fd;

    /* If this collective write is followed by an independent write,
     * it's possible to have those subsequent writes on other processes
     * race ahead and sneak in before the read-modify-write completes.
     * We carry out a collective communication at the end here so no one
     * can start independent i/o before collective I/O completes.
     *
     * need to do some gymnastics with the error codes so that if something
     * went wrong, all processes report error, but if a process has a more
     * specific error code, we can still have that process report the
     * additional information */

    vars->old_error = *error_code;
    if (*error_code != MPI_SUCCESS)
        *error_code = MPI_ERR_IO;

    /* optimization: if only one process performing i/o, we can perform
     * a less-expensive Bcast  */
#ifdef ADIOI_MPE_LOGGING
    MPE_Log_event(ADIOI_MPE_postwrite_a, 0, NULL);
#endif
    vars->error_code = *error_code;
    if (fd->hints->cb_nodes == 1) {
        *error_code = MPI_Ibcast(&vars->error_code, 1, MPI_INT,
                                 fd->hints->ranklist[0], fd->comm, &vars->req_err);
    } else {
        vars->tmp_error = *error_code;
        *error_code = MPI_Iallreduce(&vars->tmp_error, &vars->error_code, 1,
                                     MPI_INT, MPI_MAX, fd->comm, &vars->req_err);
    }

    nbc_req->data.wr.state = ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL_BCAST;
}

static void ADIOI_GEN_IwriteStridedColl_free(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_GEN_IwriteStridedColl_vars *vars = nbc_req->data.wr.wsc_vars;
    ADIO_File fd = vars->fd;
    ADIOI_Access *others_req = vars->others_req;
    int old_error = vars->old_error;

#ifdef ADIOI_MPE_LOGGING
    MPE_Log_event(ADIOI_MPE_postwrite_b, 0, NULL);
#endif
#ifdef AGGREGATION_PROFILE
    MPE_Log_event(5012, 0, NULL);
#endif

    if ((old_error != MPI_SUCCESS) && (old_error != MPI_ERR_IO))
        *error_code = old_error;

    /* free all memory allocated for collective I/O */
    ADIOI_Free(others_req[0].offsets);
    ADIOI_Free(others_req[0].mem_ptrs);
    ADIOI_Free(others_req);

    ADIOI_Free(vars->buf_idx);
    ADIOI_Free(vars->offset_list);
    ADIOI_Free(vars->st_offsets);
    ADIOI_Free(vars->fd_start);

    fd->fp_sys_posn = -1;       /* set it to null. */
#ifdef AGGREGATION_PROFILE
    MPE_Log_event(5013, 0, NULL);
#endif

    ADIOI_GEN_IwriteStridedColl_fini(nbc_req, error_code);
}

static void ADIOI_GEN_IwriteStridedColl_fini(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_GEN_IwriteStridedColl_vars *vars = nbc_req->data.wr.wsc_vars;
    MPI_Count size;

    /* This is a temporary way of filling in status. The right way is to
     * keep track of how much data was actually written during collective I/O. */
    MPI_Type_size_x(vars->datatype, &size);
    nbc_req->nbytes = size * vars->count;

    /* free the struct for parameters and variables */
    if (nbc_req->data.wr.wsc_vars) {
        ADIOI_Free(nbc_req->data.wr.wsc_vars);
        nbc_req->data.wr.wsc_vars = NULL;
    }

    /* make the request complete */
    *error_code = MPI_Grequest_complete(nbc_req->req);
    nbc_req->data.wr.state = ADIOI_IWC_STATE_COMPLETE;
}


static void ADIOI_Iexch_and_write(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_Iexch_and_write_vars *vars = nbc_req->data.wr.eaw_vars;
    ADIO_File fd = vars->fd;
    MPI_Datatype datatype = vars->datatype;
    int nprocs = vars->nprocs;
    ADIOI_Access *others_req = vars->others_req;
    MPI_Aint lb;

    /* Send data to appropriate processes and write in sizes of no more
     * than coll_bufsize.
     * The idea is to reduce the amount of extra memory required for
     * collective I/O. If all data were written all at once, which is much
     * easier, it would require temp space more than the size of user_buf,
     * which is often unacceptable. For example, to write a distributed
     * array to a file, where each local array is 8Mbytes, requiring
     * at least another 8Mbytes of temp space is unacceptable. */

    int i, j;
    ADIO_Offset st_loc = -1, end_loc = -1;
    int info_flag, coll_bufsize;
    char *value;

    *error_code = MPI_SUCCESS;  /* changed below if error */
    /* only I/O errors are currently reported */

    /* calculate the number of writes of size coll_bufsize
     * to be done by each process and the max among all processes.
     * That gives the no. of communication phases as well. */

    value = (char *) ADIOI_Malloc((MPI_MAX_INFO_VAL + 1) * sizeof(char));
    ADIOI_Info_get(fd->info, "cb_buffer_size", MPI_MAX_INFO_VAL, value, &info_flag);
    coll_bufsize = atoi(value);
    vars->coll_bufsize = coll_bufsize;
    ADIOI_Free(value);

    for (i = 0; i < nprocs; i++) {
        if (others_req[i].count) {
            st_loc = others_req[i].offsets[0];
            end_loc = others_req[i].offsets[0];
            break;
        }
    }

    for (i = 0; i < nprocs; i++)
        for (j = 0; j < others_req[i].count; j++) {
            st_loc = MPL_MIN(st_loc, others_req[i].offsets[j]);
            end_loc = MPL_MAX(end_loc, (others_req[i].offsets[j]
                                        + others_req[i].lens[j] - 1));
        }

    vars->st_loc = st_loc;
    vars->end_loc = end_loc;

    /* ntimes=ceiling_div(end_loc - st_loc + 1, coll_bufsize) */

    vars->ntimes = (int) ((end_loc - st_loc + coll_bufsize) / coll_bufsize);

    if ((st_loc == -1) && (end_loc == -1)) {
        vars->ntimes = 0;       /* this process does no writing. */
    }

    *error_code = MPI_Iallreduce(&vars->ntimes, &vars->max_ntimes, 1, MPI_INT,
                                 MPI_MAX, fd->comm, &vars->req1);

    vars->write_buf = fd->io_buf;

    vars->curr_offlen_ptr = (int *) ADIOI_Calloc(nprocs, sizeof(int));
    /* its use is explained below. calloc initializes to 0. */

    vars->count = (int *) ADIOI_Malloc(nprocs * sizeof(int));
    /* to store count of how many off-len pairs per proc are satisfied
     * in an iteration. */

    vars->partial_recv = (int *) ADIOI_Calloc(nprocs, sizeof(int));
    /* if only a portion of the last off-len pair is recd. from a process
     * in a particular iteration, the length recd. is stored here.
     * calloc initializes to 0. */

    vars->send_size = (int *) ADIOI_Malloc(nprocs * sizeof(int));
    /* total size of data to be sent to each proc. in an iteration.
     * Of size nprocs so that I can use MPI_Alltoall later. */

    vars->recv_size = (int *) ADIOI_Malloc(nprocs * sizeof(int));
    /* total size of data to be recd. from each proc. in an iteration. */

    vars->sent_to_proc = (int *) ADIOI_Calloc(nprocs, sizeof(int));
    /* amount of data sent to each proc so far. Used in
     * ADIOI_Fill_send_buffer. initialized to 0 here. */

    vars->send_buf_idx = (int *) ADIOI_Malloc(nprocs * 3 * sizeof(int));
    vars->curr_to_proc = vars->send_buf_idx + nprocs;
    vars->done_to_proc = vars->curr_to_proc + nprocs;
    /* Above three are used in ADIOI_Fill_send_buffer */

    vars->start_pos = (int *) ADIOI_Malloc(nprocs * sizeof(int));
    /* used to store the starting value of curr_offlen_ptr[i] in
     * this iteration */

    ADIOI_Datatype_iscontig(datatype, &vars->buftype_is_contig);
    if (!vars->buftype_is_contig) {
        vars->flat_buf = ADIOI_Flatten_and_find(datatype);
    }
    MPI_Type_get_extent(datatype, &lb, &vars->buftype_extent);


    /* I need to check if there are any outstanding nonblocking writes to
     * the file, which could potentially interfere with the writes taking
     * place in this collective write call. Since this is not likely to be
     * common, let me do the simplest thing possible here: Each process
     * completes all pending nonblocking operations before completing. */

    /*ADIOI_Complete_async(error_code);
     * if (*error_code != MPI_SUCCESS) return;
     * MPI_Barrier(fd->comm);
     */

    vars->done = 0;
    vars->off = st_loc;

    /* set the state to wait until MPI_Ialltoall finishes. */
    nbc_req->data.wr.state = ADIOI_IWC_STATE_IEXCH_AND_WRITE;
}

static void ADIOI_Iexch_and_write_l1_begin(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_Iexch_and_write_vars *vars = nbc_req->data.wr.eaw_vars;
    int nprocs;
    ADIOI_Access *others_req;

    int i, j;
    ADIO_Offset off, req_off;
    char *write_buf;
    int *curr_offlen_ptr, *count, req_len, *recv_size;
    int *partial_recv, *start_pos;
    ADIO_Offset size;
    static char myname[] = "ADIOI_IEXCH_AND_WRITE_L1_BEGIN";

    ADIOI_W_Iexchange_data_vars *wed_vars = NULL;

    /* loop exit condition */
    if (vars->m >= vars->ntimes) {
        ADIOI_Iexch_and_write_reset(nbc_req, error_code);
        return;
    }

    nprocs = vars->nprocs;
    others_req = vars->others_req;

    off = vars->off;
    write_buf = vars->write_buf;
    curr_offlen_ptr = vars->curr_offlen_ptr;
    count = vars->count;
    recv_size = vars->recv_size;
    partial_recv = vars->partial_recv;
    start_pos = vars->start_pos;

    /* go through all others_req and check which will be satisfied
     * by the current write */

    /* Note that MPI guarantees that displacements in filetypes are in
     * monotonically nondecreasing order and that, for writes, the
     * filetypes cannot specify overlapping regions in the file. This
     * simplifies implementation a bit compared to reads. */

    /* off = start offset in the file for the data to be written in
     * this iteration
     * size = size of data written (bytes) corresponding to off
     * req_off = off in file for a particular contiguous request
     * minus what was satisfied in previous iteration
     * req_size = size corresponding to req_off */

    /* first calculate what should be communicated */

    for (i = 0; i < nprocs; i++)
        count[i] = recv_size[i] = 0;

    size = MPL_MIN((unsigned) vars->coll_bufsize, vars->end_loc - vars->st_loc + 1 - vars->done);
    vars->size = size;

    for (i = 0; i < nprocs; i++) {
        if (others_req[i].count) {
            start_pos[i] = curr_offlen_ptr[i];
            for (j = curr_offlen_ptr[i]; j < others_req[i].count; j++) {
                if (partial_recv[i]) {
                    /* this request may have been partially
                     * satisfied in the previous iteration. */
                    req_off = others_req[i].offsets[j] + partial_recv[i];
                    req_len = others_req[i].lens[j] - partial_recv[i];
                    partial_recv[i] = 0;
                    /* modify the off-len pair to reflect this change */
                    others_req[i].offsets[j] = req_off;
                    others_req[i].lens[j] = req_len;
                } else {
                    req_off = others_req[i].offsets[j];
                    req_len = others_req[i].lens[j];
                }
                if (req_off < off + size) {
                    count[i]++;
                    ADIOI_Assert((((ADIO_Offset) (uintptr_t) write_buf) + req_off - off) ==
                                 (ADIO_Offset) (uintptr_t) (write_buf + req_off - off));
                    MPI_Get_address(write_buf + req_off - off, &(others_req[i].mem_ptrs[j]));
                    ADIOI_Assert((off + size - req_off) == (int) (off + size - req_off));
                    recv_size[i] += (int) (MPL_MIN(off + size - req_off, (unsigned) req_len));

                    if (off + size - req_off < (unsigned) req_len) {
                        partial_recv[i] = (int) (off + size - req_off);

                        /* --BEGIN ERROR HANDLING-- */
                        if ((j + 1 < others_req[i].count) &&
                            (others_req[i].offsets[j + 1] < off + size)) {
                            *error_code = MPIO_Err_create_code(MPI_SUCCESS,
                                                               MPIR_ERR_RECOVERABLE,
                                                               myname,
                                                               __LINE__,
                                                               MPI_ERR_ARG,
                                                               "Filetype specifies overlapping write regions (which is illegal according to the MPI-2 specification)",
                                                               0);
                            /* allow to continue since additional
                             * communication might have to occur
                             */
                        }
                        /* --END ERROR HANDLING-- */
                        break;
                    }
                } else
                    break;
            }
            curr_offlen_ptr[i] = j;
        }
    }

    /* create a struct for ADIOI_W_Iexchange_data() */
    wed_vars = (ADIOI_W_Iexchange_data_vars *) ADIOI_Calloc(1, sizeof(ADIOI_W_Iexchange_data_vars));
    nbc_req->data.wr.wed_vars = wed_vars;

    wed_vars->fd = vars->fd;
    wed_vars->buf = vars->buf;
    wed_vars->write_buf = vars->write_buf;
    wed_vars->flat_buf = vars->flat_buf;
    wed_vars->offset_list = vars->offset_list;
    wed_vars->len_list = vars->len_list;
    wed_vars->send_size = vars->send_size;
    wed_vars->recv_size = vars->recv_size;
    wed_vars->off = vars->off;
    wed_vars->size = vars->size;
    wed_vars->count = vars->count;
    wed_vars->start_pos = vars->start_pos;
    wed_vars->partial_recv = vars->partial_recv;
    wed_vars->sent_to_proc = vars->sent_to_proc;
    wed_vars->nprocs = vars->nprocs;
    wed_vars->myrank = vars->myrank;
    wed_vars->buftype_is_contig = vars->buftype_is_contig;
    wed_vars->contig_access_count = vars->contig_access_count;
    wed_vars->min_st_offset = vars->min_st_offset;
    wed_vars->fd_size = vars->fd_size;
    wed_vars->fd_start = vars->fd_start;
    wed_vars->fd_end = vars->fd_end;
    wed_vars->others_req = vars->others_req;
    wed_vars->send_buf_idx = vars->send_buf_idx;
    wed_vars->curr_to_proc = vars->curr_to_proc;
    wed_vars->done_to_proc = vars->done_to_proc;
    wed_vars->hole = &vars->hole;
    wed_vars->iter = vars->m;
    wed_vars->buftype_extent = vars->buftype_extent;
    wed_vars->buf_idx = vars->buf_idx;
    wed_vars->next_fn = ADIOI_Iexch_and_write_l1_body;

    ADIOI_W_Iexchange_data(nbc_req, error_code);
}

static void ADIOI_Iexch_and_write_l1_body(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_Iexch_and_write_vars *vars = nbc_req->data.wr.eaw_vars;
    ADIO_File fd = vars->fd;
    int nprocs = vars->nprocs;
    ADIO_Offset size = vars->size;
    char *write_buf = vars->write_buf;
    int *count = vars->count;
    int flag, i;

    flag = 0;
    for (i = 0; i < nprocs; i++)
        if (count[i])
            flag = 1;

    if (flag) {
        ADIOI_Assert(size == (int) size);
#if defined(ROMIO_RUN_ON_LINUX) && !defined(HAVE_AIO_LITE_H)
        MPI_Status status;
        ADIO_WriteContig(fd, write_buf, (int) size, MPI_BYTE,
                         ADIO_EXPLICIT_OFFSET, vars->off, &status, error_code);
#else
        ADIO_IwriteContig(fd, write_buf, (int) size, MPI_BYTE,
                          ADIO_EXPLICIT_OFFSET, vars->off, &vars->req3, error_code);

        nbc_req->data.wr.state = ADIOI_IWC_STATE_IEXCH_AND_WRITE_L1_BODY;
        return;
#endif
    }

    ADIOI_Iexch_and_write_l1_end(nbc_req, error_code);
}

static void ADIOI_Iexch_and_write_l1_end(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_Iexch_and_write_vars *vars = nbc_req->data.wr.eaw_vars;
    ADIO_Offset size = vars->size;

    vars->off += size;
    vars->done += size;

    /* increment m and go back to the beginning of m loop */
    vars->m++;
    ADIOI_Iexch_and_write_l1_begin(nbc_req, error_code);
}

static void ADIOI_Iexch_and_write_reset(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_Iexch_and_write_vars *vars = nbc_req->data.wr.eaw_vars;
    int nprocs = vars->nprocs;
    int *count = vars->count;
    int *recv_size = vars->recv_size;
    int i;

    for (i = 0; i < nprocs; i++)
        count[i] = recv_size[i] = 0;

    vars->m = vars->ntimes;
    ADIOI_Iexch_and_write_l2_begin(nbc_req, error_code);
}

static void ADIOI_Iexch_and_write_l2_begin(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_Iexch_and_write_vars *vars = nbc_req->data.wr.eaw_vars;
    ADIOI_W_Iexchange_data_vars *wed_vars = NULL;

    /* loop exit condition */
    if (vars->m >= vars->max_ntimes) {
        ADIOI_Iexch_and_write_fini(nbc_req, error_code);
        return;
    }

    ADIOI_Assert(vars->size == (int) vars->size);

    /* create a struct for ADIOI_W_Iexchange_data() */
    wed_vars = (ADIOI_W_Iexchange_data_vars *) ADIOI_Calloc(1, sizeof(ADIOI_W_Iexchange_data_vars));
    nbc_req->data.wr.wed_vars = wed_vars;

    wed_vars->fd = vars->fd;
    wed_vars->buf = vars->buf;
    wed_vars->write_buf = vars->write_buf;
    wed_vars->flat_buf = vars->flat_buf;
    wed_vars->offset_list = vars->offset_list;
    wed_vars->len_list = vars->len_list;
    wed_vars->send_size = vars->send_size;
    wed_vars->recv_size = vars->recv_size;
    wed_vars->off = vars->off;
    wed_vars->size = (int) vars->size;
    wed_vars->count = vars->count;
    wed_vars->start_pos = vars->start_pos;
    wed_vars->partial_recv = vars->partial_recv;
    wed_vars->sent_to_proc = vars->sent_to_proc;
    wed_vars->nprocs = vars->nprocs;
    wed_vars->myrank = vars->myrank;
    wed_vars->buftype_is_contig = vars->buftype_is_contig;
    wed_vars->contig_access_count = vars->contig_access_count;
    wed_vars->min_st_offset = vars->min_st_offset;
    wed_vars->fd_size = vars->fd_size;
    wed_vars->fd_start = vars->fd_start;
    wed_vars->fd_end = vars->fd_end;
    wed_vars->others_req = vars->others_req;
    wed_vars->send_buf_idx = vars->send_buf_idx;
    wed_vars->curr_to_proc = vars->curr_to_proc;
    wed_vars->done_to_proc = vars->done_to_proc;
    wed_vars->hole = &vars->hole;
    wed_vars->iter = vars->m;
    wed_vars->buftype_extent = vars->buftype_extent;
    wed_vars->buf_idx = vars->buf_idx;
    wed_vars->next_fn = ADIOI_Iexch_and_write_l2_end;

    /* nothing to recv, but check for send. */
    ADIOI_W_Iexchange_data(nbc_req, error_code);
}

static void ADIOI_Iexch_and_write_l2_end(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_Iexch_and_write_vars *vars = nbc_req->data.wr.eaw_vars;

    vars->m++;
    ADIOI_Iexch_and_write_l2_begin(nbc_req, error_code);
}

static void ADIOI_Iexch_and_write_fini(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_Iexch_and_write_vars *vars = nbc_req->data.wr.eaw_vars;
    void (*next_fn) (ADIOI_NBC_Request *, int *);

    ADIOI_Free(vars->curr_offlen_ptr);
    ADIOI_Free(vars->count);
    ADIOI_Free(vars->partial_recv);
    ADIOI_Free(vars->send_size);
    ADIOI_Free(vars->recv_size);
    ADIOI_Free(vars->sent_to_proc);
    ADIOI_Free(vars->start_pos);
    ADIOI_Free(vars->send_buf_idx);

    next_fn = vars->next_fn;

    /* free the struct for parameters and variables */
    ADIOI_Free(nbc_req->data.wr.eaw_vars);
    nbc_req->data.wr.eaw_vars = NULL;

    /* move to the next function */
    next_fn(nbc_req, error_code);
}


static void ADIOI_W_Iexchange_data(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_W_Iexchange_data_vars *vars = nbc_req->data.wr.wed_vars;

    /* exchange recv_size info so that each process knows how much to
     * send to whom. */

    *error_code = MPI_Ialltoall(vars->recv_size, 1, MPI_INT, vars->send_size, 1,
                                MPI_INT, vars->fd->comm, &vars->req1);

    nbc_req->data.wr.state = ADIOI_IWC_STATE_W_IEXCHANGE_DATA;
}

static void ADIOI_W_Iexchange_data_hole(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_W_Iexchange_data_vars *vars = nbc_req->data.wr.wed_vars;
    ADIO_File fd = vars->fd;
    int *recv_size = vars->recv_size;
    ADIO_Offset off = vars->off;
    int size = vars->size;
    int *count = vars->count;
    int *start_pos = vars->start_pos;
    int *partial_recv = vars->partial_recv;
    int nprocs = vars->nprocs;
    ADIOI_Access *others_req = vars->others_req;
    int *hole = vars->hole;

    int i, j, k, *tmp_len, nprocs_recv;
    MPI_Datatype *recv_types;
    int *srt_len = NULL, sum;
    ADIO_Offset *srt_off = NULL;

    /* create derived datatypes for recv */

    nprocs_recv = 0;
    for (i = 0; i < nprocs; i++)
        if (recv_size[i])
            nprocs_recv++;
    vars->nprocs_recv = nprocs_recv;

    recv_types = (MPI_Datatype *)
        ADIOI_Malloc((nprocs_recv + 1) * sizeof(MPI_Datatype));
    vars->recv_types = recv_types;
    /* +1 to avoid a 0-size malloc */

    tmp_len = (int *) ADIOI_Malloc(nprocs * sizeof(int));
    j = 0;
    for (i = 0; i < nprocs; i++) {
        if (recv_size[i]) {
            /* take care if the last off-len pair is a partial recv */
            if (partial_recv[i]) {
                k = start_pos[i] + count[i] - 1;
                tmp_len[i] = others_req[i].lens[k];
                others_req[i].lens[k] = partial_recv[i];
            }
            ADIOI_Type_create_hindexed_x(count[i],
                                         &(others_req[i].lens[start_pos[i]]),
                                         &(others_req[i].mem_ptrs[start_pos[i]]),
                                         MPI_BYTE, recv_types + j);
            /* absolute displacements; use MPI_BOTTOM in recv */
            MPI_Type_commit(recv_types + j);
            j++;
        }
    }

    /* To avoid a read-modify-write, check if there are holes in the
     * data to be written. For this, merge the (sorted) offset lists
     * others_req using a heap-merge. */

    sum = 0;
    for (i = 0; i < nprocs; i++)
        sum += count[i];
    /* valgrind-detcted optimization: if there is no work on this process we do
     * not need to search for holes */
    if (sum) {
        srt_off = (ADIO_Offset *) ADIOI_Malloc(sum * sizeof(ADIO_Offset));
        srt_len = (int *) ADIOI_Malloc(sum * sizeof(int));

        ADIOI_Heap_merge(others_req, count, srt_off, srt_len, start_pos, nprocs, nprocs_recv, sum);
    }

    /* for partial recvs, restore original lengths */
    for (i = 0; i < nprocs; i++)
        if (partial_recv[i]) {
            k = start_pos[i] + count[i] - 1;
            others_req[i].lens[k] = tmp_len[i];
        }
    ADIOI_Free(tmp_len);

    /* check if there are any holes. If yes, must do read-modify-write.
     * holes can be in three places.  'middle' is what you'd expect: the
     * processes are operating on noncontigous data.  But holes can also show
     * up at the beginning or end of the file domain (see John Bent ROMIO REQ
     * #835). Missing these holes would result in us writing more data than
     * recieved by everyone else. */

    *hole = 0;
    if (sum) {
        if (off != srt_off[0])  /* hole at the front */
            *hole = 1;
        else {  /* coalesce the sorted offset-length pairs */
            for (i = 1; i < sum; i++) {
                if (srt_off[i] <= srt_off[0] + srt_len[0]) {
                    /* ok to cast: operating on cb_buffer_size chunks */
                    int new_len = (int) srt_off[i] + srt_len[i] - (int) srt_off[0];
                    if (new_len > srt_len[0])
                        srt_len[0] = new_len;
                } else
                    break;
            }
            if (i < sum || size != srt_len[0])  /* hole in middle or end */
                *hole = 1;
        }

        ADIOI_Free(srt_off);
        ADIOI_Free(srt_len);
    }

    if (nprocs_recv) {
        if (*hole) {
            ADIO_IreadContig(fd, vars->write_buf, size, MPI_BYTE,
                             ADIO_EXPLICIT_OFFSET, off, &vars->req2, &vars->err);
            nbc_req->data.wr.state = ADIOI_IWC_STATE_W_IEXCHANGE_DATA_HOLE;
            return;
        }
    }

    ADIOI_W_Iexchange_data_send(nbc_req, error_code);
}

static void ADIOI_W_Iexchange_data_send(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_W_Iexchange_data_vars *vars = nbc_req->data.wr.wed_vars;
    ADIO_File fd = vars->fd;
    void *buf = vars->buf;
    int *send_size = vars->send_size;
    int *recv_size = vars->recv_size;
    int nprocs = vars->nprocs;
    int myrank = vars->myrank;
    int iter = vars->iter;
    MPI_Aint *buf_idx = vars->buf_idx;

    int nprocs_recv = vars->nprocs_recv;
    MPI_Datatype *recv_types = vars->recv_types;

    int i, j;
    int nprocs_send;
    char **send_buf = NULL;

    nprocs_send = 0;
    for (i = 0; i < nprocs; i++)
        if (send_size[i])
            nprocs_send++;
    vars->nprocs_send = nprocs_send;

    if (fd->atomicity) {
        /* bug fix from Wei-keng Liao and Kenin Coloma */
        vars->requests = (MPI_Request *)
            ADIOI_Malloc((nprocs_send + 1) * sizeof(MPI_Request));
        vars->send_req = vars->requests;
    } else {
        vars->requests = (MPI_Request *)
            ADIOI_Malloc((nprocs_send + nprocs_recv + 1) * sizeof(MPI_Request));
        /* +1 to avoid a 0-size malloc */

        /* post receives */
        j = 0;
        for (i = 0; i < nprocs; i++) {
            if (recv_size[i]) {
                MPI_Irecv(MPI_BOTTOM, 1, recv_types[j], i, myrank + i + 100 * iter,
                          fd->comm, vars->requests + j);
                j++;
            }
        }
        vars->send_req = vars->requests + nprocs_recv;
    }

    /* post sends. if buftype_is_contig, data can be directly sent from
     * user buf at location given by buf_idx. else use send_buf. */

#ifdef AGGREGATION_PROFILE
    MPE_Log_event(5032, 0, NULL);
#endif
    if (vars->buftype_is_contig) {
        j = 0;
        for (i = 0; i < nprocs; i++)
            if (send_size[i]) {
                MPI_Isend(((char *) buf) + buf_idx[i], send_size[i],
                          MPI_BYTE, i, myrank + i + 100 * iter, fd->comm, vars->send_req + j);
                j++;
                buf_idx[i] += send_size[i];
            }
    } else if (nprocs_send) {
        /* buftype is not contig */
        size_t msgLen = 0;
        for (i = 0; i < nprocs; i++)
            msgLen += send_size[i];
        send_buf = (char **) ADIOI_Malloc(nprocs * sizeof(char *));
        send_buf[0] = (char *) ADIOI_Malloc(msgLen * sizeof(char));
        for (i = 1; i < nprocs; i++)
            send_buf[i] = send_buf[i - 1] + send_size[i - 1];
        vars->send_buf = send_buf;

        ADIOI_Fill_send_buffer(fd, buf, vars->flat_buf, send_buf,
                               vars->offset_list, vars->len_list, send_size,
                               vars->send_req,
                               vars->sent_to_proc, nprocs, myrank,
                               vars->contig_access_count,
                               vars->min_st_offset, vars->fd_size,
                               vars->fd_start, vars->fd_end,
                               vars->send_buf_idx, vars->curr_to_proc,
                               vars->done_to_proc, iter, vars->buftype_extent);
        /* the send is done in ADIOI_Fill_send_buffer */
    }

    if (fd->atomicity) {
        vars->req3 = (MPI_Request *)
            ADIOI_Malloc((nprocs_recv + 1) * sizeof(MPI_Request));
        /* +1 to avoid a 0-size malloc */

        /* bug fix from Wei-keng Liao and Kenin Coloma */
        j = 0;
        for (i = 0; i < nprocs; i++) {
            if (recv_size[i]) {
                MPI_Irecv(MPI_BOTTOM, 1, recv_types[j], i, myrank + i + 100 * iter,
                          fd->comm, vars->req3 + j);
                j++;
            }
        }

        nbc_req->data.wr.state = ADIOI_IWC_STATE_W_IEXCHANGE_DATA_SEND;
        return;
    }

    ADIOI_W_Iexchange_data_wait(nbc_req, error_code);
}

static void ADIOI_W_Iexchange_data_wait(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_W_Iexchange_data_vars *vars = nbc_req->data.wr.wed_vars;
    ADIO_File fd = vars->fd;
    int nprocs_send = vars->nprocs_send;
    int nprocs_recv = vars->nprocs_recv;
    MPI_Datatype *recv_types = vars->recv_types;
    int i;

    for (i = 0; i < nprocs_recv; i++)
        MPI_Type_free(recv_types + i);
    ADIOI_Free(recv_types);

    i = 0;
    if (fd->atomicity) {
        /* bug fix from Wei-keng Liao and Kenin Coloma */
        MPI_Testall(nprocs_send, vars->send_req, &i, MPI_STATUSES_IGNORE);
    } else {
        MPI_Testall(nprocs_send + nprocs_recv, vars->requests, &i, MPI_STATUSES_IGNORE);
    }

    if (i) {
        ADIOI_W_Iexchange_data_fini(nbc_req, error_code);
    } else {
        nbc_req->data.wr.state = ADIOI_IWC_STATE_W_IEXCHANGE_DATA_WAIT;
    }
}

static void ADIOI_W_Iexchange_data_fini(ADIOI_NBC_Request * nbc_req, int *error_code)
{
    ADIOI_W_Iexchange_data_vars *vars = nbc_req->data.wr.wed_vars;
    void (*next_fn) (ADIOI_NBC_Request *, int *);
    ADIO_File fd = vars->fd;
    char **send_buf = vars->send_buf;

    if (fd->atomicity)
        ADIOI_Free(vars->req3);

#ifdef AGGREGATION_PROFILE
    MPE_Log_event(5033, 0, NULL);
#endif
    ADIOI_Free(vars->requests);
    if (!vars->buftype_is_contig && vars->nprocs_send) {
        ADIOI_Free(send_buf[0]);
        ADIOI_Free(send_buf);
    }

    next_fn = vars->next_fn;

    /* free the structure for parameters and variables */
    ADIOI_Free(vars);
    nbc_req->data.wr.wed_vars = NULL;

    /* move to the next function */
    next_fn(nbc_req, error_code);
}


static int ADIOI_GEN_iwc_query_fn(void *extra_state, MPI_Status * status)
{
    ADIOI_NBC_Request *nbc_req;

    nbc_req = (ADIOI_NBC_Request *) extra_state;

    MPI_Status_set_elements_x(status, MPI_BYTE, nbc_req->nbytes);

    /* can never cancel so always true */
    MPI_Status_set_cancelled(status, 0);

    /* choose not to return a value for this */
    status->MPI_SOURCE = MPI_UNDEFINED;
    /* tag has no meaning for this generalized request */
    status->MPI_TAG = MPI_UNDEFINED;

    /* this generalized request never fails */
    return MPI_SUCCESS;
}

static int ADIOI_GEN_iwc_free_fn(void *extra_state)
{
    ADIOI_NBC_Request *nbc_req;

    nbc_req = (ADIOI_NBC_Request *) extra_state;
    ADIOI_Free(nbc_req);

    return MPI_SUCCESS;
}

static int ADIOI_GEN_iwc_poll_fn(void *extra_state, MPI_Status * status)
{
    ADIOI_NBC_Request *nbc_req;
    ADIOI_GEN_IwriteStridedColl_vars *wsc_vars = NULL;
    ADIOI_Icalc_others_req_vars *cor_vars = NULL;
    ADIOI_Iexch_and_write_vars *eaw_vars = NULL;
    ADIOI_W_Iexchange_data_vars *wed_vars = NULL;
    int errcode = MPI_SUCCESS;
    int flag;

    nbc_req = (ADIOI_NBC_Request *) extra_state;

    switch (nbc_req->data.wr.state) {
        case ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL:
            wsc_vars = nbc_req->data.wr.wsc_vars;
            errcode = MPI_Testall(2, wsc_vars->req_offset, &flag, MPI_STATUSES_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                ADIOI_GEN_IwriteStridedColl_inter(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL_INDIO:
            wsc_vars = nbc_req->data.wr.wsc_vars;
            errcode = MPI_Test(&wsc_vars->req_ind_io, &flag, MPI_STATUS_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                /* call the last function */
                ADIOI_GEN_IwriteStridedColl_fini(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_GEN_IWRITESTRIDEDCOLL_BCAST:
            wsc_vars = nbc_req->data.wr.wsc_vars;
            errcode = MPI_Test(&wsc_vars->req_err, &flag, MPI_STATUS_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                errcode = wsc_vars->error_code;
                ADIOI_GEN_IwriteStridedColl_free(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_ICALC_OTHERS_REQ:
            cor_vars = nbc_req->cor_vars;
            errcode = MPI_Test(&cor_vars->req1, &flag, MPI_STATUS_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                ADIOI_Icalc_others_req_main(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_ICALC_OTHERS_REQ_MAIN:
            cor_vars = nbc_req->cor_vars;
            if (cor_vars->num_req2) {
                errcode = MPI_Testall(cor_vars->num_req2, cor_vars->req2,
                                      &flag, MPI_STATUSES_IGNORE);
                if (errcode == MPI_SUCCESS && flag) {
                    ADIOI_Icalc_others_req_fini(nbc_req, &errcode);
                }
            } else {
                ADIOI_Icalc_others_req_fini(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_IEXCH_AND_WRITE:
            eaw_vars = nbc_req->data.wr.eaw_vars;
            errcode = MPI_Test(&eaw_vars->req1, &flag, MPI_STATUS_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                eaw_vars->m = 0;
                ADIOI_Iexch_and_write_l1_begin(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_IEXCH_AND_WRITE_L1_BODY:
            eaw_vars = nbc_req->data.wr.eaw_vars;
            errcode = MPI_Test(&eaw_vars->req3, &flag, MPI_STATUS_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                ADIOI_Iexch_and_write_l1_end(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_W_IEXCHANGE_DATA:
            wed_vars = nbc_req->data.wr.wed_vars;
            errcode = MPI_Test(&wed_vars->req1, &flag, MPI_STATUS_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                ADIOI_W_Iexchange_data_hole(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_W_IEXCHANGE_DATA_HOLE:
            wed_vars = nbc_req->data.wr.wed_vars;
            errcode = MPI_Test(&wed_vars->req2, &flag, MPI_STATUSES_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                /* --BEGIN ERROR HANDLING-- */
                if (wed_vars->err != MPI_SUCCESS) {
                    errcode = MPIO_Err_create_code(wed_vars->err,
                                                   MPIR_ERR_RECOVERABLE,
                                                   "ADIOI_W_EXCHANGE_DATA",
                                                   __LINE__, MPI_ERR_IO, "**ioRMWrdwr", 0);
                    break;;
                }
                /* --END ERROR HANDLING-- */
                ADIOI_W_Iexchange_data_send(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_W_IEXCHANGE_DATA_SEND:
            wed_vars = nbc_req->data.wr.wed_vars;
            errcode = MPI_Testall(wed_vars->nprocs_recv, wed_vars->req3,
                                  &flag, MPI_STATUSES_IGNORE);
            if (errcode == MPI_SUCCESS && flag) {
                ADIOI_W_Iexchange_data_wait(nbc_req, &errcode);
            }
            break;

        case ADIOI_IWC_STATE_W_IEXCHANGE_DATA_WAIT:
            wed_vars = nbc_req->data.wr.wed_vars;
            if (wed_vars->fd->atomicity) {
                /* bug fix from Wei-keng Liao and Kenin Coloma */
                errcode = MPI_Testall(wed_vars->nprocs_send, wed_vars->send_req,
                                      &flag, MPI_STATUSES_IGNORE);
            } else {
                errcode = MPI_Testall(wed_vars->nprocs_send +
                                      wed_vars->nprocs_recv,
                                      wed_vars->requests, &flag, MPI_STATUSES_IGNORE);
            }
            if (errcode == MPI_SUCCESS && flag) {
                ADIOI_W_Iexchange_data_fini(nbc_req, &errcode);
            }
            break;

        default:
            break;
    }

    /* --BEGIN ERROR HANDLING-- */
    if (errcode != MPI_SUCCESS) {
        errcode = MPIO_Err_create_code(MPI_SUCCESS,
                                       MPIR_ERR_RECOVERABLE,
                                       "ADIOI_GEN_iwc_poll_fn", __LINE__,
                                       MPI_ERR_IO, "**mpi_grequest_complete", 0);
    }
    /* --END ERROR HANDLING-- */

    return errcode;
}

/* wait for multiple requests to complete */
static int ADIOI_GEN_iwc_wait_fn(int count, void **array_of_states,
                                 double timeout, MPI_Status * status)
{
    int i, errcode = MPI_SUCCESS;
    double starttime;
    ADIOI_NBC_Request **nbc_reqlist;

    nbc_reqlist = (ADIOI_NBC_Request **) array_of_states;

    starttime = MPI_Wtime();
    for (i = 0; i < count; i++) {
        while (nbc_reqlist[i]->data.wr.state != ADIOI_IWC_STATE_COMPLETE) {
            errcode = ADIOI_GEN_iwc_poll_fn(nbc_reqlist[i], MPI_STATUS_IGNORE);
            /* --BEGIN ERROR HANDLING-- */
            if (errcode != MPI_SUCCESS) {
                errcode = MPIO_Err_create_code(MPI_SUCCESS,
                                               MPIR_ERR_RECOVERABLE,
                                               "ADIOI_GEN_iwc_wait_fn",
                                               __LINE__, MPI_ERR_IO, "**mpi_grequest_complete", 0);
            }
            /* --END ERROR HANDLING-- */

            if ((timeout > 0) && (timeout < (MPI_Wtime() - starttime)))
                goto fn_exit;

            /* If the progress engine is blocked, we have to yield for another
             * thread to be able to unblock the progress engine. */
            /* NOTE: we're outside a critical section (safety ensured by standard),
             * we only need yield in case of user threads */
            ROMIO_THREAD_CS_YIELD();
        }
    }

  fn_exit:
    return errcode;
}

#endif /* HAVE_MPI_GREQUEST_EXTENSIONS */