Logo AND Algorithmique Numérique Distribuée

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