/* * Copyright (C) 2022 NHR@FAU, University Erlangen-Nuremberg. * All rights reserved. This file is part of nusif-solver. * Use of this source code is governed by a MIT style * license that can be found in the LICENSE file. */ #include #include #include #include "comm.h" // subroutines local to this module static int sizeOfRank(int rank, int size, int N) { return N / size + ((N % size > rank) ? 1 : 0); } static void setupCommunication(Comm* c, int direction, int layer) { MPI_Datatype type; size_t dblsize = sizeof(double); int imaxLocal = c->imaxLocal; int jmaxLocal = c->jmaxLocal; int sizes[NDIMS]; int subSizes[NDIMS]; int starts[NDIMS]; int offset = 0; sizes[IDIM] = imaxLocal + 2; sizes[JDIM] = jmaxLocal + 2; if (layer == HALO) { offset = 1; } switch (direction) { case LEFT: subSizes[IDIM] = 1; subSizes[JDIM] = jmaxLocal; starts[IDIM] = 1 - offset; starts[JDIM] = 1; break; case RIGHT: subSizes[IDIM] = 1; subSizes[JDIM] = jmaxLocal; starts[IDIM] = imaxLocal + offset; starts[JDIM] = 1; break; case BOTTOM: subSizes[IDIM] = imaxLocal; subSizes[JDIM] = 1; starts[IDIM] = 1; starts[JDIM] = 1 - offset; break; case TOP: subSizes[IDIM] = imaxLocal; subSizes[JDIM] = 1; starts[IDIM] = 1; starts[JDIM] = jmaxLocal + offset; break; } MPI_Type_create_subarray(NDIMS, sizes, subSizes, starts, MPI_ORDER_C, MPI_DOUBLE, &type); MPI_Type_commit(&type); if (layer == HALO) { c->rbufferTypes[direction] = type; } else if (layer == BULK) { c->sbufferTypes[direction] = type; } } static void assembleResult(Comm* c, double* src, double* dst, int imaxLocal[], int jmaxLocal[], int offset[], int jmax, int imax) { MPI_Request* requests; int numRequests = 1; if (c->rank == 0) { numRequests = c->size + 1; } else { numRequests = 1; } requests = (MPI_Request*)malloc(numRequests * sizeof(MPI_Request)); /* all ranks send their bulk array */ MPI_Datatype bulkType; int oldSizes[NDIMS] = { c->jmaxLocal + 2, c->imaxLocal + 2 }; int newSizes[NDIMS] = { c->jmaxLocal, c->imaxLocal }; int starts[NDIMS] = { 1, 1 }; MPI_Type_create_subarray(NDIMS, oldSizes, newSizes, starts, MPI_ORDER_C, MPI_DOUBLE, &bulkType); MPI_Type_commit(&bulkType); MPI_Isend(src, 1, bulkType, 0, 0, c->comm, &requests[0]); /* rank 0 assembles the subdomains */ if (c->rank == 0) { for (int i = 0; i < c->size; i++) { MPI_Datatype domainType; int oldSizes[NDIMS] = { jmax, imax }; int newSizes[NDIMS] = { jmaxLocal[i], imaxLocal[i] }; int starts[NDIMS] = { offset[i * NDIMS + JDIM], offset[i * NDIMS + IDIM] }; MPI_Type_create_subarray(NDIMS, oldSizes, newSizes, starts, MPI_ORDER_C, MPI_DOUBLE, &domainType); MPI_Type_commit(&domainType); MPI_Irecv(dst, 1, domainType, i, 0, c->comm, &requests[i + 1]); } } MPI_Waitall(numRequests, requests, MPI_STATUSES_IGNORE); } static int sum(int* sizes, int position) { int sum = 0; for (int i = 0; i < position; i++) { sum += sizes[i]; } return sum; } // exported subroutines void commReduction(double* v, int op) { if (op == MAX) { MPI_Allreduce(MPI_IN_PLACE, v, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD); } else if (op == SUM) { MPI_Allreduce(MPI_IN_PLACE, v, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); } } int commIsBoundary(Comm* c, int direction) { switch (direction) { case LEFT: return c->coords[IDIM] == 0; break; case RIGHT: return c->coords[IDIM] == (c->dims[IDIM] - 1); break; case BOTTOM: return c->coords[JDIM] == 0; break; case TOP: return c->coords[JDIM] == (c->dims[JDIM] - 1); break; } return 0; } void commExchange(Comm* c, double* grid) { int counts[NDIRS] = { 1, 1, 1, 1 }; MPI_Aint displs[NDIRS] = { 0, 0, 0, 0 }; MPI_Neighbor_alltoallw(grid, counts, displs, c->sbufferTypes, grid, counts, displs, c->rbufferTypes, c->comm); } void commShift(Comm* c, double* f, double* g) { MPI_Request requests[4] = { MPI_REQUEST_NULL, MPI_REQUEST_NULL, MPI_REQUEST_NULL, MPI_REQUEST_NULL }; /* shift G */ /* receive ghost cells from bottom neighbor */ MPI_Irecv(g, 1, c->rbufferTypes[BOTTOM], c->neighbours[BOTTOM], 0, c->comm, &requests[0]); /* send ghost cells to top neighbor */ MPI_Isend(g, 1, c->sbufferTypes[TOP], c->neighbours[TOP], 0, c->comm, &requests[1]); /* shift F */ /* receive ghost cells from left neighbor */ MPI_Irecv(f, 1, c->rbufferTypes[LEFT], c->neighbours[LEFT], 1, c->comm, &requests[2]); /* send ghost cells to right neighbor */ MPI_Isend(f, 1, c->sbufferTypes[RIGHT], c->neighbours[RIGHT], 1, c->comm, &requests[3]); MPI_Waitall(4, requests, MPI_STATUSES_IGNORE); } void commCollectResult(Comm* c, double* ug, double* vg, double* pg, double* u, double* v, double* p, int jmax, int imax) { int offset[c->size * NDIMS]; int imaxLocal[c->size]; int jmaxLocal[c->size]; MPI_Gather(&c->imaxLocal, 1, MPI_INT, imaxLocal, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Gather(&c->jmaxLocal, 1, MPI_INT, jmaxLocal, 1, MPI_INT, 0, MPI_COMM_WORLD); if (c->rank == 0) { for (int i = 0; i < c->size; i++) { int coords[NDIMS]; MPI_Cart_coords(c->comm, i, NDIMS, coords); offset[i * NDIMS + IDIM] = sum(imaxLocal, coords[IDIM]); offset[i * NDIMS + JDIM] = sum(jmaxLocal, coords[JDIM]); printf("Rank: %d, Coords(j,i): %d %d, Size(j,i): %d %d " "Offset(j,i): %d %d\n", i, coords[JDIM], coords[IDIM], jmaxLocal[i], imaxLocal[i], offset[i * NDIMS + JDIM], offset[i * NDIMS + IDIM]); } } /* collect P */ assembleResult(c, p, pg, imaxLocal, jmaxLocal, offset, jmax, imax); /* collect U */ assembleResult(c, u, ug, imaxLocal, jmaxLocal, offset, jmax, imax); /* collect V */ assembleResult(c, v, vg, imaxLocal, jmaxLocal, offset, jmax, imax); } void commPrintConfig(Comm* c) { fflush(stdout); MPI_Barrier(MPI_COMM_WORLD); if (commIsMaster(c)) { printf("Communication setup:\n"); } for (int i = 0; i < c->size; i++) { if (i == c->rank) { printf("\tRank %d of %d\n", c->rank, c->size); printf("\tNeighbours (bottom, top, left, right): %d %d, %d, %d\n", c->neighbours[BOTTOM], c->neighbours[TOP], c->neighbours[LEFT], c->neighbours[RIGHT]); printf("\tCoordinates (j,i) %d %d\n", c->coords[JDIM], c->coords[IDIM]); printf("\tLocal domain size (j,i) %dx%d\n", c->jmaxLocal, c->imaxLocal); fflush(stdout); } } MPI_Barrier(MPI_COMM_WORLD); } void commInit(Comm* c, int jmax, int imax) { /* setup communication */ MPI_Comm_rank(MPI_COMM_WORLD, &(c->rank)); MPI_Comm_size(MPI_COMM_WORLD, &(c->size)); int dims[NDIMS] = { 0, 0 }; int periods[NDIMS] = { 0, 0 }; MPI_Dims_create(c->size, NDIMS, dims); MPI_Cart_create(MPI_COMM_WORLD, NDIMS, dims, periods, 0, &c->comm); MPI_Cart_shift(c->comm, IDIM, 1, &c->neighbours[LEFT], &c->neighbours[RIGHT]); MPI_Cart_shift(c->comm, JDIM, 1, &c->neighbours[BOTTOM], &c->neighbours[TOP]); MPI_Cart_get(c->comm, NDIMS, c->dims, periods, c->coords); c->imaxLocal = sizeOfRank(c->rank, dims[IDIM], imax); c->jmaxLocal = sizeOfRank(c->rank, dims[JDIM], jmax); // setup buffer types for communication setupCommunication(c, LEFT, BULK); setupCommunication(c, LEFT, HALO); setupCommunication(c, RIGHT, BULK); setupCommunication(c, RIGHT, HALO); setupCommunication(c, BOTTOM, BULK); setupCommunication(c, BOTTOM, HALO); setupCommunication(c, TOP, BULK); setupCommunication(c, TOP, HALO); }