| 1 | /* matmat_mw.c : parallel multiplication of two matrices.
|
|---|
| 2 | * To execute: mpicc matmat_mw.c ; mpiexec -n 4 ./a.out N L M
|
|---|
| 3 | * Or replace "4" with however many procs you want to use.
|
|---|
| 4 | * Arguments N L M should be replaced with any integer numbers which is
|
|---|
| 5 | * no larger than the corresponding dimension decided in the "data" file.
|
|---|
| 6 | * To verify: civl verify matmat_mw.c
|
|---|
| 7 | */
|
|---|
| 8 | #include <stdlib.h>
|
|---|
| 9 | #include <stdio.h>
|
|---|
| 10 | #include <string.h>
|
|---|
| 11 | #include <mpi.h>
|
|---|
| 12 | #include <omp.h>
|
|---|
| 13 |
|
|---|
| 14 | #define comm MPI_COMM_WORLD
|
|---|
| 15 | #ifdef _CIVL
|
|---|
| 16 | #include <civlc.cvh>
|
|---|
| 17 | $input int _mpi_nprocs_lo = 2;
|
|---|
| 18 | $input int _mpi_nprocs_hi = 4;
|
|---|
| 19 | /* Dimensions of 2 matrices: a[N][L] * b[L][M] */
|
|---|
| 20 | $input int NB = 3; // upper bound of N
|
|---|
| 21 | $input int N;
|
|---|
| 22 | $assume(0 < N && N <= NB);
|
|---|
| 23 | $input int LB = 3; // upper bound of L
|
|---|
| 24 | $input int L;
|
|---|
| 25 | $assume(0 < L && L <= LB);
|
|---|
| 26 | $input int MB = 3; // upper bound of M
|
|---|
| 27 | $input int M;
|
|---|
| 28 | $assume(0 < M && M <= MB);
|
|---|
| 29 | $input double a[N][L]; // input data for matrix a
|
|---|
| 30 | $input double b[L][M]; // input data for matrix b
|
|---|
| 31 | double oracle[N][M]; // matrix stores results of a sequential run
|
|---|
| 32 | #else
|
|---|
| 33 | FILE * fp; // pointer to the data file which gives two matrices
|
|---|
| 34 | int N, L, M;
|
|---|
| 35 | #endif
|
|---|
| 36 |
|
|---|
| 37 | /* prints a matrix. In CIVL mode, it will compare the matrix with the
|
|---|
| 38 | result of the sequential run.*/
|
|---|
| 39 | void printMatrix(int numRows, int numCols, double *m) {
|
|---|
| 40 | int i, j;
|
|---|
| 41 |
|
|---|
| 42 | for (i = 0; i < numRows; i++) {
|
|---|
| 43 | for (j = 0; j < numCols; j++) {
|
|---|
| 44 | printf("%f ", m[i*numCols + j]);
|
|---|
| 45 | #ifdef _CIVL
|
|---|
| 46 | $assert(m[i*numCols + j] == oracle[i][j], "The calculated value at position [%d][%d] is %f"
|
|---|
| 47 | " but the expected one is %f", i, j, m[i*numCols+j], oracle[i][j]);
|
|---|
| 48 | #endif
|
|---|
| 49 | }
|
|---|
| 50 | printf("\n");
|
|---|
| 51 | }
|
|---|
| 52 | printf("\n");
|
|---|
| 53 | }
|
|---|
| 54 |
|
|---|
| 55 | /* Computes a vetor with length L times a matrix with dimensions [L][M] */
|
|---|
| 56 | void vecmat(double vector[L], double matrix[L][M], double result[M]) {
|
|---|
| 57 | int j, k;
|
|---|
| 58 | #pragma omp parallel for
|
|---|
| 59 | for (j = 0; j < M; j++)
|
|---|
| 60 | for (k = 0, result[j] = 0.0; k < L; k++)
|
|---|
| 61 | result[j] += vector[k]*matrix[k][j];
|
|---|
| 62 | }
|
|---|
| 63 |
|
|---|
| 64 | int main(int argc, char *argv[]) {
|
|---|
| 65 | int rank, nprocs, i, j;
|
|---|
| 66 | MPI_Status status;
|
|---|
| 67 |
|
|---|
| 68 | #ifndef _CIVL
|
|---|
| 69 | N = atoi(argv[1]);
|
|---|
| 70 | L = atoi(argv[2]);
|
|---|
| 71 | M = atoi(argv[3]);
|
|---|
| 72 | #endif
|
|---|
| 73 | MPI_Init(&argc, &argv);
|
|---|
| 74 | MPI_Comm_rank(comm, &rank);
|
|---|
| 75 | MPI_Comm_size(comm, &nprocs);
|
|---|
| 76 | if (rank == 0) {
|
|---|
| 77 | double c[N][M], tmp[M];
|
|---|
| 78 | int count;
|
|---|
| 79 | #ifndef _CIVL
|
|---|
| 80 | double a[N][L], b[L][M];
|
|---|
| 81 |
|
|---|
| 82 | fp = fopen("data", "r");
|
|---|
| 83 | for (i = 0; i < N; i++)
|
|---|
| 84 | for (j = 0; j < L; j++)
|
|---|
| 85 | fscanf(fp,"%lf", &a[i][j]);
|
|---|
| 86 | for (i = 0; i < L; i++)
|
|---|
| 87 | for (j = 0; j < M; j++)
|
|---|
| 88 | fscanf(fp,"%lf",&b[i][j]);
|
|---|
| 89 | #else
|
|---|
| 90 | $elaborate(N);
|
|---|
| 91 | $elaborate(M);
|
|---|
| 92 | $elaborate(L);
|
|---|
| 93 | for(int i=0; i < N; i++) {
|
|---|
| 94 | vecmat(a[i], b, &oracle[i][0]);
|
|---|
| 95 | }
|
|---|
| 96 |
|
|---|
| 97 | #endif
|
|---|
| 98 | MPI_Bcast(b, L*M, MPI_DOUBLE, 0, comm);
|
|---|
| 99 | for (count = 0; count < nprocs-1 && count < N; count++)
|
|---|
| 100 | MPI_Send(&a[count][0], L, MPI_DOUBLE, count+1, count+1, comm);
|
|---|
| 101 | for (i = 0; i < N; i++) {
|
|---|
| 102 | MPI_Recv(tmp, M, MPI_DOUBLE, MPI_ANY_SOURCE, MPI_ANY_TAG, comm, &status);
|
|---|
| 103 | for (j = 0; j < M; j++) c[status.MPI_TAG-1][j] = tmp[j];
|
|---|
| 104 | if (count < N) {
|
|---|
| 105 | MPI_Send(&a[count][0], L, MPI_DOUBLE, status.MPI_SOURCE, count+1, comm);
|
|---|
| 106 | count++;
|
|---|
| 107 | }
|
|---|
| 108 | }
|
|---|
| 109 | for (i = 1; i < nprocs; i++) MPI_Send(NULL, 0, MPI_INT, i, 0, comm);
|
|---|
| 110 | printMatrix(N, M, &c[0][0]);
|
|---|
| 111 | #ifndef _CIVL
|
|---|
| 112 | fclose(fp);
|
|---|
| 113 | #endif
|
|---|
| 114 | } else {
|
|---|
| 115 | double b[L][M], in[L], out[M];
|
|---|
| 116 |
|
|---|
| 117 | MPI_Bcast(b, L*M, MPI_DOUBLE, 0, comm);
|
|---|
| 118 | while (1) {
|
|---|
| 119 | MPI_Recv(in, L, MPI_DOUBLE, 0, MPI_ANY_TAG, comm, &status);
|
|---|
| 120 | if (status.MPI_TAG == 0) break;
|
|---|
| 121 | vecmat(in, b, out);
|
|---|
| 122 | MPI_Send(out, M, MPI_DOUBLE, 0, status.MPI_TAG, comm);
|
|---|
| 123 | }
|
|---|
| 124 | }
|
|---|
| 125 | MPI_Finalize();
|
|---|
| 126 | return 0;
|
|---|
| 127 | }
|
|---|
| 128 |
|
|---|