Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
make dt runnable (fix argv accesses)
[simgrid.git] / examples / smpi / NAS / dt.c
1 /*************************************************************************
2  *                                                                       * 
3  *    N  A  S   P A R A L L E L   B E N C H M A R K S  3.3     *
4  *                                     *
5  *                  D T       *
6  *                                     *
7  ************************************************************************* 
8  *                                     *
9  *   This benchmark is part of the NAS Parallel Benchmark 3.3 suite.   *
10  *                                     *
11  *   Permission to use, copy, distribute and modify this software    *
12  *   for any purpose with or without fee is hereby granted.  We      *
13  *   request, however, that all derived work reference the NAS       *
14  *   Parallel Benchmarks 3.3. This software is provided "as is"      *
15  *   without express or implied warranty.                *
16  *                                     *
17  *   Information on NPB 3.3, including the technical report, the     *
18  *   original specifications, source code, results and information     *
19  *   on how to submit new results, is available at:            *
20  *                                     *
21  *      http:  www.nas.nasa.gov/Software/NPB             *
22  *                                     *
23  *   Send comments or suggestions to  npb@nas.nasa.gov           *
24  *   Send bug reports to        npb-bugs@nas.nasa.gov        *
25  *                                     *
26  *     NAS Parallel Benchmarks Group                 *
27  *     NASA Ames Research Center                   *
28  *     Mail Stop: T27A-1                       *
29  *     Moffett Field, CA   94035-1000                *
30  *                                     *
31  *     E-mail:  npb@nas.nasa.gov                   *
32  *     Fax:   (650) 604-3957                     *
33  *                                     *
34  ************************************************************************* 
35  *                                     *
36  *   Author: M. Frumkin         *       *
37  *                                     *
38  *************************************************************************/
39
40 #include <stdlib.h>
41 #include <stdio.h>
42 #include <string.h>
43
44 #include "DGraph.h"
45 #include "smpi/mpi.h"
46 #include "nas_common.h"
47 #include "simgrid/instr.h" //TRACE_
48
49 int timer_on=0,timers_tot=64;
50 double start[64], elapsed[64];
51
52 char class;
53 int nprocs;
54 int num_samples;
55 int deviation;
56 int num_sources;
57
58 static int verify(char *bmname,double rnm2){
59   double verify_value=0.0;
60   double epsilon=1.0E-8;
61   int verified=-1;
62   if (class != 'U') {
63     if(class=='S') {
64      if(strstr(bmname,"BH")){
65        verify_value=30892725.0;
66      }else if(strstr(bmname,"WH")){
67        verify_value=67349758.0;
68      }else if(strstr(bmname,"SH")){
69        verify_value=58875767.0;
70      }else{
71        fprintf(stderr,"No such benchmark as %s.\n",bmname);
72      }
73      verified = 0;
74     }else if(class=='W') {
75       if(strstr(bmname,"BH")){
76         verify_value = 4102461.0;
77       }else if(strstr(bmname,"WH")){
78         verify_value = 204280762.0;
79       }else if(strstr(bmname,"SH")){
80         verify_value = 186944764.0;
81       }else{
82         fprintf(stderr,"No such benchmark as %s.\n",bmname);
83       }
84       verified = 0;
85     }else if(class=='A') {
86       if(strstr(bmname,"BH")){
87         verify_value = 17809491.0;
88       }else if(strstr(bmname,"WH")){
89         verify_value = 1289925229.0;
90       }else if(strstr(bmname,"SH")){
91         verify_value = 610856482.0;
92       }else{
93         fprintf(stderr,"No such benchmark as %s.\n",bmname);
94       }
95       verified = 0;
96     }else if(class=='B') {
97       if(strstr(bmname,"BH")){
98         verify_value = 4317114.0;
99       }else if(strstr(bmname,"WH")){
100         verify_value = 7877279917.0;
101       }else if(strstr(bmname,"SH")){
102         verify_value = 1836863082.0;
103       }else{
104         fprintf(stderr,"No such benchmark as %s.\n",bmname);
105         verified = 0;
106       }
107     }else if(class=='C') {
108       if(strstr(bmname,"BH")){
109         verify_value = 0.0;
110       }else if(strstr(bmname,"WH")){
111         verify_value = 0.0;
112       }else if(strstr(bmname,"SH")){
113         verify_value = 0.0;
114       }else{
115         fprintf(stderr,"No such benchmark as %s.\n",bmname);
116         verified = -1;
117       }
118     }else if(class=='D') {
119       if(strstr(bmname,"BH")){
120         verify_value = 0.0;
121       }else if(strstr(bmname,"WH")){
122         verify_value = 0.0;
123       }else if(strstr(bmname,"SH")){
124         verify_value = 0.0;
125       }else{
126         fprintf(stderr,"No such benchmark as %s.\n",bmname);
127       }
128       verified = -1;
129     }else{
130       fprintf(stderr,"No such class as %c.\n",class);
131     }
132     fprintf(stderr," %s L2 Norm = %f\n",bmname,rnm2);
133     if(verified==-1){
134       fprintf(stderr," No verification was performed.\n");
135     }else if( rnm2 - verify_value < epsilon && rnm2 - verify_value > -epsilon) {  /* abs here does not work on ALTIX */
136       verified = 1;
137       fprintf(stderr," Deviation = %f\n",(rnm2 - verify_value));
138     }else{
139       verified = 0;
140       fprintf(stderr," The correct verification value = %f\n",verify_value);
141       fprintf(stderr," Got value = %f\n",rnm2);
142     }
143   }else{
144     verified = -1;
145   }
146   return  verified;
147 }
148
149 static int ipowMod(int a,long long int n,int md){
150   int seed=1,q=a,r=1;
151   if(n<0){
152     fprintf(stderr,"ipowMod: exponent must be nonnegative exp=%lld\n",n);
153     n=-n; /* temp fix */
154     /*  return 1; */
155   }
156   if(md<=0){
157     fprintf(stderr,"ipowMod: module must be positive mod=%d",md);
158     return 1;
159   }
160   if(n==0) return 1;
161   while(n>1){
162     int n2 = n/2;
163     if (n2*2==n){
164       seed = (q*q)%md;
165       q=seed;
166       n = n2;
167     }else{
168       seed = (r*q)%md;
169       r=seed;
170       n = n-1;
171     }
172   }
173   seed = (r*q)%md;
174   return seed;
175 }
176
177 static DGraph *buildSH(const char cls){
178 /* Nodes of the graph must be topologically sorted to avoid MPI deadlock. */
179   DGraph *dg;
180   int numSources=num_sources; /* must be power of 2 */
181   int numOfLayers=0,tmpS=numSources>>1;
182   int firstLayerNode=0;
183   DGArc *ar=NULL;
184   DGNode *nd=NULL;
185   int mask=0x0,ndid=0,ndoff=0;
186   int i=0,j=0;
187   char nm[BLOCK_SIZE];
188
189   sprintf(nm,"DT_SH.%c",cls);
190   dg=newDGraph(nm);
191
192   while(tmpS>1){
193   numOfLayers++;
194   tmpS>>=1;
195   }
196   for(i=0;i<numSources;i++){
197   sprintf(nm,"Source.%d",i);
198   nd=newNode(nm);
199   AttachNode(dg,nd);
200   }
201   for(j=0;j<numOfLayers;j++){
202     mask=0x00000001<<j;
203     for(i=0;i<numSources;i++){
204       sprintf(nm,"Comparator.%d",(i+j*firstLayerNode));
205       nd=newNode(nm);
206       AttachNode(dg,nd);
207       ndoff=i&(~mask);
208       ndid=firstLayerNode+ndoff;
209       ar=newArc(dg->node[ndid],nd);
210       AttachArc(dg,ar);
211       ndoff+=mask;
212       ndid=firstLayerNode+ndoff;
213       ar=newArc(dg->node[ndid],nd);
214       AttachArc(dg,ar);
215     }
216     firstLayerNode+=numSources;
217   }
218   mask=0x00000001<<numOfLayers;
219   for(i=0;i<numSources;i++){
220     sprintf(nm,"Sink.%d",i);
221     nd=newNode(nm);
222     AttachNode(dg,nd);
223     ndoff=i&(~mask);
224     ndid=firstLayerNode+ndoff;
225     ar=newArc(dg->node[ndid],nd);
226     AttachArc(dg,ar);
227     ndoff+=mask;
228     ndid=firstLayerNode+ndoff;
229     ar=newArc(dg->node[ndid],nd);
230     AttachArc(dg,ar);
231   }
232   return dg;
233 }
234
235 static DGraph *buildWH(const char cls){
236 /*  Nodes of the graph must be topologically sorted to avoid MPI deadlock. */
237   int i=0,j=0;
238   int numSources=num_sources,maxInDeg=4;
239   int numLayerNodes=numSources,firstLayerNode=0;
240   int totComparators=0;
241   int numPrevLayerNodes=numLayerNodes;
242   int id=0,sid=0;
243   DGraph *dg;
244   DGNode *nd=NULL,*source=NULL,*tmp=NULL,*snd=NULL;
245   DGArc *ar=NULL;
246   char nm[BLOCK_SIZE];
247
248   sprintf(nm,"DT_WH.%c",cls);
249   dg=newDGraph(nm);
250
251   for(i=0;i<numSources;i++){
252     sprintf(nm,"Sink.%d",i);
253     nd=newNode(nm);
254     AttachNode(dg,nd);
255   }
256   totComparators=0;
257   numPrevLayerNodes=numLayerNodes;
258   while(numLayerNodes>maxInDeg){
259     numLayerNodes=numLayerNodes/maxInDeg;
260     if(numLayerNodes*maxInDeg<numPrevLayerNodes)numLayerNodes++;
261     for(i=0;i<numLayerNodes;i++){
262       sprintf(nm,"Comparator.%d",totComparators);
263       totComparators++;
264       nd=newNode(nm);
265       id=AttachNode(dg,nd);
266       for(j=0;j<maxInDeg;j++){
267         sid=i*maxInDeg+j;
268         if(sid>=numPrevLayerNodes) break;
269         snd=dg->node[firstLayerNode+sid];
270         ar=newArc(dg->node[id],snd);
271         AttachArc(dg,ar);
272       }
273     }
274     firstLayerNode+=numPrevLayerNodes;
275     numPrevLayerNodes=numLayerNodes;
276   }
277   source=newNode((char*)"Source");
278   AttachNode(dg,source);   
279   for(i=0;i<numPrevLayerNodes;i++){
280     nd=dg->node[firstLayerNode+i];
281     ar=newArc(source,nd);
282     AttachArc(dg,ar);
283   }
284
285   for(i=0;i<dg->numNodes/2;i++){  /* Topological sorting */
286     tmp=dg->node[i];
287     dg->node[i]=dg->node[dg->numNodes-1-i];
288     dg->node[i]->id=i;
289     dg->node[dg->numNodes-1-i]=tmp;
290     dg->node[dg->numNodes-1-i]->id=dg->numNodes-1-i;
291   }
292 return dg;
293 }
294
295 static DGraph *buildBH(const char cls){
296 /* Nodes of the graph must be topologically sorted to avoid MPI deadlock.*/
297   int i=0,j=0;
298   int numSources=num_sources,maxInDeg=4;
299   int numLayerNodes=numSources,firstLayerNode=0;
300   DGraph *dg;
301   DGNode *nd=NULL, *snd=NULL, *sink=NULL;
302   DGArc *ar=NULL;
303   int totComparators=0;
304   int numPrevLayerNodes=numLayerNodes;
305   int id=0, sid=0;
306   char nm[BLOCK_SIZE];
307
308   sprintf(nm,"DT_BH.%c",cls);
309   dg=newDGraph(nm);
310
311   for(i=0;i<numSources;i++){
312     sprintf(nm,"Source.%d",i);
313     nd=newNode(nm);
314     AttachNode(dg,nd);
315   }
316   while(numLayerNodes>maxInDeg){
317     numLayerNodes=numLayerNodes/maxInDeg;
318     if(numLayerNodes*maxInDeg<numPrevLayerNodes)numLayerNodes++;
319     for(i=0;i<numLayerNodes;i++){
320       sprintf(nm,"Comparator.%d",totComparators);
321       totComparators++;
322       nd=newNode(nm);
323       id=AttachNode(dg,nd);
324       for(j=0;j<maxInDeg;j++){
325         sid=i*maxInDeg+j;
326         if(sid>=numPrevLayerNodes) break;
327         snd=dg->node[firstLayerNode+sid];
328         ar=newArc(snd,dg->node[id]);
329         AttachArc(dg,ar);
330       }
331     }
332     firstLayerNode+=numPrevLayerNodes;
333     numPrevLayerNodes=numLayerNodes;
334   }
335   sink=newNode((char*)"Sink");
336   AttachNode(dg,sink);   
337   for(i=0;i<numPrevLayerNodes;i++){
338     nd=dg->node[firstLayerNode+i];
339     ar=newArc(nd,sink);
340     AttachArc(dg,ar);
341   }
342   return dg;
343 }
344
345 typedef struct{
346   int len;
347   double* val;
348 } Arr;
349
350 static Arr *newArr(int len){
351   Arr *arr=(Arr *)malloc(sizeof(Arr)); //Arr *arr=(Arr *)SMPI_SHARED_MALLOC(sizeof(Arr));
352   arr->len=len;
353   arr->val=(double *)malloc(len*sizeof(double)); //arr->val=(double *)SMPI_SHARED_MALLOC(len*sizeof(double));
354   return arr;
355 }
356
357 static void arrShow(Arr* a){
358   if(!a) fprintf(stderr,"-- NULL array\n");
359   else{
360     fprintf(stderr,"-- length=%d\n",a->len);
361   }
362 }
363
364 static double CheckVal(Arr *feat){
365   double csum=0.0;
366   for(int i=0;i<feat->len;i++){
367     csum+=feat->val[i]*feat->val[i]/feat->len; /* The truncation does not work since result will be 0 for large len  */
368   }
369   return csum;
370 }
371
372 static int GetFNumDPar(int* mean, int* stdev){
373   *mean=num_samples;
374   *stdev=deviation;
375   return 0;
376 }
377
378 static int GetFeatureNum(char *mbname,int id){
379   double tran=314159265.0;
380   double A=2*id+1;
381   double denom=randlc(&tran,&A);
382   char cval='S';
383   int mean=num_samples,stdev=128;
384   int rtfs=0,len=0;
385   GetFNumDPar(&mean,&stdev);
386   rtfs=ipowMod((int)(1/denom)*(int)cval,(long long int) (2*id+1),2*stdev);
387   if(rtfs<0) rtfs=-rtfs;
388   len=mean-stdev+rtfs;
389   return len;
390 }
391
392 static Arr* RandomFeatures(char *bmname,int fdim,int id){
393   int len=GetFeatureNum(bmname,id)*fdim;
394   Arr* feat=newArr(len);
395   int nxg=2,nyg=2,nzg=2,nfg=5;
396   int nx=421,ny=419,nz=1427,nf=3527;
397   long long int expon=(len*(id+1))%3141592;
398   int seedx=ipowMod(nxg,expon,nx), seedy=ipowMod(nyg,expon,ny), seedz=ipowMod(nzg,expon,nz),seedf=ipowMod(nfg,expon,nf);
399
400   if(timer_on){
401     timer_clear(id+1);
402     timer_start(id+1);
403   }
404   for(int i=0;i<len;i+=fdim){
405     seedx=(seedx*nxg)%nx;
406     seedy=(seedy*nyg)%ny;
407     seedz=(seedz*nzg)%nz;
408     seedf=(seedf*nfg)%nf;
409     feat->val[i]=seedx;
410     feat->val[i+1]=seedy;
411     feat->val[i+2]=seedz;
412     feat->val[i+3]=seedf;
413   }
414   if(timer_on){
415   timer_stop(id+1);
416   fprintf(stderr,"** RandomFeatures time in node %d = %f\n",id,timer_read(id+1));
417   }
418   return feat;
419 }
420
421 static void Resample(Arr *a,int blen){
422   long long int i=0,j=0,jlo=0,jhi=0;
423   double avval=0.0;
424   double *nval=(double *)malloc(blen*sizeof(double));
425   //double *nval=(double *)SMPI_SHARED_MALLOC(blen*sizeof(double));
426   for(i=0;i<blen;i++) nval[i]=0.0;
427   for(i=1;i<a->len-1;i++){
428     jlo=(int)(0.5*(2*i-1)*(blen/a->len));
429     jhi=(int)(0.5*(2*i+1)*(blen/a->len));
430
431     avval=a->val[i]/(jhi-jlo+1);
432     for(j=jlo;j<=jhi;j++){
433     nval[j]+=avval;
434     }
435   }
436   nval[0]=a->val[0];
437   nval[blen-1]=a->val[a->len-1];
438   free(a->val); //SMPI_SHARED_FREE(a->val);
439   a->val=nval;
440   a->len=blen;
441 }
442
443 #define fielddim 4
444 static Arr* WindowFilter(Arr *a, Arr* b,int w){
445   int i=0,j=0,k=0;
446   double rms0=0.0,rms1=0.0,rmsm1=0.0;
447   double weight=((double) (w+1))/(w+2);
448  
449   w+=1;
450   if(timer_on){
451     timer_clear(w);
452     timer_start(w);
453   }
454   if(a->len<b->len) Resample(a,b->len);
455   if(a->len>b->len) Resample(b,a->len);
456   for(i=fielddim;i<a->len-fielddim;i+=fielddim){
457     rms0=(a->val[i]-b->val[i])*(a->val[i]-b->val[i]) +(a->val[i+1]-b->val[i+1])*(a->val[i+1]-b->val[i+1])
458           +(a->val[i+2]-b->val[i+2])*(a->val[i+2]-b->val[i+2]) +(a->val[i+3]-b->val[i+3])*(a->val[i+3]-b->val[i+3]);
459     j=i+fielddim;
460     rms1=(a->val[j]-b->val[j])*(a->val[j]-b->val[j]) +(a->val[j+1]-b->val[j+1])*(a->val[j+1]-b->val[j+1])
461           +(a->val[j+2]-b->val[j+2])*(a->val[j+2]-b->val[j+2]) +(a->val[j+3]-b->val[j+3])*(a->val[j+3]-b->val[j+3]);
462     j=i-fielddim;
463     rmsm1=(a->val[j]-b->val[j])*(a->val[j]-b->val[j]) +(a->val[j+1]-b->val[j+1])*(a->val[j+1]-b->val[j+1])
464            +(a->val[j+2]-b->val[j+2])*(a->val[j+2]-b->val[j+2]) +(a->val[j+3]-b->val[j+3])*(a->val[j+3]-b->val[j+3]);
465     k=0;
466     if(rms1<rms0){
467       k=1;
468       rms0=rms1;
469     }
470     if(rmsm1<rms0) k=-1;
471     if(k==0){
472       j=i+fielddim;
473       a->val[i]=weight*b->val[i];
474       a->val[i+1]=weight*b->val[i+1];
475       a->val[i+2]=weight*b->val[i+2];
476       a->val[i+3]=weight*b->val[i+3];
477     }else if(k==1){
478       j=i+fielddim;
479       a->val[i]=weight*b->val[j];
480       a->val[i+1]=weight*b->val[j+1];
481       a->val[i+2]=weight*b->val[j+2];
482       a->val[i+3]=weight*b->val[j+3];
483     }else { /*if(k==-1)*/
484       j=i-fielddim;
485       a->val[i]=weight*b->val[j];
486       a->val[i+1]=weight*b->val[j+1];
487       a->val[i+2]=weight*b->val[j+2];
488       a->val[i+3]=weight*b->val[j+3];
489     }
490   }
491   if(timer_on){
492     timer_stop(w);
493     fprintf(stderr,"** WindowFilter time in node %d = %f\n",(w-1),timer_read(w));
494   }
495   return a;
496 }
497
498 static int SendResults(DGraph *dg,DGNode *nd,Arr *feat){
499   int i=0,tag=0;
500   DGArc *ar=NULL;
501   DGNode *head=NULL;
502   if(!feat) return 0;
503   TRACE_smpi_set_category ("SendResults");
504   for(i=0;i<nd->outDegree;i++){
505     ar=nd->outArc[i];
506     if(ar->tail!=nd) continue;
507     head=ar->head;
508     tag=ar->id;
509     if(head->address!=nd->address){
510       MPI_Send(&feat->len,1,MPI_INT,head->address,tag,MPI_COMM_WORLD);
511       MPI_Send(feat->val,feat->len,MPI_DOUBLE,head->address,tag,MPI_COMM_WORLD);
512     }
513   }
514   TRACE_smpi_set_category (NULL);
515   return 1;
516 }
517
518 static Arr* CombineStreams(DGraph *dg,DGNode *nd){
519   Arr *resfeat=newArr(num_samples*fielddim);
520   int i=0,len=0,tag=0;
521   DGArc *ar=NULL;
522   DGNode *tail=NULL;
523   MPI_Status status;
524   Arr *feat=NULL,*featp=NULL;
525
526   if(nd->inDegree==0) return NULL;
527   for(i=0;i<nd->inDegree;i++){
528     ar=nd->inArc[i];
529     if(ar->head!=nd) continue;
530     tail=ar->tail;
531     if(tail->address!=nd->address){
532       len=0;
533       tag=ar->id;
534       MPI_Recv(&len,1,MPI_INT,tail->address,tag,MPI_COMM_WORLD,&status);
535       feat=newArr(len);
536       MPI_Recv(feat->val,feat->len,MPI_DOUBLE,tail->address,tag,MPI_COMM_WORLD,&status);
537       resfeat=WindowFilter(resfeat,feat,nd->id);
538       free(feat);//SMPI_SHARED_FREE(feat);
539     }else{
540       featp=(Arr *)tail->feat;
541       feat=newArr(featp->len);
542       memcpy(feat->val,featp->val,featp->len*sizeof(double));
543       resfeat=WindowFilter(resfeat,feat,nd->id);
544       free(feat);//SMPI_SHARED_FREE(feat);
545     }
546   }
547   for(i=0;i<resfeat->len;i++) resfeat->val[i]=((int)resfeat->val[i])/nd->inDegree;
548   nd->feat=resfeat;
549   return nd->feat;
550 }
551
552 static double Reduce(Arr *a,int w){
553   double retv=0.0;
554   if(timer_on){
555     timer_clear(w);
556     timer_start(w);
557   }
558   retv=(int)(w*CheckVal(a));/* The casting needed for node and array dependent verifcation */
559   if(timer_on){
560     timer_stop(w);
561     fprintf(stderr,"** Reduce time in node %d = %f\n",(w-1),timer_read(w));
562   }
563   return retv;
564 }
565
566 static double ReduceStreams(DGraph *dg,DGNode *nd){
567   double csum=0.0;
568   int i=0,len=0,tag=0;
569   DGArc *ar=NULL;
570   DGNode *tail=NULL;
571   Arr *feat=NULL;
572   double retv=0.0;
573
574   TRACE_smpi_set_category ("ReduceStreams");
575
576   for(i=0;i<nd->inDegree;i++){
577     ar=nd->inArc[i];
578     if(ar->head!=nd) continue;
579     tail=ar->tail;
580     if(tail->address!=nd->address){
581       MPI_Status status;
582       len=0;
583       tag=ar->id;
584       MPI_Recv(&len,1,MPI_INT,tail->address,tag,MPI_COMM_WORLD,&status);
585       feat=newArr(len);
586       MPI_Recv(feat->val,feat->len,MPI_DOUBLE,tail->address,tag,MPI_COMM_WORLD,&status);
587       csum+=Reduce(feat,(nd->id+1));
588       free(feat);//SMPI_SHARED_FREE(feat);
589     }else{
590       csum+=Reduce(tail->feat,(nd->id+1));
591     }
592   }
593   if(nd->inDegree>0)csum=(((long long int)csum)/nd->inDegree);
594   retv=(nd->id+1)*csum;
595   return retv;
596 }
597
598 static int ProcessNodes(DGraph *dg,int me){
599   double chksum=0.0;
600   Arr *feat=NULL;
601   int i=0,verified=0,tag;
602   DGNode *nd=NULL;
603   double rchksum=0.0;
604   MPI_Status status;
605
606   TRACE_smpi_set_category ("ProcessNodes");
607
608   for(i=0;i<dg->numNodes;i++){
609     nd=dg->node[i];
610     if(nd->address!=me) continue;
611     if(strstr(nd->name,"Source")){
612       nd->feat=RandomFeatures(dg->name,fielddim,nd->id);
613       SendResults(dg,nd,nd->feat);
614     }else if(strstr(nd->name,"Sink")){
615       chksum=ReduceStreams(dg,nd);
616       tag=dg->numArcs+nd->id; /* make these to avoid clash with arc tags */
617       MPI_Send(&chksum,1,MPI_DOUBLE,0,tag,MPI_COMM_WORLD);
618     }else{
619       feat=CombineStreams(dg,nd);
620       SendResults(dg,nd,feat);
621     }
622   }
623
624   TRACE_smpi_set_category ("ProcessNodes");
625
626   if(me==0){ /* Report node */
627     rchksum=0.0;
628     chksum=0.0;
629     for(i=0;i<dg->numNodes;i++){
630       nd=dg->node[i];
631       if(!strstr(nd->name,"Sink")) continue;
632       tag=dg->numArcs+nd->id; /* make these to avoid clash with arc tags */
633       MPI_Recv(&rchksum,1,MPI_DOUBLE,nd->address,tag,MPI_COMM_WORLD,&status);
634       chksum+=rchksum;
635     }
636     verified=verify(dg->name,chksum);
637   }
638   return verified;
639 }
640
641 int main(int argc,char **argv ){
642   int my_rank,comm_size;
643   int i;
644   DGraph *dg=NULL;
645   int verified=0, featnum=0;
646   double bytes_sent=2.0,tot_time=0.0;
647
648   MPI_Init( &argc, &argv );
649   MPI_Comm_rank( MPI_COMM_WORLD, &my_rank );
650   MPI_Comm_size( MPI_COMM_WORLD, &comm_size );
651
652   TRACE_smpi_set_category ("begin");
653   get_info(argc, argv, &nprocs, &class);
654   check_info(DT, nprocs, class);
655
656   if    (class == 'S') { num_samples=1728; deviation=128; num_sources=4; }
657   else if (class == 'W') { num_samples=1728*8; deviation=128*2; num_sources=4*2; }
658   else if (class == 'A') { num_samples=1728*64; deviation=128*4; num_sources=4*4; }
659   else if (class == 'B') { num_samples=1728*512; deviation=128*8; num_sources=4*8; }
660   else if (class == 'C') { num_samples=1728*4096; deviation=128*16; num_sources=4*16; }
661   else if (class == 'D') { num_samples=1728*4096*8; deviation=128*32; num_sources=4*32; }
662   else {
663     printf("setparams: Internal error: invalid class type %c\n", class);
664     exit(1);
665   }
666
667   if(argc!=4 || (strncmp(argv[3],"BH",2)!=0 && strncmp(argv[3],"WH",2)!=0 && strncmp(argv[3],"SH",2)!=0)){
668     if(my_rank==0){
669     fprintf(stderr,"** Usage: mpirun -np N ../bin/dt.S GraphName\n");
670     fprintf(stderr,"** Where \n   - N is integer number of MPI processes\n");
671     fprintf(stderr,"   - S is the class S, W, or A \n");
672     fprintf(stderr,"   - GraphName is the communication graph name BH, WH, or SH.\n");
673     fprintf(stderr,"   - the number of MPI processes N should not be be less than \n");
674     fprintf(stderr,"   the number of nodes in the graph\n");
675     }
676     MPI_Finalize();
677     exit(0);
678   }
679
680   if(strncmp(argv[3],"BH",2)==0){
681     dg=buildBH(class);
682   }else if(strncmp(argv[3],"WH",2)==0){
683     dg=buildWH(class);
684   }else if(strncmp(argv[3],"SH",2)==0){
685     dg=buildSH(class);
686   }
687
688
689   if(timer_on && dg->numNodes+1>timers_tot){
690     timer_on=0;
691     if(my_rank==0)
692     fprintf(stderr,"Not enough timers. Node timeing is off. \n");
693   }
694   if(dg->numNodes>comm_size){
695     if(my_rank==0){
696     fprintf(stderr,"**  The number of MPI processes should not be less than \n");
697     fprintf(stderr,"**  the number of nodes in the graph\n");
698     fprintf(stderr,"**  Number of MPI processes = %d\n",comm_size);
699     fprintf(stderr,"**  Number nodes in the graph = %d\n",dg->numNodes);
700     }
701     MPI_Finalize();
702     exit(0);
703   }
704   for(i=0;i<dg->numNodes;i++){
705     dg->node[i]->address=i;
706   }
707   if( my_rank == 0 ){
708     printf( "\n\n NAS Parallel Benchmarks 3.3 -- DT Benchmark\n\n" );
709     graphShow(dg,0);
710     timer_clear(0);
711     timer_start(0);
712   }
713   verified=ProcessNodes(dg,my_rank);
714   TRACE_smpi_set_category ("end");
715
716   featnum=num_samples*fielddim;
717   bytes_sent=featnum*dg->numArcs;
718   bytes_sent/=1048576;
719   if(my_rank==0){
720     timer_stop(0);
721     tot_time=timer_read(0);
722     c_print_results( dg->name, class, featnum, 0, 0, dg->numNodes, 0, comm_size, tot_time, bytes_sent/tot_time,
723          "bytes transmitted", verified);
724   }
725   MPI_Finalize();
726   return 1;
727 }