Statistiques
| Révision :

root / Pi / C / MPI / Pi_MPI.c @ 184

Historique | Voir | Annoter | Télécharger (7,56 ko)

1 9 equemene
//
2 9 equemene
// Estimation of Pi using Monte Carlo exploration process
3 9 equemene
// gcc -std=c99 -O3 -o Pi Pi.c -lm
4 59 equemene
// Emmanuel Quemener <emmanuel.quemener@ens-lyon.fr>
5 9 equemene
6 56 equemene
// Needed for gethostname
7 56 equemene
#define _BSD_SOURCE
8 56 equemene
#include <sys/unistd.h>
9 56 equemene
10 9 equemene
#include <math.h>
11 9 equemene
#include <stdio.h>
12 9 equemene
#include <stdlib.h>
13 25 equemene
#include <limits.h>
14 9 equemene
#include <mpi.h>
15 126 equemene
#include <stddef.h>
16 9 equemene
17 23 equemene
#include <sys/time.h>
18 23 equemene
19 9 equemene
// Marsaglia RNG very simple implementation
20 9 equemene
#define znew  ((z=36969*(z&65535)+(z>>16))<<16)
21 9 equemene
#define wnew  ((w=18000*(w&65535)+(w>>16))&65535)
22 9 equemene
#define MWC   (znew+wnew)
23 9 equemene
#define SHR3  (jsr=(jsr=(jsr=jsr^(jsr<<17))^(jsr>>13))^(jsr<<5))
24 9 equemene
#define CONG  (jcong=69069*jcong+1234567)
25 9 equemene
#define KISS  ((MWC^CONG)+SHR3)
26 9 equemene
27 9 equemene
#define MWCfp MWC * 2.328306435454494e-10f
28 9 equemene
#define KISSfp KISS * 2.328306435454494e-10f
29 126 equemene
#define SHR3fp SHR3 * 2.328306435454494e-10f
30 126 equemene
#define CONGfp CONG * 2.328306435454494e-10f
31 9 equemene
32 9 equemene
#define ITERATIONS 1000000000
33 9 equemene
34 9 equemene
#ifdef LONG
35 25 equemene
#define LENGTH long long
36 9 equemene
#else
37 25 equemene
#define LENGTH int
38 9 equemene
#endif
39 9 equemene
40 126 equemene
typedef struct compute_node {
41 126 equemene
        LENGTH inside;
42 126 equemene
        long int useconds;
43 126 equemene
} result;
44 126 equemene
45 56 equemene
unsigned int rotl(unsigned int value, int shift) {
46 56 equemene
    return (value << shift) | (value >> (sizeof(value) * CHAR_BIT - shift));
47 56 equemene
}
48 56 equemene
49 56 equemene
unsigned int rotr(unsigned int value, int shift) {
50 56 equemene
    return (value >> shift) | (value << (sizeof(value) * CHAR_BIT - shift));
51 56 equemene
}
52 56 equemene
53 9 equemene
LENGTH MainLoopGlobal(LENGTH iterations,unsigned int seed_w,unsigned int seed_z)
54 9 equemene
{
55 126 equemene
#if defined TCONG
56 126 equemene
   unsigned int jcong=seed_z;
57 126 equemene
#elif defined TSHR3
58 126 equemene
   unsigned int jsr=seed_w;
59 126 equemene
#elif defined TMWC
60 9 equemene
   unsigned int z=seed_z;
61 9 equemene
   unsigned int w=seed_w;
62 126 equemene
#elif defined TKISS
63 126 equemene
   unsigned int jcong=seed_z;
64 126 equemene
   unsigned int jsr=seed_w;
65 126 equemene
   unsigned int z=seed_z;
66 126 equemene
   unsigned int w=seed_w;
67 126 equemene
#endif
68 9 equemene
69 126 equemene
  LENGTH total=0;
70 9 equemene
71 9 equemene
   for (LENGTH i=0;i<iterations;i++) {
72 9 equemene
73 126 equemene
#if defined TINT32
74 126 equemene
    #define THEONE 1073741824
75 126 equemene
    #if defined TCONG
76 126 equemene
        unsigned int x=CONG>>17 ;
77 126 equemene
        unsigned int y=CONG>>17 ;
78 126 equemene
    #elif defined TSHR3
79 126 equemene
        unsigned int x=SHR3>>17 ;
80 126 equemene
        unsigned int y=SHR3>>17 ;
81 126 equemene
    #elif defined TMWC
82 126 equemene
        unsigned int x=MWC>>17 ;
83 126 equemene
        unsigned int y=MWC>>17 ;
84 126 equemene
    #elif defined TKISS
85 126 equemene
        unsigned int x=KISS>>17 ;
86 126 equemene
        unsigned int y=KISS>>17 ;
87 126 equemene
    #endif
88 126 equemene
#elif defined TINT64
89 126 equemene
    #define THEONE 4611686018427387904
90 126 equemene
    #if defined TCONG
91 126 equemene
        unsigned long x=(unsigned long)(CONG>>1) ;
92 126 equemene
        unsigned long y=(unsigned long)(CONG>>1) ;
93 126 equemene
    #elif defined TSHR3
94 126 equemene
        unsigned long x=(unsigned long)(SHR3>>1) ;
95 126 equemene
        unsigned long y=(unsigned long)(SHR3>>1) ;
96 126 equemene
    #elif defined TMWC
97 126 equemene
        unsigned long x=(unsigned long)(MWC>>1) ;
98 126 equemene
        unsigned long y=(unsigned long)(MWC>>1) ;
99 126 equemene
    #elif defined TKISS
100 126 equemene
        unsigned long x=(unsigned long)(KISS>>1) ;
101 126 equemene
        unsigned long y=(unsigned long)(KISS>>1) ;
102 126 equemene
    #endif
103 126 equemene
#elif defined TFP32
104 126 equemene
    #define THEONE 1.0f
105 126 equemene
    #if defined TCONG
106 126 equemene
        float x=CONGfp ;
107 126 equemene
        float y=CONGfp ;
108 126 equemene
    #elif defined TSHR3
109 126 equemene
        float x=SHR3fp ;
110 126 equemene
        float y=SHR3fp ;
111 126 equemene
    #elif defined TMWC
112 126 equemene
        float x=MWCfp ;
113 126 equemene
        float y=MWCfp ;
114 126 equemene
    #elif defined TKISS
115 126 equemene
      float x=KISSfp ;
116 126 equemene
      float y=KISSfp ;
117 126 equemene
    #endif
118 126 equemene
#elif defined TFP64
119 126 equemene
    #define THEONE 1.0f
120 126 equemene
    #if defined TCONG
121 126 equemene
        double x=(double)CONGfp ;
122 126 equemene
        double y=(double)CONGfp ;
123 126 equemene
    #elif defined TSHR3
124 126 equemene
        double x=(double)SHR3fp ;
125 126 equemene
        double y=(double)SHR3fp ;
126 126 equemene
    #elif defined TMWC
127 126 equemene
        double x=(double)MWCfp ;
128 126 equemene
        double y=(double)MWCfp ;
129 126 equemene
    #elif defined TKISS
130 126 equemene
        double x=(double)KISSfp ;
131 126 equemene
        double y=(double)KISSfp ;
132 126 equemene
    #endif
133 126 equemene
#endif
134 9 equemene
135 9 equemene
      // Matching test
136 126 equemene
      unsigned long inside=((x*x+y*y) < THEONE) ? 1:0;
137 9 equemene
      total+=inside;
138 9 equemene
   }
139 126 equemene
140 9 equemene
   return(total);
141 9 equemene
}
142 9 equemene
143 9 equemene
int main(int argc, char *argv[]) {
144 9 equemene
145 56 equemene
  unsigned int seed_z=362436069,seed_w=52128862;
146 126 equemene
  LENGTH iterations=ITERATIONS,inside[8192],insides,part_inside,part_iterations;
147 9 equemene
  int numtasks,rank,rc,tag=1,i;
148 9 equemene
  float pi;
149 9 equemene
150 56 equemene
  char hostname[128];
151 56 equemene
152 56 equemene
  gethostname(hostname, sizeof hostname);
153 56 equemene
154 23 equemene
  struct timeval start,end;
155 23 equemene
  long int useconds;
156 23 equemene
157 9 equemene
  MPI_Status Stat;
158 80 equemene
159 9 equemene
  rc = MPI_Init(&argc,&argv);
160 9 equemene
  if (rc != MPI_SUCCESS) {
161 9 equemene
    printf ("Error starting MPI program. Terminating.\n");
162 9 equemene
    MPI_Abort(MPI_COMM_WORLD, rc);
163 9 equemene
  }
164 9 equemene
165 9 equemene
  MPI_Comm_size(MPI_COMM_WORLD,&numtasks);
166 9 equemene
  MPI_Comm_rank(MPI_COMM_WORLD,&rank);
167 9 equemene
168 126 equemene
  const int nitems=2;
169 126 equemene
  int blocklengths[2] = {1,1};
170 126 equemene
171 126 equemene
#ifdef LONG
172 126 equemene
  MPI_Datatype types[2] = {MPI_LONG, MPI_LONG};
173 126 equemene
#else
174 126 equemene
  MPI_Datatype types[2] = {MPI_INT, MPI_LONG};
175 126 equemene
#endif
176 126 equemene
177 126 equemene
  MPI_Datatype mpi_result_type;
178 126 equemene
  MPI_Aint     offsets[2];
179 126 equemene
180 126 equemene
  offsets[0] = offsetof(result, inside);
181 126 equemene
  offsets[1] = offsetof(result, useconds);
182 126 equemene
183 126 equemene
  MPI_Type_create_struct(nitems, blocklengths, offsets, types, &mpi_result_type);
184 126 equemene
  MPI_Type_commit(&mpi_result_type);
185 126 equemene
186 9 equemene
  if (rank==0) {
187 56 equemene
188 9 equemene
    if (argc > 1) {
189 25 equemene
      iterations=(LENGTH)atoll(argv[1]);
190 9 equemene
    }
191 9 equemene
    else {
192 9 equemene
      printf("\n\tPi : Estimate Pi with Monte Carlo exploration\n\n");
193 9 equemene
      printf("\t\t#1 : number of iterations (default 1 billion)\n\n");
194 9 equemene
    }
195 56 equemene
196 126 equemene
    printf ("\tSizeof int = %lld bytes.\n", (long long)sizeof(int));
197 126 equemene
    printf ("\tSizeof long = %lld bytes.\n", (long long)sizeof(long));
198 126 equemene
    printf ("\tSizeof long long = %lld bytes.\n", (long long)sizeof(long long));
199 56 equemene
200 126 equemene
    printf ("\tMax int = %u\n", INT_MAX);
201 126 equemene
    printf ("\tMax long = %ld\n", LONG_MAX);
202 126 equemene
    printf ("\tMax long long = %lld\n\n", LLONG_MAX);
203 56 equemene
204 126 equemene
    part_iterations=((iterations%numtasks) == 0) ? iterations/numtasks:iterations/numtasks+1 ;
205 56 equemene
206 9 equemene
    // Split part of code
207 9 equemene
    for (i=1;i<numtasks;i++) {
208 56 equemene
209 25 equemene
#ifdef LONG
210 80 equemene
      rc = MPI_Send(&part_iterations, 1, MPI_LONG_LONG, i, tag,
211 80 equemene
                    MPI_COMM_WORLD);
212 25 equemene
#else
213 80 equemene
      rc = MPI_Send(&part_iterations, 1, MPI_INT, i, tag,
214 80 equemene
                    MPI_COMM_WORLD);
215 80 equemene
#endif
216 9 equemene
    }
217 9 equemene
218 23 equemene
    gettimeofday(&start,(struct timezone *)0);
219 56 equemene
220 9 equemene
    insides=MainLoopGlobal(part_iterations,seed_w,seed_z);
221 56 equemene
222 23 equemene
    gettimeofday(&end,(struct timezone *)0);
223 23 equemene
    useconds=(end.tv_sec-start.tv_sec)*1000000+end.tv_usec-start.tv_usec;
224 56 equemene
225 126 equemene
    printf("\tOn %s with rank %i, find %lld inside in %lu useconds.\n",
226 56 equemene
             hostname,rank,(long long)insides,useconds);
227 56 equemene
228 9 equemene
    // Join part of code
229 56 equemene
      for (i=1;i<numtasks;i++) {
230 126 equemene
231 126 equemene
        result recv;
232 126 equemene
233 126 equemene
        rc = MPI_Recv(&recv, 1, mpi_result_type, i, tag,
234 80 equemene
                      MPI_COMM_WORLD, &Stat);
235 126 equemene
236 126 equemene
        inside[i]=recv.inside;
237 126 equemene
        useconds=recv.useconds;
238 126 equemene
239 126 equemene
        printf("\tReceive from %i, find %lld inside in %lu useconds\n",i,(long long)inside[i],useconds);
240 126 equemene
241 56 equemene
        insides+=inside[i];
242 56 equemene
      }
243 56 equemene
244 126 equemene
      pi=4.*(float)insides/(float)(part_iterations*numtasks);
245 56 equemene
246 56 equemene
      printf("\n\tPi=%.40f\n\twith error %.40f\n\twith %lld iterations\n\n",pi,
247 56 equemene
             fabs(pi-4*atan(1.))/pi,(long long)iterations);
248 56 equemene
249 9 equemene
  }
250 9 equemene
  else
251 9 equemene
    {
252 9 equemene
      // Receive information from master
253 25 equemene
#ifdef LONG
254 80 equemene
      rc = MPI_Recv(&part_iterations, 1, MPI_LONG_LONG, 0, tag,
255 80 equemene
                    MPI_COMM_WORLD, &Stat);
256 25 equemene
#else
257 80 equemene
      rc = MPI_Recv(&part_iterations, 1, MPI_INT, 0, tag,
258 80 equemene
                    MPI_COMM_WORLD, &Stat);
259 25 equemene
#endif
260 9 equemene
261 126 equemene
      /*      Not such a good idea to print information...
262 126 equemene
              printf("\tOn %s with rank %i, receive from master %lld\n",
263 126 equemene
              hostname,rank,(long long)part_iterations); */
264 56 equemene
265 23 equemene
      gettimeofday(&start,(struct timezone *)0);
266 23 equemene
267 56 equemene
      part_inside=MainLoopGlobal(part_iterations,rotr(seed_w,rank),rotl(seed_z,rank));
268 9 equemene
269 56 equemene
      gettimeofday(&end,(struct timezone *)0);
270 56 equemene
      useconds=(end.tv_sec-start.tv_sec)*1000000+end.tv_usec-start.tv_usec;
271 126 equemene
272 126 equemene
      /*
273 56 equemene
      printf("\tOn %s with rank %i, find %lld inside in %lu useconds.\n",
274 56 equemene
             hostname,rank,(long long)part_inside,useconds);
275 126 equemene
      */
276 23 equemene
277 126 equemene
      result send;
278 126 equemene
      send.inside=part_inside;
279 126 equemene
      send.useconds=useconds;
280 126 equemene
281 126 equemene
      rc = MPI_Send(&send, 1, mpi_result_type, 0, tag, MPI_COMM_WORLD);
282 126 equemene
283 9 equemene
    }
284 126 equemene
285 126 equemene
  MPI_Type_free(&mpi_result_type);
286 9 equemene
  MPI_Finalize();
287 9 equemene
}