Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright lines.
[simgrid.git] / teshsuite / smpi / coll-scatter / coll-scatter.c
1 /* Copyright (c) 2012-2021. 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 #include <stdio.h>
8 #include <mpi.h>
9
10 int main(int argc, char **argv)
11 {
12   int size;
13   int rank;
14   int success = 1;
15   int retval;
16   int sendcount = 1;            // one double to each process
17   int recvcount = 1;
18   double *sndbuf = NULL;
19   double rcvd;
20   int root = 0;                 // arbitrary choice
21
22   MPI_Init(&argc, &argv);
23   MPI_Comm_size(MPI_COMM_WORLD, &size);
24   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
25   MPI_Comm_set_errhandler(MPI_COMM_WORLD, MPI_ERRORS_RETURN);
26
27   // on root, initialize sendbuf
28   if (root == rank) {
29     sndbuf = malloc(size * sizeof(double));
30     for (int i = 0; i < size; i++) {
31       sndbuf[i] = (double) i;
32     }
33   }
34
35   retval = MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, NULL, recvcount, MPI_DOUBLE, root, MPI_COMM_WORLD);
36   if(retval!=MPI_ERR_BUFFER)
37     printf("MPI_Scatter did not return MPI_ERR_BUFFER for empty recvbuf\n");
38   retval = MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd, -1, MPI_DOUBLE, root, MPI_COMM_WORLD);
39   if(retval!=MPI_ERR_COUNT)
40     printf("MPI_Scatter did not return MPI_ERR_COUNT for -1 recvcount\n");
41   retval = MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd, recvcount, MPI_DATATYPE_NULL, root, MPI_COMM_WORLD);
42   if(retval!=MPI_ERR_TYPE)
43     printf("MPI_Scatter did not return MPI_ERR_TYPE for MPI_DATATYPE_NULL sendtype\n");
44   retval = MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd, recvcount, MPI_DOUBLE, -1, MPI_COMM_WORLD);
45   if(retval!=MPI_ERR_ROOT)
46     printf("MPI_Scatter did not return MPI_ERR_ROOT for root -1\n");
47   retval = MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd, recvcount, MPI_DOUBLE, size+1, MPI_COMM_WORLD);
48   if(retval!=MPI_ERR_ROOT)
49     printf("MPI_Scatter did not return MPI_ERR_ROOT for root > size\n");
50   retval = MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd, recvcount, MPI_DOUBLE, root, MPI_COMM_NULL);
51   if(retval!=MPI_ERR_COMM)
52     printf("MPI_Scatter did not return MPI_ERR_COMM for MPI_COMM_NULL comm\n");
53
54   retval = MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd, recvcount, MPI_DOUBLE, root, MPI_COMM_WORLD);
55   if (root == rank) {
56     free(sndbuf);
57   }
58   if (retval != MPI_SUCCESS) {
59     fprintf(stderr, "(%s:%d) MPI_Scatter() returned retval=%d\n", __FILE__, __LINE__, retval);
60     return 0;
61   }
62   // verification
63   if ((double) rank != rcvd) {
64     fprintf(stderr, "[%d] has %f instead of %d\n", rank, rcvd, rank);
65     success = 0;
66   }
67   /* test 1 */
68   if (0 == rank)
69     printf("** Small Test Result: ...\n");
70   if (!success)
71     printf("\t[%d] failed.\n", rank);
72   else
73     printf("\t[%d] ok.\n", rank);
74
75   MPI_Finalize();
76   return 0;
77 }