-
Notifications
You must be signed in to change notification settings - Fork 120
/
mpi_io.c
286 lines (249 loc) · 11.1 KB
/
mpi_io.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
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
/* File: mpi_read.c */
/* Description: This program reads data points from a file using MPI-IO */
/* that implements a simple k-means clustering algorithm */
/* 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>
#include <mpi.h>
#include "kmeans.h"
/*---< mpi_read() >----------------------------------------------------------*/
float** mpi_read(int isBinaryFile, /* flag: 0 or 1 */
char *filename, /* input file name */
int *numObjs, /* no. data objects (local) */
int *numCoords, /* no. coordinates */
MPI_Comm comm)
{
float **objects;
int i, j, len, divd, rem;
int rank, nproc;
MPI_Status status;
MPI_Comm_rank(comm, &rank);
MPI_Comm_size(comm, &nproc);
if (isBinaryFile) { /* using MPI-IO to read file concurrently */
int err;
MPI_Offset disp;
MPI_Datatype filetype;
MPI_File fh;
err = MPI_File_open(comm, filename, MPI_MODE_RDONLY, MPI_INFO_NULL, &fh);
if (err != MPI_SUCCESS) {
char errstr[MPI_MAX_ERROR_STRING];
int errlen;
MPI_Error_string(err, errstr, &errlen);
printf("Error at opening file %s (%s)\n",filename,errstr);
MPI_Finalize();
exit(1);
}
/* read numObjs & numCoords from the 1st 2 integers */
MPI_File_read(fh, numObjs, 1, MPI_INT, &status);
MPI_File_read(fh, numCoords, 1, MPI_INT, &status);
if (*numObjs <= 0 || *numCoords <= 0) {
printf("Error: file format (%s)\n",filename);
MPI_Finalize();
exit(1);
}
divd = (*numObjs) / nproc;
rem = (*numObjs) % nproc;
len = (rank < rem) ? rank*(divd+1) : rank*divd + rem;
disp = 2 * sizeof(int) + len * (*numCoords) * sizeof(float);
/* adjust numObjs to be local size */
(*numObjs) = (rank < rem) ? divd+1 : divd;
/* allocate space for data points */
objects = (float**)malloc((*numObjs) * sizeof(float*));
assert(objects != NULL);
objects[0] = (float*) malloc((*numObjs)*(*numCoords) * sizeof(float));
assert(objects[0] != NULL);
for (i=1; i<(*numObjs); i++)
objects[i] = objects[i-1] + (*numCoords);
/* define a file type for file view */
MPI_Type_contiguous((*numObjs), MPI_FLOAT, &filetype);
MPI_Type_commit(&filetype);
MPI_File_set_view(fh, disp, MPI_FLOAT, filetype, "native",
MPI_INFO_NULL);
MPI_File_read_all(fh, objects[0], (*numObjs)*(*numCoords),
MPI_FLOAT, &status);
MPI_Type_free(&filetype);
MPI_File_close(&fh);
}
else { /* ASCII format: let proc 0 read and distribute to others */
if (rank == 0) {
objects = file_read(0, filename, numObjs, numCoords);
if (objects == NULL) *numObjs = -1;
}
/* broadcast global numObjs and numCoords to the rest proc */
MPI_Bcast(numObjs, 1, MPI_INT, 0, comm);
MPI_Bcast(numCoords, 1, MPI_INT, 0, comm);
if (*numObjs == -1) {
MPI_Finalize();
exit(1);
}
divd = (*numObjs) / nproc;
rem = (*numObjs) % nproc;
if (rank == 0) {
int index = (rem > 0) ? divd+1 : divd;
/* index is the numObjs partitioned locally in proc 0 */
(*numObjs) = index;
/* distribute objects[] to other processes */
for (i=1; i<nproc; i++) {
int msg_size = (i < rem) ? (divd+1) : divd;
MPI_Send(objects[index], msg_size*(*numCoords), MPI_FLOAT,
i, i, comm);
index += msg_size;
}
/* reduce the objects[] to local size */
objects[0] = realloc(objects[0],
(*numObjs)*(*numCoords)*sizeof(float));
assert(objects[0] != NULL);
objects = realloc(objects, (*numObjs)*sizeof(float*));
assert(objects != NULL);
}
else {
/* local numObjs */
(*numObjs) = (rank < rem) ? divd+1 : divd;
/* allocate space for data points */
objects = (float**)malloc((*numObjs) *sizeof(float*));
assert(objects != NULL);
objects[0] = (float*) malloc((*numObjs)*(*numCoords)*sizeof(float));
assert(objects[0] != NULL);
for (i=1; i<(*numObjs); i++)
objects[i] = objects[i-1] + (*numCoords);
MPI_Recv(objects[0], (*numObjs)*(*numCoords), MPI_FLOAT, 0,
rank, comm, &status);
}
}
return objects;
}
/*---< mpi_write() >---------------------------------------------------------*/
int mpi_write(int isOutFileBinary, /* flag: 0 or 1 */
char *filename, /* input file name */
int numClusters, /* no. clusters */
int numObjs, /* no. data objects */
int numCoords, /* no. coordinates (local) */
float **clusters, /* [numClusters][numCoords] centers */
int *membership, /* [numObjs] */
int totalNumObjs, /* total no. data objects */
MPI_Comm comm)
{
int divd, rem, len, err;
int i, j, k, rank, nproc;
char outFileName[1024], fs_type[32], str[32], *delim;
MPI_File fh;
MPI_Status status;
MPI_Comm_rank(comm, &rank);
MPI_Comm_size(comm, &nproc);
delim = strchr(filename, ':');
if (delim != NULL) {
strncpy(fs_type, filename, delim-filename);
fs_type[delim-filename] = '\0';
/* real file name starts from delim+1 */
delim++;
}
else
delim = filename;
/* output: the coordinates of the cluster centres ----------------------*/
/* only proc 0 do this, because clusters[] are the same across all proc */
if (rank == 0) {
printf("Writing coordinates of K=%d cluster centers to file \"%s.cluster_centres\"\n",
numClusters, delim);
sprintf(outFileName, "%s.cluster_centres", filename);
err = MPI_File_open(MPI_COMM_SELF, outFileName, MPI_MODE_CREATE | MPI_MODE_WRONLY, MPI_INFO_NULL, &fh);
if (err != MPI_SUCCESS) {
char errstr[MPI_MAX_ERROR_STRING];
int errlen;
MPI_Error_string(err, errstr, &errlen);
printf("Error at opening file %s (%s)\n", outFileName,errstr);
MPI_Finalize();
exit(1);
}
if (isOutFileBinary) {
MPI_File_write(fh, &numClusters, 1, MPI_INT, &status);
MPI_File_write(fh, &numCoords, 1, MPI_INT, &status);
MPI_File_write(fh, clusters[0], numClusters*numCoords, MPI_FLOAT, &status);
}
else {
for (i=0; i<numClusters; i++) {
sprintf(str, "%d ", i);
MPI_File_write(fh, str, strlen(str), MPI_CHAR, &status);
for (j=0; j<numCoords; j++) {
sprintf(str, "%f ", clusters[i][j]);
MPI_File_write(fh, str, strlen(str), MPI_CHAR, &status);
}
MPI_File_write(fh, "\n", 1, MPI_CHAR, &status);
}
}
MPI_File_close(&fh);
}
/* output: the closest cluster centre to each of the data points --------*/
if (rank == 0)
printf("Writing membership of N=%d data objects to file \"%s.membership\"\n",
totalNumObjs, delim);
if (isOutFileBinary) {
sprintf(outFileName, "%s.membership", filename);
err = MPI_File_open(comm, outFileName, MPI_MODE_CREATE | MPI_MODE_WRONLY, MPI_INFO_NULL, &fh);
if (err != MPI_SUCCESS) {
char errstr[MPI_MAX_ERROR_STRING];
int errlen;
MPI_Error_string(err, errstr, &errlen);
printf("Error at opening file %s (%s)\n", outFileName,errstr);
MPI_Finalize();
exit(1);
}
/* write numObjs from the 1st integer */
if (rank == 0)
MPI_File_write(fh, &totalNumObjs, 1, MPI_INT, &status);
divd = totalNumObjs / nproc;
rem = totalNumObjs % nproc;
len = (rank < rem) ? rank*(divd+1) : rank*divd + rem;
MPI_Offset disp = (len + 1) * sizeof(int);
MPI_File_seek(fh, disp, MPI_SEEK_SET);
MPI_File_write(fh, membership, numObjs, MPI_INT, &status);
MPI_File_close(&fh);
}
else {
if (rank == 0) { /* gather membership[] from all processes ----------*/
int divd = totalNumObjs / nproc;
int rem = totalNumObjs % nproc;
sprintf(outFileName, "%s.membership", filename);
err = MPI_File_open(MPI_COMM_SELF, outFileName, MPI_MODE_CREATE | MPI_MODE_WRONLY, MPI_INFO_NULL, &fh);
if (err != MPI_SUCCESS) {
char errstr[MPI_MAX_ERROR_STRING];
int errlen;
MPI_Error_string(err, errstr, &errlen);
printf("Error at opening file %s (%s)\n", outFileName,errstr);
MPI_Finalize();
exit(1);
}
/* first, print out local membership[] */
for (j=0; j<numObjs; j++) {
sprintf(str, "%d %d\n", j, membership[j]);
MPI_File_write(fh, str, strlen(str), MPI_CHAR, &status);
}
k = numObjs;
for (i=1; i<nproc; i++) {
numObjs = (i < rem) ? divd+1 : divd;
MPI_Recv(membership, numObjs, MPI_INT, i, i, comm, &status);
for (j=0; j<numObjs; j++) {
sprintf(str, "%d %d\n", k++, membership[j]);
MPI_File_write(fh, str, strlen(str), MPI_CHAR, &status);
}
}
MPI_File_close(&fh);
}
else {
MPI_Send(membership, numObjs, MPI_INT, 0, rank, comm);
}
}
return 1;
}