codes-online-comm-wrkld.C 29.2 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13
/*
 * Copyright (C) 2014 University of Chicago
 * See COPYRIGHT notice in top-level directory.
 *
 */

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <mpi.h>
#include <ross.h>
#include <assert.h>
#include <deque>
14
#include <iostream>
15
#include <inttypes.h>
16
#include <fstream>
17 18
#include <boost/property_tree/json_parser.hpp>
#include <boost/property_tree/ptree.hpp>
19 20 21 22 23 24
#include "codes/codes-workload.h"
#include "codes/quickhash.h"
#include "codes/codes-jobmap.h"
#include "codes_config.h"
#include "lammps.h"
#include "nekbone_swm_user_code.h"
25
#include "nearest_neighbor_swm_user_code.h"
26
#include "all_to_one_swm_user_code.h"
27

28 29
#define ALLREDUCE_SHORT_MSG_SIZE 2048

30

31
//#define DBG_COMM 0
32 33 34 35 36 37 38

using namespace std;

static struct qhash_table *rank_tbl = NULL;
static int rank_tbl_pop = 0;
static int total_rank_cnt = 0;
ABT_thread global_prod_thread = NULL;
39
ABT_xstream self_es;
40
double cpu_freq = 1.0;
41 42 43 44 45 46 47 48 49
long num_allreduce = 0;
long num_isends = 0;
long num_irecvs = 0;
long num_barriers = 0;
long num_sends = 0;
long num_recvs = 0;
long num_sendrecv = 0;
long num_waitalls = 0;

50 51 52
//std::map<int64_t, int> send_count;
//std::map<int64_t, int> isend_count;
//std::map<int64_t, int> allreduce_count;
53 54

struct shared_context {
55
    int my_rank;
56
    uint32_t wait_id;
57 58 59 60 61
    int num_ranks;
    char workload_name[MAX_NAME_LENGTH_WKLD];
    void * swm_obj;
    ABT_thread      producer;
    std::deque<struct codes_workload_op*> fifo;
62 63 64
};

struct rank_mpi_context {
65 66 67
    struct qhash_head hash_link;
    int app_id;
    struct shared_context sctx;
68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85
};

typedef struct rank_mpi_compare {
    int app_id;
    int rank;
} rank_mpi_compare;

/*
 * peer: the receiving peer id 
 * comm_id: the communicator id being used
 * tag: tag id 
 * reqvc: virtual channel being used by the message (to be ignored)
 * rspvc: virtual channel being used by the message (to be ignored)
 * buf: the address of sender's buffer in memory
 * bytes: number of bytes to be sent 
 * reqrt and rsprt: routing types (to be ignored) */

void SWM_Send(SWM_PEER peer,
86 87 88 89 90 91 92 93 94
        SWM_COMM_ID comm_id,
        SWM_TAG tag,
        SWM_VC reqvc,
        SWM_VC rspvc,
        SWM_BUF buf,
        SWM_BYTES bytes,
        SWM_BYTES pktrspbytes,
        SWM_ROUTING_TYPE reqrt,
        SWM_ROUTING_TYPE rsprt)
95 96
{
    /* add an event in the shared queue and then yield */
97
    //    printf("\n Sending to rank %d ", comm_id);
98 99 100 101 102 103 104 105
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_SEND;
    wrkld_per_rank.u.send.tag = tag;
    wrkld_per_rank.u.send.num_bytes = bytes;
    wrkld_per_rank.u.send.dest_rank = peer;

#ifdef DBG_COMM
106
/*    if(tag != 1235 && tag != 1234) 
107 108 109 110 111 112 113 114 115 116
    {
        auto it = send_count.find(bytes);
        if(it == send_count.end())
        {
            send_count.insert(std::make_pair(bytes, 1));
        }
        else
        {
            it->second = it->second + 1;
        }
117
    }*/
118 119 120 121 122 123 124 125 126 127 128 129 130
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    wrkld_per_rank.u.send.source_rank = sctx->my_rank;
    sctx->fifo.push_back(&wrkld_per_rank);

    ABT_thread_yield_to(global_prod_thread);
131
    num_sends++;
132 133
}

134 135 136 137 138 139 140
/*
 * @param comm_id: communicator ID (For now, MPI_COMM_WORLD)
 * reqvc and rspvc: virtual channel IDs for request and response (ignore for
 * our purpose)
 * buf: buffer location for the call (ignore for our purpose)
 * reqrt and rsprt: routing types, ignore and use routing from config file instead. 
 * */
141 142 143 144 145 146 147 148 149 150 151
void SWM_Barrier(
        SWM_COMM_ID comm_id,
        SWM_VC reqvc,
        SWM_VC rspvc,
        SWM_BUF buf, 
        SWM_UNKNOWN auto1,
        SWM_UNKNOWN2 auto2,
        SWM_ROUTING_TYPE reqrt, 
        SWM_ROUTING_TYPE rsprt)
{
    /* Add an event in the shared queue and then yield */
152
#if 0
153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_DELAY;
    /* TODO: Check how to convert cycle count into delay? */
    wrkld_per_rank.u.delay.nsecs = 0.1;

#ifdef DBG_COMM
    printf("\n Barrier delay %lf ", wrkld_per_rank.u.delay.nsecs);
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    sctx->fifo.push_back(&wrkld_per_rank);

    ABT_thread_yield_to(global_prod_thread);
173 174
#endif
#ifdef DBG_COMM
175
//     printf("\n barrier ");
176 177 178 179 180
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err;
181
    int rank, size, src, dest, mask;
182 183 184 185 186 187 188 189 190 191 192 193

    err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);

    rank = sctx->my_rank;
    size = sctx->num_ranks;
    mask = 0x1;

    while(mask < size) {
194
        dest = (rank + mask) % size;
195 196
        src = (rank - mask + size) % size;

197 198
        SWM_Sendrecv(comm_id, dest, 1234, reqvc, rspvc, 0, 0, 0,
                src,  1234, 0,  reqrt, rsprt);
199 200
        mask <<= 1;
    }
201
    num_barriers++;
202 203 204
}

void SWM_Isend(SWM_PEER peer,
205 206 207 208 209 210 211 212 213 214
        SWM_COMM_ID comm_id,
        SWM_TAG tag,
        SWM_VC reqvc,
        SWM_VC rspvc,
        SWM_BUF buf,
        SWM_BYTES bytes,
        SWM_BYTES pktrspbytes,
        uint32_t * handle,
        SWM_ROUTING_TYPE reqrt,
        SWM_ROUTING_TYPE rsprt)
215 216
{
    /* add an event in the shared queue and then yield */
217
    //    printf("\n Sending to rank %d ", comm_id);
218 219 220 221 222 223 224 225
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_ISEND;
    wrkld_per_rank.u.send.tag = tag;
    wrkld_per_rank.u.send.num_bytes = bytes;
    wrkld_per_rank.u.send.dest_rank = peer;

#ifdef DBG_COMM
226
/*    if(tag != 1235 && tag != 1234) 
227 228 229 230 231 232 233 234 235 236
    {
        auto it = isend_count.find(bytes);
        if(it == isend_count.end())
        {
            isend_count.insert(std::make_pair(bytes, 1));
        }
        else
        {
            it->second = it->second + 1;
        }
237
    }*/
238 239 240 241 242 243 244 245 246 247 248 249
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    wrkld_per_rank.u.send.source_rank = sctx->my_rank;
    sctx->fifo.push_back(&wrkld_per_rank);

250 251 252 253
    *handle = sctx->wait_id;
    wrkld_per_rank.u.send.req_id = *handle;
    sctx->wait_id++;

254
    ABT_thread_yield_to(global_prod_thread);
255
    num_isends++;
256 257 258 259 260 261 262 263 264 265 266 267
}
void SWM_Recv(SWM_PEER peer,
        SWM_COMM_ID comm_id,
        SWM_TAG tag,
        SWM_BUF buf)
{
    /* Add an event in the shared queue and then yield */
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_RECV;
    wrkld_per_rank.u.recv.tag = tag;
    wrkld_per_rank.u.recv.source_rank = peer;
268
    wrkld_per_rank.u.recv.num_bytes = 0;
269 270

#ifdef DBG_COMM
271
    //printf("\n recv op tag: %d source: %d ", tag, peer);
272 273 274 275 276 277 278 279 280 281 282 283 284
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    wrkld_per_rank.u.recv.dest_rank = sctx->my_rank;
    sctx->fifo.push_back(&wrkld_per_rank);

    ABT_thread_yield_to(global_prod_thread);
285
    num_recvs++;
286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303
}

/* handle is for the request ID */
void SWM_Irecv(SWM_PEER peer,
        SWM_COMM_ID comm_id,
        SWM_TAG tag,
        SWM_BUF buf, 
        uint32_t* handle)
{
    /* Add an event in the shared queue and then yield */
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_IRECV;
    wrkld_per_rank.u.recv.tag = tag;
    wrkld_per_rank.u.recv.source_rank = peer;
    wrkld_per_rank.u.recv.num_bytes = 0;

#ifdef DBG_COMM
304
//    printf("\n irecv op tag: %d source: %d ", tag, peer);
305 306 307 308 309 310 311 312 313 314 315 316
#endif

    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    wrkld_per_rank.u.recv.dest_rank = sctx->my_rank;
    sctx->fifo.push_back(&wrkld_per_rank);
317 318 319 320
    
    *handle = sctx->wait_id;
    wrkld_per_rank.u.recv.req_id = *handle;
    sctx->wait_id++;
321 322

    ABT_thread_yield_to(global_prod_thread);
323
    num_irecvs++;
324 325 326 327
}

void SWM_Compute(long cycle_count)
{
328
    //NM: noting that cpu_frequency has been loaded in comm_online_workload_load() as GHz, e.g. cpu_freq = 2.0 means 2.0GHz
329
    if(!cpu_freq)
330
        cpu_freq = 2.0;
331 332 333
    /* Add an event in the shared queue and then yield */
    struct codes_workload_op wrkld_per_rank;

334 335 336 337
    double cpu_freq_hz = cpu_freq * 1000.0 * 1000.0 * 1000.0;
    double delay_in_seconds = cycle_count / cpu_freq_hz;
    double delay_in_ns = delay_in_seconds * 1000.0 * 1000.0 * 1000.0;

338 339
    wrkld_per_rank.op_type = CODES_WK_DELAY;
    /* TODO: Check how to convert cycle count into delay? */
340 341
    wrkld_per_rank.u.delay.nsecs = delay_in_ns;
    wrkld_per_rank.u.delay.seconds = delay_in_seconds;
342
#ifdef DBG_COMM
343
    printf("\n compute op delay: %ld ", delay_in_ns);
344 345 346 347 348 349 350 351 352 353
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    sctx->fifo.push_back(&wrkld_per_rank);
354
	
355
    ABT_thread_yield_to(global_prod_thread);
356

357 358 359 360 361 362 363 364 365 366 367 368
}

void SWM_Wait(uint32_t req_id)
{
    /* Add an event in the shared queue and then yield */
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_WAIT;
    /* TODO: Check how to convert cycle count into delay? */
    wrkld_per_rank.u.wait.req_id = req_id;

#ifdef DBG_COMM
369
//    printf("\n wait op req_id: %"PRIu32"\n", req_id);
370
//      printf("\n wait ");
371 372 373 374 375 376 377 378 379 380
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    sctx->fifo.push_back(&wrkld_per_rank);
381

382 383 384 385 386
    ABT_thread_yield_to(global_prod_thread);
}

void SWM_Waitall(int len, uint32_t * req_ids)
{
387
    num_waitalls++;
388 389 390 391 392 393 394 395 396 397 398 399
    /* Add an event in the shared queue and then yield */
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_WAITALL;
    /* TODO: Check how to convert cycle count into delay? */
    wrkld_per_rank.u.waits.count = len;
    wrkld_per_rank.u.waits.req_ids = (unsigned int*)calloc(len, sizeof(int));    

    for(int i = 0; i < len; i++)
        wrkld_per_rank.u.waits.req_ids[i] = req_ids[i];

#ifdef DBG_COMM
400
//    for(int i = 0; i < len; i++)
401
//        printf("\n wait op len %d req_id: %"PRIu32"\n", len, req_ids[i]);
402 403 404 405 406 407 408 409 410 411
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    sctx->fifo.push_back(&wrkld_per_rank);
412

413 414 415 416
    ABT_thread_yield_to(global_prod_thread);
}

void SWM_Sendrecv(
417 418 419 420 421 422 423 424 425 426 427 428 429
        SWM_COMM_ID comm_id,
        SWM_PEER sendpeer,
        SWM_TAG sendtag,
        SWM_VC sendreqvc,
        SWM_VC sendrspvc,
        SWM_BUF sendbuf,
        SWM_BYTES sendbytes,
        SWM_BYTES pktrspbytes,
        SWM_PEER recvpeer,
        SWM_TAG recvtag,
        SWM_BUF recvbuf,
        SWM_ROUTING_TYPE reqrt,
        SWM_ROUTING_TYPE rsprt)
430
{
431
    //    printf("\n Sending to %d receiving from %d ", sendpeer, recvpeer);
432 433 434 435 436 437 438 439 440 441 442 443 444
    struct codes_workload_op send_op;

    send_op.op_type = CODES_WK_SEND;
    send_op.u.send.tag = sendtag;
    send_op.u.send.num_bytes = sendbytes;
    send_op.u.send.dest_rank = sendpeer;

    /* Add an event in the shared queue and then yield */
    struct codes_workload_op recv_op;

    recv_op.op_type = CODES_WK_RECV;
    recv_op.u.recv.tag = recvtag;
    recv_op.u.recv.source_rank = recvpeer;
445
    recv_op.u.recv.num_bytes = 0;
446

447
#ifdef DBG_COMM
448
/*    if(sendtag != 1235 && sendtag != 1234) 
449 450 451 452 453 454 455 456 457 458
    {
        auto it = send_count.find(sendbytes);
        if(it == send_count.end())
        {
            send_count.insert(std::make_pair(sendbytes, 1));
        }
        else
        {
            it->second = it->second + 1;
        }
459
    }*/
460 461 462 463 464 465 466 467 468 469 470 471 472 473 474
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    recv_op.u.recv.dest_rank = sctx->my_rank;
    send_op.u.send.source_rank = sctx->my_rank;
    sctx->fifo.push_back(&send_op);
    sctx->fifo.push_back(&recv_op);

    ABT_thread_yield_to(global_prod_thread);
475
    num_sendrecv++;
476 477
}

478
/* @param count: number of bytes in Allreduce
479 480 481 482 483 484 485 486 487
 * @param respbytes: number of bytes to be sent in response (ignore for our
 * purpose)
 * $params comm_id: communicator ID (MPI_COMM_WORLD for our case)
 * @param sendreqvc: virtual channel of the sender request (ignore for our
 * purpose)
 * @param sendrspvc: virtual channel of the response request (ignore for our
 * purpose)
 * @param sendbuf and rcvbuf: buffers for send and receive calls (ignore for
 * our purpose) */
488
void SWM_Allreduce(
489
        SWM_BYTES count,
490 491 492 493 494 495 496
        SWM_BYTES respbytes,
        SWM_COMM_ID comm_id,
        SWM_VC sendreqvc,
        SWM_VC sendrspvc,
        SWM_BUF sendbuf,
        SWM_BUF rcvbuf)
{
497
#if 0
498
    /* TODO: For now, simulate a constant delay for ALlreduce*/
499
    //    printf("\n Allreduce bytes %d ", bytes);
500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520
    /* Add an event in the shared queue and then yield */
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_DELAY;
    /* TODO: Check how to convert cycle count into delay? */
    wrkld_per_rank.u.delay.nsecs = bytes + 0.1;

#ifdef DBG_COMM
    printf("\n Allreduce delay %lf ", wrkld_per_rank.u.delay.nsecs);
#endif
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    sctx->fifo.push_back(&wrkld_per_rank);

    ABT_thread_yield_to(global_prod_thread);
521 522
#endif

523 524 525 526 527 528 529 530 531 532 533
#ifdef DBG_COMM
        auto it = allreduce_count.find(count);
        if(it == allreduce_count.end())
        {
            allreduce_count.insert(std::make_pair(count, 1));
        }
        else
        {
            it->second = it->second + 1;
        }
#endif
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
    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);

    int comm_size, i, send_idx, recv_idx, last_idx, send_cnt, recv_cnt;
    int pof2, mask, rem, newrank, newdst, dst, *cnts, *disps;
    int rank = sctx->my_rank;
    comm_size = sctx->num_ranks;

    cnts = disps = NULL;

    pof2 = 1;
    while (pof2 <= comm_size) pof2 <<= 1;
    pof2 >>=1;

    rem = comm_size - pof2;

    /* In the non-power-of-two case, all even-numbered
       processes of rank < 2*rem send their data to
       (rank+1). These even-numbered processes no longer
       participate in the algorithm until the very end. The
       remaining processes form a nice power-of-two. */
    if (rank < 2*rem) {
        if (rank % 2 == 0) { /* even */
563
            SWM_Send(rank+1, comm_id, 1235, sendreqvc, sendrspvc, 0, count, 1, 0, 0);
564 565
            newrank = -1;
        } else { /* odd */
566
            SWM_Recv(rank-1, comm_id, 1235, 0);
567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588
            newrank = rank / 2;
        }
    } else {
        newrank = rank - rem;
    }

    /* If op is user-defined or count is less than pof2, use
       recursive doubling algorithm. Otherwise do a reduce-scatter
       followed by allgather. (If op is user-defined,
       derived datatypes are allowed and the user could pass basic
       datatypes on one process and derived on another as long as
       the type maps are the same. Breaking up derived
       datatypes to do the reduce-scatter is tricky, therefore
       using recursive doubling in that case.) */
    if (newrank != -1) {
        if ((count <= ALLREDUCE_SHORT_MSG_SIZE) || (count < pof2)) {

            mask = 0x1;
            while (mask < pof2) {
                newdst = newrank ^ mask;
                dst = (newdst < rem) ? newdst*2 + 1 : newdst + rem;

589 590
                SWM_Sendrecv(comm_id, dst, 1235, sendreqvc, sendrspvc, 0,
                        count, 1, dst, 1235, 0, 0, 0);
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

                mask <<= 1;
            }
        } else {
            /* do a reduce-scatter followed by allgather */
            /* for the reduce-scatter, calculate the count that
               each process receives and the displacement within
               the buffer */

            cnts = (int*)malloc(pof2*sizeof(int));
            disps = (int*)malloc(pof2*sizeof(int));

            for (i=0; i<(pof2-1); i++)
                cnts[i] = count/pof2;
            cnts[pof2-1] = count - (count/pof2)*(pof2-1);

            disps[0] = 0;
            for (i=1; i<pof2; i++)
                disps[i] = disps[i-1] + cnts[i-1];

            mask = 0x1;
            send_idx = recv_idx = 0;
            last_idx = pof2;
            while (mask < pof2) {
                newdst = newrank ^ mask;
                dst = (newdst < rem) ? newdst*2 + 1 : newdst + rem;
                send_cnt = recv_cnt = 0;
                if (newrank < newdst) {
                    send_idx = recv_idx + pof2/(mask*2);
                    for (i=send_idx; i<last_idx; i++)
                        send_cnt += cnts[i];
                    for (i=recv_idx; i<send_idx; i++)
                        recv_cnt += cnts[i];
                } else {
                    recv_idx = send_idx + pof2/(mask*2);
                    for (i=send_idx; i<recv_idx; i++)
                        send_cnt += cnts[i];
                    for (i=recv_idx; i<last_idx; i++)
                        recv_cnt += cnts[i];
                }

632 633
                SWM_Sendrecv(comm_id, dst, 1235, sendreqvc, sendrspvc, 0,
                        send_cnt, 1, dst, 1235, 0, 0, 0);
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

                send_idx = recv_idx;
                mask <<= 1;

                if(mask < pof2)
                    last_idx = recv_idx + pof2/mask;
            }

            /* now do the allgather */
            mask >>= 1;
            while (mask > 0) {
                newdst = newrank ^ mask;
                /* find real rank of dest */
                dst = (newdst < rem) ? newdst*2 + 1 : newdst + rem;

                send_cnt = recv_cnt = 0;
                if (newrank < newdst) {
                    if (mask != pof2/2)
                        last_idx = last_idx + pof2/(mask*2);

                    recv_idx = send_idx + pof2/(mask*2);
                    for (i=send_idx; i<recv_idx; i++)
                        send_cnt += cnts[i];
                    for (i=recv_idx; i<last_idx; i++)
                        recv_cnt += cnts[i];
                } else {
                    recv_idx = send_idx - pof2/(mask*2);
                    for (i=send_idx; i<last_idx; i++)
                        send_cnt += cnts[i];
                    for (i=recv_idx; i<send_idx; i++)
                        recv_cnt += cnts[i];
                }

667 668
                SWM_Sendrecv(comm_id, dst, 1235, sendreqvc, sendrspvc, 0,
                        send_cnt, 1, dst, 1235, 0, 0, 0);
669 670 671 672 673 674 675 676 677 678

                if (newrank > newdst) send_idx = recv_idx;

                mask >>= 1;
            }
        }
    }

    if(rank < 2*rem) {
        if(rank % 2) {/* odd */
679
            SWM_Send(rank-1, comm_id, 1235, sendreqvc, sendrspvc, 0, count, 1, 0, 0);
680
        } else {
681
            SWM_Recv(rank+1, comm_id, 1235, 0);
682 683 684 685 686
        }
    }

    if(cnts) free(cnts);
    if(disps) free(disps);
687 688

    num_allreduce++;
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
}

void SWM_Allreduce(
        SWM_BYTES bytes,
        SWM_BYTES respbytes,
        SWM_COMM_ID comm_id,
        SWM_VC sendreqvc,
        SWM_VC sendrspvc,
        SWM_BUF sendbuf,
        SWM_BUF rcvbuf,
        SWM_UNKNOWN auto1,
        SWM_UNKNOWN2 auto2,
        SWM_ROUTING_TYPE reqrt,
        SWM_ROUTING_TYPE rsprt)
{
    SWM_Allreduce(bytes, respbytes, comm_id, sendreqvc, sendrspvc, sendbuf, rcvbuf);
}

void SWM_Finalize()
{
    /* Add an event in the shared queue and then yield */
    struct codes_workload_op wrkld_per_rank;

    wrkld_per_rank.op_type = CODES_WK_END;

    /* Retreive the shared context state */
    ABT_thread prod;
    void * arg;
    int err = ABT_thread_self(&prod);
    assert(err == ABT_SUCCESS);
    err =  ABT_thread_get_arg(prod, &arg);
    assert(err == ABT_SUCCESS);
    struct shared_context * sctx = static_cast<shared_context*>(arg);
    sctx->fifo.push_back(&wrkld_per_rank);

724
#ifdef DBG_COMM 
725
/*    auto it = allreduce_count.begin();
726 727 728 729 730 731 732 733 734 735 736 737 738 739 740
    for(; it != allreduce_count.end(); it++)
    {
        cout << "\n Allreduce " << it->first << " " << it->second;
    }
    
    it = send_count.begin();
    for(; it != send_count.end(); it++)
    {
        cout << "\n Send " << it->first << " " << it->second;
    }
    
    it = isend_count.begin();
    for(; it != isend_count.end(); it++)
    {
        cout << "\n isend " << it->first << " " << it->second;
741
    }*/
742
#endif
743 744
//#ifdef DBG_COMM
//    printf("\n finalize workload for rank %d ", sctx->my_rank);
745
//    printf("\n finalize workload for rank %d num_sends %d num_recvs %d num_isends %d num_irecvs %d num_allreduce %d num_barrier %d num_waitalls %d", sctx->my_rank, num_sends, num_recvs, num_isends, num_irecvs, num_allreduce, num_barriers, num_waitalls);
746
//#endif
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
    ABT_thread_yield_to(global_prod_thread);
}

static int hash_rank_compare(void *key, struct qhash_head *link)
{
    rank_mpi_compare *in = (rank_mpi_compare*)key;
    rank_mpi_context *tmp;

    tmp = qhash_entry(link, rank_mpi_context, hash_link);
    if (tmp->sctx.my_rank == in->rank && tmp->app_id == in->app_id)
        return 1;
    return 0;
}
static void workload_caller(void * arg)
{
    shared_context* sctx = static_cast<shared_context*>(arg);

    if(strcmp(sctx->workload_name, "lammps") == 0)
    {
        LAMMPS_SWM * lammps_swm = static_cast<LAMMPS_SWM*>(sctx->swm_obj);
        lammps_swm->call();
    }
    else if(strcmp(sctx->workload_name, "nekbone") == 0) 
    {
        NEKBONESWMUserCode * nekbone_swm = static_cast<NEKBONESWMUserCode*>(sctx->swm_obj);
        nekbone_swm->call();
    }
774 775 776 777 778
    else if(strcmp(sctx->workload_name, "nearest_neighbor") == 0)
    {
       NearestNeighborSWMUserCode * nn_swm = static_cast<NearestNeighborSWMUserCode*>(sctx->swm_obj);
       nn_swm->call();
    }
779 780 781 782 783
    else if(strcmp(sctx->workload_name, "incast") == 0 || strcmp(sctx->workload_name, "incast1") == 0 || strcmp(sctx->workload_name, "incast2") == 0)
    {
       AllToOneSWMUserCode * incast_swm = static_cast<AllToOneSWMUserCode*>(sctx->swm_obj);
       incast_swm->call();
    }
784 785 786 787 788 789
}
static int comm_online_workload_load(const char * params, int app_id, int rank)
{
    /* LOAD parameters from JSON file*/
    online_comm_params * o_params = (online_comm_params*)params;
    int nprocs = o_params->nprocs;
790

791 792
    rank_mpi_context *my_ctx = new rank_mpi_context;
    //my_ctx = (rank_mpi_context*)caloc(1, sizeof(rank_mpi_context));  
793 794
    assert(my_ctx); 
    my_ctx->sctx.my_rank = rank; 
795
    my_ctx->sctx.num_ranks = nprocs;
796
    my_ctx->sctx.wait_id = 0;
797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816
    my_ctx->app_id = app_id;

    void** generic_ptrs;
    int array_len = 1;
    generic_ptrs = (void**)calloc(array_len,  sizeof(void*));
    generic_ptrs[0] = (void*)&rank;

    strcpy(my_ctx->sctx.workload_name, o_params->workload_name);
    boost::property_tree::ptree root;
    string path;
    path.append(SWM_DATAROOTDIR);

    if(strcmp(o_params->workload_name, "lammps") == 0)
    {
        path.append("/lammps_workload.json");
    }
    else if(strcmp(o_params->workload_name, "nekbone") == 0)
    {
        path.append("/workload.json"); 
    }
817 818 819 820
    else if(strcmp(o_params->workload_name, "nearest_neighbor") == 0)
    {
        path.append("/skeleton.json"); 
    }
821 822 823 824 825 826 827 828 829 830 831 832
    else if(strcmp(o_params->workload_name, "incast") == 0)
    {
        path.append("/incast.json"); 
    }
    else if(strcmp(o_params->workload_name, "incast1") == 0)
    {
        path.append("/incast1.json"); 
    }
    else if(strcmp(o_params->workload_name, "incast2") == 0)
    {
        path.append("/incast2.json"); 
    }
833 834 835
    else
        tw_error(TW_LOC, "\n Undefined workload type %s ", o_params->workload_name);

836
    try {
837
        std::ifstream jsonFile(path.c_str());
838 839
        boost::property_tree::json_parser::read_json(jsonFile, root);
        uint32_t process_cnt = root.get<uint32_t>("jobs.size", 1);
840
        cpu_freq = root.get<double>("jobs.cfg.cpu_freq") / 1e9; 
841 842 843 844 845 846
    }
    catch(std::exception & e)
    {
        printf("%s \n", e.what());
        return -1;
    }
847 848 849 850 851 852 853
    if(strcmp(o_params->workload_name, "lammps") == 0)
    {
        LAMMPS_SWM * lammps_swm = new LAMMPS_SWM(root, generic_ptrs);
        my_ctx->sctx.swm_obj = (void*)lammps_swm;
    }
    else if(strcmp(o_params->workload_name, "nekbone") == 0)
    {
854 855
        NEKBONESWMUserCode * nekbone_swm = new NEKBONESWMUserCode(root, generic_ptrs);
        my_ctx->sctx.swm_obj = (void*)nekbone_swm;
856
    }
857 858 859 860 861
    else if(strcmp(o_params->workload_name, "nearest_neighbor") == 0)
    {
        NearestNeighborSWMUserCode * nn_swm = new NearestNeighborSWMUserCode(root, generic_ptrs);
        my_ctx->sctx.swm_obj = (void*)nn_swm;
    }
862 863 864 865 866
    else if(strcmp(o_params->workload_name, "incast") == 0 || strcmp(o_params->workload_name, "incast1") == 0 || strcmp(o_params->workload_name, "incast2") == 0)
    {
        AllToOneSWMUserCode * incast_swm = new AllToOneSWMUserCode(root, generic_ptrs);
        my_ctx->sctx.swm_obj = (void*)incast_swm;
    }
867 868

    if(global_prod_thread == NULL)
869 870
    {
        ABT_xstream_self(&self_es);
871
        ABT_thread_self(&global_prod_thread);
872
    }
873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894
    ABT_thread_create_on_xstream(self_es, 
            &workload_caller, (void*)&(my_ctx->sctx),
            ABT_THREAD_ATTR_NULL, &(my_ctx->sctx.producer));

    rank_mpi_compare cmp;
    cmp.app_id = app_id;
    cmp.rank = rank;

    if(!rank_tbl)
    {
        rank_tbl = qhash_init(hash_rank_compare, quickhash_64bit_hash, nprocs);
        if(!rank_tbl)
            return -1;
    }
    qhash_add(rank_tbl, &cmp, &(my_ctx->hash_link));
    rank_tbl_pop++;

    return 0;
}

static void comm_online_workload_get_next(int app_id, int rank, struct codes_workload_op * op)
{
895 896 897
    /* At this point, we will use the "call" function. The send/receive/wait
     * definitions will be replaced by our own function definitions that will do a
     * yield to argobots if an event is not available. */
898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918
    /* if shared queue is empty then yield */

    rank_mpi_context * temp_data;
    struct qhash_head * hash_link = NULL;
    rank_mpi_compare cmp;
    cmp.rank = rank;
    cmp.app_id = app_id;
    hash_link = qhash_search(rank_tbl, &cmp);
    if(!hash_link)
    {
        printf("\n not found for rank id %d , %d", rank, app_id);
        op->op_type = CODES_WK_END;
        return;
    }
    temp_data = qhash_entry(hash_link, rank_mpi_context, hash_link);
    assert(temp_data);
    while(temp_data->sctx.fifo.empty())
    {
        ABT_thread_yield_to(temp_data->sctx.producer); 
    }
    struct codes_workload_op * front_op = temp_data->sctx.fifo.front();
919
    assert(front_op);
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
    *op = *front_op;
    temp_data->sctx.fifo.pop_front();
    return;
}
static int comm_online_workload_get_rank_cnt(const char *params, int app_id)
{
    online_comm_params * o_params = (online_comm_params*)params;
    int nprocs = o_params->nprocs;
    return nprocs;
}

static int comm_online_workload_finalize(const char* params, int app_id, int rank)
{
    rank_mpi_context * temp_data;
    struct qhash_head * hash_link = NULL;
    rank_mpi_compare cmp;
    cmp.rank = rank;
    cmp.app_id = app_id;
    hash_link = qhash_search(rank_tbl, &cmp);
    if(!hash_link)
    {
        printf("\n not found for rank id %d , %d", rank, app_id);
        return -1;
    }
    temp_data = qhash_entry(hash_link, rank_mpi_context, hash_link);
    assert(temp_data);

    ABT_thread_join(temp_data->sctx.producer);    
    ABT_thread_free(&(temp_data->sctx.producer));
    return 0;
}
951
extern "C" {
952
/* workload method name and function pointers for the CODES workload API */
953
struct codes_workload_method online_comm_workload_method =
954
{
955
    //.method_name =
956
    (char*)"online_comm_workload",
957 958 959 960 961 962 963 964 965 966 967 968
    //.codes_workload_read_config = 
    NULL,
    //.codes_workload_load = 
    comm_online_workload_load,
    //.codes_workload_get_next = 
    comm_online_workload_get_next,
    // .codes_workload_get_next_rc2 = 
    NULL,
    // .codes_workload_get_rank_cnt
    comm_online_workload_get_rank_cnt,
    // .codes_workload_finalize = 
    comm_online_workload_finalize
969
};
970
} // closing brace for extern "C"
971