MOAB: Mesh Oriented datABase  (version 5.4.1)
mhdf_parallel.c
Go to the documentation of this file.
00001 /* Test mhdf library on top of Parallel HDF5.
00002  *
00003  * Call mhdf API similar to how WriteHDF5Parallel does,
00004  * but avoid any of our own parallel communiation.
00005  *
00006  * Output will be composed of:
00007  *  a 1-D array of hexes, one per processor, where the first
00008  *    process writes all the nodes for its hex and every other
00009  *    process writes the right 4 nodes of its hex.
00010  *  a set containing all of the hexes
00011  *  a set containing all of the nodes
00012  *  a set per processor containing all of the entities written by that proc
00013  *  a tag specifying the process ID that wrote each entity.
00014  */
00015 
00016 #include "mhdf.h"
00017 
00018 #include <time.h>
00019 #include <stdlib.h>
00020 #include <stdio.h>
00021 #include <assert.h>
00022 
00023 #include <H5Ppublic.h>
00024 #include <H5Tpublic.h>
00025 #include <H5Epublic.h>
00026 #include <H5FDmpi.h>
00027 #include <H5FDmpio.h>
00028 
00029 #ifdef H5_HAVE_PARALLEL
00030 
00031 const char filename[]      = "mhdf_ll.h5m";
00032 const char proc_tag_name[] = "proc_id";
00033 const char elem_handle[]   = "Hex";
00034 
00035 int RANK;
00036 int NUM_PROC;
00037 
00038 #define CHECK( A ) check( &( A ), __LINE__ )
00039 void check( mhdf_Status* status, int line )
00040 {
00041     if( mhdf_isError( status ) )
00042     {
00043         fprintf( stderr, "%d: ERROR at line %d: \"%s\"\n", RANK, line, mhdf_message( status ) );
00044         abort();
00045     }
00046 }
00047 
00048 typedef int handle_id_t;
00049 #define id_type H5T_NATIVE_INT
00050 
00051 /* create file layout in serial */
00052 void create_file()
00053 {
00054     const char* elem_types[1] = { elem_handle };
00055     const int num_elem_types  = sizeof( elem_types ) / sizeof( elem_types[0] );
00056     const int num_nodes       = 4 + 4 * NUM_PROC;
00057     const int num_hexes       = NUM_PROC;
00058     const int num_sets        = 2 + NUM_PROC;
00059     const int set_data_size   = num_hexes + num_nodes + num_hexes + num_nodes;
00060     const int num_tag_vals    = num_hexes + num_nodes + num_sets - 2;
00061     const char* history[]     = { __FILE__, NULL };
00062     const int history_size    = sizeof( history ) / sizeof( history[0] );
00063     int default_tag_val       = NUM_PROC;
00064     hid_t data, tagdata[2];
00065     long junk;
00066 
00067     mhdf_Status status;
00068     mhdf_FileHandle file;
00069 
00070     time_t t;
00071     t          = time( NULL );
00072     history[1] = ctime( &t );
00073 
00074     file = mhdf_createFile( filename, 1, elem_types, num_elem_types, id_type, &status );
00075     CHECK( status );
00076 
00077     mhdf_writeHistory( file, history, history_size, &status );
00078     CHECK( status );
00079 
00080     data = mhdf_createNodeCoords( file, 3, num_nodes, &junk, &status );
00081     CHECK( status );
00082     mhdf_closeData( file, data, &status );
00083     CHECK( status );
00084 
00085     mhdf_addElement( file, elem_types[0], 0, &status );
00086     CHECK( status );
00087     data = mhdf_createConnectivity( file, elem_handle, 8, num_hexes, &junk, &status );
00088     CHECK( status );
00089     mhdf_closeData( file, data, &status );
00090     CHECK( status );
00091 
00092     data = mhdf_createSetMeta( file, num_sets, &junk, &status );
00093     CHECK( status );
00094     mhdf_closeData( file, data, &status );
00095     CHECK( status );
00096 
00097     data = mhdf_createSetData( file, set_data_size, &status );
00098     CHECK( status );
00099     mhdf_closeData( file, data, &status );
00100     CHECK( status );
00101 
00102     mhdf_createTag( file, proc_tag_name, mhdf_INTEGER, 1, mhdf_DENSE_TYPE, &default_tag_val, &default_tag_val, 0, 0,
00103                     &status );
00104     CHECK( status );
00105 
00106     mhdf_createSparseTagData( file, proc_tag_name, num_tag_vals, tagdata, &status );
00107     CHECK( status );
00108     mhdf_closeData( file, tagdata[0], &status );
00109     CHECK( status );
00110     mhdf_closeData( file, tagdata[1], &status );
00111     CHECK( status );
00112 
00113     mhdf_closeFile( file, &status );
00114     CHECK( status );
00115 }
00116 
00117 void write_file_data()
00118 {
00119     const int total_num_nodes = 4 + 4 * NUM_PROC;
00120     const int total_num_hexes = NUM_PROC;
00121     long first_node, first_elem, first_set, count, ntag;
00122     unsigned long ucount;
00123     mhdf_index_t set_desc[4] = { 0, -1, -1, 0 };
00124     hid_t handle, handles[2];
00125     mhdf_Status status;
00126     mhdf_FileHandle file;
00127     int num_node, offset, dim;
00128     double coords[3 * 8] = { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 1.0, 0.0, 0.0, 1.0,
00129                              0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 1.0, 0.0, 0.0, 1.0 };
00130     handle_id_t list[10];
00131     int i, tagdata[10];
00132     for( i = 0; i < 10; i++ )
00133         tagdata[i] = RANK;
00134 
00135     /* open file */
00136     handle = H5Pcreate( H5P_FILE_ACCESS );
00137     H5Pset_fapl_mpio( handle, MPI_COMM_WORLD, MPI_INFO_NULL );
00138     file = mhdf_openFileWithOpt( filename, 1, &ucount, id_type, handle, &status );
00139     CHECK( status );
00140     H5Pclose( handle );
00141 
00142     /* write node coordinates */
00143     if( 0 == RANK )
00144     {
00145         num_node   = 8;
00146         offset     = 0;
00147         coords[12] = coords[15] = coords[18] = coords[21] = 1.0;
00148     }
00149     else
00150     {
00151         num_node  = 4;
00152         offset    = 4 + 4 * RANK;
00153         coords[0] = coords[3] = coords[6] = coords[9] = 1.0 + RANK;
00154     }
00155     handle = mhdf_openNodeCoords( file, &count, &dim, &first_node, &status );
00156     CHECK( status );
00157     assert( count == total_num_nodes );
00158     assert( dim == 3 );
00159     mhdf_writeNodeCoords( handle, offset, num_node, coords, &status );
00160     CHECK( status );
00161     mhdf_closeData( file, handle, &status );
00162     CHECK( status );
00163 
00164     /* write hex connectivity */
00165     for( i = 0; i < 8; ++i )
00166         list[i] = 4 * RANK + i + first_node;
00167     handle = mhdf_openConnectivity( file, elem_handle, &dim, &count, &first_elem, &status );
00168     CHECK( status );
00169     assert( count == total_num_hexes );
00170     assert( 8 == dim );
00171     mhdf_writeConnectivity( handle, RANK, 1, id_type, list, &status );
00172     CHECK( status );
00173     mhdf_closeData( file, handle, &status );
00174     CHECK( status );
00175 
00176     /* write set descriptions */
00177     handle = mhdf_openSetMeta( file, &count, &first_set, &status );
00178     CHECK( status );
00179     assert( count == 2 + NUM_PROC );
00180 
00181     /* write set descriptions for per-proc sets */
00182     set_desc[0] = total_num_nodes + total_num_hexes + 9 + 5 * RANK - 1;
00183     mhdf_writeSetMeta( handle, 2 + RANK, 1, MHDF_INDEX_TYPE, set_desc, &status );
00184     CHECK( status );
00185 
00186     /* write set descriptions for multi-proc sets */
00187     if( 0 == RANK )
00188     {
00189         set_desc[0] = total_num_nodes - 1;
00190         mhdf_writeSetMeta( handle, 0, 1, MHDF_INDEX_TYPE, set_desc, &status );
00191         CHECK( status );
00192         set_desc[0] += total_num_hexes;
00193         mhdf_writeSetMeta( handle, 1, 1, MHDF_INDEX_TYPE, set_desc, &status );
00194         CHECK( status );
00195     }
00196 
00197     mhdf_closeData( file, handle, &status );
00198     CHECK( status );
00199 
00200     /* write set contents */
00201 
00202     handle = mhdf_openSetData( file, &count, &status );
00203     CHECK( status );
00204     assert( 2 * total_num_nodes + 2 * total_num_hexes == count );
00205 
00206     /* write per-proc sets */
00207     if( 0 == RANK )
00208     {
00209         count  = 9;
00210         offset = total_num_nodes + total_num_hexes;
00211         for( i = 0; i < 8; ++i )
00212             list[i] = first_node + i;
00213         list[8] = first_elem;
00214     }
00215     else
00216     {
00217         count  = 5;
00218         offset = total_num_nodes + total_num_hexes + 4 + 5 * RANK;
00219         for( i = 0; i < 4; ++i )
00220             list[i] = 4 + 4 * RANK + first_node + i;
00221         list[4] = RANK + first_elem;
00222     }
00223     mhdf_writeSetData( handle, offset, count, id_type, list, &status );
00224     CHECK( status );
00225 
00226     /* write multi-proc sets */
00227     /* write nodes */
00228     offset = ( 0 == RANK ) ? 0 : 4 + 4 * RANK;
00229     mhdf_writeSetData( handle, offset, count - 1, id_type, list, &status );
00230     CHECK( status );
00231     /* write hex */
00232     offset = total_num_nodes + RANK;
00233     mhdf_writeSetData( handle, offset, 1, id_type, list + count - 1, &status );
00234     CHECK( status );
00235 
00236     mhdf_closeData( file, handle, &status );
00237     CHECK( status );
00238 
00239     /* write tag data */
00240     offset = ( 0 == RANK ) ? 0 : 4 + 4 * RANK;
00241     offset += 2 * RANK;
00242     list[count++] = first_set + 2 + RANK;
00243     mhdf_openSparseTagData( file, proc_tag_name, &ntag, &ntag, handles, &status );
00244     CHECK( status );
00245     assert( ntag == total_num_nodes + total_num_hexes + NUM_PROC );
00246     mhdf_writeSparseTagEntities( handles[0], offset, count, id_type, list, &status );
00247     CHECK( status );
00248     mhdf_writeTagValues( handles[1], offset, count, H5T_NATIVE_INT, tagdata, &status );
00249     CHECK( status );
00250     mhdf_closeData( file, handles[0], &status );
00251     CHECK( status );
00252     mhdf_closeData( file, handles[1], &status );
00253     CHECK( status );
00254 
00255     /* done */
00256     mhdf_closeFile( file, &status );
00257     CHECK( status );
00258 }
00259 
00260 /* Define a dummy error handler to register with HDF5. 
00261  * This function doesn't do anything except pass the
00262  * error on to the default handler that would have been
00263  * called anyway.  It's only purpose is to provide a 
00264  * spot to set a break point so we can figure out where
00265  * (in our code) that we made an invalid call into HDF5
00266  */
00267 #if defined( H5E_auto_t_vers ) && H5E_auto_t_vers > 1
00268 herr_t ( *default_handler )( hid_t, void* );
00269 static herr_t handle_hdf5_error( hid_t stack, void* data )
00270 #else
00271 herr_t ( *default_handler )( void* );
00272 static herr_t handle_hdf5_error( void* data )
00273 #endif
00274 {
00275     herr_t result = 0;
00276     if( default_handler )
00277 #if defined( H5E_auto_t_vers ) && H5E_auto_t_vers > 1
00278         result = (default_handler)( stack, data );
00279 #else
00280         result = (default_handler)( data );
00281 #endif
00282     assert( 0 );
00283     return result;
00284 }
00285 
00286 #endif /* #ifdef H5_HAVE_PARALLEL */
00287 
00288 int main( int argc, char* argv[] )
00289 {
00290 #ifdef H5_HAVE_PARALLEL
00291     int rval;
00292     void* data;
00293     herr_t err;
00294 
00295 #if defined( H5Eget_auto_vers ) && H5Eget_auto_vers > 1
00296     err = H5Eget_auto( H5E_DEFAULT, &default_handler, &data );
00297 #else
00298     err = H5Eget_auto( &default_handler, &data );
00299 #endif
00300     if( err >= 0 )
00301     {
00302 #if defined( H5Eset_auto_vers ) && H5Eset_auto_vers > 1
00303         H5Eset_auto( H5E_DEFAULT, &handle_hdf5_error, data );
00304 #else
00305         H5Eset_auto( &handle_hdf5_error, data );
00306 #endif
00307     }
00308 
00309     rval = MPI_Init( &argc, &argv );
00310     if( rval ) return rval;
00311     rval = MPI_Comm_rank( MPI_COMM_WORLD, &RANK );
00312     if( rval ) return rval;
00313     rval = MPI_Comm_size( MPI_COMM_WORLD, &NUM_PROC );
00314     if( rval ) return rval;
00315 
00316     if( RANK == 0 ) create_file();
00317 
00318     /* Wait for rank 0 to finish creating the file, otherwise rank 1 may find it to be invalid */
00319     rval = MPI_Barrier( MPI_COMM_WORLD );
00320     if( rval ) return rval;
00321 
00322     write_file_data();
00323 
00324     MPI_Finalize();
00325 #endif
00326     return 0;
00327 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines