1 /* Test mhdf library on top of Parallel HDF5.
2  *
3  * Call mhdf API similar to how WriteHDF5Parallel does,
4  * but avoid any of our own parallel communiation.
5  *
6  * Output will be composed of:
7  *  a 1-D array of hexes, one per processor, where the first
8  *    process writes all the nodes for its hex and every other
9  *    process writes the right 4 nodes of its hex.
10  *  a set containing all of the hexes
11  *  a set containing all of the nodes
12  *  a set per processor containing all of the entities written by that proc
13  *  a tag specifying the process ID that wrote each entity.
14  */
15 
16 #include "mhdf.h"
17 
18 #include <time.h>
19 #include <stdlib.h>
20 #include <stdio.h>
21 #include <assert.h>
22 
23 #include <H5Ppublic.h>
24 #include <H5Tpublic.h>
25 #include <H5Epublic.h>
26 #include <H5FDmpi.h>
27 #include <H5FDmpio.h>
28 
29 #ifdef H5_HAVE_PARALLEL
30 
31 const char filename[] = "mhdf_ll.h5m";
32 const char proc_tag_name[] = "proc_id";
33 const char elem_handle[] = "Hex";
34 
35 int RANK;
36 int NUM_PROC;
37 
38 #define CHECK( A ) check( &(A), __LINE__ )
check(mhdf_Status * status,int line)39 void check( mhdf_Status* status, int line )
40 {
41   if (mhdf_isError( status )) {
42     fprintf(stderr,"%d: ERROR at line %d: \"%s\"\n", RANK, line, mhdf_message( status ) );
43     abort();
44   }
45 }
46 
47 typedef int handle_id_t;
48 #define id_type H5T_NATIVE_INT
49 
50 /* create file layout in serial */
create_file()51 void create_file()
52 {
53   const char* elem_types[1] = {elem_handle};
54   const int num_elem_types = sizeof(elem_types)/sizeof(elem_types[0]);
55   const int num_nodes = 4+4*NUM_PROC;
56   const int num_hexes = NUM_PROC;
57   const int num_sets = 2+NUM_PROC;
58   const int set_data_size = num_hexes + num_nodes + num_hexes + num_nodes;
59   const int num_tag_vals = num_hexes + num_nodes + num_sets - 2;
60   const char* history[] = { __FILE__, NULL };
61   const int history_size = sizeof(history) / sizeof(history[0]);
62   int default_tag_val = NUM_PROC;
63   hid_t data, tagdata[2];
64   long junk;
65 
66   mhdf_Status status;
67   mhdf_FileHandle file;
68 
69   time_t t;
70   t = time(NULL);
71   history[1] = ctime(&t);
72 
73   file = mhdf_createFile( filename, 1, elem_types, num_elem_types, id_type, &status );
74   CHECK(status);
75 
76   mhdf_writeHistory( file, history, history_size, &status );
77   CHECK(status);
78 
79   data = mhdf_createNodeCoords( file, 3, num_nodes, &junk, &status );
80   CHECK(status);
81   mhdf_closeData( file, data, &status );
82   CHECK(status);
83 
84   mhdf_addElement( file, elem_types[0], 0, &status );
85   CHECK(status);
86   data = mhdf_createConnectivity( file, elem_handle, 8, num_hexes, &junk, &status );
87   CHECK(status);
88   mhdf_closeData( file, data, &status );
89   CHECK(status);
90 
91   data = mhdf_createSetMeta( file, num_sets, &junk, &status );
92   CHECK(status);
93   mhdf_closeData( file, data, &status );
94   CHECK(status);
95 
96   data = mhdf_createSetData( file, set_data_size, &status );
97   CHECK(status);
98   mhdf_closeData( file, data, &status );
99   CHECK(status);
100 
101   mhdf_createTag( file, proc_tag_name, mhdf_INTEGER, 1, mhdf_DENSE_TYPE,
102                   &default_tag_val, &default_tag_val, 0, 0, &status );
103   CHECK(status);
104 
105   mhdf_createSparseTagData( file, proc_tag_name, num_tag_vals, tagdata, &status );
106   CHECK(status);
107   mhdf_closeData( file, tagdata[0], &status );
108   CHECK(status);
109   mhdf_closeData( file, tagdata[1], &status );
110   CHECK(status);
111 
112   mhdf_closeFile( file, &status );
113   CHECK(status);
114 }
115 
write_file_data()116 void write_file_data()
117 {
118   const int total_num_nodes = 4+4*NUM_PROC;
119   const int total_num_hexes = NUM_PROC;
120   long first_node, first_elem, first_set, count, ntag;
121   unsigned long ucount;
122   mhdf_index_t set_desc[4] = { 0, -1, -1, 0 };
123   hid_t handle, handles[2];
124   mhdf_Status status;
125   mhdf_FileHandle file;
126   int num_node, offset, dim;
127   double coords[3*8] = { 0.0, 0.0, 0.0,
128                          0.0, 1.0, 0.0,
129                          0.0, 1.0, 1.0,
130                          0.0, 0.0, 1.0,
131                          0.0, 0.0, 0.0,
132                          0.0, 1.0, 0.0,
133                          0.0, 1.0, 1.0,
134                          0.0, 0.0, 1.0 };
135   handle_id_t list[10];
136   int i, tagdata[10];
137   for (i = 0; i < 10; i++) tagdata[i] = RANK;
138 
139 
140     /* open file */
141   handle = H5Pcreate( H5P_FILE_ACCESS );
142   H5Pset_fapl_mpio( handle, MPI_COMM_WORLD, MPI_INFO_NULL );
143   file = mhdf_openFileWithOpt( filename, 1, &ucount, id_type, handle, &status );
144   CHECK(status);
145   H5Pclose( handle );
146 
147     /* write node coordinates */
148   if (0 == RANK) {
149     num_node = 8;
150     offset = 0;
151     coords[12] = coords[15] = coords[18] = coords[21] = 1.0;
152   }
153   else {
154     num_node = 4;
155     offset = 4 + 4*RANK;
156     coords[0] = coords[3] = coords[6] = coords[9] = 1.0 + RANK;
157   }
158   handle = mhdf_openNodeCoords( file, &count, &dim, &first_node, &status );
159   CHECK(status);
160   assert( count == total_num_nodes );
161   assert( dim == 3 );
162   mhdf_writeNodeCoords( handle, offset, num_node, coords, &status );
163   CHECK(status);
164   mhdf_closeData( file, handle, &status );
165   CHECK(status);
166 
167     /* write hex connectivity */
168   for (i = 0; i < 8; ++i)
169     list[i] = 4*RANK + i + first_node;
170   handle = mhdf_openConnectivity( file, elem_handle, &dim, &count, &first_elem, &status );
171   CHECK(status);
172   assert( count == total_num_hexes );
173   assert( 8 == dim );
174   mhdf_writeConnectivity( handle, RANK, 1, id_type, list, &status );
175   CHECK(status);
176   mhdf_closeData( file, handle, &status );
177   CHECK(status);
178 
179     /* write set descriptions */
180   handle = mhdf_openSetMeta( file, &count, &first_set, &status );
181   CHECK(status);
182   assert( count == 2+NUM_PROC );
183 
184     /* write set descriptions for per-proc sets */
185   set_desc[0] = total_num_nodes + total_num_hexes + 9 + 5*RANK - 1;
186   mhdf_writeSetMeta( handle, 2+RANK, 1, MHDF_INDEX_TYPE, set_desc, &status );
187   CHECK(status);
188 
189     /* write set descriptions for multi-proc sets */
190   if (0 == RANK) {
191     set_desc[0] = total_num_nodes - 1;
192     mhdf_writeSetMeta( handle, 0, 1, MHDF_INDEX_TYPE, set_desc, &status );
193     CHECK(status);
194     set_desc[0] += total_num_hexes;
195     mhdf_writeSetMeta( handle, 1, 1, MHDF_INDEX_TYPE, set_desc, &status );
196     CHECK(status);
197   }
198 
199   mhdf_closeData( file, handle, &status );
200   CHECK(status);
201 
202     /* write set contents */
203 
204   handle = mhdf_openSetData( file, &count, &status );
205   CHECK(status);
206   assert( 2*total_num_nodes + 2*total_num_hexes == count );
207 
208     /* write per-proc sets */
209   if (0 == RANK) {
210     count = 9;
211     offset = total_num_nodes + total_num_hexes;
212     for (i = 0; i < 8; ++i)
213       list[i] = first_node + i;
214     list[8] = first_elem;
215   }
216   else {
217     count = 5;
218     offset = total_num_nodes + total_num_hexes + 4 + 5*RANK;
219     for (i = 0; i < 4; ++i)
220       list[i] = 4 + 4*RANK + first_node + i;
221     list[4] = RANK + first_elem;
222   }
223   mhdf_writeSetData( handle, offset, count, id_type, list, &status );
224   CHECK(status);
225 
226     /* write multi-proc sets */
227     /* write nodes */
228   offset = (0 == RANK) ? 0 : 4 + 4*RANK;
229   mhdf_writeSetData( handle, offset, count-1, id_type, list, &status );
230   CHECK(status);
231     /* write hex */
232   offset = total_num_nodes + RANK;
233   mhdf_writeSetData( handle, offset, 1, id_type, list + count - 1, &status );
234   CHECK(status);
235 
236   mhdf_closeData( file, handle, &status );
237   CHECK(status);
238 
239     /* write tag data */
240   offset = (0 == RANK) ? 0 : 4 + 4*RANK;
241   offset += 2*RANK;
242   list[count++] = first_set + 2 + RANK;
243   mhdf_openSparseTagData( file, proc_tag_name, &ntag, &ntag, handles, &status );
244   CHECK(status);
245   assert( ntag == total_num_nodes + total_num_hexes + NUM_PROC );
246   mhdf_writeSparseTagEntities( handles[0], offset, count, id_type, list, &status );
247   CHECK(status);
248   mhdf_writeTagValues( handles[1], offset, count, H5T_NATIVE_INT, tagdata, &status );
249   CHECK(status);
250   mhdf_closeData( file, handles[0], &status );
251   CHECK(status);
252   mhdf_closeData( file, handles[1], &status );
253   CHECK(status);
254 
255     /* done */
256   mhdf_closeFile( file, &status );
257   CHECK(status);
258 }
259 
260 /* Define a dummy error handler to register with HDF5.
261  * This function doesn't do anything except pass the
262  * error on to the default handler that would have been
263  * called anyway.  It's only purpose is to provide a
264  * spot to set a break point so we can figure out where
265  * (in our code) that we made an invalid call into HDF5
266  */
267 #if defined(H5E_auto_t_vers) && H5E_auto_t_vers > 1
268 herr_t (*default_handler)( hid_t, void* );
handle_hdf5_error(hid_t stack,void * data)269 static herr_t handle_hdf5_error( hid_t stack, void* data )
270 #else
271 herr_t (*default_handler)( void* );
272 static herr_t handle_hdf5_error( void* data )
273 #endif
274 {
275   herr_t result = 0;
276   if (default_handler)
277 #if defined(H5E_auto_t_vers) && H5E_auto_t_vers > 1
278     result = (default_handler)(stack,data);
279 #else
280     result = (default_handler)(data);
281 #endif
282   assert(0);
283   return result;
284 }
285 
286 #endif /* #ifdef H5_HAVE_PARALLEL */
287 
main(int argc,char * argv[])288 int main( int argc, char* argv[] )
289 {
290 #ifdef H5_HAVE_PARALLEL
291   int rval;
292   void* data;
293   herr_t err;
294 
295 #if defined(H5Eget_auto_vers) && H5Eget_auto_vers > 1
296   err = H5Eget_auto( H5E_DEFAULT, &default_handler, &data );
297 #else
298   err = H5Eget_auto( &default_handler, &data );
299 #endif
300   if (err >= 0) {
301 #if defined(H5Eset_auto_vers) && H5Eset_auto_vers > 1
302     H5Eset_auto( H5E_DEFAULT, &handle_hdf5_error, data );
303 #else
304     H5Eset_auto( &handle_hdf5_error, data );
305 #endif
306   }
307 
308   rval = MPI_Init( &argc, &argv );
309   if (rval)
310     return rval;
311   rval = MPI_Comm_rank( MPI_COMM_WORLD, &RANK );
312   if (rval)
313     return rval;
314   rval = MPI_Comm_size( MPI_COMM_WORLD, &NUM_PROC );
315   if (rval)
316     return rval;
317 
318   if (RANK == 0)
319     create_file();
320 
321   /* Wait for rank 0 to finish creating the file, otherwise rank 1 may find it to be invalid */
322   rval = MPI_Barrier(MPI_COMM_WORLD);
323   if (rval)
324     return rval;
325 
326   write_file_data();
327 
328   MPI_Finalize();
329 #endif
330   return 0;
331 }
332