Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
test more extensively error returns for collectives.
[simgrid.git] / teshsuite / smpi / coll-reduce-scatter / coll-reduce-scatter.c
1 /* Copyright (c) 2013-2019. The SimGrid Team.
2  * All rights reserved.                           */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 /*
8  * Test of reduce scatter.
9  * Each processor contributes its rank + the index to the reduction,  then receives the ith sum
10  * Can be called with any number of processors.
11  */
12
13 #include "mpi.h"
14 #include <stdio.h>
15 #include <stdlib.h>
16
17 int main( int argc, char **argv )
18 {
19   int err = 0;
20   int toterr;
21   int size;
22   int rank;
23   int i;
24   MPI_Comm comm;
25
26   MPI_Init( &argc, &argv );
27   comm = MPI_COMM_WORLD;
28
29   MPI_Comm_size( comm, &size );
30   MPI_Comm_rank( comm, &rank );
31   int* sendbuf = (int *) malloc( size * sizeof(int) );
32   for (i=0; i<size; i++)
33     sendbuf[i] = rank + i;
34   int* recvcounts = (int*) malloc (size * sizeof(int));
35   int* recvbuf  = (int*) malloc (size * sizeof(int));
36   for (i=0; i<size; i++)
37     recvcounts[i] = 1;
38   int retval;
39
40   retval = MPI_Reduce_scatter(NULL, recvbuf, recvcounts, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
41   if(retval!=MPI_ERR_BUFFER)
42     printf("MPI_Reduce_scatter did not return MPI_ERR_BUFFER for empty sendbuf\n");
43   retval = MPI_Reduce_scatter(sendbuf, NULL, recvcounts, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
44   if(retval!=MPI_ERR_BUFFER)
45     printf("MPI_Reduce_scatter did not return MPI_ERR_BUFFER for empty recvbuf\n");
46   retval = MPI_Reduce_scatter(sendbuf, recvbuf, NULL, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
47   if(retval!=MPI_ERR_ARG)
48     printf("MPI_Reduce_scatter did not return MPI_ERR_ARG for NULL recvcounts\n");
49   retval = MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, MPI_DATATYPE_NULL, MPI_SUM, MPI_COMM_WORLD);
50   if(retval!=MPI_ERR_TYPE)
51     printf("MPI_Reduce_scatter did not return MPI_ERR_TYPE for MPI_DATATYPE_NULL type\n");
52   retval = MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, MPI_DOUBLE, MPI_OP_NULL, MPI_COMM_WORLD);
53   if(retval!=MPI_ERR_OP)
54     printf("MPI_Reduce_scatter did not return MPI_ERR_OP for MPI_OP_NULL op\n");
55   retval = MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, MPI_DOUBLE, MPI_SUM, MPI_COMM_NULL);
56   if(retval!=MPI_ERR_COMM)
57     printf("MPI_Reduce_scatter did not return MPI_ERR_COMM for MPI_COMM_NULL comm\n");
58
59   MPI_Reduce_scatter( sendbuf, recvbuf, recvcounts, MPI_INT, MPI_SUM, comm );
60   int sumval = size * rank + ((size - 1) * size)/2;
61   /* recvbuf should be size * (rank + i) */
62   if (recvbuf[0] != sumval) {
63     err++;
64     fprintf( stdout, "Did not get expected value for reduce scatter\n" );
65     fprintf( stdout, "[%d] Got %d expected %d\n", rank, recvbuf[0], sumval );
66   }
67
68   MPI_Allreduce( &err, &toterr, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD );
69   if (rank == 0 && toterr == 0) {
70     printf( " No Errors\n" );
71   }
72   free(sendbuf);
73   free(recvcounts);
74   free(recvbuf);
75
76   MPI_Finalize();
77
78   return toterr;
79 }