Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
617bfaa190902c9b94743a736e85d34f7d5f76c1
[simgrid.git] / examples / smpi / MM / Summa.c
1 /*!
2  * Classical Block Matrix Multiplication example
3  *
4  * Authors: Quintin Jean-Noël
5  */
6 #include "Matrix_init.h"
7 #include "Summa.h"
8 #include "timer.h"
9 #include "xbt/log.h"
10  XBT_LOG_NEW_DEFAULT_CATEGORY(MM_Summa,
11                              "Messages specific for this msg example");
12
13 inline double Summa(
14                      double *a, double *b, double *c,
15                      size_t lda, size_t ldb, size_t ldc,
16                      size_t m, size_t k_a, size_t k_b, size_t n,
17                      size_t Block_size, size_t start, size_t end,
18                      size_t row, size_t col, size_t size_row, size_t size_col,
19                      double *a_local, double *b_local,
20                      MPI_Datatype Block_a, MPI_Datatype Block_a_local,
21                      MPI_Datatype Block_b,
22                      MPI_Comm row_comm, MPI_Comm col_comm, int subs)
23 {
24   double *B_a     , *B_b     ; //matrix blocks
25   size_t err;
26   double alpha = 1, beta = 1;  //C := alpha * a * b + beta * c
27   size_t B_proc_col, B_proc_row; // Number of bloc(row or col) on one processor
28   B_proc_col =  k_b / Block_size;  // Number of block on one processor
29   B_proc_row = k_a / Block_size; // Number of block on one processor
30
31   //size_t lda = k_a, ldb = n, ldc = n;
32   size_t lda_local = lda;
33   size_t ldb_local = ldb;
34
35
36   double time, computation_time = 0, communication_time = 0;
37   struct timespec start_time, end_time; //time mesure
38   struct timespec start_time_intern, end_time_intern; //time mesure
39
40
41
42
43   get_time(&start_time);
44
45   /*-------------Distributed Matrix Multiplication algorithm-----------------*/
46   size_t iter;
47   for( iter = start; iter < end; iter++ ){
48     size_t pivot_row, pivot_col, pos_a, pos_b;
49 #ifdef CYCLIC
50     // pivot row on processor layer
51     pivot_row = (iter % size_col);
52     pivot_col = (iter % size_row);
53     //position of the block
54     if(subs == 1){
55       pos_a = (size_t)((iter - start) / size_row) * Block_size;
56       pos_b = (size_t)((iter - start) / size_col) * ldb * Block_size;
57     }else{
58       pos_a = (size_t)(iter / size_row) * Block_size;
59       pos_b = (size_t)(iter / size_col) * ldb * Block_size;
60     }
61 #else
62     // pivot row on processor layer
63     pivot_row = (size_t)(iter / B_proc_col) % size_col;
64     pivot_col = (size_t)(iter / B_proc_row) % size_row;
65     //position of the block
66     if(subs == 1){
67       pos_a = ((iter - start) % B_proc_row) * Block_size;
68       pos_b = ((iter - start) % B_proc_col) * ldb * Block_size;
69     }else{
70       pos_a = (iter % B_proc_row) * Block_size;
71       pos_b = (iter % B_proc_col) * ldb * Block_size;
72     }
73 #endif
74     XBT_DEBUG( "pivot: %zu, iter: %zu, B_proc_col: %zu, "
75                 "size_col:%zu, size_row: %zu\n",
76                 pivot_row, iter, B_proc_row,size_col,size_row);
77     MPI_Barrier(row_comm);
78     MPI_Barrier(col_comm);
79
80     get_time(&start_time_intern);
81     //Broadcast the row
82     if(size_row > 1){
83       MPI_Datatype * Block;
84       if( pivot_col != col ){
85         B_a = a_local;
86         lda_local = Block_size;
87         XBT_DEBUG("recieve B_a %zu,%zu \n",m , Block_size);
88         Block = &Block_a_local;
89       }else{
90         B_a = a + pos_a;
91         lda_local = lda;
92         XBT_DEBUG("sent B_a %zu,%zu \n",m , Block_size);
93         Block = &Block_a;
94       }
95       err = MPI_Bcast(B_a, 1, *Block, pivot_col, row_comm);
96       if (err != MPI_SUCCESS) {
97         perror("Error Bcast A\n");
98         return -1;
99       }
100     }else{
101       B_a = a + pos_a;
102       XBT_DEBUG("position of B_a: %zu \n", pos_a);
103     }
104
105     //Broadcast the col
106     if(size_col > 1){
107       if( pivot_row == row ){
108         B_b = b + pos_b;
109         XBT_DEBUG("sent B_b Block_size: %zu, pos:%zu \n",
110                     ldb, pos_b);
111       }else{
112         B_b = b_local;
113         XBT_DEBUG("recieve B_b %zu,%zu \n", Block_size,n);
114       }
115       err = MPI_Bcast(B_b, 1, Block_b, pivot_row, col_comm );
116       if (err != MPI_SUCCESS) {
117         perror("Error Bcast B\n");
118         MPI_Finalize();
119         exit(-1);
120       }
121     }else{
122       B_b = b + pos_b;
123       XBT_DEBUG("position of B_b: %zu \n", pos_b);
124     }
125     get_time(&end_time_intern);
126     communication_time += get_timediff(&start_time_intern,&end_time_intern);
127
128     MPI_Barrier(row_comm);
129     MPI_Barrier(col_comm);
130     get_time(&start_time_intern);
131     XBT_DEBUG("execute Gemm number: %zu\n", iter);
132     //We have recieved a line of block and a colomn
133    //              cblas_dgemm( CblasRowMajor, CblasNoTrans, CblasNoTrans,
134    //               m, n, Block_size, alpha, B_a, lda_local, B_b, ldb_local,
135    //               beta, c, ldc );
136     int i, j, k;
137     for(i = 0; i < m; i++)
138       for(j = 0; j < n; j++)
139         for(k = 0; k < Block_size; k++)
140           c[i*ldc+j] += B_a[j*lda_local+k]*B_b[k*ldb_local+j];
141
142     get_time(&end_time_intern);
143     computation_time += get_timediff(&start_time_intern,&end_time_intern);
144
145   }
146   MPI_Barrier(row_comm);
147   MPI_Barrier(col_comm);
148
149   get_time(&end_time);
150   time = get_timediff(&start_time,&end_time);
151   printf("communication time: %le nanoseconds, "
152          "computation time: %le nanoseconds\n",
153          communication_time, computation_time);
154
155
156   return time;
157 }