1 /*************************************************************************
3 * N A S P A R A L L E L B E N C H M A R K S 3.3 *
7 *************************************************************************
9 * This benchmark is part of the NAS Parallel Benchmark 3.3 suite. *
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. *
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: *
21 * http: www.nas.nasa.gov/Software/NPB *
23 * Send comments or suggestions to npb@nas.nasa.gov *
24 * Send bug reports to npb-bugs@nas.nasa.gov *
26 * NAS Parallel Benchmarks Group *
27 * NASA Ames Research Center *
29 * Moffett Field, CA 94035-1000 *
31 * E-mail: npb@nas.nasa.gov *
32 * Fax: (650) 604-3957 *
34 *************************************************************************
36 * Author: M. Frumkin * *
38 *************************************************************************/
46 #include "nas_common.h"
47 #include "simgrid/instr.h" //TRACE_
49 int timer_on=0,timers_tot=64;
50 double start[64], elapsed[64];
58 static int verify(char *bmname,double rnm2){
59 double verify_value=0.0;
60 double epsilon=1.0E-8;
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;
71 fprintf(stderr,"No such benchmark as %s.\n",bmname);
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;
82 fprintf(stderr,"No such benchmark as %s.\n",bmname);
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;
93 fprintf(stderr,"No such benchmark as %s.\n",bmname);
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;
104 fprintf(stderr,"No such benchmark as %s.\n",bmname);
107 }else if(class=='C') {
108 if(strstr(bmname,"BH")){
110 }else if(strstr(bmname,"WH")){
112 }else if(strstr(bmname,"SH")){
115 fprintf(stderr,"No such benchmark as %s.\n",bmname);
118 }else if(class=='D') {
119 if(strstr(bmname,"BH")){
121 }else if(strstr(bmname,"WH")){
123 }else if(strstr(bmname,"SH")){
126 fprintf(stderr,"No such benchmark as %s.\n",bmname);
130 fprintf(stderr,"No such class as %c.\n",class);
132 fprintf(stderr," %s L2 Norm = %f\n",bmname,rnm2);
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 */
137 fprintf(stderr," Deviation = %f\n",(rnm2 - verify_value));
140 fprintf(stderr," The correct verification value = %f\n",verify_value);
141 fprintf(stderr," Got value = %f\n",rnm2);
149 static int ipowMod(int a,long long int n,int md){
152 fprintf(stderr,"ipowMod: exponent must be nonnegative exp=%lld\n",n);
157 fprintf(stderr,"ipowMod: module must be positive mod=%d",md);
177 static DGraph *buildSH(const char cls){
178 /* Nodes of the graph must be topologically sorted to avoid MPI deadlock. */
180 int numSources=num_sources; /* must be power of 2 */
181 int numOfLayers=0,tmpS=numSources>>1;
182 int firstLayerNode=0;
185 int mask=0x0,ndid=0,ndoff=0;
189 sprintf(nm,"DT_SH.%c",cls);
196 for(i=0;i<numSources;i++){
197 sprintf(nm,"Source.%d",i);
201 for(j=0;j<numOfLayers;j++){
203 for(i=0;i<numSources;i++){
204 sprintf(nm,"Comparator.%d",(i+j*firstLayerNode));
208 ndid=firstLayerNode+ndoff;
209 ar=newArc(dg->node[ndid],nd);
212 ndid=firstLayerNode+ndoff;
213 ar=newArc(dg->node[ndid],nd);
216 firstLayerNode+=numSources;
218 mask=0x00000001<<numOfLayers;
219 for(i=0;i<numSources;i++){
220 sprintf(nm,"Sink.%d",i);
224 ndid=firstLayerNode+ndoff;
225 ar=newArc(dg->node[ndid],nd);
228 ndid=firstLayerNode+ndoff;
229 ar=newArc(dg->node[ndid],nd);
235 static DGraph *buildWH(const char cls){
236 /* Nodes of the graph must be topologically sorted to avoid MPI deadlock. */
238 int numSources=num_sources,maxInDeg=4;
239 int numLayerNodes=numSources,firstLayerNode=0;
240 int totComparators=0;
241 int numPrevLayerNodes=numLayerNodes;
244 DGNode *nd=NULL,*source=NULL,*tmp=NULL,*snd=NULL;
248 sprintf(nm,"DT_WH.%c",cls);
251 for(i=0;i<numSources;i++){
252 sprintf(nm,"Sink.%d",i);
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);
265 id=AttachNode(dg,nd);
266 for(j=0;j<maxInDeg;j++){
268 if(sid>=numPrevLayerNodes) break;
269 snd=dg->node[firstLayerNode+sid];
270 ar=newArc(dg->node[id],snd);
274 firstLayerNode+=numPrevLayerNodes;
275 numPrevLayerNodes=numLayerNodes;
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);
285 for(i=0;i<dg->numNodes/2;i++){ /* Topological sorting */
287 dg->node[i]=dg->node[dg->numNodes-1-i];
289 dg->node[dg->numNodes-1-i]=tmp;
290 dg->node[dg->numNodes-1-i]->id=dg->numNodes-1-i;
295 static DGraph *buildBH(const char cls){
296 /* Nodes of the graph must be topologically sorted to avoid MPI deadlock.*/
298 int numSources=num_sources,maxInDeg=4;
299 int numLayerNodes=numSources,firstLayerNode=0;
301 DGNode *nd=NULL, *snd=NULL, *sink=NULL;
303 int totComparators=0;
304 int numPrevLayerNodes=numLayerNodes;
308 sprintf(nm,"DT_BH.%c",cls);
311 for(i=0;i<numSources;i++){
312 sprintf(nm,"Source.%d",i);
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);
323 id=AttachNode(dg,nd);
324 for(j=0;j<maxInDeg;j++){
326 if(sid>=numPrevLayerNodes) break;
327 snd=dg->node[firstLayerNode+sid];
328 ar=newArc(snd,dg->node[id]);
332 firstLayerNode+=numPrevLayerNodes;
333 numPrevLayerNodes=numLayerNodes;
335 sink=newNode((char*)"Sink");
337 for(i=0;i<numPrevLayerNodes;i++){
338 nd=dg->node[firstLayerNode+i];
350 static Arr *newArr(int len){
351 Arr *arr=(Arr *)malloc(sizeof(Arr)); //Arr *arr=(Arr *)SMPI_SHARED_MALLOC(sizeof(Arr));
353 arr->val=(double *)malloc(len*sizeof(double)); //arr->val=(double *)SMPI_SHARED_MALLOC(len*sizeof(double));
357 static void arrShow(Arr* a){
358 if(!a) fprintf(stderr,"-- NULL array\n");
360 fprintf(stderr,"-- length=%d\n",a->len);
364 static double CheckVal(Arr *feat){
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 */
372 static int GetFNumDPar(int* mean, int* stdev){
378 static int GetFeatureNum(char *mbname,int id){
379 double tran=314159265.0;
381 double denom=randlc(&tran,&A);
383 int mean=num_samples,stdev=128;
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;
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);
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;
410 feat->val[i+1]=seedy;
411 feat->val[i+2]=seedz;
412 feat->val[i+3]=seedf;
416 fprintf(stderr,"** RandomFeatures time in node %d = %f\n",id,timer_read(id+1));
421 static void Resample(Arr *a,int blen){
422 long long int i=0,j=0,jlo=0,jhi=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));
431 avval=a->val[i]/(jhi-jlo+1);
432 for(j=jlo;j<=jhi;j++){
437 nval[blen-1]=a->val[a->len-1];
438 free(a->val); //SMPI_SHARED_FREE(a->val);
444 static Arr* WindowFilter(Arr *a, Arr* b,int w){
446 double rms0=0.0,rms1=0.0,rmsm1=0.0;
447 double weight=((double) (w+1))/(w+2);
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]);
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]);
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]);
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];
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)*/
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];
493 fprintf(stderr,"** WindowFilter time in node %d = %f\n",(w-1),timer_read(w));
498 static int SendResults(DGraph *dg,DGNode *nd,Arr *feat){
503 TRACE_smpi_set_category ("SendResults");
504 for(i=0;i<nd->outDegree;i++){
506 if(ar->tail!=nd) continue;
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);
514 TRACE_smpi_set_category (NULL);
518 static Arr* CombineStreams(DGraph *dg,DGNode *nd){
519 Arr *resfeat=newArr(num_samples*fielddim);
524 Arr *feat=NULL,*featp=NULL;
526 if(nd->inDegree==0) return NULL;
527 for(i=0;i<nd->inDegree;i++){
529 if(ar->head!=nd) continue;
531 if(tail->address!=nd->address){
534 MPI_Recv(&len,1,MPI_INT,tail->address,tag,MPI_COMM_WORLD,&status);
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);
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);
547 for(i=0;i<resfeat->len;i++) resfeat->val[i]=((int)resfeat->val[i])/nd->inDegree;
552 static double Reduce(Arr *a,int w){
558 retv=(int)(w*CheckVal(a));/* The casting needed for node and array dependent verifcation */
561 fprintf(stderr,"** Reduce time in node %d = %f\n",(w-1),timer_read(w));
566 static double ReduceStreams(DGraph *dg,DGNode *nd){
574 TRACE_smpi_set_category ("ReduceStreams");
576 for(i=0;i<nd->inDegree;i++){
578 if(ar->head!=nd) continue;
580 if(tail->address!=nd->address){
584 MPI_Recv(&len,1,MPI_INT,tail->address,tag,MPI_COMM_WORLD,&status);
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);
590 csum+=Reduce(tail->feat,(nd->id+1));
593 if(nd->inDegree>0)csum=(((long long int)csum)/nd->inDegree);
594 retv=(nd->id+1)*csum;
598 static int ProcessNodes(DGraph *dg,int me){
601 int i=0,verified=0,tag;
606 TRACE_smpi_set_category ("ProcessNodes");
608 for(i=0;i<dg->numNodes;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);
619 feat=CombineStreams(dg,nd);
620 SendResults(dg,nd,feat);
624 TRACE_smpi_set_category ("ProcessNodes");
626 if(me==0){ /* Report node */
629 for(i=0;i<dg->numNodes;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);
636 verified=verify(dg->name,chksum);
641 int main(int argc,char **argv ){
642 int my_rank,comm_size;
645 int verified=0, featnum=0;
646 double bytes_sent=2.0,tot_time=0.0;
648 MPI_Init( &argc, &argv );
649 MPI_Comm_rank( MPI_COMM_WORLD, &my_rank );
650 MPI_Comm_size( MPI_COMM_WORLD, &comm_size );
652 TRACE_smpi_set_category ("begin");
653 get_info(argc, argv, &nprocs, &class);
654 check_info(DT, nprocs, class);
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; }
663 printf("setparams: Internal error: invalid class type %c\n", class);
667 if(argc!=4 || (strncmp(argv[3],"BH",2)!=0 && strncmp(argv[3],"WH",2)!=0 && strncmp(argv[3],"SH",2)!=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");
680 if(strncmp(argv[3],"BH",2)==0){
682 }else if(strncmp(argv[3],"WH",2)==0){
684 }else if(strncmp(argv[3],"SH",2)==0){
688 if(timer_on && dg->numNodes+1>timers_tot){
691 fprintf(stderr,"Not enough timers. Node timeing is off. \n");
693 if(dg->numNodes>comm_size){
695 fprintf(stderr,"** The number of MPI processes should not be less than \n");
696 fprintf(stderr,"** the number of nodes in the graph\n");
697 fprintf(stderr,"** Number of MPI processes = %d\n",comm_size);
698 fprintf(stderr,"** Number nodes in the graph = %d\n",dg->numNodes);
703 for(i=0;i<dg->numNodes;i++){
704 dg->node[i]->address=i;
707 printf( "\n\n NAS Parallel Benchmarks 3.3 -- DT Benchmark\n\n" );
712 verified=ProcessNodes(dg,my_rank);
713 TRACE_smpi_set_category ("end");
715 featnum=num_samples*fielddim;
716 bytes_sent=featnum*dg->numArcs;
720 tot_time=timer_read(0);
721 c_print_results( dg->name, class, featnum, 0, 0, dg->numNodes, 0, comm_size, tot_time, bytes_sent/tot_time,
722 "bytes transmitted", verified);