Logo AND Algorithmique Numérique Distribuée

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