Logo AND Algorithmique Numérique Distribuée

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