Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
9f909f28a9caa3c40ebc37e87230fb211db4383b
[simgrid.git] / teshsuite / smpi / mpich3-test / datatype / tfree.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *
4  *  (C) 2003 by Argonne National Laboratory.
5  *      See COPYRIGHT in top-level directory.
6  */
7 #include "mpi.h"
8 #include <stdio.h>
9 #include <stdlib.h>
10 #include "mpitest.h"
11
12 /*
13 static char MTEST_Descrip[] = "Test that freed datatypes have reference count semantics";
14 */
15
16 #define VEC_NELM 128
17 #define VEC_STRIDE 8
18
19 int main( int argc, char *argv[] )
20 {
21     int errs = 0;
22     int rank, size, source, dest, i;
23     MPI_Comm      comm;
24     MPI_Status    status;
25     MPI_Request   req;
26     MPI_Datatype  strideType;
27     MPI_Datatype  tmpType[1024];
28     int           *buf = 0;
29
30     MTest_Init( &argc, &argv );
31
32     comm = MPI_COMM_WORLD;
33
34     MPI_Comm_rank( comm, &rank );
35     MPI_Comm_size( comm, &size );
36
37     if (size < 2) {
38         fprintf( stderr, "This test requires at least two processes." );
39         MPI_Abort( MPI_COMM_WORLD, 1 );
40         exit(1);
41     }
42
43     source  = 0;
44     dest    = size - 1;
45
46     /* 
47        The idea here is to create a simple but non-contig datatype,
48        perform an irecv with it, free it, and then create 
49        many new datatypes.  While not a complete test, if the datatype
50        was freed and the space was reused, this test may detect 
51        that error 
52        A similar test for sends might work by sending a large enough message
53        to force the use of rendezvous send. 
54     */
55     MPI_Type_vector( VEC_NELM, 1, VEC_STRIDE, MPI_INT, &strideType );
56     MPI_Type_commit( &strideType );
57
58     if (rank == dest) {
59         buf = (int *)malloc( VEC_NELM * VEC_STRIDE * sizeof(int) );
60         for (i=0; i<VEC_NELM*VEC_STRIDE; i++) buf[i] = -i;
61         MPI_Irecv( buf, 1, strideType, source, 0, comm, &req );
62         MPI_Type_free( &strideType );
63
64         for (i=0; i<1024; i++) {
65             MPI_Type_vector( VEC_NELM, 1, 1, MPI_INT, &tmpType[i] );
66             MPI_Type_commit( &tmpType[i] );
67         }
68
69         MPI_Sendrecv( 0, 0, MPI_INT, source, 1, 
70                       0, 0, MPI_INT, source, 1, comm, &status );
71
72         MPI_Wait( &req, &status );
73         for (i=0; i<VEC_NELM; i++) {
74             if (buf[VEC_STRIDE*i] != i) {
75                 errs++;
76                 if (errs < 10) {
77                     printf( "buf[%d] = %d, expected %d\n", VEC_STRIDE*i, 
78                             buf[VEC_STRIDE*i], i );
79                 }
80             }
81         }
82         for (i=0; i<1024; i++) {
83             MPI_Type_free( &tmpType[i] );
84         }
85         free( buf );
86     }
87     else if (rank == source) {
88         buf = (int *)malloc( VEC_NELM * sizeof(int) );
89         for (i=0; i<VEC_NELM; i++) buf[i] = i;
90         /* Synchronize with the receiver */
91         MPI_Sendrecv( 0, 0, MPI_INT, dest, 1, 
92                       0, 0, MPI_INT, dest, 1, comm, &status );
93         MPI_Send( buf, VEC_NELM, MPI_INT, dest, 0, comm );
94         free( buf );
95     }
96
97     /* Clean up the strideType */
98     if (rank != dest) {
99         MPI_Type_free( &strideType );
100     }
101
102
103     MTest_Finalize( errs );
104     MPI_Finalize();
105     return 0;
106 }