Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
our coll-* tests include bad calls for coverage, they need MPI_ERRORS_RETURN
[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   MPI_Comm_set_errhandler(MPI_COMM_WORLD, MPI_ERRORS_RETURN);
32
33   int* sendbuf = (int *) malloc( size * sizeof(int) );
34   for (i=0; i<size; i++)
35     sendbuf[i] = rank + i;
36   int* recvcounts = (int*) malloc (size * sizeof(int));
37   int* recvbuf  = (int*) malloc (size * sizeof(int));
38   for (i=0; i<size; i++)
39     recvcounts[i] = 1;
40   int retval;
41
42   retval = MPI_Reduce_scatter(NULL, recvbuf, recvcounts, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
43   if(retval!=MPI_ERR_BUFFER)
44     printf("MPI_Reduce_scatter did not return MPI_ERR_BUFFER for empty sendbuf\n");
45   retval = MPI_Reduce_scatter(sendbuf, NULL, recvcounts, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
46   if(retval!=MPI_ERR_BUFFER)
47     printf("MPI_Reduce_scatter did not return MPI_ERR_BUFFER for empty recvbuf\n");
48   retval = MPI_Reduce_scatter(sendbuf, recvbuf, NULL, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
49   if(retval!=MPI_ERR_ARG)
50     printf("MPI_Reduce_scatter did not return MPI_ERR_ARG for NULL recvcounts\n");
51   retval = MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, MPI_DATATYPE_NULL, MPI_SUM, MPI_COMM_WORLD);
52   if(retval!=MPI_ERR_TYPE)
53     printf("MPI_Reduce_scatter did not return MPI_ERR_TYPE for MPI_DATATYPE_NULL type\n");
54   retval = MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, MPI_DOUBLE, MPI_OP_NULL, MPI_COMM_WORLD);
55   if(retval!=MPI_ERR_OP)
56     printf("MPI_Reduce_scatter did not return MPI_ERR_OP for MPI_OP_NULL op\n");
57   retval = MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, MPI_DOUBLE, MPI_SUM, MPI_COMM_NULL);
58   if(retval!=MPI_ERR_COMM)
59     printf("MPI_Reduce_scatter did not return MPI_ERR_COMM for MPI_COMM_NULL comm\n");
60
61   MPI_Reduce_scatter( sendbuf, recvbuf, recvcounts, MPI_INT, MPI_SUM, comm );
62   int sumval = size * rank + ((size - 1) * size)/2;
63   /* recvbuf should be size * (rank + i) */
64   if (recvbuf[0] != sumval) {
65     err++;
66     fprintf( stdout, "Did not get expected value for reduce scatter\n" );
67     fprintf( stdout, "[%d] Got %d expected %d\n", rank, recvbuf[0], sumval );
68   }
69
70   MPI_Allreduce( &err, &toterr, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD );
71   if (rank == 0 && toterr == 0) {
72     printf( " No Errors\n" );
73   }
74   free(sendbuf);
75   free(recvcounts);
76   free(recvbuf);
77
78   MPI_Finalize();
79
80   return toterr;
81 }