Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
first commit to add the mpich-test suite to smpi tesh suite. Obviously all tests...
[simgrid.git] / teshsuite / smpi / mpich-test / pt2pt / nbtest.c
1 #include <stdio.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include "mpi.h"
5
6 #if defined(NEEDS_STDLIB_PROTOTYPES)
7 #include "protofix.h"
8 #endif
9
10 /*
11    Test to make sure that nonblocking routines actually work
12    In this example, we assume that we do not know the message
13    sizes ahead of time.
14
15    Just like nblock, but with the probe test.
16 */
17
18 int main( int argc, char **argv )
19 {
20     int count, tag, nsend, myid, np, rcnt, scnt, i, j, *send_buf;
21     int length, finished;
22     int baselen = 1;
23     int **recv_buf;
24     MPI_Status status, rtn_status;
25     MPI_Request *rsend, *rrecv;
26
27     MPI_Init( &argc, &argv );
28     MPI_Comm_rank( MPI_COMM_WORLD, &myid );
29     MPI_Comm_size( MPI_COMM_WORLD, &np );
30 /*
31   MPE_Errors_call_dbx_in_xterm( (argv)[0], (char *)0 ); 
32   MPE_Signals_call_debugger();
33   */
34     if (argc > 2 && argv[1] && strcmp( argv[1], "-first" ) == 0) 
35         baselen = atoi(argv[2]);
36
37 /* malloc buffers */
38     nsend = 3 * np;
39     rsend = (MPI_Request *) malloc ( nsend * sizeof(MPI_Request) );
40     rrecv = (MPI_Request *) malloc ( nsend * sizeof(MPI_Request) );
41     recv_buf = (int **) malloc ( nsend * sizeof(int *) );
42     if (!rsend || !rrecv) {
43         fprintf( stderr, "Failed to allocate space for requests\n" );
44         MPI_Abort( MPI_COMM_WORLD, 1 );
45     }
46
47     for (count = baselen; count < 10000; count *= 2) {
48         /* We'll send/recv from everyone */
49         scnt = 0;
50         rcnt = 0;
51
52         /* do sends */
53         send_buf   = (int *)malloc( count * sizeof(int) );
54         for (j=0; j<3; j++) {
55             tag = j;
56             for (i=0; i<np; i++) {
57                 if (i != myid) 
58                     MPI_Isend( send_buf, count, MPI_INT, i, tag, 
59                                MPI_COMM_WORLD, &rsend[scnt++] );
60             }
61             /* Check sends, one could free memory here if they are done */
62             for (i=0; i<scnt; i++) {
63                 MPI_Test( &rsend[i], &finished, &status );
64             }
65         }
66
67         /* do recvs */
68         for (j=0; j<3; j++) {
69             tag = j;
70             for (i=0; i<np; i++) {
71                 if (i != myid)  {
72                     MPI_Probe(MPI_ANY_SOURCE,tag,MPI_COMM_WORLD,&status);
73                     MPI_Get_count(&status,MPI_INT,&length); 
74                     /* printf("[%d] length = %d\n",myid,length); 
75                        fflush(stdout); */
76                     recv_buf[rcnt] = (int *)malloc(length * sizeof(int));
77                     MPI_Recv(recv_buf[rcnt],length,MPI_INT,status.MPI_SOURCE, 
78                              status.MPI_TAG,MPI_COMM_WORLD,&rtn_status);
79                     rcnt++;
80                 }
81             }
82         }
83
84         /* Wait on sends */
85         for (i=0; i<scnt; i++) {
86             MPI_Wait( &rsend[i], &status );
87         }
88
89         /* free buffers */
90         for (i=0; i<rcnt; i++) free(recv_buf[i]);
91         free( send_buf );
92         
93         MPI_Barrier( MPI_COMM_WORLD );
94         if (myid == 0 && (count % 64) == 0) {
95             printf( "All processes completed for count = %ld ints of data\n", 
96                     (long)count ); fflush(stdout);
97         }
98     }
99
100     MPI_Finalize();
101     return 0;
102 }
103