Logo AND Algorithmique Numérique Distribuée

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