-
Notifications
You must be signed in to change notification settings - Fork 120
/
mpi_main.c
214 lines (179 loc) · 8.25 KB
/
mpi_main.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
/* File: mpi_main.c (an MPI version) */
/* Description: This program shows an example on how to call a subroutine */
/* that implements a simple k-means clustering algorithm */
/* based on Euclid distance. */
/* Input file format: */
/* ascii file: each line contains 1 data object */
/* binary file: first 4-byte integer is the number of data */
/* objects and 2nd integer is the no. of features (or */
/* coordinates) of each object */
/* */
/* Author: Wei-keng Liao */
/* ECE Department Northwestern University */
/* email: [email protected] */
/* Copyright, 2005, Wei-keng Liao */
/* */
/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
#include <stdio.h>
#include <stdlib.h>
#include <string.h> /* strtok() */
#include <sys/types.h> /* open() */
#include <sys/stat.h>
#include <fcntl.h>
#include <unistd.h> /* getopt() */
#include <mpi.h>
int _debug;
#include "kmeans.h"
int mpi_kmeans(float**, int, int, int, float, int*, float**, MPI_Comm);
float** mpi_read(int, char*, int*, int*, MPI_Comm);
int mpi_write(int, char*, int, int, int, float**, int*, int, MPI_Comm);
/*---< usage() >------------------------------------------------------------*/
static void usage(char *argv0, float threshold) {
char *help =
"Usage: %s [switches] -i filename -n num_clusters\n"
" -i filename : file containing data to be clustered\n"
" -b : input file is in binary format (default no)\n"
" -r : output file in binary format (default no)\n"
" -n num_clusters: number of clusters (K must > 1)\n"
" -t threshold : threshold value (default %.4f)\n"
" -o : output timing results (default no)\n"
" -d : enable debug mode\n";
fprintf(stderr, help, argv0, threshold);
}
/*---< main() >-------------------------------------------------------------*/
int main(int argc, char **argv) {
int opt;
extern char *optarg;
extern int optind;
int i, j;
int isInFileBinary, isOutFileBinary;
int is_output_timing, is_print_usage;
int numClusters, numCoords, numObjs, totalNumObjs;
int *membership; /* [numObjs] */
char *filename;
float **objects; /* [numObjs][numCoords] data objects */
float **clusters; /* [numClusters][numCoords] cluster center */
float threshold;
double timing, io_timing, clustering_timing;
int rank, nproc, mpi_namelen;
char mpi_name[MPI_MAX_PROCESSOR_NAME];
MPI_Status status;
MPI_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &nproc);
MPI_Get_processor_name(mpi_name,&mpi_namelen);
/* some default values */
_debug = 0;
threshold = 0.001;
numClusters = 0;
isInFileBinary = 0;
isOutFileBinary = 0;
is_output_timing = 0;
is_print_usage = 0;
filename = NULL;
while ( (opt=getopt(argc,argv,"p:i:n:t:abdorh"))!= EOF) {
switch (opt) {
case 'i': filename=optarg;
break;
case 'b': isInFileBinary = 1;
break;
case 'r': isOutFileBinary = 1;
break;
case 't': threshold=atof(optarg);
break;
case 'n': numClusters = atoi(optarg);
break;
case 'o': is_output_timing = 1;
break;
case 'd': _debug = 1;
break;
case 'h': is_print_usage = 1;
break;
default: is_print_usage = 1;
break;
}
}
if (filename == 0 || numClusters <= 1 || is_print_usage == 1) {
if (rank == 0) usage(argv[0], threshold);
MPI_Finalize();
exit(1);
}
if (_debug) printf("Proc %d of %d running on %s\n", rank, nproc, mpi_name);
MPI_Barrier(MPI_COMM_WORLD);
io_timing = MPI_Wtime();
/* read data points from file ------------------------------------------*/
objects = mpi_read(isInFileBinary, filename, &numObjs, &numCoords,
MPI_COMM_WORLD);
if (_debug) { /* print the first 4 objects' coordinates */
int num = (numObjs < 4) ? numObjs : 4;
for (i=0; i<num; i++) {
char strline[1024], strfloat[16];
sprintf(strline,"%d: objects[%d]= ",rank,i);
for (j=0; j<numCoords; j++) {
sprintf(strfloat,"%10f",objects[i][j]);
strcat(strline, strfloat);
}
strcat(strline, "\n");
printf("%s",strline);
}
}
timing = MPI_Wtime();
io_timing = timing - io_timing;
clustering_timing = timing;
/* allocate a 2D space for clusters[] (coordinates of cluster centers)
this array should be the same across all processes */
clusters = (float**) malloc(numClusters * sizeof(float*));
assert(clusters != NULL);
clusters[0] = (float*) malloc(numClusters * numCoords * sizeof(float));
assert(clusters[0] != NULL);
for (i=1; i<numClusters; i++)
clusters[i] = clusters[i-1] + numCoords;
MPI_Allreduce(&numObjs, &totalNumObjs, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
/* pick first numClusters elements in feature[] as initial cluster centers*/
if (rank == 0) {
for (i=0; i<numClusters; i++)
for (j=0; j<numCoords; j++)
clusters[i][j] = objects[i][j];
}
MPI_Bcast(clusters[0], numClusters*numCoords, MPI_FLOAT, 0, MPI_COMM_WORLD);
/* membership: the cluster id for each data object */
membership = (int*) malloc(numObjs * sizeof(int));
assert(membership != NULL);
/* start the core computation -------------------------------------------*/
mpi_kmeans(objects, numCoords, numObjs, numClusters, threshold, membership,
clusters, MPI_COMM_WORLD);
free(objects[0]);
free(objects);
timing = MPI_Wtime();
clustering_timing = timing - clustering_timing;
/* output: the coordinates of the cluster centres ----------------------*/
mpi_write(isOutFileBinary, filename, numClusters, numObjs, numCoords,
clusters, membership, totalNumObjs, MPI_COMM_WORLD);
free(membership);
free(clusters[0]);
free(clusters);
/*---- output performance numbers ---------------------------------------*/
if (is_output_timing) {
double max_io_timing, max_clustering_timing;
io_timing += MPI_Wtime() - timing;
/* get the max timing measured among all processes */
MPI_Reduce(&io_timing, &max_io_timing, 1, MPI_DOUBLE,
MPI_MAX, 0, MPI_COMM_WORLD);
MPI_Reduce(&clustering_timing, &max_clustering_timing, 1, MPI_DOUBLE,
MPI_MAX, 0, MPI_COMM_WORLD);
if (rank == 0) {
printf("\nPerforming **** Simple Kmeans (MPI) ****\n");
printf("Num of processes = %d\n", nproc);
printf("Input file: %s\n", filename);
printf("numObjs = %d\n", totalNumObjs);
printf("numCoords = %d\n", numCoords);
printf("numClusters = %d\n", numClusters);
printf("threshold = %.4f\n", threshold);
printf("I/O time = %10.4f sec\n", max_io_timing);
printf("Computation timing = %10.4f sec\n", max_clustering_timing);
}
}
MPI_Finalize();
return(0);
}