Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
remove few tests which may never finish, and change one that used too much stack...
[simgrid.git] / teshsuite / smpi / mpich-test / pt2pt / bsendtest.c
1 /*
2  * Program to test that the "no overtaking messages" semantics
3  * of point to point communications in MPI is satisfied. 
4  * A long message is sent using MPI_BSend and received using MPI_Recv,
5  * followed by lots of short messages sent the same way.
6  *
7  *                              Patrick Bridges
8  *                              bridges@mcs.anl.gov
9  *                              patrick@CS.MsState.Edu
10  */
11
12 #include <stdio.h>
13 /* Needed for malloc declaration */
14 #include <stdlib.h>
15 #include "test.h"
16 #include "mpi.h"
17
18 #define SIZE 10000
19
20 static int src  = 0;
21 static int dest = 1;
22
23 /* Which tests to perform (not yet implemented) */
24 /* static int Do_Buffer = 1; */
25 /* static int Do_Standard = 1; */
26
27 /* Prototypes for picky compilers */
28 void Generate_Data ( double *, int );
29 void Normal_Test_Recv ( double *, int );
30 void Buffered_Test_Send ( double *, int );
31 void Buffered_Test_Ibsend ( double *, int );
32 int Check_Data ( double *, int );
33 void Clear_Buffer ( double *, int );
34
35 void Generate_Data(buffer, buff_size)
36 double *buffer;
37 int buff_size;
38 {
39     int i;
40
41     for (i = 0; i < buff_size; i++)
42         buffer[i] = (double)i+1;
43 }
44
45 void Normal_Test_Recv(buffer, buff_size)
46 double *buffer;
47 int buff_size;
48 {
49     int i, j;
50     MPI_Status Stat;
51     double     *b;
52
53     b = buffer;
54     for (j = 0; j < 2; j++) {
55         /* Receive a long message */
56         MPI_Recv(b, (buff_size/2 - 10), MPI_DOUBLE, src, 
57                  2000, MPI_COMM_WORLD, &Stat);
58         b += buff_size/2 - 10;
59         /* Followed by 10 short ones */
60         for (i = 0; i < 10; i++) {
61             MPI_Recv(b++, 1, MPI_DOUBLE, src, 2000, MPI_COMM_WORLD, &Stat);
62         }
63     }
64 }
65
66 void Buffered_Test_Send(buffer, buff_size)
67 double *buffer;
68 int buff_size;
69 {
70     int i, j;
71     void *bbuffer;
72     int  size;
73
74     for (j = 0; j < 2; j++) {
75         /* send a long message */
76         MPI_Bsend(buffer, (buff_size/2 - 10), MPI_DOUBLE, dest, 2000, 
77                  MPI_COMM_WORLD);
78         buffer += buff_size/2 - 10;
79         /* Followed by 10 short ones */
80         for (i = 0; i < 10; i++)
81             MPI_Bsend(buffer++, 1, MPI_DOUBLE, 
82                       dest, 2000, MPI_COMM_WORLD);
83         /* Force this set of Bsends to complete */
84         MPI_Buffer_detach( &bbuffer, &size );
85         MPI_Buffer_attach( bbuffer, size );
86     }
87 }
88
89 void Buffered_Test_Ibsend(buffer, buff_size)
90 double *buffer;
91 int buff_size;
92 {
93     int i, j;
94     void *bbuffer;
95     int  size;
96     int  cnt;
97     MPI_Request req[20];
98     MPI_Status  statuses[20];
99
100     for (j = 0; j < 2; j++) {
101         /* send a long message */
102         cnt = 0;
103         MPI_Ibsend(buffer, (buff_size/2 - 10), MPI_DOUBLE, dest, 2000, 
104                  MPI_COMM_WORLD, &req[cnt++]);
105         buffer += buff_size/2 - 10;
106         /* Followed by 10 short ones */
107         for (i = 0; i < 10; i++)
108             MPI_Ibsend(buffer++, 1, MPI_DOUBLE, 
109                       dest, 2000, MPI_COMM_WORLD, &req[cnt++]);
110         /* Wait for these to finish (should finish immediately) */
111         MPI_Waitall( cnt, req, statuses );
112
113         /* Force this set of Bsends to complete; this may take longer than
114            the Waitall */
115         MPI_Buffer_detach( &bbuffer, &size );
116         MPI_Buffer_attach( bbuffer, size );
117     }
118 }
119
120 int Check_Data(buffer, buff_size)
121 double *buffer;
122 int buff_size;
123 {
124     int i;
125     int err = 0;
126
127     for (i = 0; i < buff_size; i++)
128         if (buffer[i] != (i + 1)) {
129             err++;
130             fprintf( stderr, "Value at %d is %f, should be %f\n", i, 
131                     buffer[i], (double)(i+1) );
132             if (err > 10) return 1;
133             }
134     return err;
135 }
136
137 void Clear_Buffer(buffer, buff_size)
138 double *buffer;
139 int buff_size;
140 {
141     int i;
142     for (i = 0; i < buff_size; i++)
143         buffer[i] = -1;
144 }
145
146
147 int main(int argc, char **argv)
148 {
149     int rank; /* My Rank (0 or 1) */
150     double buffer[SIZE], *tmpbuffer, *tmpbuf;
151     int tsize, bsize;
152     char *Current_Test = NULL;
153
154     MPI_Init(&argc, &argv);
155     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
156
157     if (rank == src) { 
158         Generate_Data(buffer, SIZE);
159         MPI_Pack_size( SIZE, MPI_DOUBLE, MPI_COMM_WORLD, &bsize );
160         tmpbuffer = (double *) malloc( bsize + 22*MPI_BSEND_OVERHEAD );
161         if (!tmpbuffer) {
162             fprintf( stderr, "Could not allocate bsend buffer of size %d\n",
163                      bsize );
164             MPI_Abort( MPI_COMM_WORLD, 1 );
165             }
166         MPI_Buffer_attach( tmpbuffer, bsize + 22*MPI_BSEND_OVERHEAD );
167         Buffered_Test_Send(buffer, SIZE);
168         Buffered_Test_Ibsend(buffer, SIZE);
169         MPI_Buffer_detach( &tmpbuf, &tsize );
170         Test_Waitforall( );
171         MPI_Finalize();
172
173     } else if (rank == dest) {
174         Test_Init("bsendtest", rank);
175         /* Test 3 */
176         Current_Test = (char*)"Overtaking Test (Buffered Send -> Normal Receive)";
177         Clear_Buffer(buffer, SIZE);
178         /* For Bsend */
179         Normal_Test_Recv(buffer, SIZE);
180         if (Check_Data(buffer, SIZE))
181             Test_Failed(Current_Test);
182         else
183             Test_Passed(Current_Test);
184
185         /* For Ibsend */
186         Current_Test = (char*)"Overtaking Test (Buffered Isend -> Normal Receive)";
187         Clear_Buffer(buffer, SIZE);
188         Normal_Test_Recv(buffer, SIZE);
189         if (Check_Data(buffer, SIZE))
190             Test_Failed(Current_Test);
191         else
192             Test_Passed(Current_Test);
193
194         Test_Waitforall( );
195         {
196             int rval = Summarize_Test_Results(); /* Returns number of tests;
197                                                     that failed */
198             Test_Finalize();
199             MPI_Finalize();
200             return rval;
201         }
202     } else {
203         fprintf(stderr, "*** This program uses exactly 2 processes! ***\n");
204         MPI_Abort( MPI_COMM_WORLD, 1 );
205     }
206
207     return 0;
208 }
209
210
211