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  * Parallel tests for file image operations
16  */
17 
18 #include "testphdf5.h"
19 
20 /* file_image_daisy_chain_test
21  *
22  * Process zero:
23  *
24  *	1) Creates a core file with an integer vector data set of
25  * 	   length n (= mpi_size),
26  *
27  *	2) Initializes the vector to zero in * location 0, and to -1
28  *         everywhere else.
29  *
30  *	3) Flushes the core file, and gets an image of it.  Closes
31  *	   the core file.
32  *
33  *	4) Sends the image to process 1.
34  *
35  *	5) Awaits receipt on a file image from process n-1.
36  *
37  *	6) opens the image received from process n-1, verifies that
38  *	   it contains a vector of length equal to mpi_size, and
39  *	   that the vector contains (0, 1, 2, ... n-1)
40  *
41  *	7) closes the core file and exits.
42  *
43  * Process i (0 < i < n)
44  *
45  *	1) Await receipt of file image from process (i - 1).
46  *
47  *	2) Open the image with the core file driver, verify that i
48  *	   contains a vector v of length, and that v[j] = j for
49  *	   0 <= j < i, and that v[j] == -1 for i <= j < n
50  *
51  *	3) Set v[i] = i in the core file.
52  *
53  *	4) Flush the core file and send it to process (i + 1) % n.
54  *
55  *	5) close the core file and exit.
56  *
57  * Test fails on a hang (if an image is not received), or on invalid data.
58  *
59  *                                               JRM -- 11/28/11
60  */
61 void
file_image_daisy_chain_test(void)62 file_image_daisy_chain_test(void)
63 {
64     char file_name[1024] = "\0";
65     int mpi_size, mpi_rank;
66     int mpi_result;
67     int i;
68     int space_ndims;
69     MPI_Status rcvstat;
70     int * vector_ptr = NULL;
71     hid_t fapl_id = -1;
72     hid_t file_id;			/* file IDs */
73     hid_t dset_id = -1;
74     hid_t dset_type_id = -1;
75     hid_t space_id = -1;
76     herr_t err;
77     hsize_t dims[1];
78     void * image_ptr = NULL;
79     ssize_t bytes_read;
80     ssize_t image_len;
81     hbool_t vector_ok = TRUE;
82     htri_t tri_result;
83 
84 
85     /* set up MPI parameters */
86     MPI_Comm_size(MPI_COMM_WORLD, &mpi_size);
87     MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);
88 
89     /* setup file name */
90     HDsnprintf(file_name, 1024, "file_image_daisy_chain_test_%05d.h5",
91                (int)mpi_rank);
92 
93     if(mpi_rank == 0) {
94 
95 	/* 1) Creates a core file with an integer vector data set
96          *    of length mpi_size,
97          */
98 	fapl_id = H5Pcreate(H5P_FILE_ACCESS);
99         VRFY((fapl_id >= 0), "creating fapl");
100 
101 	err = H5Pset_fapl_core(fapl_id, (size_t)(64 *1024), FALSE);
102         VRFY((err >= 0), "setting core file driver in fapl.");
103 
104         file_id = H5Fcreate(file_name, 0, H5P_DEFAULT, fapl_id);
105         VRFY((file_id >= 0), "created core file");
106 
107         dims[0] = (hsize_t)mpi_size;
108 	space_id = H5Screate_simple(1, dims, dims);
109         VRFY((space_id >= 0), "created data space");
110 
111         dset_id = H5Dcreate2(file_id, "v", H5T_NATIVE_INT, space_id,
112                              H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
113         VRFY((dset_id >= 0), "created data set");
114 
115 
116 	/* 2) Initialize the vector to zero in location 0, and
117          *    to -1 everywhere else.
118          */
119 
120 	vector_ptr = (int *)HDmalloc((size_t)(mpi_size) * sizeof(int));
121         VRFY((vector_ptr != NULL), "allocated in memory representation of vector");
122 
123         vector_ptr[0] = 0;
124         for(i = 1; i < mpi_size; i++)
125             vector_ptr[i] = -1;
126 
127 	err = H5Dwrite(dset_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL,
128                        H5P_DEFAULT, (void *)vector_ptr);
129         VRFY((err >= 0), "wrote initial data to vector.");
130 
131         HDfree(vector_ptr);
132         vector_ptr = NULL;
133 
134 
135         /* 3) Flush the core file, and get an image of it.  Close
136          *    the core file.
137          */
138 	err = H5Fflush(file_id, H5F_SCOPE_GLOBAL);
139         VRFY((err >= 0), "flushed core file.");
140 
141         image_len = H5Fget_file_image(file_id, NULL, (size_t)0);
142         VRFY((image_len > 0), "got image file size");
143 
144         image_ptr = (void *)HDmalloc((size_t)image_len);
145         VRFY(image_ptr != NULL, "allocated file image buffer.");
146 
147         bytes_read = H5Fget_file_image(file_id, image_ptr, (size_t)image_len);
148         VRFY(bytes_read == image_len, "wrote file into image buffer");
149 
150         err = H5Sclose(space_id);
151 	VRFY((err >= 0), "closed data space.");
152 
153 	err = H5Dclose(dset_id);
154 	VRFY((err >= 0), "closed data set.");
155 
156 	err = H5Fclose(file_id);
157 	VRFY((err >= 0), "closed core file(1).");
158 
159 	err = H5Pclose(fapl_id);
160 	VRFY((err >= 0), "closed fapl(1).");
161 
162 
163         /* 4) Send the image to process 1. */
164 
165         mpi_result = MPI_Ssend((void *)(&image_len), (int)sizeof(ssize_t),
166 			       MPI_BYTE, 1, 0, MPI_COMM_WORLD);
167 	VRFY((mpi_result == MPI_SUCCESS), "sent image size to process 1");
168 
169         mpi_result = MPI_Ssend((void *)image_ptr, (int)image_len,
170 			       MPI_BYTE, 1, 0, MPI_COMM_WORLD);
171 	VRFY((mpi_result == MPI_SUCCESS), "sent image to process 1");
172 
173         HDfree(image_ptr);
174         image_ptr = NULL;
175         image_len = 0;
176 
177 
178 	/* 5) Await receipt on a file image from process n-1. */
179 
180 	mpi_result = MPI_Recv((void *)(&image_len), (int)sizeof(ssize_t),
181                               MPI_BYTE, mpi_size - 1, 0, MPI_COMM_WORLD,
182                               &rcvstat);
183 	VRFY((mpi_result == MPI_SUCCESS), "received image len from process n-1");
184 
185         image_ptr = (void *)HDmalloc((size_t)image_len);
186         VRFY(image_ptr != NULL, "allocated file image receive buffer.");
187 
188 	mpi_result = MPI_Recv((void *)image_ptr, (int)image_len,
189                               MPI_BYTE, mpi_size - 1, 0, MPI_COMM_WORLD,
190                               &rcvstat);
191 	VRFY((mpi_result == MPI_SUCCESS), \
192              "received file image from process n-1");
193 
194 	/* 6) open the image received from process n-1, verify that
195          *    it contains a vector of length equal to mpi_size, and
196 	 *    that the vector contains (0, 1, 2, ... n-1).
197          */
198 	fapl_id = H5Pcreate(H5P_FILE_ACCESS);
199         VRFY((fapl_id >= 0), "creating fapl");
200 
201 	err = H5Pset_fapl_core(fapl_id, (size_t)(64 *1024), FALSE);
202         VRFY((err >= 0), "setting core file driver in fapl.");
203 
204 	err = H5Pset_file_image(fapl_id, image_ptr, (size_t)image_len);
205         VRFY((err >= 0), "set file image in fapl.");
206 
207         file_id = H5Fopen(file_name, H5F_ACC_RDWR, fapl_id);
208         VRFY((file_id >= 0), "opened received file image file");
209 
210 	dset_id = H5Dopen2(file_id, "v", H5P_DEFAULT);
211         VRFY((dset_id >= 0), "opened data set");
212 
213 	dset_type_id = H5Dget_type(dset_id);
214         VRFY((dset_type_id >= 0), "obtained data set type");
215 
216 	tri_result = H5Tequal(dset_type_id, H5T_NATIVE_INT);
217         VRFY((tri_result == TRUE), "verified data set type");
218 
219 	space_id = H5Dget_space(dset_id);
220         VRFY((space_id >= 0), "opened data space");
221 
222 	space_ndims = H5Sget_simple_extent_ndims(space_id);
223 	VRFY((space_ndims == 1), "verified data space num dims(1)");
224 
225 	space_ndims = H5Sget_simple_extent_dims(space_id, dims, NULL);
226 	VRFY((space_ndims == 1), "verified data space num dims(2)");
227 	VRFY((dims[0] == (hsize_t)mpi_size), "verified data space dims");
228 
229 	vector_ptr = (int *)HDmalloc((size_t)(mpi_size) * sizeof(int));
230         VRFY((vector_ptr != NULL), "allocated in memory rep of vector");
231 
232 	err = H5Dread(dset_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL,
233                       H5P_DEFAULT, (void *)vector_ptr);
234         VRFY((err >= 0), "read received vector.");
235 
236 	vector_ok = TRUE;
237 	for(i = 0; i < mpi_size; i++)
238             if(vector_ptr[i] != i)
239                 vector_ok = FALSE;
240         VRFY((vector_ok), "verified received vector.");
241 
242         HDfree(vector_ptr);
243         vector_ptr = NULL;
244 
245 	/* 7) closes the core file and exit. */
246 
247         err = H5Sclose(space_id);
248 	VRFY((err >= 0), "closed data space.");
249 
250 	err = H5Dclose(dset_id);
251 	VRFY((err >= 0), "closed data set.");
252 
253 	err = H5Fclose(file_id);
254 	VRFY((err >= 0), "closed core file(1).");
255 
256 	err = H5Pclose(fapl_id);
257 	VRFY((err >= 0), "closed fapl(1).");
258 
259         HDfree(image_ptr);
260         image_ptr = NULL;
261         image_len = 0;
262     } else {
263         /* 1) Await receipt of file image from process (i - 1). */
264 
265 	mpi_result = MPI_Recv((void *)(&image_len), (int)sizeof(ssize_t),
266                               MPI_BYTE, mpi_rank - 1, 0, MPI_COMM_WORLD,
267                               &rcvstat);
268 	VRFY((mpi_result == MPI_SUCCESS), \
269              "received image size from process mpi_rank-1");
270 
271         image_ptr = (void *)HDmalloc((size_t)image_len);
272         VRFY(image_ptr != NULL, "allocated file image receive buffer.");
273 
274 	mpi_result = MPI_Recv((void *)image_ptr, (int)image_len,
275                               MPI_BYTE, mpi_rank - 1, 0, MPI_COMM_WORLD,
276                               &rcvstat);
277 	VRFY((mpi_result == MPI_SUCCESS), \
278              "received file image from process mpi_rank-1");
279 
280 	/* 2) Open the image with the core file driver, verify that it
281 	 *    contains a vector v of length, and that v[j] = j for
282 	 *    0 <= j < i, and that v[j] == -1 for i <= j < n
283 	 */
284 	fapl_id = H5Pcreate(H5P_FILE_ACCESS);
285         VRFY((fapl_id >= 0), "creating fapl");
286 
287 	err = H5Pset_fapl_core(fapl_id, (size_t)(64 * 1024), FALSE);
288         VRFY((err >= 0), "setting core file driver in fapl.");
289 
290 	err = H5Pset_file_image(fapl_id, image_ptr, (size_t)image_len);
291         VRFY((err >= 0), "set file image in fapl.");
292 
293         file_id = H5Fopen(file_name, H5F_ACC_RDWR, fapl_id);
294 	H5Eprint2(H5P_DEFAULT, stderr);
295         VRFY((file_id >= 0), "opened received file image file");
296 
297 	dset_id = H5Dopen2(file_id, "v", H5P_DEFAULT);
298         VRFY((dset_id >= 0), "opened data set");
299 
300 	dset_type_id = H5Dget_type(dset_id);
301         VRFY((dset_type_id >= 0), "obtained data set type");
302 
303 	tri_result = H5Tequal(dset_type_id, H5T_NATIVE_INT);
304         VRFY((tri_result == TRUE), "verified data set type");
305 
306 	space_id = H5Dget_space(dset_id);
307         VRFY((space_id >= 0), "opened data space");
308 
309 	space_ndims = H5Sget_simple_extent_ndims(space_id);
310 	VRFY((space_ndims == 1), "verified data space num dims(1)");
311 
312 	space_ndims = H5Sget_simple_extent_dims(space_id, dims, NULL);
313 	VRFY((space_ndims == 1), "verified data space num dims(2)");
314 	VRFY((dims[0] == (hsize_t)mpi_size), "verified data space dims");
315 
316 	vector_ptr = (int *)HDmalloc((size_t)(mpi_size) * sizeof(int));
317         VRFY((vector_ptr != NULL), "allocated in memory rep of vector");
318 
319 	err = H5Dread(dset_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL,
320                        H5P_DEFAULT, (void *)vector_ptr);
321         VRFY((err >= 0), "read received vector.");
322 
323 	vector_ok = TRUE;
324 	for(i = 0; i < mpi_size; i++){
325             if(i < mpi_rank) {
326                 if(vector_ptr[i] != i)
327                     vector_ok = FALSE;
328 	    } else {
329                 if(vector_ptr[i] != -1)
330                     vector_ok = FALSE;
331 	    }
332         }
333         VRFY((vector_ok), "verified received vector.");
334 
335 
336 	/* 3) Set v[i] = i in the core file. */
337 
338 	vector_ptr[mpi_rank] = mpi_rank;
339 
340 	err = H5Dwrite(dset_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL,
341                        H5P_DEFAULT, (void *)vector_ptr);
342         VRFY((err >= 0), "wrote modified data to vector.");
343 
344         HDfree(vector_ptr);
345         vector_ptr = NULL;
346 
347 
348 	/* 4) Flush the core file and send it to process (mpi_rank + 1) % n. */
349 
350 	err = H5Fflush(file_id, H5F_SCOPE_GLOBAL);
351         VRFY((err >= 0), "flushed core file.");
352 
353         image_len = H5Fget_file_image(file_id, NULL, (size_t)0);
354         VRFY((image_len > 0), "got (possibly modified) image file len");
355 
356         image_ptr = (void *)HDrealloc((void *)image_ptr, (size_t)image_len);
357         VRFY(image_ptr != NULL, "re-allocated file image buffer.");
358 
359         bytes_read = H5Fget_file_image(file_id, image_ptr, (size_t)image_len);
360         VRFY(bytes_read == image_len, "wrote file into image buffer");
361 
362         mpi_result = MPI_Ssend((void *)(&image_len), (int)sizeof(ssize_t),
363 			       MPI_BYTE, (mpi_rank + 1) % mpi_size, 0,
364                                MPI_COMM_WORLD);
365 	VRFY((mpi_result == MPI_SUCCESS), \
366              "sent image size to process (mpi_rank + 1) % mpi_size");
367 
368         mpi_result = MPI_Ssend((void *)image_ptr, (int)image_len,
369 			       MPI_BYTE, (mpi_rank + 1) % mpi_size, 0,
370                                MPI_COMM_WORLD);
371 	VRFY((mpi_result == MPI_SUCCESS), \
372               "sent image to process (mpi_rank + 1) % mpi_size");
373 
374         HDfree(image_ptr);
375         image_ptr = NULL;
376         image_len = 0;
377 
378 	/* 5) close the core file and exit. */
379 
380         err = H5Sclose(space_id);
381 	VRFY((err >= 0), "closed data space.");
382 
383 	err = H5Dclose(dset_id);
384 	VRFY((err >= 0), "closed data set.");
385 
386 	err = H5Fclose(file_id);
387 	VRFY((err >= 0), "closed core file(1).");
388 
389 	err = H5Pclose(fapl_id);
390 	VRFY((err >= 0), "closed fapl(1).");
391     }
392 
393     return;
394 
395 } /* file_image_daisy_chain_test() */
396 
397