00001
00002
00003
00004
00005
00006 #include "mpi.h"
00007 #include <stdio.h>
00008 #include <string.h>
00009 #include <stdlib.h>
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021 int main(int argc, char **argv)
00022 {
00023 MPI_Datatype newtype;
00024 int i, ndims, array_of_gsizes[3], array_of_distribs[3];
00025 int order, nprocs, len, *buf, bufcount, mynod;
00026 int array_of_dargs[3], array_of_psizes[3];
00027 MPI_File fh;
00028 MPI_Status status;
00029 double stim, write_tim, new_write_tim, write_bw;
00030 double read_tim, new_read_tim, read_bw;
00031 char *filename;
00032
00033 MPI_Init(&argc,&argv);
00034 MPI_Comm_rank(MPI_COMM_WORLD, &mynod);
00035 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
00036
00037
00038
00039 if (!mynod) {
00040 i = 1;
00041 while ((i < argc) && strcmp("-fname", *argv)) {
00042 i++;
00043 argv++;
00044 }
00045 if (i >= argc) {
00046 fprintf(stderr, "\n*# Usage: coll_perf -fname filename\n\n");
00047 MPI_Abort(MPI_COMM_WORLD, 1);
00048 }
00049 argv++;
00050 len = strlen(*argv);
00051 filename = (char *) malloc(len+1);
00052 strcpy(filename, *argv);
00053 MPI_Bcast(&len, 1, MPI_INT, 0, MPI_COMM_WORLD);
00054 MPI_Bcast(filename, len+1, MPI_CHAR, 0, MPI_COMM_WORLD);
00055 }
00056 else {
00057 MPI_Bcast(&len, 1, MPI_INT, 0, MPI_COMM_WORLD);
00058 filename = (char *) malloc(len+1);
00059 MPI_Bcast(filename, len+1, MPI_CHAR, 0, MPI_COMM_WORLD);
00060 }
00061
00062
00063 ndims = 3;
00064 order = MPI_ORDER_C;
00065
00066 array_of_gsizes[0] = 128;
00067 array_of_gsizes[1] = 128;
00068 array_of_gsizes[2] = 128;
00069
00070 array_of_distribs[0] = MPI_DISTRIBUTE_BLOCK;
00071 array_of_distribs[1] = MPI_DISTRIBUTE_BLOCK;
00072 array_of_distribs[2] = MPI_DISTRIBUTE_BLOCK;
00073
00074 array_of_dargs[0] = MPI_DISTRIBUTE_DFLT_DARG;
00075 array_of_dargs[1] = MPI_DISTRIBUTE_DFLT_DARG;
00076 array_of_dargs[2] = MPI_DISTRIBUTE_DFLT_DARG;
00077
00078 for (i=0; i<ndims; i++) array_of_psizes[i] = 0;
00079 MPI_Dims_create(nprocs, ndims, array_of_psizes);
00080
00081 MPI_Type_create_darray(nprocs, mynod, ndims, array_of_gsizes,
00082 array_of_distribs, array_of_dargs,
00083 array_of_psizes, order, MPI_INT, &newtype);
00084 MPI_Type_commit(&newtype);
00085
00086 MPI_Type_size(newtype, &bufcount);
00087 bufcount = bufcount/sizeof(int);
00088 buf = (int *) malloc(bufcount * sizeof(int));
00089
00090
00091
00092
00093 MPI_File_open(MPI_COMM_WORLD, filename, MPI_MODE_CREATE | MPI_MODE_RDWR,
00094 MPI_INFO_NULL, &fh);
00095 MPI_File_set_view(fh, 0, MPI_INT, newtype, "native", MPI_INFO_NULL);
00096 MPI_File_write_all(fh, buf, bufcount, MPI_INT, &status);
00097 MPI_File_seek(fh, 0, MPI_SEEK_SET);
00098 MPI_File_read_all(fh, buf, bufcount, MPI_INT, &status);
00099 MPI_File_close(&fh);
00100
00101 MPI_Barrier(MPI_COMM_WORLD);
00102
00103
00104 MPI_File_open(MPI_COMM_WORLD, filename, MPI_MODE_CREATE | MPI_MODE_RDWR,
00105 MPI_INFO_NULL, &fh);
00106 MPI_File_set_view(fh, 0, MPI_INT, newtype, "native", MPI_INFO_NULL);
00107
00108 MPI_Barrier(MPI_COMM_WORLD);
00109 stim = MPI_Wtime();
00110 MPI_File_write_all(fh, buf, bufcount, MPI_INT, &status);
00111 write_tim = MPI_Wtime() - stim;
00112 MPI_File_close(&fh);
00113
00114 MPI_Allreduce(&write_tim, &new_write_tim, 1, MPI_DOUBLE, MPI_MAX,
00115 MPI_COMM_WORLD);
00116
00117 if (mynod == 0) {
00118 write_bw = (array_of_gsizes[0]*array_of_gsizes[1]*array_of_gsizes[2]*sizeof(int))/(new_write_tim*1024.0*1024.0);
00119 fprintf(stderr, "Global array size %d x %d x %d integers\n", array_of_gsizes[0], array_of_gsizes[1], array_of_gsizes[2]);
00120 fprintf(stderr, "Collective write time = %f sec, Collective write bandwidth = %f Mbytes/sec\n", new_write_tim, write_bw);
00121 }
00122
00123 MPI_Barrier(MPI_COMM_WORLD);
00124
00125
00126 MPI_File_open(MPI_COMM_WORLD, filename, MPI_MODE_CREATE | MPI_MODE_RDWR,
00127 MPI_INFO_NULL, &fh);
00128 MPI_File_set_view(fh, 0, MPI_INT, newtype, "native", MPI_INFO_NULL);
00129
00130 MPI_Barrier(MPI_COMM_WORLD);
00131 stim = MPI_Wtime();
00132 MPI_File_read_all(fh, buf, bufcount, MPI_INT, &status);
00133 read_tim = MPI_Wtime() - stim;
00134 MPI_File_close(&fh);
00135
00136 MPI_Allreduce(&read_tim, &new_read_tim, 1, MPI_DOUBLE, MPI_MAX,
00137 MPI_COMM_WORLD);
00138
00139 if (mynod == 0) {
00140 read_bw = (array_of_gsizes[0]*array_of_gsizes[1]*array_of_gsizes[2]*sizeof(int))/(new_read_tim*1024.0*1024.0);
00141 fprintf(stderr, "Collective read time = %f sec, Collective read bandwidth = %f Mbytes/sec\n", new_read_tim, read_bw);
00142 }
00143
00144 MPI_Type_free(&newtype);
00145 free(buf);
00146 free(filename);
00147
00148 MPI_Finalize();
00149 return 0;
00150 }