Skip to content
Snippets Groups Projects
dispatch.c 27 KiB
Newer Older
/**
 * 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;
        }
baptiste.coudray's avatar
baptiste.coudray committed

        int dimensions[2] = {nb_rows, nb_columns};
        chunk_info_init(&dc->chunks_info[i], dc->type, dimensions, y, x, i == dc->my_rank);

        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);
}

baptiste.coudray's avatar
baptiste.coudray committed
extern struct dispatch_context *dispatch_context_new(const 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) {
baptiste.coudray's avatar
baptiste.coudray committed
    printf("[dispatch_context] my_rank = %d, world_size = %d, network_dimensions = [%d][%d], n_dimensions = %d, data_dimensions = [%d][%d]\n",
           dc->my_rank, dc->world_size, dc->network_dimensions[0], dc->network_dimensions[1], dc->n_dimensions,
           dc->data_dimensions[0], dc->data_dimensions[1]);
}

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]);
baptiste.coudray's avatar
baptiste.coudray committed
    int dimensions[2] = {1, thickness_x};
baptiste.coudray's avatar
baptiste.coudray committed
    futhark_entry_get_envelope_1d(fc, &fut_west, &fut_east, fut_chunk_data, dimensions[1] * dc->type);

    envelope_t inner_envelope = (envelope_t) {0};
baptiste.coudray's avatar
baptiste.coudray committed
    // West
    {
        int start_x = dc->chunk_info->x;
        chunk_info_init(&inner_envelope.west, dc->type, dimensions, 0, start_x, true);
        futhark_values_u8_1d(fc, fut_west, inner_envelope.west.data);
        futhark_free_u8_1d(fc, fut_west);
    }
    // East
    {
        int start_x = dc->chunk_info->x + dc->chunk_info->dimensions[1] - thickness_x;
        chunk_info_init(&inner_envelope.east, dc->type, dimensions, 0, start_x, true);
        futhark_values_u8_1d(fc, fut_east, inner_envelope.east.data);
        futhark_free_u8_1d(fc, fut_east);
    }
baptiste.coudray's avatar
baptiste.coudray committed
    futhark_free_u8_1d(fc, fut_chunk_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);

baptiste.coudray's avatar
baptiste.coudray committed
    struct futhark_u8_2d *fut_north_west;
    struct futhark_u8_2d *fut_north;
    struct futhark_u8_2d *fut_north_east;
baptiste.coudray's avatar
baptiste.coudray committed
    struct futhark_u8_2d *fut_west;
    struct futhark_u8_2d *fut_east;
    struct futhark_u8_2d *fut_south_west;
baptiste.coudray's avatar
baptiste.coudray committed
    struct futhark_u8_2d *fut_south;
    struct futhark_u8_2d *fut_south_east;

    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);

    envelope_t inner_envelope = (envelope_t) {0};

baptiste.coudray's avatar
baptiste.coudray committed
    // North-West
    {
        int dimensions[2] = {thickness_y, thickness_x};
        int start_y = dc->chunk_info->y;
        int start_x = dc->chunk_info->x;
        chunk_info_init(&inner_envelope.north_west, dc->type, dimensions, start_y, start_x, true);
        futhark_values_u8_2d(fc, fut_north_west, inner_envelope.north_west.data);
        futhark_free_u8_2d(fc, fut_north_west);
    }
    // North
    {
        int dimensions[2] = {thickness_y, dc->chunk_info->dimensions[1]};
        int start_y = dc->chunk_info->y;
        int start_x = dc->chunk_info->x;
        chunk_info_init(&inner_envelope.north, dc->type, dimensions, start_y, start_x, true);
        futhark_values_u8_2d(fc, fut_north, inner_envelope.north.data);
        futhark_free_u8_2d(fc, fut_north);
    }
    // North-East
    {
        int dimensions[2] = {thickness_y, thickness_x};
        int start_y = dc->chunk_info->y;
        int start_x = dc->chunk_info->x + dc->chunk_info->dimensions[1] - thickness_x;
        chunk_info_init(&inner_envelope.north_east, dc->type, dimensions, start_y, start_x, true);
        futhark_values_u8_2d(fc, fut_north_east, inner_envelope.north_east.data);
        futhark_free_u8_2d(fc, fut_north_east);
    }
    // West
    {
        int dimensions[2] = {dc->chunk_info->dimensions[0], thickness_x};
        int start_y = dc->chunk_info->y;
        int start_x = dc->chunk_info->x;
        chunk_info_init(&inner_envelope.west, dc->type, dimensions, start_y, start_x, true);
        futhark_values_u8_2d(fc, fut_west, inner_envelope.west.data);
        futhark_free_u8_2d(fc, fut_west);
    }
    // East
    {
        int dimensions[2] = {dc->chunk_info->dimensions[0], thickness_x};
        int start_y = dc->chunk_info->y;
        int start_x = dc->chunk_info->x + dc->chunk_info->dimensions[1] - thickness_x;
        chunk_info_init(&inner_envelope.east, dc->type, dimensions, start_y, start_x, true);
        futhark_values_u8_2d(fc, fut_east, inner_envelope.east.data);
        futhark_free_u8_2d(fc, fut_east);
    }
    // South-West
    {
        int dimensions[2] = {thickness_y, thickness_x};
        int start_y = dc->chunk_info->y + dc->chunk_info->dimensions[0] - thickness_y;
        int start_x = dc->chunk_info->x;
        chunk_info_init(&inner_envelope.south_west, dc->type, dimensions, start_y, start_x, true);
        futhark_values_u8_2d(fc, fut_south_west, inner_envelope.south_west.data);
        futhark_free_u8_2d(fc, fut_south_west);
    }
    // South
    {
        int dimensions[2] = {thickness_y, dc->chunk_info->dimensions[1]};
        int start_y = dc->chunk_info->y + dc->chunk_info->dimensions[0] - thickness_y;
        int start_x = dc->chunk_info->x;
        chunk_info_init(&inner_envelope.south, dc->type, dimensions, start_y, start_x, true);
        futhark_values_u8_2d(fc, fut_south, inner_envelope.south.data);
        futhark_free_u8_2d(fc, fut_south);
    }
    // South-East
    {
        int dimensions[2] = {thickness_y, thickness_x};
        int start_y = dc->chunk_info->y + dc->chunk_info->dimensions[0] - thickness_y;
        int start_x = dc->chunk_info->x + dc->chunk_info->dimensions[1] - thickness_x;
        chunk_info_init(&inner_envelope.south_east, dc->type, dimensions, start_y, start_x, true);
        futhark_values_u8_2d(fc, fut_south_east, inner_envelope.south_east.data);
        futhark_free_u8_2d(fc, fut_south_east);
    }
    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]);

baptiste.coudray's avatar
baptiste.coudray committed
        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++]);

baptiste.coudray's avatar
baptiste.coudray committed
        int dimensions[2] = {1, min(thickness, dc->chunks_info[dest_source].dimensions[1])};
        int start_x = dc->chunks_info[dest_source].x + dc->chunks_info[dest_source].dimensions[1] - dimensions[1];
        chunk_info_init(&outer_envelope.west, dc->type, dimensions, 0, start_x, true);

        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++]);
    }

    // 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]);

        void *inner_envelope_east = ((uint8_t *) inner_envelope->east.data)
baptiste.coudray's avatar
baptiste.coudray committed
                                    + ((inner_envelope->east.dimensions[1] - send_count) * 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++]);

baptiste.coudray's avatar
baptiste.coudray committed
        int dimensions[2] = {1, min(thickness, dc->chunks_info[dest_source].dimensions[1])};
        int start_x = dc->chunks_info[dest_source].x + dimensions[1];
        chunk_info_init(&outer_envelope.east, dc->type, dimensions, 0, start_x, true);

        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]);
    }

    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];

        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 */
baptiste.coudray's avatar
baptiste.coudray committed
        int dimensions[2] = {
                min(thickness, dc->chunks_info[dest_source].dimensions[0]),
                dc->chunks_info[dest_source].dimensions[1]
        };
        int start_y = dc->chunks_info[dest_source].y + dc->chunks_info[dest_source].dimensions[0] - dimensions[0];
        int start_x = dc->chunks_info[dest_source].x;
        chunk_info_init(&outer_envelope.north, dc->type, dimensions, start_y, start_x, true);

        MPI_Irecv(outer_envelope.north.data, outer_envelope.north.count * dc->type, MPI_UINT8_T, dest_source_y,
baptiste.coudray's avatar
baptiste.coudray committed
                  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];

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

        /* Neighbour send west column, which correspond to east envelope */
baptiste.coudray's avatar
baptiste.coudray committed
        int dimensions[2] = {
                dc->chunks_info[dest_source].dimensions[0],
                min(thickness, dc->chunks_info[dest_source].dimensions[1])
        };
        int start_y = dc->chunks_info[dest_source].y;
        int start_x = dc->chunks_info[dest_source].x + dc->chunks_info[dest_source].dimensions[1] - dimensions[1];
        chunk_info_init(&outer_envelope.east, dc->type, dimensions, start_y, start_x, true);
        MPI_Irecv(outer_envelope.east.data, outer_envelope.east.count * dc->type, MPI_UINT8_T, dest_source_x,
baptiste.coudray's avatar
baptiste.coudray committed
                  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];

        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++]);

baptiste.coudray's avatar
baptiste.coudray committed
        /* Neighbour send north row, which correspond to south envelope */
        int dimensions[2] = {
                min(thickness, dc->chunks_info[dest_source].dimensions[0]),
                dc->chunks_info[dest_source].dimensions[1]
        };
        int start_y = dc->chunks_info[dest_source].y;
        int start_x = dc->chunks_info[dest_source].x;
        chunk_info_init(&outer_envelope.south, dc->type, dimensions, start_y, start_x, true);
        MPI_Irecv(outer_envelope.south.data, outer_envelope.south.count * dc->type, MPI_UINT8_T, dest_source_y,
baptiste.coudray's avatar
baptiste.coudray committed
                  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];
baptiste.coudray's avatar
baptiste.coudray committed
        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++]);
baptiste.coudray's avatar
baptiste.coudray committed
        /* Neighbour send west column, which correspond to east envelope */
        int dimensions[2] = {
                dc->chunks_info[dest_source].dimensions[0],
                min(thickness, dc->chunks_info[dest_source].dimensions[1])
        };
baptiste.coudray's avatar
baptiste.coudray committed
        int start_y = dc->chunks_info[dest_source].y;
        int start_x = dc->chunks_info[dest_source].x;
        chunk_info_init(&outer_envelope.west, dc->type, dimensions, start_y, start_x, true);
        MPI_Irecv(outer_envelope.west.data, outer_envelope.west.count * dc->type, MPI_UINT8_T, dest_source_x,
baptiste.coudray's avatar
baptiste.coudray committed
                  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;
        MPI_Isend(inner_envelope->north_east.data, send_count * dc->type, MPI_UINT8_T, dest_source,
baptiste.coudray's avatar
baptiste.coudray committed
                  NORTH_EAST_CELLS_TAG, MPI_COMM_WORLD, &requests[i_request++]);

        /* Neighbour send south-west cell, which correspond to north-east cell */
baptiste.coudray's avatar
baptiste.coudray committed
        int dimensions[2] = {thickness, thickness};
        int start_y = dc->chunks_info[dest_source].y + dc->chunks_info[dest_source].dimensions[0] + dimensions[0];
        int start_x = dc->chunks_info[dest_source].x;
        chunk_info_init(&outer_envelope.north_east, dc->type, dimensions, start_y, start_x, true);
        MPI_Irecv(outer_envelope.north_east.data, outer_envelope.north_east.count * dc->type, MPI_UINT8_T, dest_source,
baptiste.coudray's avatar
baptiste.coudray committed
                  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];
baptiste.coudray's avatar
baptiste.coudray committed
        int dest_source = INDEX_2D_TO_1D(dest_source_y, dest_source_x, dc->network_dimensions[1]);

        int send_count = thickness * thickness;
baptiste.coudray's avatar
baptiste.coudray committed
        MPI_Isend(inner_envelope->south_east.data, send_count * dc->type, MPI_UINT8_T, dest_source,
                  SOUTH_EAST_CELLS_TAG, MPI_COMM_WORLD, &requests[i_request++]);

        /* Neighbour send north-west cell, which correspond to south-east cell */
baptiste.coudray's avatar
baptiste.coudray committed
        int dimensions[2] = {thickness, thickness};
        int start_y = dc->chunks_info[dest_source].y;
        int start_x = dc->chunks_info[dest_source].x;
        chunk_info_init(&outer_envelope.south_east, dc->type, dimensions, start_y, start_x, true);
        MPI_Irecv(outer_envelope.south_east.data, outer_envelope.south_east.count * dc->type, MPI_UINT8_T, dest_source,
                  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;
        MPI_Isend(inner_envelope->south_west.data, send_count * dc->type, MPI_UINT8_T, dest_source,
baptiste.coudray's avatar
baptiste.coudray committed
                  SOUTH_WEST_CELLS_TAG, MPI_COMM_WORLD, &requests[i_request++]);

        /* Neighbour send north-east cell, which correspond to south-west cell */
baptiste.coudray's avatar
baptiste.coudray committed
        int dimensions[2] = {thickness, thickness};
        int start_y = dc->chunks_info[dest_source].y;
        int start_x = dc->chunks_info[dest_source].x + dc->chunk_info[dest_source].dimensions[1] - thickness;
        chunk_info_init(&outer_envelope.south_west, dc->type, dimensions, start_y, start_x, true);
        MPI_Irecv(outer_envelope.south_west.data, outer_envelope.south_west.count * dc->type, MPI_UINT8_T, dest_source,
baptiste.coudray's avatar
baptiste.coudray committed
                  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;
        MPI_Isend(inner_envelope->north_west.data, send_count * dc->type, MPI_INT8_T, dest_source,
baptiste.coudray's avatar
baptiste.coudray committed
                  NORTH_WEST_CELLS_TAG, MPI_COMM_WORLD, &requests[i_request++]);

        /* Neighbour send south-east cell, which correspond to north-west cell */
baptiste.coudray's avatar
baptiste.coudray committed
        int dimensions[2] = {thickness, thickness};
        int start_y = dc->chunks_info[dest_source].y + dc->chunk_info[dest_source].dimensions[0] - thickness;
        int start_x = dc->chunks_info[dest_source].x + dc->chunk_info[dest_source].dimensions[1] - thickness;
        chunk_info_init(&outer_envelope.north_west, dc->type, dimensions, start_y, start_x, true);
        MPI_Irecv(outer_envelope.north_west.data, outer_envelope.north_west.count * dc->type, MPI_INT8_T, dest_source,
baptiste.coudray's avatar
baptiste.coudray committed
                  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) {
baptiste.coudray's avatar
baptiste.coudray committed
    chunk_info_free(dc->chunk_info);
    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) {
baptiste.coudray's avatar
baptiste.coudray committed
    chunk_info_free(&envelope->north);
    chunk_info_free(&envelope->north_east);
    chunk_info_free(&envelope->east);
    chunk_info_free(&envelope->south_east);
    chunk_info_free(&envelope->south_west);
    chunk_info_free(&envelope->west);
    chunk_info_free(&envelope->north_west);