Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
further refactoring
[simgrid.git] / teshsuite / smpi / coll-scatter / coll-scatter.c
1 /* Copyright (c) 2012-2014. 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  * MESSAGE PASSING INTERFACE TEST CASE SUITE
9  *
10  * Copyright IBM Corp. 1995
11  *
12  * IBM Corp. hereby grants a non-exclusive license to use, copy, modify, and
13  *distribute this software for any purpose and without fee provided that the
14  *above copyright notice and the following paragraphs appear in all copies.
15
16  * IBM Corp. makes no representation that the test cases comprising this
17  * suite are correct or are an accurate representation of any standard.
18
19  * In no event shall IBM be liable to any party for direct, indirect, special
20  * incidental, or consequential damage arising out of the use of this software
21  * even if IBM Corp. has been advised of the possibility of such damage.
22
23  * IBM CORP. SPECIFICALLY DISCLAIMS ANY WARRANTIES INCLUDING, BUT NOT LIMITED
24  * TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
25  * PURPOSE.  THE SOFTWARE PROVIDED HEREUNDER IS ON AN "AS IS" BASIS AND IBM
26  * CORP. HAS NO OBLIGATION TO PROVIDE MAINTENANCE, SUPPORT, UPDATES,
27  * ENHANCEMENTS, OR MODIFICATIONS.
28  **/
29 #include <stdio.h>
30 #include <mpi.h>
31
32 static int ibm_test(int rank, int size)
33 {
34
35 #define MAXLEN  10000
36
37   int success = 1;
38   int root = 0;
39   int i, j, k;
40   int *out;
41   int *in;
42
43   out = malloc(MAXLEN * 64 * sizeof(int));
44   in = malloc(MAXLEN * sizeof(int));
45
46   for (j = 1; j <= MAXLEN; j *= 10) {
47     root = (root + 1) % size;
48     if (rank == root)
49       for (i = 0; i < j * size; i++)
50         out[i] = i;
51
52     MPI_Scatter(out, j, MPI_INT, in, j, MPI_INT, root, MPI_COMM_WORLD);
53
54     for (k = 0; k < j; k++) {
55       if (in[k] != k + rank * j) {
56         fprintf(stderr, "task %d bad answer (%d) at index %d k of %d (should be %d)",rank, in[k], k, j, (k + rank * j));
57         return (0);
58       }
59     }
60   }
61   free(out);
62   free(in);
63   MPI_Barrier(MPI_COMM_WORLD);
64   return (success);
65 }
66
67 /** small test: the root sends a single distinct double to other processes */
68 static int small_test(int rank, int size)
69 {
70   int success = 1;
71   int retval;
72   int sendcount = 1;            // one double to each process
73   int recvcount = 1;
74   int i;
75   double *sndbuf = NULL;
76   double rcvd;
77   int root = 0;                 // arbitrary choice
78
79   // on root, initialize sendbuf
80   if (root == rank) {
81     sndbuf = malloc(size * sizeof(double));
82     for (i = 0; i < size; i++) {
83       sndbuf[i] = (double) i;
84     }
85   }
86
87   retval = MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd, recvcount, MPI_DOUBLE, root, MPI_COMM_WORLD);
88   if (root == rank) {
89     free(sndbuf);
90     }
91   if (retval != MPI_SUCCESS) {
92     fprintf(stderr, "(%s:%d) MPI_Scatter() returned retval=%d\n", __FILE__, __LINE__, retval);
93     return 0;
94   }
95   // verification
96   if ((double) rank != rcvd) {
97     fprintf(stderr, "[%d] has %f instead of %d\n", rank, rcvd, rank);
98     success = 0;
99   }
100   return (success);
101 }
102
103 int main(int argc, char **argv)
104 {
105   int size, rank;
106
107   MPI_Init(&argc, &argv);
108   MPI_Comm_size(MPI_COMM_WORLD, &size);
109   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
110
111   /* test 1 */
112   if (0 == rank)
113     printf("** Small Test Result: ...\n");
114   if (!small_test(rank, size))
115     printf("\t[%d] failed.\n", rank);
116   else
117     printf("\t[%d] ok.\n", rank);
118
119   MPI_Barrier(MPI_COMM_WORLD);
120
121   /* test 2 */
122   if (0 == rank)
123     printf("** IBM Test Result: ...\n");
124   if (!ibm_test(rank, size))
125     printf("\t[%d] failed.\n", rank);
126   else
127     printf("\t[%d] ok.\n", rank);
128
129   MPI_Finalize();
130   return 0;
131 }