mpi line with matrix exchange

This commit is contained in:
Fabio Salvini 2016-11-13 13:11:24 +01:00
parent df121930ad
commit 14541623a4
3 changed files with 73 additions and 28 deletions

View File

@ -3,14 +3,16 @@
*/
#include <stdio.h>
#include <string.h>
#include <math.h>
#include <mpi.h>
#include "../config/config.h"
#include "../utils/utils.h"
#define TAG_BORDER 0
#define TAG_MATRIX 1
double *compute_jacobi(int n, double init_value, double threshold, borders b, int *rows, int *iterations);
double *compute_jacobi(int n, double init_value, double threshold, borders b, int *iterations);
int rank;
int numprocs;
@ -25,7 +27,6 @@ int main(int argc, char* argv[]) {
double *x;
double startwtime = 0.0, endwtime;
int iterations;
int rows;
MPI_Init(&argc, &argv);
@ -62,12 +63,13 @@ int main(int argc, char* argv[]) {
startwtime = MPI_Wtime();
}
x = compute_jacobi(n, init_value, threshold, b, &rows, &iterations);
x = compute_jacobi(n, init_value, threshold, b, &iterations);
if (rank == 0) {
endwtime = MPI_Wtime();
printf("Wall clock time: %fs\n", endwtime - startwtime);
printf("Iterations: %d\n", iterations);
/* print_sa_matrix(x, n + 2, n + 2); */
}
destroy_sa_matrix(x);
@ -77,28 +79,31 @@ int main(int argc, char* argv[]) {
return 0;
}
double *compute_jacobi(int n, double init_value, double threshold, borders b, int *rows, int *iterations) {
double *compute_jacobi(int n, double init_value, double threshold, borders b, int *iterations) {
double *complete_x;
double *x;
double max_diff, global_max_diff, new_x;
int i, j;
int rows, rows_to_transmit;
int receive_pos;
MPI_Status status;
if (rank == 0) {
*rows = n - (n / numprocs) * (numprocs - 1);
rows = n - (n / numprocs) * (numprocs - 1);
} else {
*rows = n / numprocs;
rows = n / numprocs;
}
LOG(printf("[Process %d/%d] rows: %d\n", rank, numprocs, *rows));
LOG(printf("[Process %d/%d] rows: %d\n", rank, numprocs, rows));
/* LOG(printf("[Process %d/%d] initializing matrix\n", rank, numprocs)); */
/* Initialize the matrix */
x = create_sa_matrix(*rows + 2, n + 2);
for (i = 0; i < *rows + 2; i++) {
x = create_sa_matrix(rows + 2, n + 2);
for (i = 0; i < rows + 2; i++) {
for (j = 1; j <= n; j++) {
x[sa_index(n + 2, i, j)] = init_value;
}
}
/* Initialize boundary regions */
for (i = 0; i < *rows + 2; i++) {
for (i = 0; i < rows + 2; i++) {
x[sa_index(n + 2, i, 0)] = b.west;
x[sa_index(n + 2, i, n + 1)] = b.east;
}
@ -109,7 +114,7 @@ double *compute_jacobi(int n, double init_value, double threshold, borders b, in
}
if (rank == numprocs - 1){
for (i = 1; i < n + 1; i++) {
x[sa_index(n + 2, *rows + 1, i)] = b.south;
x[sa_index(n + 2, rows + 1, i)] = b.south;
}
}
/* LOG(printf("[Process %d/%d] matrix initialized\n", rank, numprocs)); */
@ -118,7 +123,7 @@ double *compute_jacobi(int n, double init_value, double threshold, borders b, in
do {
max_diff = 0;
global_max_diff = 0;
for (i = 1; i <= *rows; i++) {
for (i = 1; i <= rows; i++) {
for (j = 1; j <= n; j++) {
new_x = 0.25 * (x[sa_index(n + 2, i - 1, j)] + x[sa_index(n + 2, i, j + 1)] + x[sa_index(n + 2, i + 1, j)] + x[sa_index(n + 2, i, j - 1)]);
max_diff = (double) fmax(max_diff, fabs(new_x - x[sa_index(n + 2, i, j)]));
@ -128,8 +133,8 @@ double *compute_jacobi(int n, double init_value, double threshold, borders b, in
if (rank % 2 == 0) {
if (rank != numprocs - 1) {
// Send and receive south border
MPI_Send(&x[sa_index(n + 2, *rows, 0)], n + 2, MPI_DOUBLE, rank + 1, TAG_BORDER, MPI_COMM_WORLD);
MPI_Recv(&x[sa_index(n + 2, *rows + 1, 0)], n + 2, MPI_DOUBLE, rank + 1, TAG_BORDER, MPI_COMM_WORLD, &status);
MPI_Send(&x[sa_index(n + 2, rows, 0)], n + 2, MPI_DOUBLE, rank + 1, TAG_BORDER, MPI_COMM_WORLD);
MPI_Recv(&x[sa_index(n + 2, rows + 1, 0)], n + 2, MPI_DOUBLE, rank + 1, TAG_BORDER, MPI_COMM_WORLD, &status);
}
if (rank != 0) {
// Send and receive north border
@ -142,8 +147,8 @@ double *compute_jacobi(int n, double init_value, double threshold, borders b, in
MPI_Send(&x[sa_index(n + 2, 1, 0)], n + 2, MPI_DOUBLE, rank - 1, TAG_BORDER, MPI_COMM_WORLD);
if (rank != numprocs - 1) {
// Receive and send south border
MPI_Recv(&x[sa_index(n + 2, *rows + 1, 0)], n + 2, MPI_DOUBLE, rank + 1, TAG_BORDER, MPI_COMM_WORLD, &status);
MPI_Send(&x[sa_index(n + 2, *rows, 0)], n + 2, MPI_DOUBLE, rank + 1, TAG_BORDER, MPI_COMM_WORLD);
MPI_Recv(&x[sa_index(n + 2, rows + 1, 0)], n + 2, MPI_DOUBLE, rank + 1, TAG_BORDER, MPI_COMM_WORLD, &status);
MPI_Send(&x[sa_index(n + 2, rows, 0)], n + 2, MPI_DOUBLE, rank + 1, TAG_BORDER, MPI_COMM_WORLD);
}
}
/* LOG(printf("[Process %d/%d] max_diff: %f\n", rank, numprocs, max_diff)); */
@ -151,5 +156,27 @@ double *compute_jacobi(int n, double init_value, double threshold, borders b, in
/* LOG(printf("[Process %d/%d] global_max_diff: %f\n", rank, numprocs, global_max_diff)); */
(*iterations)++;
} while (global_max_diff > threshold);
return x;
if (rank == 0) {
complete_x = create_sa_matrix(n + 2, n + 2);
memcpy(complete_x, x, (rows + 1) * (n + 2) * sizeof(double)); /* rows + 2 if rank == numprocs - 1 */
rows_to_transmit = n / numprocs;
receive_pos = rows + 1;
for (i = 1; i < numprocs; i++) {
if (i == numprocs - 1) {
rows_to_transmit++;
}
MPI_Recv(&complete_x[sa_index(n + 2, receive_pos, 0)], rows_to_transmit * (n + 2), MPI_DOUBLE, i, TAG_MATRIX, MPI_COMM_WORLD, &status);
receive_pos += n / numprocs;
}
} else {
complete_x = NULL;
rows_to_transmit = rows;
if (rank == numprocs - 1) {
rows_to_transmit++;
}
MPI_Send(&x[sa_index(n + 2, 1, 0)], rows_to_transmit * (n + 2), MPI_DOUBLE, 0, TAG_MATRIX, MPI_COMM_WORLD);
}
return complete_x;
}

View File

@ -1,9 +1,6 @@
#include <stdio.h>
#include <stdlib.h>
/*
* Create a matrix stored in a single array.
*/
double *create_sa_matrix(int rows, int cols) {
double *x;
@ -11,21 +8,26 @@ double *create_sa_matrix(int rows, int cols) {
return x;
}
/*
* Destroy a single array matrix.
*/
void destroy_sa_matrix(double *x) {
free(x);
}
/*
* Get the index for the single array matrix
* that correspond to the given row and column.
*/
int sa_index(int cols, int r, int c) {
return r * cols + c;
}
void print_sa_matrix(double *x, int rows, int cols) {
int i, j;
for (i = 0; i < rows; i++) {
for (j = 0; j < cols; j++) {
printf("%f\t", x[sa_index(cols, i, j)]);
}
printf("\n");
}
fflush(stdout);
}
double **create_matrix(int rows, int cols) {
int i;
double **x;

View File

@ -17,11 +17,27 @@ typedef struct borders {
} borders;
/*
* Create a matrix stored in a single array.
*/
double *create_sa_matrix(int rows, int cols);
/*
* Destroy a single array matrix.
*/
void destroy_sa_matrix(double *x);
/*
* Get the index for the single array matrix
* that correspond to the given row and column.
*/
int sa_index(int cols, int r, int c);
/*
* Print a single array matrix.
*/
void print_sa_matrix(double *x, int rows, int cols);
double **create_matrix(int rows, int cols);
void destroy_matrix(double **x, int rows);
void print_matrix(double **x, int rows, int cols);