1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
2  * Copyright by The HDF Group.                                               *
3  * Copyright by the Board of Trustees of the University of Illinois.         *
4  * All rights reserved.                                                      *
5  *                                                                           *
6  * This file is part of HDF5.  The full HDF5 copyright notice, including     *
7  * terms governing use, modification, and redistribution, is contained in    *
8  * the COPYING file, which can be found at the root of the source code       *
9  * distribution tree, or in https://support.hdfgroup.org/ftp/HDF5/releases.  *
10  * If you do not have access to either file, you may request a copy from     *
11  * help@hdfgroup.org.                                                        *
12  * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
13 
14 /*
15  * Programmer:	Raymond Lu <slu@ncsa.uiuc.edu>
16  *		Friday, Oct 3, 2004
17  *
18  * Purpose:	Tests performance of metadata
19  */
20 
21 #include "h5test.h"
22 
23 #ifdef H5_HAVE_PARALLEL
24 #define MAINPROCESS	(!mpi_rank)	/* define process 0 as main process */
25 #endif /*H5_HAVE_PARALLEL*/
26 
27 /* File_Access_type bits */
28 #define FACC_DEFAULT	0x0	/* serial as default */
29 #define FACC_MPIO	0x1	/* MPIO */
30 
31 /* Which test to run */
32 int RUN_TEST = 0x0;     /* all tests as default */
33 int TEST_1   = 0x1;     /* Test 1 */
34 int TEST_2   = 0x2;     /* Test 2 */
35 int TEST_3   = 0x4;     /* Test 3 */
36 
37 
38 const char *FILENAME[] = {
39     "meta_perf_1",
40     "meta_perf_2",
41     "meta_perf_3",
42     NULL
43 };
44 
45 /* Default values for performance. Can be changed through command line options */
46 int 	NUM_DSETS = 16;
47 int 	NUM_ATTRS = 8;
48 int 	BATCH_ATTRS = 2;
49 hbool_t flush_dset = FALSE;
50 hbool_t flush_attr = FALSE;
51 int 	nerrors = 0;			/* errors count */
52 hid_t	fapl;
53 
54 /* Data space IDs */
55 hid_t	space;
56 hid_t	small_space;
57 
58 /* Performance data */
59 typedef struct p_time {
60     double total;
61     double avg;
62     double max;
63     double min;
64     double start;
65     char   func[32];
66 } p_time;
67 
68 /*Test file access type for parallel.  MPIO as default */
69 int facc_type = FACC_DEFAULT;
70 
71 double  retrieve_time(void);
72 void    perf(p_time *perf_t, double start_t, double end_t);
73 void    print_perf(p_time, p_time, p_time);
74 
75 
76 /*-------------------------------------------------------------------------
77  * Function:	parse_options
78  *
79   Purpose:	Parse command line options
80  *
81  * Programmer:	Raymond Lu
82  *		Friday, Oct 3, 2003
83  *
84  * Modifications:
85  *
86  *-------------------------------------------------------------------------
87  */
88 static int
parse_options(int argc,char ** argv)89 parse_options(int argc, char **argv)
90 {
91     int t;
92 
93     /* Use default values */
94     if(argc==1)
95 	return(0);
96 
97     while (--argc){
98 	if (**(++argv) != '-'){
99 	    break;
100 	}else{
101 	    switch(*(*argv+1)){
102                 case 'h':   /* Help page */
103                             return(1);
104 
105 		case 'd':   /* Number of datasets */
106                             NUM_DSETS = atoi((*argv+1)+1);
107 			    if (NUM_DSETS < 0){
108 				nerrors++;
109 				return(1);
110 			    }
111 			    break;
112 
113 		case 'a':   /* Number of attributes per dataset */
114                             NUM_ATTRS = atoi((*argv+1)+1);
115 			    if (NUM_ATTRS < 0){
116 				nerrors++;
117 				return(1);
118 			    }
119 			    break;
120 
121 		case 'n':   /* Number of attributes to be created in batch */
122                             BATCH_ATTRS = atoi((*argv+1)+1);
123 			    if (BATCH_ATTRS < 0){
124 				nerrors++;
125 				return(1);
126 			    }
127 			    break;
128 
129 		case 'm':   /* Use the MPI-IO driver */
130 			    facc_type = FACC_MPIO;
131 			    break;
132 
133                 case 'f':   /* Call H5Fflush for each dataset or attribute */
134                             if(!strcmp("a", (*argv+2)))
135                                 flush_attr = TRUE;
136                             else if(!strcmp("d", (*argv+2)))
137                                 flush_dset = TRUE;
138                             else {
139                                 nerrors++;
140                                 return(1);
141                             }
142                             break;
143 
144 		case 't':   /* Which test to run */
145                             t = atoi((*argv+1)+1);
146 			    if (t < 1 || t > 3){
147 				nerrors++;
148 				return(1);
149 			    }
150                             if(t == 1)
151                                 RUN_TEST |= TEST_1;
152                             else if(t == 2)
153                                 RUN_TEST |= TEST_2;
154                             else
155                                 RUN_TEST |= TEST_3;
156 
157 			    break;
158 
159  		default:    nerrors++;
160 			    return(1);
161 	    }
162 	}
163     } /*while*/
164 
165     /* Check valid values */
166 #ifndef H5_HAVE_PARALLEL
167     if(facc_type == FACC_MPIO)
168     {
169         nerrors++;
170         return(1);
171     }
172 #endif /*H5_HAVE_PARALLEL*/
173 
174     if(NUM_ATTRS && !BATCH_ATTRS)
175         NUM_ATTRS = 0;
176 
177     if(!NUM_ATTRS && BATCH_ATTRS)
178         BATCH_ATTRS = 0;
179 
180     if(!NUM_DSETS) {
181         nerrors++;
182         return(1);
183     }
184 
185     if(NUM_ATTRS && BATCH_ATTRS) {
186         if(BATCH_ATTRS > NUM_ATTRS || NUM_ATTRS % BATCH_ATTRS) {
187 	    nerrors++;
188             return(1);
189         }
190     }
191 
192     return(0);
193 }
194 
195 
196 /*-------------------------------------------------------------------------
197  * Function:	usage
198  *
199   Purpose:	Prints help page
200  *
201  * Programmer:	Raymond Lu
202  *		Friday, Oct 3, 2003
203  *
204  * Modifications:
205  *
206  *-------------------------------------------------------------------------
207  */
208 static void
usage(void)209 usage(void)
210 {
211     printf("Usage: perf_meta [-h] [-m] [-d<num_datasets>]"
212            "[-a<num_attributes>]\n"
213            "\t[-n<batch_attributes>] [-f<option>] [-t<test>]\n");
214     printf("\t-h"
215 	"\t\t\thelp page.\n");
216     printf("\t-m"
217 	"\t\t\tset MPIO as the file driver when parallel HDF5\n"
218         "\t\t\t\tis enabled.  -m must be specified\n"
219         "\t\t\t\twhen running parallel program.\n");
220     printf("\t-d<num_datasets>"
221 	"\tset number of datasets for meta data \n"
222         "\t\t\t\tperformance test\n");
223     printf("\t-a<num_attributes>"
224         "\tset number of attributes per dataset for meta \n"
225         "\t\t\t\tdata performance test.\n");
226     printf("\t-n<batch_attributes>"
227 	"\tset batch number of attributes for dataset \n"
228         "\t\t\t\tfor meta data performance test.\n");
229     printf("\t-f<option>"
230 	"\t\tflush data to disk after closing a dataset \n"
231         "\t\t\t\tor attribute.  Valid options are \"d\" for \n"
232         "\t\t\t\tdataset, \"a\" for attribute.  Disabled is \n"
233         "\t\t\t\tthe default.\n");
234     printf("\t-t<tests>"
235 	"\t\trun specific test.  Give only one number each \n"
236         "\t\t\t\ttime. i.e. \"-t1 -t3\" will run test 1 and 3. \n"
237         "\t\t\t\tDefault is all three tests.  The 3 tests are: \n\n"
238         "\t\t\t\t1. Create <num_attributes> attributes for each \n"
239         "\t\t\t\t   of <num_datasets> existing datasets.\n"
240         "\t\t\t\t2. Create <num_attributes> attributes for each \n"
241         "\t\t\t\t   of <num_datasets> new datasets.\n"
242         "\t\t\t\t3. Create <batch_attributes> attributes for \n"
243         "\t\t\t\t   each of <num_dataset> new datasets for \n"
244         "\t\t\t\t   <num_attributes>/<batch_attributes> times.\n");
245 }
246 
247 
248 /*-------------------------------------------------------------------------
249  * Function:	create_dspace
250  *
251  * Purpose:	Attempts to create data space.
252  *
253  * Return:	Success:	0
254  *
255  *		Failure:	-1
256  *
257  * Programmer:	Raymond Lu
258  *		Friday, Oct 3, 2003
259  *
260  * Modifications:
261  *
262  *-------------------------------------------------------------------------
263  */
264 static herr_t
create_dspace(void)265 create_dspace(void)
266 {
267     hsize_t	dims[2];
268     hsize_t	small_dims[2];
269 
270     /* Create the data space */
271     dims[0] = 256;
272     dims[1] = 512;
273     if((space = H5Screate_simple(2, dims, NULL)) < 0)
274 	    goto error;
275 
276     /* Create a small data space for attributes */
277     small_dims[0] = 16;
278     small_dims[1] = 8;
279     if((small_space = H5Screate_simple(2, small_dims, NULL)) < 0)
280 	    goto error;
281 
282     return 0;
283 
284 error:
285     return -1;
286 }
287 
288 
289 /*-------------------------------------------------------------------------
290  * Function:	create_dsets
291  *
292  * Purpose:	Attempts to create some datasets.
293  *
294  * Return:	Success:	0
295  *
296  *		Failure:	-1
297  *
298  * Programmer:	Raymond Lu
299  *		Friday, Oct 3, 2003
300  *
301  * Modifications:
302  *
303  *-------------------------------------------------------------------------
304  */
305 static herr_t
create_dsets(hid_t file)306 create_dsets(hid_t file)
307 {
308     hid_t	dataset;
309     char	dset_name[32];
310     int		i;
311 
312     /*
313      * Create a dataset using the default dataset creation properties.
314      */
315     for(i = 0; i < NUM_DSETS; i++) {
316 	HDsprintf(dset_name, "dataset %d", i);
317     	if((dataset = H5Dcreate2(file, dset_name, H5T_NATIVE_DOUBLE, space,
318                 H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)) < 0)
319             goto error;
320 
321     	if(H5Dclose(dataset) < 0)
322             goto error;
323     } /* end for */
324 
325     return 0;
326 
327 error:
328     return -1;
329 
330 }
331 
332 
333 /*-------------------------------------------------------------------------
334  * Function:	create_attrs_1
335  *
336  * Purpose:	Attempts to create all attributes for each existing dataset.
337  *
338  * Return:	Success:	0
339  *
340  *		Failure:	-1
341  *
342  * Programmer:	Raymond Lu
343  *		Friday, Oct 3, 2003
344  *
345  * Modifications:
346  *
347  *-------------------------------------------------------------------------
348  */
349 static herr_t
create_attrs_1(void)350 create_attrs_1(void)
351 {
352     hid_t	file, dataset, attr;
353     char	filename[128];
354     char	dset_name[64];
355     char	attr_name[128];
356     int		i, j;
357     p_time      attr_t  = {0, 0, 0, 1000000, 0, ""};
358     p_time      open_t  = {0, 0, 0, 1000000, 0, "H5Dopen2"};
359     p_time      close_t = {0, 0, 0, 1000000, 0, ""};
360 
361 #ifdef H5_HAVE_PARALLEL
362     /* need the rank for printing data */
363     int         mpi_rank;
364     if(facc_type == FACC_MPIO)
365         MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
366 #endif /*H5_HAVE_PARALLEL*/
367 
368     h5_fixname(FILENAME[0], fapl, filename, sizeof filename);
369 
370     if ((file=H5Fcreate(filename, H5F_ACC_TRUNC, H5P_DEFAULT,
371 	fapl)) < 0)
372 	goto error;
373 
374     if(create_dsets(file) < 0)
375 	goto error;
376 
377     /*
378      * Create all(user specifies the number) attributes for each dataset
379      */
380     for(i = 0; i < NUM_DSETS; i++) {
381 	HDsprintf(dset_name, "dataset %d", i);
382         open_t.start = retrieve_time();
383 	if((dataset = H5Dopen2(file, dset_name, H5P_DEFAULT)) < 0)
384 		goto error;
385 	perf(&open_t, open_t.start, retrieve_time());
386 
387 	for(j = 0; j < NUM_ATTRS; j++) {
388             HDsprintf(attr_name, "all attrs for each dset %d", j);
389             attr_t.start = retrieve_time();
390             if((attr = H5Acreate2(dataset, attr_name, H5T_NATIVE_DOUBLE,
391                     small_space, H5P_DEFAULT, H5P_DEFAULT)) < 0)
392                 goto error;
393             if(H5Aclose(attr) < 0)
394                 goto error;
395             perf(&attr_t, attr_t.start, retrieve_time());
396             if(flush_attr && H5Fflush(file, H5F_SCOPE_LOCAL) < 0)
397                 goto error;
398     	} /* end for */
399 
400 	close_t.start = retrieve_time();
401     	if(H5Dclose(dataset) < 0)
402             goto error;
403 	perf(&close_t, close_t.start, retrieve_time());
404         if(flush_dset && H5Fflush(file,  H5F_SCOPE_LOCAL) < 0)
405             goto error;
406     } /* end for */
407 
408     if(facc_type == FACC_MPIO) {
409 #ifdef H5_HAVE_PARALLEL
410         MPI_Barrier(MPI_COMM_WORLD);
411 #endif /*H5_HAVE_PARALLEL*/
412     }
413 
414 #ifdef H5_HAVE_PARALLEL
415     if (facc_type == FACC_DEFAULT || (facc_type != FACC_DEFAULT && MAINPROCESS)) /* only process 0 reports */
416 #endif /*H5_HAVE_PARALLEL*/
417     {
418         /* Calculate the average time */
419         open_t.avg  = open_t.total / NUM_DSETS;
420         close_t.avg = close_t.total / NUM_DSETS;
421         if(NUM_ATTRS)
422             attr_t.avg  = attr_t.total / (NUM_ATTRS*NUM_DSETS);
423 
424         /* Print out the performance result */
425         HDfprintf(stderr, "1.  Create %d attributes for each of %d existing datasets\n",
426             NUM_ATTRS, NUM_DSETS);
427         print_perf(open_t, close_t, attr_t);
428     }
429 
430     if (H5Fclose(file) < 0) goto error;
431 
432     return 0;
433 
434 error:
435     return -1;
436 }
437 
438 
439 /*-------------------------------------------------------------------------
440  * Function:	create_attrs_2
441  *
442  * Purpose:	Attempts to create all attributes for each new dataset.
443  *
444  * Return:	Success:	0
445  *
446  *		Failure:	-1
447  *
448  * Programmer:	Raymond Lu
449  *		Friday, Oct 3, 2003
450  *
451  * Modifications:
452  *
453  *-------------------------------------------------------------------------
454  */
455 static herr_t
create_attrs_2(void)456 create_attrs_2(void)
457 {
458     hid_t	file, dataset, attr;
459     char	filename[128];
460     char	dset_name[64];
461     char	attr_name[128];
462     int		i, j;
463     p_time      attr_t  = {0, 0, 0, 1000000, 0, ""};
464     p_time      create_t  = {0, 0, 0, 1000000, 0, "H5Dcreate2"};
465     p_time      close_t = {0, 0, 0, 1000000, 0, ""};
466 
467 #ifdef H5_HAVE_PARALLEL
468     /* need the rank for printing data */
469     int         mpi_rank;
470     if(facc_type == FACC_MPIO)
471         MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
472 #endif /*H5_HAVE_PARALLEL*/
473 
474     h5_fixname(FILENAME[1], fapl, filename, sizeof filename);
475 
476     if ((file=H5Fcreate(filename, H5F_ACC_TRUNC, H5P_DEFAULT, fapl)) < 0)
477 	goto error;
478 
479     /*
480      * Create all(user specifies the number) attributes for each new dataset
481      */
482     for(i = 0; i < NUM_DSETS; i++) {
483 	HDsprintf(dset_name, "dataset %d", i);
484         create_t.start = retrieve_time();
485    	if((dataset = H5Dcreate2(file, dset_name, H5T_NATIVE_DOUBLE,
486                 space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)) < 0)
487             goto error;
488 	perf(&create_t, create_t.start, retrieve_time());
489 
490 	for(j = 0; j < NUM_ATTRS; j++) {
491             HDsprintf(attr_name, "all attrs for each dset %d", j);
492             attr_t.start = retrieve_time();
493             if((attr = H5Acreate2(dataset, attr_name, H5T_NATIVE_DOUBLE,
494                     small_space, H5P_DEFAULT, H5P_DEFAULT)) < 0)
495                 goto error;
496             if(H5Aclose(attr) < 0)
497                 goto error;
498             perf(&attr_t, attr_t.start, retrieve_time());
499             if(flush_attr && H5Fflush(file,  H5F_SCOPE_LOCAL) < 0)
500                 goto error;
501 	} /* end for */
502 
503 	close_t.start = retrieve_time();
504     	if(H5Dclose(dataset) < 0)
505             goto error;
506 	perf(&close_t, close_t.start, retrieve_time());
507         if(flush_dset && H5Fflush(file,  H5F_SCOPE_LOCAL) < 0)
508             goto error;
509     } /* end for */
510 
511 #ifdef H5_HAVE_PARALLEL
512     if(facc_type == FACC_MPIO)
513         MPI_Barrier(MPI_COMM_WORLD);
514 #endif /*H5_HAVE_PARALLEL*/
515 
516 #ifdef H5_HAVE_PARALLEL
517     /* only process 0 reports if parallel */
518     if (facc_type == FACC_DEFAULT || (facc_type != FACC_DEFAULT && MAINPROCESS))
519 #endif /*H5_HAVE_PARALLEL*/
520     {
521         /* Calculate the average time */
522         create_t.avg = create_t.total / NUM_DSETS;
523         close_t.avg  = close_t.total / NUM_DSETS;
524         if(NUM_ATTRS)
525             attr_t.avg = attr_t.total / (NUM_ATTRS*NUM_DSETS);
526 
527         /* Print out the performance result */
528         HDfprintf(stderr, "2.  Create %d attributes for each of %d new datasets\n",
529             NUM_ATTRS, NUM_DSETS);
530         print_perf(create_t, close_t, attr_t);
531     }
532 
533     if (H5Fclose(file) < 0) goto error;
534 
535     return 0;
536 
537 error:
538     return -1;
539 }
540 
541 
542 /*-------------------------------------------------------------------------
543  * Function:	create_attrs_3
544  *
545  * Purpose:	Attempts to create some attributes for each dataset in a
546  * 		loop.
547  *
548  * Return:	Success:	0
549  *
550  *		Failure:	-1
551  *
552  * Programmer:	Raymond Lu
553  *		Friday, Oct 3, 2003
554  *
555  * Modifications:
556  *
557  *-------------------------------------------------------------------------
558  */
559 static herr_t
create_attrs_3(void)560 create_attrs_3(void)
561 {
562     hid_t	file, dataset, attr;
563     char	filename[128];
564     char	dset_name[64];
565     char	attr_name[128];
566     int		loop_num;
567     int		i, j, k;
568     p_time      attr_t  = {0, 0, 0, 1000000, 0, ""};
569     p_time      open_t  = {0, 0, 0, 1000000, 0, "H5Dopen2"};
570     p_time      close_t = {0, 0, 0, 1000000, 0, ""};
571 
572 #ifdef H5_HAVE_PARALLEL
573     /* need the rank for printing data */
574     int         mpi_rank;
575     if(facc_type == FACC_MPIO)
576         MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
577 #endif /*H5_HAVE_PARALLEL*/
578 
579     h5_fixname(FILENAME[2], fapl, filename, sizeof filename);
580 
581     if ((file=H5Fcreate(filename, H5F_ACC_TRUNC, H5P_DEFAULT,
582 	fapl)) < 0)
583 	goto error;
584 
585     if(create_dsets(file) < 0)
586 	goto error;
587 
588     /*
589      * Create some(user specifies the number) attributes for each dataset
590      * in a loop
591      */
592     loop_num = NUM_ATTRS/BATCH_ATTRS;
593 
594     for(i = 0; i < loop_num; i++) {
595     	for(j = 0; j < NUM_DSETS; j++) {
596             HDsprintf(dset_name, "dataset %d", j);
597             open_t.start = retrieve_time();
598             if((dataset = H5Dopen2(file, dset_name, H5P_DEFAULT)) < 0)
599                 goto error;
600             perf(&open_t, open_t.start, retrieve_time());
601 
602             for(k = 0; k < BATCH_ATTRS; k++) {
603                 HDsprintf(attr_name, "some attrs for each dset %d %d", i, k);
604                 attr_t.start = retrieve_time();
605                 if((attr = H5Acreate2(dataset, attr_name, H5T_NATIVE_DOUBLE,
606                         small_space, H5P_DEFAULT, H5P_DEFAULT)) < 0)
607                     goto error;
608                 if(H5Aclose(attr) < 0)
609                     goto error;
610                 perf(&attr_t, attr_t.start, retrieve_time());
611                 if(flush_attr && H5Fflush(file,  H5F_SCOPE_LOCAL) < 0)
612                     goto error;
613             } /* end for */
614 
615             close_t.start = retrieve_time();
616             if(H5Dclose(dataset) < 0)
617                 goto error;
618             perf(&close_t, close_t.start, retrieve_time());
619             if(flush_dset && H5Fflush(file,  H5F_SCOPE_LOCAL) < 0)
620                 goto error;
621     	} /* end for */
622     } /* end for */
623 
624 #ifdef H5_HAVE_PARALLEL
625     if(facc_type == FACC_MPIO)
626         MPI_Barrier(MPI_COMM_WORLD);
627 #endif /*H5_HAVE_PARALLEL*/
628 
629 #ifdef H5_HAVE_PARALLEL
630     /* only process 0 reports if parallel */
631     if (facc_type == FACC_DEFAULT || (facc_type != FACC_DEFAULT && MAINPROCESS))
632 #endif /*H5_HAVE_PARALLEL*/
633     {
634         /* Calculate the average time */
635         open_t.avg = open_t.total / (loop_num*NUM_DSETS);
636         close_t.avg = close_t.total / (loop_num*NUM_DSETS);
637         attr_t.avg = attr_t.total / (NUM_ATTRS*NUM_DSETS);
638 
639         /* Print out the performance result */
640         HDfprintf(stderr, "3.  Create %d attributes for each of %d existing datasets for %d times\n",
641             BATCH_ATTRS, NUM_DSETS, loop_num);
642         print_perf(open_t, close_t, attr_t);
643     }
644 
645     if (H5Fclose(file) < 0) goto error;
646 
647     return 0;
648 
649 error:
650     return -1;
651 }
652 
653 
654 /*-------------------------------------------------------------------------
655  * Function:	retrieve_time
656  *
657  * Purpose:     Returns time in seconds, in a double number.
658  *
659  * Programmer:	Raymond Lu
660  *		Friday, Oct 3, 2003
661  *
662  * Modifications:
663  *
664  *-------------------------------------------------------------------------
665  */
retrieve_time(void)666 double retrieve_time(void)
667 {
668 #ifdef H5_HAVE_PARALLEL
669     if(facc_type == FACC_DEFAULT) {
670 #endif /*H5_HAVE_PARALLEL*/
671         struct timeval t;
672         HDgettimeofday(&t, NULL);
673         return ((double)t.tv_sec + (double)t.tv_usec / 1000000);
674 #ifdef H5_HAVE_PARALLEL
675     } else {
676         return MPI_Wtime();
677     }
678 #endif /*H5_HAVE_PARALLEL*/
679 }
680 
681 
682 /*-------------------------------------------------------------------------
683  * Function:	perf
684  *
685  * Purpose:	Calculate total time, maximal and minimal time of
686  *              performance.
687  *
688  * Programmer:	Raymond Lu
689  *		Friday, Oct 3, 2003
690  *
691  * Modifications:
692  *
693  *-------------------------------------------------------------------------
694  */
perf(p_time * perf_t,double start_t,double end_t)695 void perf(p_time *perf_t, double start_t, double end_t)
696 {
697 	double t = end_t - start_t;
698 
699 #ifdef H5_HAVE_PARALLEL
700     if(facc_type == FACC_MPIO) {
701         double reduced_t;
702         double t_max, t_min;
703         int    mpi_size, mpi_rank;
704 
705         MPI_Comm_size(MPI_COMM_WORLD, &mpi_size);
706         MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
707         MPI_Barrier(MPI_COMM_WORLD);
708 
709         MPI_Reduce(&t, &reduced_t, 1, MPI_DOUBLE, MPI_SUM, 0,
710                 MPI_COMM_WORLD);
711         reduced_t /= mpi_size;
712 
713         MPI_Reduce(&t, &t_max, 1, MPI_DOUBLE, MPI_MAX, 0,
714                 MPI_COMM_WORLD);
715         MPI_Reduce(&t, &t_min, 1, MPI_DOUBLE, MPI_MIN, 0,
716                 MPI_COMM_WORLD);
717 
718         if (MAINPROCESS) {
719             perf_t->total += reduced_t;
720 
721 	    if(t_max > perf_t->max)
722 		perf_t->max = t_max;
723 	    if(t_min < perf_t->min)
724 		perf_t->min = t_min;
725         }
726     } else
727 #endif /*H5_HAVE_PARALLEL*/
728     {
729 	perf_t->total += t;
730 
731 	if(t > perf_t->max)
732 		perf_t->max = t;
733 	if(t < perf_t->min)
734 		perf_t->min = t;
735     }
736 }
737 
738 
739 /*-------------------------------------------------------------------------
740  * Function:	print_perf
741  *
742  * Purpose:	Print out performance data.
743  *
744  * Programmer:	Raymond Lu
745  *		Friday, Oct 3, 2003
746  *
747  * Modifications:
748  *
749  *-------------------------------------------------------------------------
750  */
print_perf(p_time open_t,p_time close_t,p_time attr_t)751 void print_perf(p_time open_t, p_time close_t, p_time attr_t)
752 {
753     HDfprintf(stderr, "\t%s:\t\tavg=%.6fs;\tmax=%.6fs;\tmin=%.6fs\n",
754 		open_t.func, open_t.avg, open_t.max, open_t.min);
755     HDfprintf(stderr, "\tH5Dclose:\t\tavg=%.6fs;\tmax=%.6fs;\tmin=%.6fs\n",
756 		close_t.avg, close_t.max, close_t.min);
757     if(NUM_ATTRS)
758         HDfprintf(stderr, "\tH5A(create & close):\tavg=%.6fs;\tmax=%.6fs;\tmin=%.6fs\n",
759 		attr_t.avg, attr_t.max, attr_t.min);
760 }
761 
762 
763 /*-------------------------------------------------------------------------
764  * Function:	main
765  *
766  * Purpose:	Tests
767  *
768  * Return:	Success:	exit(0)
769  *
770  *		Failure:	exit(1)
771  *
772  * Programmer:	Raymond Lu
773  *		Friday, Oct 3, 2003
774  *
775  * Modifications:
776  *
777  *-------------------------------------------------------------------------
778  */
779 int
main(int argc,char ** argv)780 main(int argc, char **argv)
781 {
782 #ifdef H5_HAVE_PARALLEL
783     int mpi_size, mpi_rank;				/* mpi variables */
784 #endif /*H5_HAVE_PARALLEL*/
785 
786     if(parse_options(argc, argv) != 0) {
787        usage();
788        return 0;
789     }
790 
791 #ifdef H5_HAVE_PARALLEL
792     if(facc_type == FACC_MPIO) {
793         MPI_Init(&argc, &argv);
794         MPI_Comm_size(MPI_COMM_WORLD, &mpi_size);
795         MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
796     }
797 #endif /*H5_HAVE_PARALLEL*/
798 
799 #ifdef H5_HAVE_PARALLEL
800     if (facc_type == FACC_DEFAULT || (facc_type != FACC_DEFAULT && MAINPROCESS))
801 #endif /*H5_HAVE_PARALLEL*/
802         HDfprintf(stderr, "\t\tPerformance result of metadata for datasets and attributes\n\n");
803 
804     fapl = H5Pcreate (H5P_FILE_ACCESS);
805 #ifdef H5_HAVE_PARALLEL
806     if(facc_type == FACC_MPIO)
807         H5Pset_fapl_mpio(fapl, MPI_COMM_WORLD, MPI_INFO_NULL);
808 #endif /*H5_HAVE_PARALLEL*/
809 
810     nerrors += create_dspace() < 0 	?1:0;
811 
812     if((RUN_TEST & TEST_1) || !RUN_TEST)
813         nerrors += create_attrs_1() < 0 	?1:0;
814     if((RUN_TEST & TEST_2) || !RUN_TEST)
815         nerrors += create_attrs_2() < 0 	?1:0;
816     if(((RUN_TEST & TEST_3) || !RUN_TEST) && BATCH_ATTRS && NUM_ATTRS)
817         nerrors += create_attrs_3() < 0 	?1:0;
818 
819     if (H5Sclose(space) < 0) goto error;
820     if (H5Sclose(small_space) < 0) goto error;
821 
822     h5_clean_files(FILENAME, fapl);
823 
824 #ifdef H5_HAVE_PARALLEL
825     if(facc_type == FACC_MPIO)
826         /* MPI_Finalize must be called AFTER H5close which may use MPI calls */
827         MPI_Finalize();
828 #endif /*H5_HAVE_PARALLEL*/
829 
830     if (nerrors) goto error;
831 #ifdef H5_HAVE_PARALLEL
832     if (facc_type != FACC_DEFAULT && MAINPROCESS)
833 #endif /*H5_HAVE_PARALLEL*/
834         printf("All metadata performance tests passed.\n");
835 
836     return 0;
837 
838  error:
839     nerrors = MAX(1, nerrors);
840 #ifdef H5_HAVE_PARALLEL
841     if (facc_type != FACC_DEFAULT && MAINPROCESS)
842 #endif /*H5_HAVE_PARALLEL*/
843         printf("***** %d PERFORMANCE TEST%s FAILED! *****\n",
844 	   nerrors, 1 == nerrors ? "" : "S");
845 
846     return 1;
847 }
848 
849