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' || class == 'D') {
110 fprintf(stderr,"No such class as %c.\n",class);
112 fprintf(stderr," %s L2 Norm = %f\n",bmname,rnm2);
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 */
117 fprintf(stderr," Deviation = %f\n",(rnm2 - verify_value));
120 fprintf(stderr," The correct verification value = %f\n",verify_value);
121 fprintf(stderr," Got value = %f\n",rnm2);
129 static int ipowMod(int a,long long int n,int md){
133 fprintf(stderr,"ipowMod: exponent must be nonnegative exp=%lld\n",n);
134 exp=-n; /* temp fix */
138 fprintf(stderr,"ipowMod: module must be positive mod=%d",md);
159 static DGraph *buildSH(const char cls){
160 /* Nodes of the graph must be topologically sorted to avoid MPI deadlock. */
162 unsigned int numSources=num_sources; /* must be power of 2 */
164 unsigned int tmpS=numSources>>1;
165 int firstLayerNode=0;
168 unsigned int mask=0x0;
170 unsigned int i=0,j=0;
173 snprintf(nm,BLOCK_SIZE,"DT_SH.%c",cls);
180 for(i=0;i<numSources;i++){
181 snprintf(nm,BLOCK_SIZE,"Source.%d",i);
185 for(j=0;j<numOfLayers;j++){
187 for(i=0;i<numSources;i++){
188 snprintf(nm,BLOCK_SIZE,"Comparator.%d",(i+j*firstLayerNode));
192 ndid=firstLayerNode+ndoff;
193 ar=newArc(dg->node[ndid],nd);
196 ndid=firstLayerNode+ndoff;
197 ar=newArc(dg->node[ndid],nd);
200 firstLayerNode+=numSources;
202 mask=0x00000001<<numOfLayers;
203 for(i=0;i<numSources;i++){
204 snprintf(nm,BLOCK_SIZE,"Sink.%d",i);
208 ndid=firstLayerNode+ndoff;
209 ar=newArc(dg->node[ndid],nd);
212 ndid=firstLayerNode+ndoff;
213 ar=newArc(dg->node[ndid],nd);
219 static DGraph *buildWH(const char cls){
220 /* Nodes of the graph must be topologically sorted to avoid MPI deadlock. */
222 int numSources=num_sources,maxInDeg=4;
223 int numLayerNodes=numSources,firstLayerNode=0;
224 int totComparators=0;
225 int numPrevLayerNodes=numLayerNodes;
228 DGNode *nd=NULL,*source=NULL,*tmp=NULL,*snd=NULL;
232 snprintf(nm,BLOCK_SIZE,"DT_WH.%c",cls);
235 for(i=0;i<numSources;i++){
236 snprintf(nm,BLOCK_SIZE,"Sink.%d",i);
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 snprintf(nm,BLOCK_SIZE,"Comparator.%d",totComparators);
249 id=AttachNode(dg,nd);
250 for(j=0;j<maxInDeg;j++){
252 if(sid>=numPrevLayerNodes) break;
253 snd=dg->node[firstLayerNode+sid];
254 ar=newArc(dg->node[id],snd);
258 firstLayerNode+=numPrevLayerNodes;
259 numPrevLayerNodes=numLayerNodes;
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);
269 for(i=0;i<dg->numNodes/2;i++){ /* Topological sorting */
271 dg->node[i]=dg->node[dg->numNodes-1-i];
273 dg->node[dg->numNodes-1-i]=tmp;
274 dg->node[dg->numNodes-1-i]->id=dg->numNodes-1-i;
279 static DGraph *buildBH(const char cls){
280 /* Nodes of the graph must be topologically sorted to avoid MPI deadlock.*/
282 int numSources=num_sources,maxInDeg=4;
283 int numLayerNodes=numSources,firstLayerNode=0;
285 DGNode *nd=NULL, *snd=NULL, *sink=NULL;
287 int totComparators=0;
288 int numPrevLayerNodes=numLayerNodes;
292 snprintf(nm,BLOCK_SIZE,"DT_BH.%c",cls);
295 for(i=0;i<numSources;i++){
296 snprintf(nm,BLOCK_SIZE,"Source.%d",i);
300 while(numLayerNodes>maxInDeg){
301 numLayerNodes=numLayerNodes/maxInDeg;
302 if(numLayerNodes*maxInDeg<numPrevLayerNodes)numLayerNodes++;
303 for(i=0;i<numLayerNodes;i++){
304 snprintf(nm,BLOCK_SIZE,"Comparator.%d",totComparators);
307 id=AttachNode(dg,nd);
308 for(j=0;j<maxInDeg;j++){
310 if(sid>=numPrevLayerNodes) break;
311 snd=dg->node[firstLayerNode+sid];
312 ar=newArc(snd,dg->node[id]);
316 firstLayerNode+=numPrevLayerNodes;
317 numPrevLayerNodes=numLayerNodes;
319 sink=newNode((char*)"Sink");
321 for(i=0;i<numPrevLayerNodes;i++){
322 nd=dg->node[firstLayerNode+i];
334 static Arr *newArr(int len){
335 Arr *arr=(Arr *)malloc(sizeof(Arr)); //Arr *arr=(Arr *)SMPI_SHARED_MALLOC(sizeof(Arr));
337 arr->val=(double *)malloc(len*sizeof(double)); //arr->val=(double *)SMPI_SHARED_MALLOC(len*sizeof(double));
341 static double CheckVal(Arr *feat){
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 */
349 static int GetFNumDPar(int* mean, int* stdev){
355 static int GetFeatureNum(char *mbname,int id){
356 double tran=314159265.0;
358 double denom=randlc(&tran,&A);
360 int mean=num_samples,stdev=128;
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;
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);
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;
387 feat->val[i+1]=seedy;
388 feat->val[i+2]=seedz;
389 feat->val[i+3]=seedf;
393 fprintf(stderr,"** RandomFeatures time in node %d = %f\n",id,timer_read(id+1));
398 static void Resample(Arr *a,int blen){
399 long long int i=0,j=0,jlo=0,jhi=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));
408 avval=a->val[i]/(jhi-jlo+1);
409 for(j=jlo;j<=jhi;j++){
414 nval[blen-1]=a->val[a->len-1];
415 free(a->val); //SMPI_SHARED_FREE(a->val);
421 static Arr* WindowFilter(Arr *a, Arr* b,int w){
423 double rms0=0.0,rms1=0.0,rmsm1=0.0;
424 double weight=((double) (w+1))/(w+2);
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]);
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]);
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]);
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];
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)*/
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];
470 fprintf(stderr,"** WindowFilter time in node %d = %f\n",(w-1),timer_read(w));
475 static int SendResults(DGraph *dg,DGNode *nd,Arr *feat){
480 TRACE_smpi_set_category ("SendResults");
481 for(i=0;i<nd->outDegree;i++){
483 if(ar->tail!=nd) continue;
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);
491 TRACE_smpi_set_category (NULL);
495 static Arr* CombineStreams(DGraph *dg,DGNode *nd){
496 Arr *resfeat=newArr(num_samples*fielddim);
501 Arr *feat=NULL,*featp=NULL;
503 if(nd->inDegree==0) return NULL;
504 for(i=0;i<nd->inDegree;i++){
506 if(ar->head!=nd) continue;
508 if(tail->address!=nd->address){
511 MPI_Recv(&len,1,MPI_INT,tail->address,tag,MPI_COMM_WORLD,&status);
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);
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);
524 for(i=0;i<resfeat->len;i++) resfeat->val[i]=((int)resfeat->val[i])/nd->inDegree;
529 static double Reduce(Arr *a,int w){
535 retv=(int)(w*CheckVal(a));/* The casting needed for node and array dependent verifcation */
538 fprintf(stderr,"** Reduce time in node %d = %f\n",(w-1),timer_read(w));
543 static double ReduceStreams(DGraph *dg,DGNode *nd){
551 TRACE_smpi_set_category ("ReduceStreams");
553 for(i=0;i<nd->inDegree;i++){
555 if(ar->head!=nd) continue;
557 if(tail->address!=nd->address){
561 MPI_Recv(&len,1,MPI_INT,tail->address,tag,MPI_COMM_WORLD,&status);
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);
567 csum+=Reduce(tail->feat,(nd->id+1));
570 if(nd->inDegree>0)csum=(((long long int)csum)/nd->inDegree);
571 retv=(nd->id+1)*csum;
575 static int ProcessNodes(DGraph *dg,int me){
578 int i=0,verified=0,tag;
583 TRACE_smpi_set_category ("ProcessNodes");
585 for(i=0;i<dg->numNodes;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);
596 feat=CombineStreams(dg,nd);
597 SendResults(dg,nd,feat);
601 TRACE_smpi_set_category ("ProcessNodes");
603 if(me==0){ /* Report node */
606 for(i=0;i<dg->numNodes;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);
613 verified=verify(dg->name,chksum);
618 int main(int argc,char **argv ){
619 int my_rank,comm_size;
622 int verified=0, featnum=0;
623 double bytes_sent=2.0,tot_time=0.0;
625 MPI_Init( &argc, &argv );
626 MPI_Comm_rank( MPI_COMM_WORLD, &my_rank );
627 MPI_Comm_size( MPI_COMM_WORLD, &comm_size );
629 TRACE_smpi_set_category ("begin");
630 get_info(argc, argv, &nprocs, &class);
631 check_info(DT, nprocs, class);
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; }
640 printf("setparams: Internal error: invalid class type %c\n", class);
644 if(argc!=4 || (strncmp(argv[3],"BH",2)!=0 && strncmp(argv[3],"WH",2)!=0 && strncmp(argv[3],"SH",2)!=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");
657 if(strncmp(argv[3],"BH",2)==0){
659 }else if(strncmp(argv[3],"WH",2)==0){
661 }else if(strncmp(argv[3],"SH",2)==0){
665 if(timer_on && dg->numNodes+1>timers_tot){
668 fprintf(stderr,"Not enough timers. Node timeing is off. \n");
670 if(dg->numNodes>comm_size){
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);
680 for(i=0;i<dg->numNodes;i++){
681 dg->node[i]->address=i;
684 printf( "\n\n NAS Parallel Benchmarks 3.3 -- DT Benchmark\n\n" );
689 verified=ProcessNodes(dg,my_rank);
690 TRACE_smpi_set_category ("end");
692 featnum=num_samples*fielddim;
693 bytes_sent=featnum*dg->numArcs;
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);