Skip to content
Snippets Groups Projects
dispatch.c 24.5 KiB
Newer Older
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
/**
 * Author: Baptiste Coudray
 * School: HEPIA
 * Class: ITI-3
 * Year: 2020-2021
 */

#include <stdio.h>
#include <stdint.h>
#include <mpi.h>
#include <math.h>
#include <stdbool.h>
#include <stdlib.h>
#include <string.h>
#include "dispatch.h"
#include "envelope.h"

#define min(a, b) (((a) <= (b)) ? (a) : (b))
#define max(a, b) (((a) >= (b)) ? (a) : (b))

#define INDEX_2D_TO_1D(y, x, nb_columns) ((y) * (nb_columns) + (x))

#define NORTH_ROW_TAG 0
#define EAST_COLUMN_TAG 1
#define SOUTH_ROW_TAG 2
#define WEST_COLUMN_TAG 3

#define NORTH_EAST_CELLS_TAG 4
#define SOUTH_EAST_CELLS_TAG 5
#define SOUTH_WEST_CELLS_TAG 6
#define NORTH_WEST_CELLS_TAG 7

#define CHUNK_DATA_TAG 8
#define ROOT_RANK 0

struct dispatch_context {
    int my_rank;
    int world_size;
    int my_cart_rank;
    int coordinates[2];
    MPI_Comm communicators[3]; /* cart_comm, row_comm, column_comm */
    int network_dimensions[2];
    int data_dimensions[2];
    int n_dimensions;
    chunk_info_t *chunk_info;
    int type;
    chunk_info_t *chunks_info;
};

static void get_world_size(struct dispatch_context *dc) {
    MPI_Comm_size(MPI_COMM_WORLD, &dc->world_size);
}

static void get_my_rank(struct dispatch_context *dc) {
    MPI_Comm_rank(MPI_COMM_WORLD, &dc->my_rank);
}

static void find_network_dimensions(struct dispatch_context *dc) {
    /* 1D */
    // | R0 <-> R1 <-> R2 <-> ... <-> Rn |
    if (dc->data_dimensions[0] == 1) {
        dc->network_dimensions[0] = 1;
        dc->network_dimensions[1] = dc->world_size;
    } else {
        /* 2D */
        int limit = (int) sqrt(dc->world_size);
        bool first_pass = true;
        for (int i = 1; i <= limit; ++i) {
            if (dc->world_size % i == 0) {
                int new_grid_ny = i;
                int new_grid_nx = dc->world_size / i;
                int current_difference = abs(dc->network_dimensions[0] - dc->network_dimensions[1]);
                int new_difference = abs(new_grid_ny - new_grid_nx);
                if (first_pass || current_difference > new_difference) {
                    dc->network_dimensions[0] = new_grid_ny;
                    dc->network_dimensions[1] = new_grid_nx;
                }
                first_pass = false;
            }
        }
    }
}

static void create_network_communicators(struct dispatch_context *dc) {
    int periods[2] = {true, true}; // Cyclic on row-column-depth

    MPI_Cart_create(MPI_COMM_WORLD, 2, dc->network_dimensions, periods, 1, &dc->communicators[0]);

    /* Create row communicator */
    int remain_dims[2] = {false, true};
    MPI_Cart_sub(dc->communicators[0], remain_dims, &dc->communicators[1]);

    /* Create column communicator */
    remain_dims[0] = true; // row
    remain_dims[1] = false; // column
    MPI_Cart_sub(dc->communicators[0], remain_dims, &dc->communicators[2]);

    MPI_Comm_rank(dc->communicators[0], &dc->my_cart_rank);

    MPI_Cart_coords(dc->communicators[0], dc->my_cart_rank, 2, dc->coordinates);
}

static void divide_data(struct dispatch_context *dc) {
    dc->chunks_info = calloc((size_t) dc->world_size, sizeof(chunk_info_t));

    int nb_rows_per_process = dc->data_dimensions[0] / dc->network_dimensions[0];
    int remaining_rows = dc->data_dimensions[0] % dc->network_dimensions[0];

    int nb_columns_per_process = dc->data_dimensions[1] / dc->network_dimensions[1];
    int remaining_columns = dc->data_dimensions[1] % dc->network_dimensions[1];

    for (int i = 0, y = 0, x = 0; i < dc->world_size; ++i) {
        int nb_rows = nb_rows_per_process;
        if (remaining_rows > 0) {
            ++nb_rows;
        }
        int nb_columns = nb_columns_per_process;
        if (remaining_columns > 0) {
            ++nb_columns;
            --remaining_columns;
        }
        dc->chunks_info[i].y = y;
        dc->chunks_info[i].x = x;
        dc->chunks_info[i].dimensions[0] = nb_rows;
        dc->chunks_info[i].dimensions[1] = nb_columns;
        dc->chunks_info[i].count = nb_rows * nb_columns;

        x += nb_columns;
        if (x >= dc->data_dimensions[1]) {
            x = 0;
            y += dc->chunks_info[max(i - 1, 0)].dimensions[0];
            remaining_columns = dc->data_dimensions[1] % dc->network_dimensions[1];
            --remaining_rows;
        }
    }
    dc->chunk_info = &dc->chunks_info[dc->my_rank];
}

static void create_chunk_data(struct dispatch_context *dc, int type) {
    dc->type = type;
    dc->chunk_info->data = calloc((size_t) dc->chunk_info->count, (size_t) dc->type);
}

extern struct dispatch_context *dispatch_context_new(int *dimensions, int type, int n_dimensions) {
    struct dispatch_context *dispatch_context = calloc(1, sizeof(struct dispatch_context));
    get_world_size(dispatch_context);
    get_my_rank(dispatch_context);
    dispatch_context->n_dimensions = n_dimensions;

    switch (n_dimensions) {
        case 1:
            dispatch_context->data_dimensions[0] = 1;
            dispatch_context->data_dimensions[1] = dimensions[0];
            break;
        case 2:
            dispatch_context->data_dimensions[0] = dimensions[0];
            dispatch_context->data_dimensions[1] = dimensions[1];
            break;
        case 3:
            return NULL;
        default:
            fprintf(stderr, "Invalid dimensions size.");
            MPI_Abort(MPI_COMM_WORLD, 1);
            break;
    }

    find_network_dimensions(dispatch_context);
    create_network_communicators(dispatch_context);
    divide_data(dispatch_context);
    create_chunk_data(dispatch_context, type);
    return dispatch_context;
}

extern void dispatch_context_print(struct dispatch_context *dc) {
    printf("R%d: Nd = [%d][%d], CDd = [%d][%d], Indexes = [%d][%d]\n",
           dc->my_rank, dc->network_dimensions[0], dc->network_dimensions[1],
           dc->chunk_info->dimensions[0], dc->chunk_info->dimensions[1],
           dc->chunk_info->y, dc->chunk_info->x
    );
}

static envelope_t get_inner_envelope_1d(struct dispatch_context *dc, struct futhark_context *fc, int thickness) {
    struct futhark_u8_1d *fut_chunk_data = futhark_new_u8_1d(fc, dc->chunk_info->data,
                                                             dc->chunk_info->dimensions[1] * dc->type);
    struct futhark_u8_1d *fut_west;
    struct futhark_u8_1d *fut_east;
    int thickness_x = min(thickness, dc->chunk_info->dimensions[1]);

    futhark_entry_get_envelope_1d(fc, &fut_west, &fut_east, fut_chunk_data, thickness_x * dc->type);

    envelope_t inner_envelope = (envelope_t) {0};
    inner_envelope.west.data = calloc((size_t) thickness_x, (size_t) dc->type);
    inner_envelope.west.dimensions[0] = 1;
    inner_envelope.west.dimensions[1] = thickness_x;

    inner_envelope.east.data = calloc((size_t) thickness_x, (size_t) dc->type);
    inner_envelope.east.dimensions[0] = 1;
    inner_envelope.east.dimensions[1] = thickness_x;

    futhark_values_u8_1d(fc, fut_west, inner_envelope.west.data);
    futhark_values_u8_1d(fc, fut_east, inner_envelope.east.data);

    return inner_envelope;
}

static envelope_t get_inner_envelope_2d(struct dispatch_context *dc, struct futhark_context *fc, int thickness) {
    struct futhark_u8_2d *fut_chunk_data = futhark_new_u8_2d(fc, dc->chunk_info->data,
                                                             dc->chunk_info->dimensions[0],
                                                             dc->chunk_info->dimensions[1] * dc->type);

    struct futhark_u8_2d *fut_north;
    struct futhark_u8_2d *fut_north_east;
    struct futhark_u8_2d *fut_east;
    struct futhark_u8_2d *fut_south_east;
    struct futhark_u8_2d *fut_south;
    struct futhark_u8_2d *fut_south_west;
    struct futhark_u8_2d *fut_west;
    struct futhark_u8_2d *fut_north_west;

    int thickness_y = min(thickness, dc->chunk_info->dimensions[0]);
    int thickness_x = min(thickness, dc->chunk_info->dimensions[1]);

    futhark_entry_get_envelope_2d(fc, &fut_north, &fut_north_east, &fut_east, &fut_south_east, &fut_south,
                                  &fut_south_west, &fut_west, &fut_north_west, fut_chunk_data, thickness_y,
                                  thickness_x * dc->type);

    size_t s_thickness_y = (size_t) thickness_x;
    size_t s_thickness_x = (size_t) thickness_y;
    size_t s_type = (size_t) dc->type;

    envelope_t inner_envelope = (envelope_t) {0};
    inner_envelope.north.data = calloc(s_thickness_y * (size_t) dc->chunk_info->dimensions[1], s_type);
    inner_envelope.east.data = calloc((size_t) dc->chunk_info->dimensions[0] * s_thickness_x, s_type);
    inner_envelope.south.data = calloc(s_thickness_y * (size_t) dc->chunk_info->dimensions[1], s_type);
    inner_envelope.west.data = calloc((size_t) dc->chunk_info->dimensions[0] * s_thickness_x, s_type);
    inner_envelope.north_east.data = calloc(s_thickness_y * s_thickness_x, s_type);
    inner_envelope.south_east.data = calloc(s_thickness_y * s_thickness_x, s_type);
    inner_envelope.south_west.data = calloc(s_thickness_y * s_thickness_x, s_type);
    inner_envelope.north_west.data = calloc(s_thickness_y * s_thickness_x, s_type);

    futhark_values_u8_2d(fc, fut_north, inner_envelope.north.data);
    futhark_values_u8_2d(fc, fut_east, inner_envelope.east.data);
    futhark_values_u8_2d(fc, fut_south, inner_envelope.south.data);
    futhark_values_u8_2d(fc, fut_west, inner_envelope.west.data);
    futhark_values_u8_2d(fc, fut_north_east, inner_envelope.north_east.data);
    futhark_values_u8_2d(fc, fut_south_east, inner_envelope.south_east.data);
    futhark_values_u8_2d(fc, fut_south_west, inner_envelope.south_west.data);
    futhark_values_u8_2d(fc, fut_north_west, inner_envelope.north_west.data);

    return inner_envelope;
}

static envelope_t get_outer_envelope_1d(struct dispatch_context *dc, int thickness, envelope_t *inner_envelope) {
    int coordinate_x = dc->coordinates[1];
    MPI_Request requests[4];
    int i_request = 0;
    envelope_t outer_envelope = (envelope_t) {0};

    // West-part
    {
        int dest_source_x = (coordinate_x - 1 >= 0) ? coordinate_x - 1 : dc->network_dimensions[1] - 1;
        int dest_source = INDEX_2D_TO_1D(0, dest_source_x, dc->network_dimensions[1]);

        int send_count = min(thickness, dc->chunk_info->dimensions[1]);
        int recv_count = min(thickness, dc->chunks_info[dest_source].dimensions[1]);

        outer_envelope.west.data = calloc((size_t) recv_count, (size_t) dc->type);

        void *inner_envelope_west = ((uint8_t *) inner_envelope->west.data)
                                    + INDEX_2D_TO_1D(0, 0, inner_envelope->west.dimensions[1])
                                      * dc->type;
        MPI_Isend(inner_envelope_west, send_count * dc->type, MPI_UINT8_T, dest_source_x, WEST_COLUMN_TAG,
                  dc->communicators[1], &requests[i_request++]);

        void *outer_envelope_west = ((uint8_t *) outer_envelope.west.data)
                                    + INDEX_2D_TO_1D(0, 0, recv_count)
                                      * dc->type;
        MPI_Irecv(outer_envelope_west, recv_count * dc->type, MPI_UINT8_T, dest_source_x, EAST_COLUMN_TAG,
                  dc->communicators[1], &requests[i_request++]);
    }

    // East-part
    {
        int dest_source_x = (coordinate_x + 1) % dc->network_dimensions[1];
        int dest_source = INDEX_2D_TO_1D(0, dest_source_x, dc->network_dimensions[1]);

        int send_count = min(thickness, dc->chunk_info->dimensions[1]);
        int recv_count = min(thickness, dc->chunks_info[dest_source].dimensions[1]);

        outer_envelope.east.data = calloc((size_t) recv_count, (size_t) dc->type);

        void *inner_envelope_east = ((uint8_t *) inner_envelope->east.data)
                                    + INDEX_2D_TO_1D(0, inner_envelope->east.dimensions[1] - send_count,
                                                     inner_envelope->east.dimensions[1])
                                      * dc->type;

        MPI_Isend(inner_envelope_east, send_count * dc->type, MPI_UINT8_T, dest_source_x, EAST_COLUMN_TAG,
                  dc->communicators[1], &requests[i_request++]);

        void *outer_envelope_east = ((uint8_t *) outer_envelope.east.data)
                                    + INDEX_2D_TO_1D(0, 0, dc->chunk_info->dimensions[1])
                                      * dc->type;
        MPI_Irecv(outer_envelope_east, recv_count * dc->type, MPI_UINT8_T, dest_source_x, WEST_COLUMN_TAG,
                  dc->communicators[1], &requests[i_request]);
    }

    MPI_Waitall(4, requests, MPI_STATUSES_IGNORE);
    return outer_envelope;
}

static envelope_t get_outer_envelope_2d(struct dispatch_context *dc, int thickness, envelope_t *inner_envelope) {
    int coordinate_y = dc->coordinates[0];
    int coordinate_x = dc->coordinates[1];
    MPI_Request requests[16] = {0};
    int i_request = 0;

    envelope_t outer_envelope = (envelope_t) {0};

    // North
    {
        int dest_source_y = (coordinate_y - 1) >= 0 ? (coordinate_y - 1) : (dc->network_dimensions[0] - 1);
        int dest_source = INDEX_2D_TO_1D(dest_source_y, coordinate_x, dc->network_dimensions[1]);

        int send_count = min(thickness, dc->chunk_info->dimensions[0]) * dc->chunk_info->dimensions[1];
        int dimensions[2] = {
                min(thickness, dc->chunks_info[dest_source].dimensions[0]),
                dc->chunks_info[dest_source].dimensions[1]
        };
        chunk_info_init(&outer_envelope.north, dc->type, dimensions, 0, 0);

        MPI_Isend(inner_envelope->north.data, send_count * dc->type, MPI_UINT8_T, dest_source_y, NORTH_ROW_TAG,
                  dc->communicators[2], &requests[i_request++]);

        /* Neighbour send south row, which correspond to north envelope */
        MPI_Irecv(outer_envelope.north.data, outer_envelope.north.count * dc->type, MPI_UINT8_T, dest_source_y,
                  SOUTH_ROW_TAG,
                  dc->communicators[2], &requests[i_request++]);
    }

    // East
    {
        int dest_source_x = (coordinate_x + 1) % dc->network_dimensions[1];
        int dest_source = INDEX_2D_TO_1D(coordinate_y, dest_source_x, dc->network_dimensions[1]);

        int send_count = min(thickness, dc->chunk_info->dimensions[1]) * dc->chunk_info->dimensions[0];
        int dimensions[2] = {
                dc->chunks_info[dest_source].dimensions[0],
                min(thickness, dc->chunks_info[dest_source].dimensions[1])
        };
        chunk_info_init(&outer_envelope.east, dc->type, dimensions, 0, 0);

        MPI_Isend(inner_envelope->east.data, send_count * dc->type, MPI_UINT8_T, dest_source_x, EAST_COLUMN_TAG,
                  dc->communicators[1],
                  &requests[i_request++]);

        /* Neighbour send west column, which correspond to east envelope */
        MPI_Irecv(outer_envelope.east.data, outer_envelope.east.count * dc->type, MPI_UINT8_T, dest_source_x,
                  WEST_COLUMN_TAG,
                  dc->communicators[1],
                  &requests[i_request++]);
    }

    // South
    {
        int dest_source_y = (coordinate_y + 1) % dc->network_dimensions[0];
        int dest_source = INDEX_2D_TO_1D(dest_source_y, coordinate_x, dc->network_dimensions[1]);

        int send_count = min(thickness, dc->chunk_info->dimensions[0]) * dc->chunk_info->dimensions[1];
        int dimensions[2] = {
                min(thickness, dc->chunks_info[dest_source].dimensions[0]),
                dc->chunks_info[dest_source].dimensions[1]
        };
        chunk_info_init(&outer_envelope.south, dc->type, dimensions, 0, 0);

        MPI_Isend(inner_envelope->south.data, send_count * dc->type, MPI_UINT8_T, dest_source_y, SOUTH_ROW_TAG,
                  dc->communicators[2], &requests[i_request++]);

        /* Neighbour send south row, which correspond to north envelope */
        MPI_Irecv(outer_envelope.south.data, outer_envelope.south.count * dc->type, MPI_UINT8_T, dest_source_y,
                  NORTH_ROW_TAG,
                  dc->communicators[2], &requests[i_request++]);
    }

    // West
    {
        int dest_source_x = (coordinate_x - 1) >= 0 ? coordinate_x - 1 : dc->network_dimensions[1] - 1;
        int dest_source = INDEX_2D_TO_1D(coordinate_y, dest_source_x, dc->network_dimensions[1]);

        int send_count = min(thickness, dc->chunk_info->dimensions[1]) * dc->chunk_info->dimensions[0];

        int dimensions[2] = {
                dc->chunks_info[dest_source].dimensions[0],
                min(thickness, dc->chunks_info[dest_source].dimensions[1])
        };
        chunk_info_init(&outer_envelope.west, dc->type, dimensions, 0, 0);

        MPI_Isend(inner_envelope->west.data, send_count * dc->type, MPI_UINT8_T, dest_source_x, WEST_COLUMN_TAG,
                  dc->communicators[1],
                  &requests[i_request++]);

        /* Neighbour send west column, which correspond to east envelope */
        MPI_Irecv(outer_envelope.west.data, outer_envelope.west.count * dc->type, MPI_UINT8_T, dest_source_x,
                  EAST_COLUMN_TAG,
                  dc->communicators[1],
                  &requests[i_request++]);
    }

    // North-East
    {
        int dest_source_y = (coordinate_y - 1) >= 0 ? coordinate_y - 1 : dc->network_dimensions[0] - 1;
        int dest_source_x = (coordinate_x + 1) % dc->network_dimensions[1];
        int dest_source = INDEX_2D_TO_1D(dest_source_y, dest_source_x, dc->network_dimensions[1]);

        int send_count = thickness * thickness;
        int dimensions[2] = {
                thickness,
                thickness,
        };
        chunk_info_init(&outer_envelope.north_east, dc->type, dimensions, 0, 0);

        MPI_Isend(inner_envelope->north_east.data, send_count * dc->type, MPI_UINT8_T, dest_source,
                  NORTH_EAST_CELLS_TAG, MPI_COMM_WORLD,
                  &requests[i_request++]);
        /* Neighbour send south-west cell, which correspond to north-east cell */
        MPI_Irecv(outer_envelope.north_east.data, outer_envelope.north_east.count * dc->type, MPI_UINT8_T, dest_source,
                  SOUTH_WEST_CELLS_TAG, MPI_COMM_WORLD,
                  &requests[i_request++]);
    }

    // South-East
    {
        int dest_source_y = (coordinate_y + 1) % dc->network_dimensions[0];
        int dest_source_x = (coordinate_x + 1) % dc->network_dimensions[1];
        int destSource = INDEX_2D_TO_1D(dest_source_y, dest_source_x, dc->network_dimensions[1]);

        int send_count = thickness * thickness;
        int dimensions[2] = {
                thickness,
                thickness,
        };
        chunk_info_init(&outer_envelope.south_east, dc->type, dimensions, 0, 0);

        MPI_Isend(inner_envelope->south_east.data, send_count * dc->type, MPI_UINT8_T, destSource,
                  SOUTH_EAST_CELLS_TAG, MPI_COMM_WORLD,
                  &requests[i_request++]);
        /* Neighbour send north-west cell, which correspond to south-east cell */
        MPI_Irecv(outer_envelope.south_east.data, outer_envelope.south_east.count * dc->type, MPI_UINT8_T, destSource,
                  NORTH_WEST_CELLS_TAG, MPI_COMM_WORLD,
                  &requests[i_request++]);
    }

    // South-West
    {
        int dest_source_y = (coordinate_y + 1) % dc->network_dimensions[0];
        int dest_source_x = (coordinate_x - 1) >= 0 ? coordinate_x - 1 : dc->network_dimensions[1] - 1;
        int dest_source = INDEX_2D_TO_1D(dest_source_y, dest_source_x, dc->network_dimensions[1]);

        int send_count = thickness * thickness;
        int dimensions[2] = {
                thickness,
                thickness,
        };
        chunk_info_init(&outer_envelope.south_west, dc->type, dimensions, 0, 0);

        MPI_Isend(inner_envelope->south_west.data, send_count * dc->type, MPI_UINT8_T, dest_source,
                  SOUTH_WEST_CELLS_TAG, MPI_COMM_WORLD,
                  &requests[i_request++]);
        /* Neighbour send north-east cell, which correspond to south-west cell */
        MPI_Irecv(outer_envelope.south_west.data, outer_envelope.south_west.count * dc->type, MPI_UINT8_T, dest_source,
                  NORTH_EAST_CELLS_TAG, MPI_COMM_WORLD,
                  &requests[i_request++]);
    }

    // North-West
    {
        int dest_source_y = (coordinate_y - 1) >= 0 ? coordinate_y - 1 : dc->network_dimensions[0] - 1;
        int dest_source_x = (coordinate_x - 1) >= 0 ? coordinate_x - 1 : dc->network_dimensions[1] - 1;
        int dest_source = INDEX_2D_TO_1D(dest_source_y, dest_source_x, dc->network_dimensions[1]);

        int send_count = thickness * thickness;
        int dimensions[2] = {
                thickness,
                thickness,
        };
        chunk_info_init(&outer_envelope.north_west, dc->type, dimensions, 0, 0);

        MPI_Isend(inner_envelope->north_west.data, send_count * dc->type, MPI_INT8_T, dest_source,
                  NORTH_WEST_CELLS_TAG, MPI_COMM_WORLD,
                  &requests[i_request++]);
        /* Neighbour send south-east cell, which correspond to north-west cell */
        MPI_Irecv(outer_envelope.north_west.data, outer_envelope.north_west.count * dc->type, MPI_INT8_T, dest_source,
                  SOUTH_EAST_CELLS_TAG, MPI_COMM_WORLD,
                  &requests[i_request++]);
    }

    MPI_Waitall(i_request, requests, MPI_STATUSES_IGNORE);
    return outer_envelope;
}

envelope_t get_inner_envelope(struct dispatch_context *dc, struct futhark_context *fc, int thickness) {
    envelope_t inner_envelope = {0};
    switch (dc->n_dimensions) {
        case 1:
            inner_envelope = get_inner_envelope_1d(dc, fc, thickness);
            break;
        case 2:
            inner_envelope = get_inner_envelope_2d(dc, fc, thickness);
            break;
        case 3:
            break;
        default:
            fprintf(stderr, "Invalid dimensions size.");
            MPI_Abort(MPI_COMM_WORLD, 1);
            break;
    }
    return inner_envelope;
}

envelope_t get_outer_envelope(struct dispatch_context *dc, struct futhark_context *fc, int thickness) {
    envelope_t inner_envelope = get_inner_envelope(dc, fc, thickness);
    envelope_t outer_envelope = {0};
    switch (dc->n_dimensions) {
        case 1:
            outer_envelope = get_outer_envelope_1d(dc, thickness, &inner_envelope);
            break;
        case 2:
            outer_envelope = get_outer_envelope_2d(dc, thickness, &inner_envelope);
            break;
        case 3:
            break;
        default:
            fprintf(stderr, "Invalid dimensions size.");
            MPI_Abort(MPI_COMM_WORLD, 1);
            break;
    }
    futhark_context_sync(fc);
    envelope_free(&inner_envelope);
    return outer_envelope;
}

static void chunk_data_to_data(struct dispatch_context *dc, void *chunk_data, void *data, int rank) {
    int y = dc->chunks_info[rank].y;
    int x = dc->chunks_info[rank].x;
    uint8_t *data8 = (uint8_t *) data;
    uint8_t *chunk_data8 = (uint8_t *) chunk_data;

    for (int i = 0; i < dc->chunks_info[rank].dimensions[0]; ++i) {
        memcpy(
                data8 + INDEX_2D_TO_1D(y + i, x, dc->data_dimensions[1]),
                chunk_data8 + INDEX_2D_TO_1D(i, 0, dc->chunks_info[rank].dimensions[1]),
                (size_t) dc->chunks_info[rank].dimensions[1] * (size_t) dc->type
        );
//        data8[INDEX_2D_TO_1D(y + i, x + j, n)] = chunkBoard[INDEX_2D_TO_1D(i, j, chunkN)];
    }
}

void *get_chunk_data(struct dispatch_context *dc) {
    return dc->chunk_info->data;
}

void *get_data(struct dispatch_context *dc) {
    void *data = NULL;
    if (dc->my_rank == ROOT_RANK) {
        data = calloc((size_t) dc->data_dimensions[0] * (size_t) dc->data_dimensions[1], (size_t) dc->type);
        chunk_data_to_data(dc, dc->chunk_info->data, data, dc->my_rank);
        for (int i = 0; i < dc->world_size; ++i) {
            if (i != dc->my_rank) {
                dc->chunks_info[i].data = calloc((size_t) dc->chunks_info[i].count, (size_t) dc->type);
                MPI_Recv(dc->chunks_info[i].data, dc->chunks_info[i].count * dc->type, MPI_UINT8_T, i, CHUNK_DATA_TAG,
                         MPI_COMM_WORLD, MPI_STATUS_IGNORE);
                chunk_data_to_data(dc, dc->chunks_info[i].data, data, i);
                free(dc->chunks_info[i].data);
            }
        }
    } else {
        MPI_Send(dc->chunk_info->data, dc->chunk_info->count * dc->type, MPI_UINT8_T, ROOT_RANK, CHUNK_DATA_TAG,
                 MPI_COMM_WORLD);
    }
    return data;
}

chunk_info_t get_chunk_info(struct dispatch_context *dc) {
    return *dc->chunk_info;
}

void dispatch_context_free(struct dispatch_context *dc) {
    free(dc->chunk_info->data);
    free(dc->chunks_info);
    MPI_Comm_free(&dc->communicators[2]);
    MPI_Comm_free(&dc->communicators[1]);
    MPI_Comm_free(&dc->communicators[0]);
}

void envelope_free(envelope_t *envelope) {
    free(envelope->north.data);
    free(envelope->north_east.data);
    free(envelope->east.data);
    free(envelope->south_east.data);
    free(envelope->south_west.data);
    free(envelope->west.data);
    free(envelope->north_west.data);
}