MOAB: Mesh Oriented datABase  (version 5.4.1)
gcrm_par.cpp
Go to the documentation of this file.
00001 #include "TestUtil.hpp"
00002 #include "moab/Core.hpp"
00003 #include "moab/ParallelComm.hpp"
00004 #include "moab/ProgOptions.hpp"
00005 #include "MBParallelConventions.h"
00006 #include "moab/ReadUtilIface.hpp"
00007 #include "MBTagConventions.hpp"
00008 
00009 #include <sstream>
00010 
00011 using namespace moab;
00012 
00013 std::string example = TestDir + "unittest/io/gcrm_r3.nc";
00014 
00015 void test_read_onevar_trivial();
00016 #if defined( MOAB_HAVE_MPI ) && defined( MOAB_HAVE_ZOLTAN )
00017 void test_read_onevar_rcbzoltan();
00018 #endif
00019 
00020 void test_read_mesh_parallel_trivial();
00021 #if defined( MOAB_HAVE_MPI ) && defined( MOAB_HAVE_ZOLTAN )
00022 void test_read_mesh_parallel_rcbzoltan();
00023 #endif
00024 
00025 void test_gather_onevar_on_rank0();
00026 void test_gather_onevar_on_rank1();
00027 
00028 void test_multiple_loads_of_same_file();
00029 
00030 // Helper functions
00031 void read_one_cell_var( bool rcbzoltan );
00032 void read_mesh_parallel( bool rcbzoltan );
00033 void gather_one_cell_var( int gather_set_rank );
00034 void multiple_loads_of_same_file();
00035 
00036 std::string read_options;
00037 const double eps = 1e-6;
00038 const int layers = 3;
00039 
00040 int main( int argc, char* argv[] )
00041 {
00042     MPI_Init( &argc, &argv );
00043     int result = 0;
00044 
00045     result += RUN_TEST( test_read_onevar_trivial );
00046 #if defined( MOAB_HAVE_MPI ) && defined( MOAB_HAVE_ZOLTAN )
00047     result += RUN_TEST( test_read_onevar_rcbzoltan );
00048 #endif
00049 
00050     result += RUN_TEST( test_read_mesh_parallel_trivial );
00051 #if defined( MOAB_HAVE_MPI ) && defined( MOAB_HAVE_ZOLTAN )
00052     result += RUN_TEST( test_read_mesh_parallel_rcbzoltan );
00053 #endif
00054 
00055     result += RUN_TEST( test_gather_onevar_on_rank0 );
00056     result += RUN_TEST( test_gather_onevar_on_rank1 );
00057 
00058     result += RUN_TEST( test_multiple_loads_of_same_file );
00059 
00060     MPI_Finalize();
00061     return result;
00062 }
00063 
00064 void test_read_onevar_trivial()
00065 {
00066     read_one_cell_var( false );
00067 }
00068 
00069 void test_read_onevar_rcbzoltan()
00070 {
00071     read_one_cell_var( true );
00072 }
00073 
00074 void test_read_mesh_parallel_trivial()
00075 {
00076     read_mesh_parallel( false );
00077 }
00078 
00079 void test_read_mesh_parallel_rcbzoltan()
00080 {
00081     read_mesh_parallel( true );
00082 }
00083 
00084 void test_gather_onevar_on_rank0()
00085 {
00086     gather_one_cell_var( 0 );
00087 }
00088 
00089 void test_gather_onevar_on_rank1()
00090 {
00091     gather_one_cell_var( 1 );
00092 }
00093 
00094 void test_multiple_loads_of_same_file()
00095 {
00096     multiple_loads_of_same_file();
00097 }
00098 
00099 // Helper functions
00100 void read_one_cell_var( bool rcbzoltan )
00101 {
00102     Core moab;
00103     Interface& mb = moab;
00104 
00105     read_options = "PARALLEL=READ_PART;PARTITION_METHOD=TRIVIAL;NO_EDGES;VARIABLE=vorticity";
00106     if( rcbzoltan )
00107         read_options = "PARALLEL=READ_PART;PARTITION_METHOD=RCBZOLTAN;NO_EDGES;VARIABLE=vorticity;DEBUG_IO=1";
00108 
00109     ErrorCode rval = mb.load_file( example.c_str(), NULL, read_options.c_str() );CHECK_ERR( rval );
00110 
00111     // Get local edges
00112     Range local_edges;
00113     rval = mb.get_entities_by_type( 0, MBEDGE, local_edges );CHECK_ERR( rval );
00114     CHECK_EQUAL( (size_t)0, local_edges.size() );
00115 
00116     // Get local cells
00117     Range local_cells;
00118     rval = mb.get_entities_by_type( 0, MBPOLYGON, local_cells );CHECK_ERR( rval );
00119     // No mixed elements
00120     CHECK_EQUAL( (size_t)1, local_cells.psize() );
00121 
00122     Tag gid_tag = mb.globalId_tag();
00123 
00124     std::vector< int > gids( local_cells.size() );
00125     rval = mb.tag_get_data( gid_tag, local_cells, &gids[0] );
00126     Range local_cell_gids;
00127     std::copy( gids.rbegin(), gids.rend(), range_inserter( local_cell_gids ) );
00128 
00129     ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 );
00130     int procs           = pcomm->proc_config().proc_size();
00131     int rank            = pcomm->proc_config().proc_rank();
00132 
00133     // Make check runs this test on two processors
00134     if( 2 == procs )
00135     {
00136         // Check tag for cell variable vorticity at timestep 0
00137         Tag vorticity_tag0;
00138         rval = mb.tag_get_handle( "vorticity0", layers, MB_TYPE_DOUBLE, vorticity_tag0 );CHECK_ERR( rval );
00139 
00140         // Check tag for cell variable vorticity at timestep 1
00141         Tag vorticity_tag1;
00142         rval = mb.tag_get_handle( "vorticity1", layers, MB_TYPE_DOUBLE, vorticity_tag1 );CHECK_ERR( rval );
00143 
00144         // Get vorticity0 and vorticity1 tag values on 3 local cells
00145         double vorticity0_val[3 * layers];
00146         double vorticity1_val[3 * layers];
00147 
00148         if( rcbzoltan )
00149         {
00150             CHECK_EQUAL( (size_t)14, local_cell_gids.psize() );
00151 
00152             if( 0 == rank )
00153             {
00154                 CHECK_EQUAL( (size_t)321, local_cells.size() );
00155                 CHECK_EQUAL( (size_t)321, local_cell_gids.size() );
00156 
00157                 CHECK_EQUAL( 3, (int)local_cell_gids[0] );
00158                 CHECK_EQUAL( 162, (int)local_cell_gids[159] );
00159                 CHECK_EQUAL( 640, (int)local_cell_gids[318] );
00160 
00161                 EntityHandle cell_ents[] = { local_cells[0], local_cells[159], local_cells[318] };
00162                 rval                     = mb.tag_get_data( vorticity_tag0, cell_ents, 3, vorticity0_val );CHECK_ERR( rval );
00163 
00164                 // Timestep 0
00165                 // Layer 0
00166                 CHECK_REAL_EQUAL( -0.725999, vorticity0_val[0 * layers], eps );
00167                 CHECK_REAL_EQUAL( -1.814997, vorticity0_val[1 * layers], eps );
00168                 CHECK_REAL_EQUAL( 0.131688, vorticity0_val[2 * layers], eps );
00169                 // Layer 1
00170                 CHECK_REAL_EQUAL( -0.725989, vorticity0_val[0 * layers + 1], eps );
00171                 CHECK_REAL_EQUAL( -1.814972, vorticity0_val[1 * layers + 1], eps );
00172                 CHECK_REAL_EQUAL( 0.131686, vorticity0_val[2 * layers + 1], eps );
00173 
00174                 rval = mb.tag_get_data( vorticity_tag1, cell_ents, 3, vorticity1_val );CHECK_ERR( rval );
00175 
00176                 // Timestep 1
00177                 // Layer 0
00178                 CHECK_REAL_EQUAL( -0.706871, vorticity1_val[0 * layers], eps );
00179                 CHECK_REAL_EQUAL( -1.767178, vorticity1_val[1 * layers], eps );
00180                 CHECK_REAL_EQUAL( 0.128218, vorticity1_val[2 * layers], eps );
00181                 // Layer 1
00182                 CHECK_REAL_EQUAL( -0.706861, vorticity1_val[0 * layers + 1], eps );
00183                 CHECK_REAL_EQUAL( -1.767153, vorticity1_val[1 * layers + 1], eps );
00184                 CHECK_REAL_EQUAL( 0.128216, vorticity1_val[2 * layers + 1], eps );
00185             }
00186             else if( 1 == rank )
00187             {
00188                 CHECK_EQUAL( (size_t)321, local_cells.size() );
00189                 CHECK_EQUAL( (size_t)321, local_cell_gids.size() );
00190 
00191                 CHECK_EQUAL( 1, (int)local_cell_gids[0] );
00192                 CHECK_EQUAL( 366, (int)local_cell_gids[161] );
00193                 CHECK_EQUAL( 1, (int)local_cell_gids[322] );
00194 
00195                 EntityHandle cell_ents[] = { local_cells[0], local_cells[161], local_cells[322] };
00196                 rval                     = mb.tag_get_data( vorticity_tag0, cell_ents, 3, vorticity0_val );CHECK_ERR( rval );
00197 
00198                 // Timestep 0
00199                 // Layer 0
00200                 CHECK_REAL_EQUAL( 3.629994, vorticity0_val[0 * layers], eps );
00201                 CHECK_REAL_EQUAL( -1.512985, vorticity0_val[1 * layers], eps );
00202                 CHECK_REAL_EQUAL( 3.629994, vorticity0_val[2 * layers], eps );
00203                 // Layer 1
00204                 CHECK_REAL_EQUAL( 3.629944, vorticity0_val[0 * layers + 1], eps );
00205                 CHECK_REAL_EQUAL( -1.512964, vorticity0_val[1 * layers + 1], eps );
00206                 CHECK_REAL_EQUAL( 3.629944, vorticity0_val[2 * layers + 1], eps );
00207 
00208                 rval = mb.tag_get_data( vorticity_tag1, cell_ents, 3, vorticity1_val );CHECK_ERR( rval );
00209 
00210                 // Timestep 1
00211                 // Layer 0
00212                 CHECK_REAL_EQUAL( 3.534355, vorticity1_val[0 * layers], eps );
00213                 CHECK_REAL_EQUAL( -1.473122, vorticity1_val[1 * layers], eps );
00214                 CHECK_REAL_EQUAL( 3.534355, vorticity1_val[2 * layers], eps );
00215                 // Layer 1
00216                 CHECK_REAL_EQUAL( 3.534306, vorticity1_val[0 * layers + 1], eps );
00217                 CHECK_REAL_EQUAL( -1.473102, vorticity1_val[1 * layers + 1], eps );
00218                 CHECK_REAL_EQUAL( 3.534306, vorticity1_val[2 * layers + 1], eps );
00219             }
00220         }
00221         else
00222         {
00223             CHECK_EQUAL( (size_t)321, local_cells.size() );
00224             CHECK_EQUAL( (size_t)321, local_cell_gids.size() );
00225             CHECK_EQUAL( (size_t)1, local_cell_gids.psize() );
00226 
00227             EntityHandle cell_ents[] = { local_cells[0], local_cells[160], local_cells[320] };
00228             rval                     = mb.tag_get_data( vorticity_tag0, cell_ents, 3, vorticity0_val );CHECK_ERR( rval );
00229 
00230             rval = mb.tag_get_data( vorticity_tag1, cell_ents, 3, vorticity1_val );CHECK_ERR( rval );
00231 
00232             if( 0 == rank )
00233             {
00234                 CHECK_EQUAL( 1, (int)local_cell_gids[0] );
00235                 CHECK_EQUAL( 161, (int)local_cell_gids[160] );
00236                 CHECK_EQUAL( 321, (int)local_cell_gids[320] );
00237 
00238                 // Timestep 0
00239                 // Layer 0
00240                 CHECK_REAL_EQUAL( 3.629994, vorticity0_val[0 * layers], eps );
00241                 CHECK_REAL_EQUAL( -1.708188, vorticity0_val[1 * layers], eps );
00242                 CHECK_REAL_EQUAL( 0.131688, vorticity0_val[2 * layers], eps );
00243                 // Layer 1
00244                 CHECK_REAL_EQUAL( 3.629944, vorticity0_val[0 * layers + 1], eps );
00245                 CHECK_REAL_EQUAL( -1.708164, vorticity0_val[1 * layers + 1], eps );
00246                 CHECK_REAL_EQUAL( 0.131686, vorticity0_val[2 * layers + 1], eps );
00247 
00248                 // Timestep 1
00249                 // Layer 0
00250                 CHECK_REAL_EQUAL( 3.534355, vorticity1_val[0 * layers], eps );
00251                 CHECK_REAL_EQUAL( -1.663182, vorticity1_val[1 * layers], eps );
00252                 CHECK_REAL_EQUAL( 0.128218, vorticity1_val[2 * layers], eps );
00253                 // Layer 1
00254                 CHECK_REAL_EQUAL( 3.534306, vorticity1_val[0 * layers + 1], eps );
00255                 CHECK_REAL_EQUAL( -1.663160, vorticity1_val[1 * layers + 1], eps );
00256                 CHECK_REAL_EQUAL( 0.128216, vorticity1_val[2 * layers + 1], eps );
00257             }
00258             else if( 1 == rank )
00259             {
00260                 CHECK_EQUAL( 322, (int)local_cell_gids[0] );
00261                 CHECK_EQUAL( 482, (int)local_cell_gids[160] );
00262                 CHECK_EQUAL( 642, (int)local_cell_gids[320] );
00263 
00264                 // Timestep 0
00265                 // Layer 0
00266                 CHECK_REAL_EQUAL( -0.554888, vorticity0_val[0 * layers], eps );
00267                 CHECK_REAL_EQUAL( 2.434397, vorticity0_val[1 * layers], eps );
00268                 CHECK_REAL_EQUAL( -0.554888, vorticity0_val[2 * layers], eps );
00269                 // Layer 1
00270                 CHECK_REAL_EQUAL( -0.554881, vorticity0_val[0 * layers + 1], eps );
00271                 CHECK_REAL_EQUAL( 2.434363, vorticity0_val[1 * layers + 1], eps );
00272                 CHECK_REAL_EQUAL( -0.554881, vorticity0_val[2 * layers + 1], eps );
00273 
00274                 // Timestep 1
00275                 // Layer 0
00276                 CHECK_REAL_EQUAL( -0.540269, vorticity1_val[0 * layers], eps );
00277                 CHECK_REAL_EQUAL( 2.370258, vorticity1_val[1 * layers], eps );
00278                 CHECK_REAL_EQUAL( -0.540269, vorticity1_val[2 * layers], eps );
00279                 // Layer 1
00280                 CHECK_REAL_EQUAL( -0.540262, vorticity1_val[0 * layers + 1], eps );
00281                 CHECK_REAL_EQUAL( 2.370226, vorticity1_val[1 * layers + 1], eps );
00282                 CHECK_REAL_EQUAL( -0.540262, vorticity1_val[2 * layers + 1], eps );
00283             }
00284         }
00285     }
00286 }
00287 
00288 void read_mesh_parallel( bool rcbzoltan )
00289 {
00290     Core moab;
00291     Interface& mb = moab;
00292 
00293     read_options = "PARALLEL=READ_PART;PARTITION_METHOD=TRIVIAL;PARALLEL_RESOLVE_SHARED_ENTS;VARIABLE=";
00294     if( rcbzoltan )
00295         read_options = "PARALLEL=READ_PART;PARTITION_METHOD=RCBZOLTAN;PARALLEL_RESOLVE_SHARED_ENTS;VARIABLE=";
00296 
00297     ErrorCode rval = mb.load_file( example.c_str(), NULL, read_options.c_str() );CHECK_ERR( rval );
00298 
00299     ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 );
00300     int procs           = pcomm->proc_config().proc_size();
00301     int rank            = pcomm->proc_config().proc_rank();
00302 
00303     rval = pcomm->check_all_shared_handles();CHECK_ERR( rval );
00304 
00305     // Get local vertices
00306     Range local_verts;
00307     rval = mb.get_entities_by_type( 0, MBVERTEX, local_verts );CHECK_ERR( rval );
00308 
00309     int verts_num = local_verts.size();
00310     if( 2 == procs )
00311     {
00312         if( rcbzoltan )
00313         {
00314             if( 0 == rank )
00315                 CHECK_EQUAL( 688, verts_num );
00316             else if( 1 == rank )
00317                 CHECK_EQUAL( 687, verts_num );  // Not owned vertices included
00318         }
00319         else
00320         {
00321             if( 0 == rank )
00322                 CHECK_EQUAL( 687, verts_num );
00323             else if( 1 == rank )
00324                 CHECK_EQUAL( 688, verts_num );  // Not owned vertices included
00325         }
00326     }
00327 
00328     rval = pcomm->filter_pstatus( local_verts, PSTATUS_NOT_OWNED, PSTATUS_NOT );CHECK_ERR( rval );
00329 
00330     verts_num = local_verts.size();
00331     if( 2 == procs )
00332     {
00333         if( rcbzoltan )
00334         {
00335             if( 0 == rank )
00336                 CHECK_EQUAL( 688, verts_num );
00337             else if( 1 == rank )
00338                 CHECK_EQUAL( 592, verts_num );  // Not owned vertices excluded
00339         }
00340         else
00341         {
00342             if( 0 == rank )
00343                 CHECK_EQUAL( 687, verts_num );
00344             else if( 1 == rank )
00345                 CHECK_EQUAL( 593, verts_num );  // Not owned vertices excluded
00346         }
00347     }
00348 
00349     // Get local edges
00350     Range local_edges;
00351     rval = mb.get_entities_by_type( 0, MBEDGE, local_edges );CHECK_ERR( rval );
00352 
00353     int edges_num = local_edges.size();
00354     if( 2 == procs )
00355     {
00356         if( rcbzoltan )
00357         {
00358             if( 0 == rank )
00359                 CHECK_EQUAL( 1008, edges_num );
00360             else if( 1 == rank )
00361                 CHECK_EQUAL( 1007, edges_num );  // Not owned edges included
00362         }
00363         else
00364         {
00365             if( 0 == rank )
00366                 CHECK_EQUAL( 1007, edges_num );
00367             else if( 1 == rank )
00368                 CHECK_EQUAL( 1008, edges_num );  // Not owned edges included
00369         }
00370     }
00371 
00372     rval = pcomm->filter_pstatus( local_edges, PSTATUS_NOT_OWNED, PSTATUS_NOT );CHECK_ERR( rval );
00373 
00374     edges_num = local_edges.size();
00375     if( 2 == procs )
00376     {
00377         if( rcbzoltan )
00378         {
00379             if( 0 == rank )
00380                 CHECK_EQUAL( 1008, edges_num );
00381             else if( 1 == rank )
00382                 CHECK_EQUAL( 912, edges_num );  // Not owned edges excluded
00383         }
00384         else
00385         {
00386             if( 0 == rank )
00387                 CHECK_EQUAL( 1007, edges_num );
00388             else if( 1 == rank )
00389                 CHECK_EQUAL( 913, edges_num );  // Not owned edges excluded
00390         }
00391     }
00392 
00393     // Get local cells
00394     Range local_cells;
00395     rval = mb.get_entities_by_type( 0, MBPOLYGON, local_cells );CHECK_ERR( rval );
00396     // No mixed elements
00397     CHECK_EQUAL( (size_t)1, local_cells.psize() );
00398 
00399     int cells_num = local_cells.size();
00400     if( 2 == procs )
00401     {
00402         if( rcbzoltan )
00403         {
00404             if( 0 == rank )
00405                 CHECK_EQUAL( 321, cells_num );
00406             else
00407                 CHECK_EQUAL( 321, cells_num );
00408         }
00409         else
00410             CHECK_EQUAL( 321, cells_num );
00411     }
00412 
00413     rval = pcomm->filter_pstatus( local_cells, PSTATUS_NOT_OWNED, PSTATUS_NOT );CHECK_ERR( rval );
00414 
00415     cells_num = local_cells.size();
00416     if( 2 == procs )
00417     {
00418         if( rcbzoltan )
00419         {
00420             if( 0 == rank )
00421                 CHECK_EQUAL( 321, cells_num );
00422             else
00423                 CHECK_EQUAL( 321, cells_num );
00424         }
00425         else
00426             CHECK_EQUAL( 321, cells_num );
00427     }
00428 
00429     std::cout << "proc: " << rank << " verts:" << verts_num << "\n";
00430 
00431     int total_verts_num;
00432     MPI_Reduce( &verts_num, &total_verts_num, 1, MPI_INT, MPI_SUM, 0, pcomm->proc_config().proc_comm() );
00433     if( 0 == rank )
00434     {
00435         std::cout << "total vertices: " << total_verts_num << "\n";
00436         CHECK_EQUAL( 1280, total_verts_num );
00437     }
00438 
00439     std::cout << "proc: " << rank << " edges:" << edges_num << "\n";
00440 
00441     int total_edges_num;
00442     MPI_Reduce( &edges_num, &total_edges_num, 1, MPI_INT, MPI_SUM, 0, pcomm->proc_config().proc_comm() );
00443     if( 0 == rank )
00444     {
00445         std::cout << "total edges: " << total_edges_num << "\n";
00446         CHECK_EQUAL( 1920, total_edges_num );
00447     }
00448 
00449     std::cout << "proc: " << rank << " cells:" << cells_num << "\n";
00450 
00451     int total_cells_num;
00452     MPI_Reduce( &cells_num, &total_cells_num, 1, MPI_INT, MPI_SUM, 0, pcomm->proc_config().proc_comm() );
00453     if( 0 == rank )
00454     {
00455         std::cout << "total cells: " << total_cells_num << "\n";
00456         CHECK_EQUAL( 642, total_cells_num );
00457     }
00458 
00459 #ifdef MOAB_HAVE_HDF5_PARALLEL
00460     std::string write_options( "PARALLEL=WRITE_PART;" );
00461 
00462     std::string output_file = "test_gcrm";
00463     if( rcbzoltan ) output_file += "_rcbzoltan";
00464     output_file += ".h5m";
00465 
00466     mb.write_file( output_file.c_str(), NULL, write_options.c_str() );
00467 #endif
00468 }
00469 
00470 void gather_one_cell_var( int gather_set_rank )
00471 {
00472     Core moab;
00473     Interface& mb = moab;
00474 
00475     EntityHandle file_set;
00476     ErrorCode rval = mb.create_meshset( MESHSET_SET, file_set );CHECK_ERR( rval );
00477 
00478     read_options = "PARALLEL=READ_PART;PARTITION_METHOD=TRIVIAL;PARALLEL_RESOLVE_SHARED_ENTS";
00479     std::ostringstream gather_set_option;
00480     gather_set_option << ";GATHER_SET=" << gather_set_rank;
00481     read_options += gather_set_option.str();
00482 
00483     rval = mb.load_file( example.c_str(), &file_set, read_options.c_str() );CHECK_ERR( rval );
00484 
00485     ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 );
00486     int procs           = pcomm->proc_config().proc_size();
00487     int rank            = pcomm->proc_config().proc_rank();
00488 
00489     // Make sure gather_set_rank is valid
00490     if( gather_set_rank < 0 || gather_set_rank >= procs ) return;
00491 
00492     Range cells, cells_owned;
00493     rval = mb.get_entities_by_type( file_set, MBPOLYGON, cells );CHECK_ERR( rval );
00494 
00495     // Get local owned cells
00496     rval = pcomm->filter_pstatus( cells, PSTATUS_NOT_OWNED, PSTATUS_NOT, -1, &cells_owned );CHECK_ERR( rval );
00497 
00498     EntityHandle gather_set = 0;
00499     if( gather_set_rank == rank )
00500     {
00501         // Get gather set
00502         ReadUtilIface* readUtilIface;
00503         mb.query_interface( readUtilIface );
00504         rval = readUtilIface->get_gather_set( gather_set );CHECK_ERR( rval );
00505         assert( gather_set != 0 );
00506     }
00507 
00508     Tag vorticity_tag0, gid_tag;
00509     rval = mb.tag_get_handle( "vorticity0", layers, MB_TYPE_DOUBLE, vorticity_tag0, MB_TAG_DENSE );CHECK_ERR( rval );
00510 
00511     gid_tag = mb.globalId_tag();
00512 
00513     pcomm->gather_data( cells_owned, vorticity_tag0, gid_tag, gather_set, gather_set_rank );
00514 
00515     if( gather_set_rank == rank )
00516     {
00517         // Get gather set cells
00518         Range gather_set_cells;
00519         rval = mb.get_entities_by_type( gather_set, MBPOLYGON, gather_set_cells );CHECK_ERR( rval );
00520         CHECK_EQUAL( (size_t)642, gather_set_cells.size() );
00521         CHECK_EQUAL( (size_t)1, gather_set_cells.psize() );
00522 
00523         // Check vorticity0 tag values on 4 gather set cells: first cell, two median cells, and last
00524         // cell
00525         EntityHandle cell_ents[] = { gather_set_cells[0], gather_set_cells[320], gather_set_cells[321],
00526                                      gather_set_cells[641] };
00527         double vorticity0_val[4 * layers];
00528         rval = mb.tag_get_data( vorticity_tag0, &cell_ents[0], 4, vorticity0_val );CHECK_ERR( rval );
00529 
00530         // Only check first two layers
00531         // Layer 0
00532         CHECK_REAL_EQUAL( 3.629994, vorticity0_val[0 * layers], eps );
00533         CHECK_REAL_EQUAL( 0.131688, vorticity0_val[1 * layers], eps );
00534         CHECK_REAL_EQUAL( -0.554888, vorticity0_val[2 * layers], eps );
00535         CHECK_REAL_EQUAL( -0.554888, vorticity0_val[3 * layers], eps );
00536         // Layer 1
00537         CHECK_REAL_EQUAL( 3.629944, vorticity0_val[0 * layers + 1], eps );
00538         CHECK_REAL_EQUAL( 0.131686, vorticity0_val[1 * layers + 1], eps );
00539         CHECK_REAL_EQUAL( -0.554881, vorticity0_val[2 * layers + 1], eps );
00540         CHECK_REAL_EQUAL( -0.554881, vorticity0_val[3 * layers + 1], eps );
00541     }
00542 }
00543 
00544 void multiple_loads_of_same_file()
00545 {
00546     Core moab;
00547     Interface& mb = moab;
00548 
00549     // Need a file set for nomesh to work right
00550     EntityHandle file_set;
00551     ErrorCode rval;
00552     rval = mb.create_meshset( MESHSET_SET, file_set );CHECK_ERR( rval );
00553 
00554     // Read first only header information, no mesh, no variable
00555     read_options = "PARALLEL=READ_PART;PARTITION;NOMESH;VARIABLE=;PARTITION_METHOD=TRIVIAL";
00556 
00557     rval = mb.load_file( example.c_str(), &file_set, read_options.c_str() );CHECK_ERR( rval );
00558 
00559     // Create mesh, no variable
00560     read_options = "PARALLEL=READ_PART;PARTITION;PARALLEL_RESOLVE_SHARED_ENTS;PARTITION_METHOD="
00561                    "TRIVIAL;VARIABLE=";
00562 
00563     rval = mb.load_file( example.c_str(), &file_set, read_options.c_str() );CHECK_ERR( rval );
00564 
00565     // Read variable vorticity at timestep 0, no mesh
00566     read_options = "PARALLEL=READ_PART;PARTITION;PARTITION_METHOD=TRIVIAL;NOMESH;VARIABLE="
00567                    "vorticity;TIMESTEP=0";
00568 
00569     rval = mb.load_file( example.c_str(), &file_set, read_options.c_str() );CHECK_ERR( rval );
00570 
00571     Range local_verts;
00572     rval = mb.get_entities_by_type( file_set, MBVERTEX, local_verts );CHECK_ERR( rval );
00573 
00574     Range local_edges;
00575     rval = mb.get_entities_by_type( file_set, MBEDGE, local_edges );CHECK_ERR( rval );
00576 
00577     Range local_cells;
00578     rval = mb.get_entities_by_type( file_set, MBPOLYGON, local_cells );CHECK_ERR( rval );
00579     // No mixed elements
00580     CHECK_EQUAL( (size_t)1, local_cells.psize() );
00581 
00582     ParallelComm* pcomm = ParallelComm::get_pcomm( &mb, 0 );
00583     int procs           = pcomm->proc_config().proc_size();
00584     int rank            = pcomm->proc_config().proc_rank();
00585 
00586     // Make check runs this test on two processors
00587     if( 2 == procs )
00588     {
00589         CHECK_EQUAL( (size_t)321, local_cells.size() );
00590 
00591         // Check tag for cell variable vorticity at timestep 0
00592         Tag vorticity_tag0;
00593         rval = mb.tag_get_handle( "vorticity0", layers, MB_TYPE_DOUBLE, vorticity_tag0 );CHECK_ERR( rval );
00594 
00595         // Get vorticity0 tag values on 3 local cells
00596         double vorticity0_val[3 * layers];
00597         EntityHandle cell_ents[] = { local_cells[0], local_cells[160], local_cells[320] };
00598         rval                     = mb.tag_get_data( vorticity_tag0, cell_ents, 3, vorticity0_val );CHECK_ERR( rval );
00599 
00600         if( 0 == rank )
00601         {
00602             CHECK_EQUAL( (size_t)687, local_verts.size() );
00603             CHECK_EQUAL( (size_t)1007, local_edges.size() );
00604 
00605             // Layer 0
00606             CHECK_REAL_EQUAL( 3.629994, vorticity0_val[0 * layers], eps );
00607             CHECK_REAL_EQUAL( -1.708188, vorticity0_val[1 * layers], eps );
00608             CHECK_REAL_EQUAL( 0.131688, vorticity0_val[2 * layers], eps );
00609             // Layer 1
00610             CHECK_REAL_EQUAL( 3.629944, vorticity0_val[0 * layers + 1], eps );
00611             CHECK_REAL_EQUAL( -1.708164, vorticity0_val[1 * layers + 1], eps );
00612             CHECK_REAL_EQUAL( 0.131686, vorticity0_val[2 * layers + 1], eps );
00613         }
00614         else if( 1 == rank )
00615         {
00616             CHECK_EQUAL( (size_t)688, local_verts.size() );
00617             CHECK_EQUAL( (size_t)1008, local_edges.size() );
00618 
00619             // Layer 0
00620             CHECK_REAL_EQUAL( -0.554888, vorticity0_val[0 * layers], eps );
00621             CHECK_REAL_EQUAL( 2.434397, vorticity0_val[1 * layers], eps );
00622             CHECK_REAL_EQUAL( -0.554888, vorticity0_val[2 * layers], eps );
00623             // Layer 1
00624             CHECK_REAL_EQUAL( -0.554881, vorticity0_val[0 * layers + 1], eps );
00625             CHECK_REAL_EQUAL( 2.434363, vorticity0_val[1 * layers + 1], eps );
00626             CHECK_REAL_EQUAL( -0.554881, vorticity0_val[2 * layers + 1], eps );
00627         }
00628     }
00629 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines