MOAB: Mesh Oriented datABase
(version 5.4.1)
|
00001 #include "iMeshP.h" 00002 #include "moab_mpi.h" 00003 #include <iostream> 00004 #include <algorithm> 00005 #include <vector> 00006 #include <sstream> 00007 #include <cassert> 00008 #include <cmath> 00009 #include <map> 00010 #include <cstring> 00011 #include <cstdio> // remove() 00012 00013 #if !defined( _MSC_VER ) && !defined( __MINGW32__ ) 00014 #include <unistd.h> 00015 #endif 00016 00017 #define STRINGIFY_( X ) #X 00018 #define STRINGIFY( X ) STRINGIFY_( X ) 00019 const char* const FILENAME = "iMeshP_test_file"; 00020 00021 /************************************************************************** 00022 Error Checking 00023 **************************************************************************/ 00024 00025 #define CHKERR \ 00026 do \ 00027 { \ 00028 if( ierr ) \ 00029 { \ 00030 std::cerr << "Error code " << ierr << " at " << __FILE__ << ":" << __LINE__ << std::endl; \ 00031 return ierr; \ 00032 } \ 00033 } while( false ) 00034 00035 #define PCHECK \ 00036 do \ 00037 { \ 00038 ierr = is_any_proc_error( ierr ); \ 00039 CHKERR; \ 00040 } while( false ) 00041 00042 // Use my_rank instead of rank to avoid shadowed declaration 00043 #define ASSERT( A ) \ 00044 do \ 00045 { \ 00046 if( is_any_proc_error( !( A ) ) ) \ 00047 { \ 00048 int my_rank = 0; \ 00049 MPI_Comm_rank( MPI_COMM_WORLD, &my_rank ); \ 00050 if( 0 == my_rank ) \ 00051 std::cerr << "Failed assertion: " #A << std::endl << " at " __FILE__ ":" << __LINE__ << std::endl; \ 00052 return 1; \ 00053 } \ 00054 } while( false ) 00055 00056 // Test if is_my_error is non-zero on any processor in MPI_COMM_WORLD 00057 int is_any_proc_error( int is_my_error ) 00058 { 00059 int result; 00060 int err = MPI_Allreduce( &is_my_error, &result, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD ); 00061 return err || result; 00062 } 00063 00064 /************************************************************************** 00065 Test Declarations 00066 **************************************************************************/ 00067 00068 class PartMap; 00069 00070 /**\brief Consistency check for parallel load 00071 * 00072 * All other tests depend on this one. 00073 */ 00074 int test_load( iMesh_Instance, iMeshP_PartitionHandle prtn, PartMap& map, int comm_size ); 00075 00076 /**\brief Test partition query methods 00077 * 00078 * Test: 00079 * - iMeshP_getPartitionComm 00080 * - iMeshP_getNumPartitions 00081 * - iMeshP_getPartitions 00082 */ 00083 int test_get_partitions( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00084 00085 /**\brief Test part quyery methods 00086 * 00087 * Test: 00088 * - iMeshP_getNumGlobalParts 00089 * - iMeshP_getNumLocalParts 00090 * - iMeshP_getLocalParts 00091 */ 00092 int test_get_parts( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00093 00094 /**\brief Test query by entity type 00095 * 00096 * Test: 00097 * - iMeshP_getNumOfTypeAll 00098 * - iMeshP_getNumOfType 00099 * - iMeshP_getEntities 00100 * - 00101 */ 00102 int test_get_by_type( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00103 00104 /**\brief Test query by entity topology 00105 * 00106 * Test: 00107 * - iMeshP_getNumOfTopoAll 00108 * - iMeshP_getNumOfTopo 00109 * - iMeshP_getEntities 00110 * - 00111 */ 00112 int test_get_by_topo( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00113 00114 /**\brief Test mapping from part id to part handle 00115 * 00116 * Test: 00117 * - iMeshP_getPartIdFromPartHandle 00118 * - iMeshP_getPartIdsFromPartHandlesArr 00119 * - iMeshP_getPartHandleFromPartId 00120 * - iMeshP_getPartHandlesFromPartsIdsArr 00121 * - iMeshP_getRankOfPart 00122 * - iMeshP_getRankOfPartArr 00123 */ 00124 int test_part_id_handle( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00125 00126 /**\brief Test get part rank 00127 * 00128 * Tests: 00129 * - iMeshP_getRankOfPart 00130 * - iMeshP_getRankOfPartArr 00131 */ 00132 int test_part_rank( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00133 00134 /**\brief Test querying of part neighbors 00135 * 00136 * Test: 00137 * - iMeshP_getNumPartNbors 00138 * - iMeshP_getNumPartNborsArr 00139 * - iMeshP_getPartNbors 00140 * - iMeshP_getPartNborsArr 00141 */ 00142 int test_get_neighbors( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00143 00144 /**\brief Test querying of part boundary entities 00145 * 00146 * Test: 00147 * - iMeshP_getNumPartBdryEnts 00148 * - iMeshP_getPartBdryEnts 00149 */ 00150 int test_get_part_boundary( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00151 00152 /**\brief Test querying of part boundary entities 00153 * 00154 * Test: 00155 * - iMeshP_initPartBdryEntIter 00156 * - iMeshP_initPartBdryEntArrIter 00157 */ 00158 int test_part_boundary_iter( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00159 00160 /**\brief Test adjacent entity query 00161 * 00162 * Test: 00163 * - iMeshP_getAdjEntities 00164 */ 00165 int test_get_adjacencies( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00166 00167 /**\brief Test entity iterators 00168 * 00169 * Test: 00170 * - iMeshP_initEntIter 00171 * - iMeshP_initEntArrIter 00172 */ 00173 int test_entity_iterator( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00174 00175 /**\brief Test entity owner queries 00176 * 00177 * Test: 00178 * - iMeshP_getEntOwnerPart 00179 * - iMeshP_getEntOwnerPartArr 00180 * - iMeshP_isEntOwner 00181 * - iMeshP_isEntOwnerArr 00182 */ 00183 int test_entity_owner( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00184 00185 /**\brief Test entity status 00186 * 00187 * Test: 00188 * - iMeshP_getEntStatus 00189 * - iMeshP_getEntStatusArr 00190 */ 00191 int test_entity_status( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00192 00193 /**\brief Test information about entity copies for interface entities 00194 * 00195 * Test: 00196 * - iMeshP_getNumCopies 00197 * - iMeshP_getCopyParts 00198 */ 00199 int test_entity_copy_parts( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00200 00201 /**\brief Test information about entity copies for interface entities 00202 * 00203 * Test: 00204 * - iMeshP_getCopies 00205 * - iMeshP_getCopyOnPart 00206 * - iMeshP_getOwnerCopy 00207 */ 00208 int test_entity_copies( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00209 00210 /**\brief Test creation of ghost entities 00211 * 00212 * Test: 00213 * - iMeshP_createGhostEntsAll 00214 */ 00215 int test_create_ghost_ents( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00216 00217 /**\brief Test commuinication of tag data 00218 * 00219 * Test: 00220 * - iMeshP_pushTags 00221 * - iMeshP_pushTagsEnt 00222 */ 00223 int test_push_tag_data_iface( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00224 int test_push_tag_data_ghost( iMesh_Instance, iMeshP_PartitionHandle prtn, const PartMap& ); 00225 00226 int test_exchange_ents( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ); 00227 00228 /************************************************************************** 00229 Helper Funcions 00230 **************************************************************************/ 00231 00232 class PartMap 00233 { 00234 public: 00235 int num_parts() const 00236 { 00237 return sortedPartList.size(); 00238 } 00239 00240 iMeshP_Part part_id_from_local_id( int local_id ) const 00241 { 00242 return sortedPartList[idx_from_local_id( local_id )]; 00243 } 00244 00245 int local_id_from_part_id( iMeshP_Part part ) const 00246 { 00247 return partLocalIds[idx_from_part_id( part )]; 00248 } 00249 00250 int rank_from_part_id( iMeshP_Part part ) const 00251 { 00252 return partRanks[idx_from_part_id( part )]; 00253 } 00254 00255 int rank_from_local_id( int id ) const 00256 { 00257 return partRanks[idx_from_local_id( id )]; 00258 } 00259 00260 int count_from_rank( int rank ) const 00261 { 00262 return std::count( partRanks.begin(), partRanks.end(), rank ); 00263 } 00264 00265 void part_id_from_rank( int rank, std::vector< iMeshP_Part >& parts ) const; 00266 00267 void local_id_from_rank( int rank, std::vector< int >& ids ) const; 00268 00269 const std::vector< iMeshP_Part >& get_parts() const 00270 { 00271 return sortedPartList; 00272 } 00273 00274 const std::vector< int >& get_ranks() const 00275 { 00276 return partRanks; 00277 } 00278 00279 int build_map( iMesh_Instance imesh, iMeshP_PartitionHandle partition, int num_expected_parts ); 00280 00281 static int part_from_coords( iMesh_Instance imesh, iMeshP_PartHandle part, int& id_out ); 00282 00283 private: 00284 inline int idx_from_part_id( iMeshP_Part id ) const 00285 { 00286 return std::lower_bound( sortedPartList.begin(), sortedPartList.end(), id ) - sortedPartList.begin(); 00287 } 00288 inline int idx_from_local_id( int id ) const 00289 { 00290 return localIdReverseMap[id]; 00291 } 00292 00293 std::vector< iMeshP_Part > sortedPartList; 00294 std::vector< int > partRanks; 00295 std::vector< int > partLocalIds; 00296 std::vector< int > localIdReverseMap; 00297 }; 00298 00299 /**\brief Create mesh for use in parallel tests */ 00300 int create_mesh( const char* filename, int num_parts ); 00301 00302 int create_mesh_in_memory( int rank, int size, iMesh_Instance imesh, iMeshP_PartitionHandle& prtn, PartMap& map ); 00303 00304 /**\brief get unique identifier for each vertex */ 00305 int vertex_tag( iMesh_Instance imesh, iBase_EntityHandle vertex, int& tag ); 00306 00307 int get_local_parts( iMesh_Instance instance, 00308 iMeshP_PartitionHandle prtn, 00309 std::vector< iMeshP_PartHandle >& handles, 00310 std::vector< iMeshP_Part >* ids = 0 ) 00311 { 00312 iMeshP_PartHandle* arr = 0; 00313 int ierr, alloc = 0, size; 00314 iMeshP_getLocalParts( instance, prtn, &arr, &alloc, &size, &ierr );CHKERR; 00315 handles.resize( size ); 00316 std::copy( arr, arr + size, handles.begin() ); 00317 free( arr ); 00318 if( !ids ) return iBase_SUCCESS; 00319 00320 ids->resize( size ); 00321 alloc = size; 00322 iMeshP_Part* ptr = &( *ids )[0]; 00323 iMeshP_getPartIdsFromPartHandlesArr( instance, prtn, &handles[0], handles.size(), &ptr, &alloc, &size, &ierr );CHKERR; 00324 assert( size == (int)ids->size() ); 00325 assert( ptr == &( *ids )[0] ); 00326 return iBase_SUCCESS; 00327 } 00328 00329 static int get_entities( iMesh_Instance imesh, 00330 iBase_EntitySetHandle set, 00331 iBase_EntityType type, 00332 iMesh_EntityTopology topo, 00333 std::vector< iBase_EntityHandle >& entities ) 00334 { 00335 iBase_EntityHandle* array = 0; 00336 int junk = 0, size = 0, err; 00337 iMesh_getEntities( imesh, set, type, topo, &array, &junk, &size, &err ); 00338 if( !err ) 00339 { 00340 entities.clear(); 00341 entities.resize( size ); 00342 std::copy( array, array + size, entities.begin() ); 00343 free( array ); 00344 } 00345 return err; 00346 } 00347 00348 static int get_part_quads_and_verts( iMesh_Instance imesh, 00349 iMeshP_PartHandle part, 00350 std::vector< iBase_EntityHandle >& elems, 00351 std::vector< iBase_EntityHandle >& verts ) 00352 { 00353 int ierr = get_entities( imesh, part, iBase_FACE, iMesh_QUADRILATERAL, elems );CHKERR; 00354 00355 verts.resize( 4 * elems.size() ); 00356 std::vector< int > junk( elems.size() + 1 ); 00357 int junk1 = verts.size(), count, junk2 = junk.size(), junk3; 00358 iBase_EntityHandle* junk4 = &verts[0]; 00359 int* junk5 = &junk[0]; 00360 iMesh_getEntArrAdj( imesh, &elems[0], elems.size(), iBase_VERTEX, &junk4, &junk1, &count, &junk5, &junk2, &junk3, 00361 &ierr );CHKERR; 00362 assert( junk1 == (int)verts.size() ); 00363 assert( count == (int)( 4 * elems.size() ) ); 00364 assert( junk2 == (int)junk.size() ); 00365 assert( junk4 == &verts[0] ); 00366 assert( junk5 == &junk[0] ); 00367 std::sort( verts.begin(), verts.end() ); 00368 verts.erase( std::unique( verts.begin(), verts.end() ), verts.end() ); 00369 return iBase_SUCCESS; 00370 } 00371 00372 static int get_coords( iMesh_Instance imesh, const iBase_EntityHandle* verts, int num_verts, double* coords ) 00373 { 00374 double* junk1 = coords; 00375 int junk2 = 3 * num_verts; 00376 int junk3; 00377 int ierr; 00378 iMesh_getVtxArrCoords( imesh, verts, num_verts, iBase_INTERLEAVED, &junk1, &junk2, &junk3, &ierr ); 00379 if( iBase_SUCCESS != ierr ) return ierr; 00380 assert( junk1 == coords ); 00381 assert( junk2 == 3 * num_verts ); 00382 assert( junk3 == 3 * num_verts ); 00383 return iBase_SUCCESS; 00384 } 00385 00386 /************************************************************************** 00387 Main Method 00388 **************************************************************************/ 00389 00390 #define RUN_TEST( A ) run_test( &( A ), #A ) 00391 00392 int run_test( int ( *func )( iMesh_Instance, iMeshP_PartitionHandle, const PartMap& ), const char* func_name ) 00393 { 00394 int rank, size, ierr; 00395 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 00396 MPI_Comm_size( MPI_COMM_WORLD, &size ); 00397 iMesh_Instance imesh; 00398 iMesh_newMesh( 0, &imesh, &ierr, 0 ); 00399 PCHECK; 00400 00401 iMeshP_PartitionHandle prtn; 00402 iMeshP_createPartitionAll( imesh, MPI_COMM_WORLD, &prtn, &ierr ); 00403 PCHECK; 00404 00405 PartMap map; 00406 00407 #ifdef MOAB_HAVE_HDF5 00408 if( rank == 0 ) 00409 { 00410 ierr = create_mesh( FILENAME, size ); 00411 } 00412 MPI_Bcast( &ierr, 1, MPI_INT, 0, MPI_COMM_WORLD ); 00413 if( ierr ) 00414 { 00415 if( rank == 0 ) 00416 { 00417 std::cerr << "Failed to create input test file on root processor. Aborting." << std::endl; 00418 } 00419 abort(); 00420 } 00421 00422 ierr = test_load( imesh, prtn, map, size ); 00423 if( ierr ) 00424 { 00425 if( rank == 0 ) 00426 { 00427 std::cerr << "Failed to load input mesh." << std::endl 00428 << "Cannot run further tests." << std::endl 00429 << "ABORTING" << std::endl; 00430 } 00431 abort(); 00432 } 00433 #else 00434 // so we have MPI and no HDF5; in order to run the test we need to create the 00435 // model in memory, and then call sync to resolve shared ents, as if it was read 00436 ierr = create_mesh_in_memory( rank, size, imesh, prtn, map ); 00437 MPI_Bcast( &ierr, 1, MPI_INT, 0, MPI_COMM_WORLD ); 00438 if( ierr ) 00439 { 00440 if( rank == 0 ) 00441 { 00442 std::cerr << "Failed to create mesh. Aborting." << std::endl; 00443 } 00444 abort(); 00445 } 00446 00447 #endif 00448 int result = ( *func )( imesh, prtn, map ); 00449 int is_err = is_any_proc_error( result ); 00450 if( rank == 0 ) 00451 { 00452 if( is_err ) 00453 std::cout << func_name << " : FAILED!!" << std::endl; 00454 else 00455 std::cout << func_name << " : success" << std::endl; 00456 } 00457 00458 iMesh_dtor( imesh, &ierr );CHKERR; 00459 return is_err; 00460 } 00461 00462 int main( int argc, char* argv[] ) 00463 { 00464 MPI_Init( &argc, &argv ); 00465 int size, rank; 00466 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 00467 MPI_Comm_size( MPI_COMM_WORLD, &size ); 00468 00469 if( argc > 2 && !strcmp( argv[1], "-p" ) ) 00470 { 00471 #if !defined( _MSC_VER ) && !defined( __MINGW32__ ) 00472 std::cout << "Processor " << rank << " of " << size << " with PID " << getpid() << std::endl; 00473 std::cout.flush(); 00474 #endif 00475 // loop forever on requested processor, giving the user time 00476 // to attach a debugger. Once the debugger in attached, user 00477 // can change 'pause'. E.g. on gdb do "set var pause = 0" 00478 if( atoi( argv[2] ) == rank ) 00479 { 00480 volatile int pause = 1; 00481 while( pause ) 00482 ; 00483 } 00484 MPI_Barrier( MPI_COMM_WORLD ); 00485 } 00486 00487 int num_errors = 0; 00488 num_errors += RUN_TEST( test_get_partitions ); 00489 num_errors += RUN_TEST( test_get_parts ); 00490 num_errors += RUN_TEST( test_get_by_type ); 00491 num_errors += RUN_TEST( test_get_by_topo ); 00492 num_errors += RUN_TEST( test_part_id_handle ); 00493 num_errors += RUN_TEST( test_part_rank ); 00494 num_errors += RUN_TEST( test_get_neighbors ); 00495 num_errors += RUN_TEST( test_get_part_boundary ); 00496 num_errors += RUN_TEST( test_part_boundary_iter ); 00497 num_errors += RUN_TEST( test_get_adjacencies ); 00498 num_errors += RUN_TEST( test_entity_iterator ); 00499 num_errors += RUN_TEST( test_entity_owner ); 00500 num_errors += RUN_TEST( test_entity_status ); 00501 num_errors += RUN_TEST( test_entity_copy_parts ); 00502 num_errors += RUN_TEST( test_entity_copies ); 00503 num_errors += RUN_TEST( test_push_tag_data_iface ); 00504 num_errors += RUN_TEST( test_push_tag_data_ghost ); 00505 num_errors += RUN_TEST( test_create_ghost_ents ); 00506 num_errors += RUN_TEST( test_exchange_ents ); 00507 00508 // wait until all procs are done before writing summary data 00509 std::cout.flush(); 00510 MPI_Barrier( MPI_COMM_WORLD ); 00511 00512 #ifdef MOAB_HAVE_HDF5 00513 // clean up output file 00514 if( rank == 0 ) remove( FILENAME ); 00515 #endif 00516 00517 if( rank == 0 ) 00518 { 00519 if( !num_errors ) 00520 std::cout << "All tests passed" << std::endl; 00521 else 00522 std::cout << num_errors << " TESTS FAILED!" << std::endl; 00523 } 00524 00525 MPI_Finalize(); 00526 return num_errors; 00527 } 00528 00529 // Create a mesh 00530 // 00531 // 00532 // Groups of four quads will be arranged into parts as follows: 00533 // +------+------+------+------+------+----- 00534 // | | | 00535 // | | | 00536 // + Part 0 + Part 2 + Part 4 00537 // | | | 00538 // | | | 00539 // +------+------+------+------+------+----- 00540 // | | | 00541 // | | | 00542 // + Part 1 + Part 3 + Part 5 00543 // | | | 00544 // | | | 00545 // +------+------+------+------+------+----- 00546 // 00547 // Vertices will be enumerated as follows: 00548 // 1------6-----11-----16-----21-----26----- 00549 // | | | 00550 // | | | 00551 // 2 7 12 17 22 27 00552 // | | | 00553 // | | | 00554 // 3------8-----13-----18-----23-----28----- 00555 // | | | 00556 // | | | 00557 // 4 9 14 19 24 29 00558 // | | | 00559 // | | | 00560 // 5-----10-----15-----20-----25-----30----- 00561 // 00562 // Element IDs will be [4*rank+1,4*rank+5] 00563 template < int size > 00564 struct EHARR 00565 { 00566 iBase_EntityHandle h[size]; 00567 iBase_EntityHandle& operator[]( int i ) 00568 { 00569 return h[i]; 00570 } 00571 operator iBase_EntityHandle*() 00572 { 00573 return h; 00574 } 00575 }; 00576 int create_mesh( const char* filename, int num_parts ) 00577 { 00578 const char* tagname = "GLOBAL_ID"; 00579 int ierr; 00580 00581 iMesh_Instance imesh; 00582 iMesh_newMesh( 0, &imesh, &ierr, 0 );CHKERR; 00583 00584 const int num_full_cols = 2 * ( num_parts / 2 ); 00585 const int need_half_cols = num_parts % 2; 00586 const int num_cols = num_full_cols + 2 * need_half_cols; 00587 const int num_vtx = 5 + 5 * num_cols - 4 * ( num_parts % 2 ); 00588 std::vector< EHARR< 5 > > vertices( num_cols + 1 ); 00589 std::vector< EHARR< 4 > > elements( num_cols ); 00590 std::vector< int > vertex_ids( num_vtx ); 00591 std::vector< iBase_EntityHandle > vertex_list( num_vtx ); 00592 for( int i = 0; i < num_vtx; ++i ) 00593 vertex_ids[i] = i + 1; 00594 00595 // create vertices 00596 int vl_pos = 0; 00597 for( int i = 0; i <= num_cols; ++i ) 00598 { 00599 double coords[15] = { static_cast< double >( i ), 0, 0, static_cast< double >( i ), 1, 0, 00600 static_cast< double >( i ), 2, 0, static_cast< double >( i ), 3, 0, 00601 static_cast< double >( i ), 4, 0 }; 00602 iBase_EntityHandle* ptr = vertices[i]; 00603 const int n = ( num_full_cols == num_cols || i <= num_full_cols ) ? 5 : 3; 00604 int junk1 = n, junk2 = n; 00605 iMesh_createVtxArr( imesh, n, iBase_INTERLEAVED, coords, 3 * n, &ptr, &junk1, &junk2, &ierr );CHKERR; 00606 assert( ptr == vertices[i] ); 00607 assert( junk1 == n ); 00608 assert( junk2 == n ); 00609 for( int j = 0; j < n; ++j ) 00610 vertex_list[vl_pos++] = vertices[i][j]; 00611 } 00612 00613 // create elements 00614 for( int i = 0; i < num_cols; ++i ) 00615 { 00616 iBase_EntityHandle conn[16]; 00617 for( int j = 0; j < 4; ++j ) 00618 { 00619 conn[4 * j] = vertices[i][j]; 00620 conn[4 * j + 1] = vertices[i][j + 1]; 00621 conn[4 * j + 2] = vertices[i + 1][j + 1]; 00622 conn[4 * j + 3] = vertices[i + 1][j]; 00623 } 00624 iBase_EntityHandle* ptr = elements[i]; 00625 const int n = ( i < num_full_cols ) ? 4 : 2; 00626 int junk1 = n, junk2 = n, junk3 = n, junk4 = n; 00627 int stat[4]; 00628 int* ptr2 = stat; 00629 iMesh_createEntArr( imesh, iMesh_QUADRILATERAL, conn, 4 * n, &ptr, &junk1, &junk2, &ptr2, &junk3, &junk4, 00630 &ierr );CHKERR; 00631 assert( ptr == elements[i] ); 00632 assert( junk1 == n ); 00633 assert( junk2 == n ); 00634 assert( ptr2 == stat ); 00635 assert( junk3 == n ); 00636 assert( junk4 == n ); 00637 } 00638 00639 // create partition 00640 iMeshP_PartitionHandle partition; 00641 iMeshP_createPartitionAll( imesh, MPI_COMM_SELF, &partition, &ierr );CHKERR; 00642 for( int i = 0; i < num_parts; ++i ) 00643 { 00644 iMeshP_PartHandle part; 00645 iMeshP_createPart( imesh, partition, &part, &ierr );CHKERR; 00646 iBase_EntityHandle quads[] = { elements[2 * ( i / 2 )][2 * ( i % 2 )], 00647 elements[2 * ( i / 2 ) + 1][2 * ( i % 2 )], 00648 elements[2 * ( i / 2 )][2 * ( i % 2 ) + 1], 00649 elements[2 * ( i / 2 ) + 1][2 * ( i % 2 ) + 1] }; 00650 iMesh_addEntArrToSet( imesh, quads, 4, part, &ierr );CHKERR; 00651 } 00652 00653 // assign global ids to vertices 00654 iBase_TagHandle id_tag = 0; 00655 iMesh_getTagHandle( imesh, tagname, &id_tag, &ierr, strlen( tagname ) ); 00656 if( iBase_SUCCESS == ierr ) 00657 { 00658 int tag_size, tag_type; 00659 iMesh_getTagSizeValues( imesh, id_tag, &tag_size, &ierr );CHKERR; 00660 if( tag_size != 1 ) return iBase_TAG_ALREADY_EXISTS; 00661 iMesh_getTagType( imesh, id_tag, &tag_type, &ierr );CHKERR; 00662 if( tag_type != iBase_INTEGER ) return iBase_TAG_ALREADY_EXISTS; 00663 } 00664 else 00665 { 00666 iMesh_createTag( imesh, tagname, 1, iBase_INTEGER, &id_tag, &ierr, strlen( tagname ) );CHKERR; 00667 } 00668 iMesh_setIntArrData( imesh, &vertex_list[0], num_vtx, id_tag, &vertex_ids[0], num_vtx, &ierr );CHKERR; 00669 00670 // write file 00671 iBase_EntitySetHandle root_set; 00672 iMesh_getRootSet( imesh, &root_set, &ierr ); 00673 iMeshP_saveAll( imesh, partition, root_set, filename, 0, &ierr, strlen( filename ), 0 );CHKERR; 00674 00675 iMesh_dtor( imesh, &ierr );CHKERR; 00676 00677 return 0; 00678 } 00679 int create_mesh_in_memory( int rank, 00680 int num_parts, 00681 iMesh_Instance imesh, 00682 iMeshP_PartitionHandle& partition, 00683 PartMap& map ) 00684 { 00685 const char* tagname = "GLOBAL_ID"; 00686 int ierr; 00687 00688 const int num_cols = 2; 00689 const int num_vtx = 9; 00690 // we are on the top or botton row 00691 int bottom = rank % 2; // 0 2 00692 // 1 3 00693 std::vector< EHARR< 3 > > vertices( 3 ); 00694 std::vector< EHARR< 2 > > elements( 2 ); // 4 elements per process 00695 std::vector< int > vertex_ids( num_vtx ); 00696 std::vector< iBase_EntityHandle > vertex_list( num_vtx ); 00697 int start = 1 + 2 * bottom + 10 * ( rank / 2 ); 00698 00699 for( int i = 0; i < 3; ++i ) 00700 for( int j = 0; j < 3; ++j ) 00701 { 00702 vertex_ids[i + 3 * j] = start + i + 5 * j; 00703 } 00704 00705 // create vertices 00706 int vl_pos = 0; 00707 int startI = 2 * ( rank / 2 ); // so it will be 0, 0, 2, 2, ...) 00708 for( int i = 0; i <= 2; ++i ) 00709 { 00710 double coords[9] = { 00711 static_cast< double >( i + startI ), 2. * bottom, 0, 00712 static_cast< double >( i + startI ), 1 + 2. * bottom, 0, 00713 static_cast< double >( i + startI ), 2 + 2. * bottom, 0, 00714 }; 00715 iBase_EntityHandle* ptr = vertices[i]; 00716 const int n = 3; 00717 int junk1 = n, junk2 = n; 00718 iMesh_createVtxArr( imesh, n, iBase_INTERLEAVED, coords, 3 * n, &ptr, &junk1, &junk2, &ierr );CHKERR; 00719 assert( ptr == vertices[i] ); 00720 assert( junk1 == n ); 00721 assert( junk2 == n ); 00722 for( int j = 0; j < n; ++j ) 00723 vertex_list[vl_pos++] = vertices[i][j]; 00724 } 00725 00726 // create elements 00727 for( int i = 0; i < num_cols; ++i ) 00728 { 00729 iBase_EntityHandle conn[8]; 00730 for( int j = 0; j < 2; ++j ) 00731 { 00732 conn[4 * j] = vertices[i][j]; 00733 conn[4 * j + 1] = vertices[i][j + 1]; 00734 conn[4 * j + 2] = vertices[i + 1][j + 1]; 00735 conn[4 * j + 3] = vertices[i + 1][j]; 00736 } 00737 iBase_EntityHandle* ptr = elements[i]; 00738 const int n = 2; 00739 int junk1 = n, junk2 = n, junk3 = n, junk4 = n; 00740 int stat[4]; 00741 int* ptr2 = stat; 00742 iMesh_createEntArr( imesh, iMesh_QUADRILATERAL, conn, 4 * n, &ptr, &junk1, &junk2, &ptr2, &junk3, &junk4, 00743 &ierr );CHKERR; 00744 assert( ptr == elements[i] ); 00745 assert( junk1 == n ); 00746 assert( junk2 == n ); 00747 assert( ptr2 == stat ); 00748 assert( junk3 == n ); 00749 assert( junk4 == n ); 00750 } 00751 00752 // create partition 00753 iMeshP_createPartitionAll( imesh, MPI_COMM_WORLD, &partition, &ierr );CHKERR; 00754 00755 iMeshP_PartHandle part; 00756 iMeshP_createPart( imesh, partition, &part, &ierr );CHKERR; 00757 iBase_EntityHandle quads[] = { elements[0][0], elements[0][1], elements[1][0], elements[1][1] }; 00758 iMesh_addEntArrToSet( imesh, quads, 4, part, &ierr );CHKERR; 00759 00760 // assign global ids to vertices 00761 iBase_TagHandle id_tag = 0; 00762 iMesh_getTagHandle( imesh, tagname, &id_tag, &ierr, strlen( tagname ) ); 00763 if( iBase_SUCCESS == ierr ) 00764 { 00765 int tag_size, tag_type; 00766 iMesh_getTagSizeValues( imesh, id_tag, &tag_size, &ierr );CHKERR; 00767 if( tag_size != 1 ) return iBase_TAG_ALREADY_EXISTS; 00768 iMesh_getTagType( imesh, id_tag, &tag_type, &ierr );CHKERR; 00769 if( tag_type != iBase_INTEGER ) return iBase_TAG_ALREADY_EXISTS; 00770 } 00771 else 00772 { 00773 iMesh_createTag( imesh, tagname, 1, iBase_INTEGER, &id_tag, &ierr, strlen( tagname ) );CHKERR; 00774 } 00775 iMesh_setIntArrData( imesh, &vertex_list[0], num_vtx, id_tag, &vertex_ids[0], num_vtx, &ierr );CHKERR; 00776 00777 // some mesh sync 00778 iMeshP_syncPartitionAll( imesh, partition, &ierr );CHKERR; 00779 iMeshP_syncMeshAll( imesh, partition, &ierr );CHKERR; 00780 00781 ierr = map.build_map( imesh, partition, num_parts );CHKERR; 00782 return 0; 00783 } 00784 // generate unique for each vertex from coordinates. 00785 // Assume integer coordinate values with x in [0,inf] and y in [0,4] 00786 // as generated by create_mean(..). 00787 int vertex_tag( iMesh_Instance imesh, iBase_EntityHandle vertex, int& tag ) 00788 { 00789 int ierr; 00790 double x, y, z; 00791 iMesh_getVtxCoord( imesh, vertex, &x, &y, &z, &ierr );CHKERR; 00792 00793 int xc = (int)round( x ); 00794 int yc = (int)round( y ); 00795 tag = 5 * xc + yc + 1; 00796 return ierr; 00797 } 00798 00799 /************************************************************************** 00800 Test Implementations 00801 **************************************************************************/ 00802 00803 int test_load( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, PartMap& map, int proc_size ) 00804 { 00805 int ierr; 00806 00807 iBase_EntitySetHandle root_set; 00808 iMesh_getRootSet( imesh, &root_set, &ierr ); 00809 const char* opt = "moab:PARTITION=PARALLEL_PARTITION"; 00810 iMeshP_loadAll( imesh, prtn, root_set, FILENAME, opt, &ierr, strlen( FILENAME ), strlen( opt ) ); 00811 PCHECK; 00812 00813 ierr = map.build_map( imesh, prtn, proc_size );CHKERR; 00814 return iBase_SUCCESS; 00815 } 00816 00817 /**\brief Test partition query methods 00818 * 00819 * Test: 00820 * - iMeshP_getPartitionComm 00821 * - iMeshP_getNumPartitions 00822 * - iMeshP_getPartitions 00823 */ 00824 int test_get_partitions( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& ) 00825 { 00826 int ierr; 00827 00828 // test iMeshP_getPartitionCom 00829 MPI_Comm comm = MPI_COMM_SELF; 00830 iMeshP_getPartitionComm( imesh, prtn, &comm, &ierr ); 00831 PCHECK; 00832 ASSERT( comm == MPI_COMM_WORLD ); 00833 00834 // test iMeshP_getPartitions 00835 iMeshP_PartitionHandle* array = 0; 00836 int alloc = 0, size = -1; 00837 iMeshP_getPartitions( imesh, &array, &alloc, &size, &ierr ); 00838 PCHECK; 00839 ASSERT( array != 0 ); 00840 ASSERT( alloc == size ); 00841 ASSERT( size > 0 ); 00842 int idx = std::find( array, array + size, prtn ) - array; 00843 free( array ); 00844 ASSERT( idx < size ); 00845 00846 // test iMesP_getNumPartitions 00847 int size2 = -1; 00848 iMeshP_getNumPartitions( imesh, &size2, &ierr ); 00849 PCHECK; 00850 ASSERT( size2 == size ); 00851 return 0; 00852 } 00853 00854 /**\brief Test part quyery methods 00855 * 00856 * Test: 00857 * - iMeshP_getNumGlobalParts 00858 * - iMeshP_getNumLocalParts 00859 * - iMeshP_getLocalParts 00860 */ 00861 int test_get_parts( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 00862 { 00863 int size, rank, ierr; 00864 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 00865 MPI_Comm_size( MPI_COMM_WORLD, &size ); 00866 00867 int num_part_g; 00868 iMeshP_getNumGlobalParts( imesh, prtn, &num_part_g, &ierr ); 00869 PCHECK; 00870 ASSERT( num_part_g == map.num_parts() ); 00871 00872 int num_part_l; 00873 iMeshP_getNumLocalParts( imesh, prtn, &num_part_l, &ierr ); 00874 PCHECK; 00875 ASSERT( num_part_l == map.count_from_rank( rank ) ); 00876 00877 std::vector< iMeshP_PartHandle > parts( num_part_l ); 00878 iMeshP_PartHandle* ptr = &parts[0]; 00879 int junk1 = num_part_l, count = -1; 00880 iMeshP_getLocalParts( imesh, prtn, &ptr, &junk1, &count, &ierr ); 00881 PCHECK; 00882 assert( ptr == &parts[0] ); 00883 assert( junk1 == num_part_l ); 00884 ASSERT( count == num_part_l ); 00885 00886 return iBase_SUCCESS; 00887 } 00888 00889 static int test_get_by_type_topo_all( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, bool test_type, int num_parts ) 00890 { 00891 // calculate number of quads and vertices in entire mesh 00892 // from number of parts (see create_mesh(..) function.) 00893 const int expected_global_quad_count = 4 * num_parts; 00894 const int num_col = 2 * ( num_parts / 2 + num_parts % 2 ); 00895 const int expected_global_vtx_count = num_parts == 1 ? 9 : num_parts % 2 ? 1 + 5 * num_col : 5 + 5 * num_col; 00896 00897 // test getNumOf*All for root set 00898 int ierr, count; 00899 iBase_EntitySetHandle root; 00900 iMesh_getRootSet( imesh, &root, &ierr ); 00901 if( test_type ) 00902 iMeshP_getNumOfTypeAll( imesh, prtn, root, iBase_VERTEX, &count, &ierr ); 00903 else 00904 iMeshP_getNumOfTopoAll( imesh, prtn, root, iMesh_POINT, &count, &ierr ); 00905 PCHECK; 00906 ASSERT( count == expected_global_vtx_count ); 00907 if( test_type ) 00908 iMeshP_getNumOfTypeAll( imesh, prtn, root, iBase_FACE, &count, &ierr ); 00909 else 00910 iMeshP_getNumOfTopoAll( imesh, prtn, root, iMesh_QUADRILATERAL, &count, &ierr ); 00911 PCHECK; 00912 ASSERT( count == expected_global_quad_count ); 00913 00914 // create an entity set containing half of the quads 00915 std::vector< iBase_EntityHandle > all_quads, half_quads; 00916 ierr = get_entities( imesh, root, iBase_FACE, iMesh_QUADRILATERAL, all_quads ); 00917 assert( 0 == all_quads.size() % 2 ); 00918 half_quads.resize( all_quads.size() / 2 ); 00919 for( size_t i = 0; i < all_quads.size() / 2; ++i ) 00920 half_quads[i] = all_quads[2 * i]; 00921 iBase_EntitySetHandle set; 00922 iMesh_createEntSet( imesh, 1, &set, &ierr );CHKERR; 00923 iMesh_addEntArrToSet( imesh, &half_quads[0], half_quads.size(), set, &ierr );CHKERR; 00924 00925 // test getNumOf*All with defined set 00926 if( test_type ) 00927 iMeshP_getNumOfTypeAll( imesh, prtn, set, iBase_VERTEX, &count, &ierr ); 00928 else 00929 iMeshP_getNumOfTopoAll( imesh, prtn, set, iMesh_POINT, &count, &ierr ); 00930 PCHECK; 00931 ASSERT( count == 0 ); 00932 if( test_type ) 00933 iMeshP_getNumOfTypeAll( imesh, prtn, set, iBase_FACE, &count, &ierr ); 00934 else 00935 iMeshP_getNumOfTopoAll( imesh, prtn, set, iMesh_QUADRILATERAL, &count, &ierr ); 00936 PCHECK; 00937 ASSERT( count == expected_global_quad_count / 2 ); 00938 00939 return 0; 00940 } 00941 00942 static int test_get_by_type_topo_local( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, bool test_type ) 00943 { 00944 int ierr; 00945 iBase_EntitySetHandle root; 00946 iMesh_getRootSet( imesh, &root, &ierr ); 00947 00948 // select a single part 00949 std::vector< iMeshP_PartHandle > parts; 00950 ierr = get_local_parts( imesh, prtn, parts );CHKERR; 00951 iMeshP_PartHandle part = parts.front(); 00952 00953 // get the entities contained in the part 00954 std::vector< iBase_EntityHandle > part_quads, part_all; 00955 ierr = get_entities( imesh, part, iBase_FACE, iMesh_QUADRILATERAL, part_quads );CHKERR; 00956 ierr = get_entities( imesh, part, iBase_ALL_TYPES, iMesh_ALL_TOPOLOGIES, part_all );CHKERR; 00957 00958 // compare local counts (using root set) 00959 00960 int count; 00961 if( test_type ) 00962 iMeshP_getNumOfType( imesh, prtn, part, root, iBase_FACE, &count, &ierr ); 00963 else 00964 iMeshP_getNumOfTopo( imesh, prtn, part, root, iMesh_QUADRILATERAL, &count, &ierr );CHKERR; 00965 ASSERT( count == (int)part_quads.size() ); 00966 00967 if( test_type ) 00968 iMeshP_getNumOfType( imesh, prtn, part, root, iBase_ALL_TYPES, &count, &ierr ); 00969 else 00970 iMeshP_getNumOfTopo( imesh, prtn, part, root, iMesh_ALL_TOPOLOGIES, &count, &ierr );CHKERR; 00971 ASSERT( count == (int)part_all.size() ); 00972 00973 // compare local contents (using root set) 00974 00975 iBase_EntityHandle* ptr = 0; 00976 int num_ent, junk1 = 0; 00977 iMeshP_getEntities( imesh, prtn, part, root, test_type ? iBase_FACE : iBase_ALL_TYPES, 00978 test_type ? iMesh_ALL_TOPOLOGIES : iMesh_QUADRILATERAL, &ptr, &junk1, &num_ent, &ierr );CHKERR; 00979 std::vector< iBase_EntityHandle > act_quads( ptr, ptr + num_ent ); 00980 free( ptr ); 00981 junk1 = num_ent = 0; 00982 ptr = 0; 00983 iMeshP_getEntities( imesh, prtn, part, root, iBase_ALL_TYPES, iMesh_ALL_TOPOLOGIES, &ptr, &junk1, &num_ent, &ierr );CHKERR; 00984 std::vector< iBase_EntityHandle > act_all( ptr, ptr + num_ent ); 00985 free( ptr ); 00986 std::sort( part_quads.begin(), part_quads.end() ); 00987 std::sort( part_all.begin(), part_all.end() ); 00988 std::sort( act_quads.begin(), act_quads.end() ); 00989 std::sort( act_all.begin(), act_all.end() ); 00990 ASSERT( part_quads == act_quads ); 00991 ASSERT( part_all == act_all ); 00992 00993 // create an entity set containing half of the quads from the part 00994 std::vector< iBase_EntityHandle > half_quads( part_quads.size() / 2 ); 00995 for( size_t i = 0; i < half_quads.size(); ++i ) 00996 half_quads[i] = part_quads[2 * i]; 00997 iBase_EntitySetHandle set; 00998 iMesh_createEntSet( imesh, 1, &set, &ierr );CHKERR; 00999 iMesh_addEntArrToSet( imesh, &half_quads[0], half_quads.size(), set, &ierr );CHKERR; 01000 01001 // check if there exists any quads not in the part that we 01002 // can add to the set 01003 std::vector< iBase_EntityHandle > all_quads, other_quads; 01004 ierr = get_entities( imesh, root, iBase_FACE, iMesh_QUADRILATERAL, all_quads );CHKERR; 01005 std::sort( all_quads.begin(), all_quads.end() ); 01006 std::sort( part_quads.begin(), part_quads.end() ); 01007 std::set_difference( all_quads.begin(), all_quads.end(), part_quads.begin(), part_quads.end(), 01008 std::back_inserter( other_quads ) ); 01009 iMesh_addEntArrToSet( imesh, &other_quads[0], other_quads.size(), set, &ierr );CHKERR; 01010 01011 // compare local counts (using non-root set) 01012 01013 if( test_type ) 01014 iMeshP_getNumOfType( imesh, prtn, part, set, iBase_FACE, &count, &ierr ); 01015 else 01016 iMeshP_getNumOfTopo( imesh, prtn, part, set, iMesh_QUADRILATERAL, &count, &ierr );CHKERR; 01017 ASSERT( count == (int)half_quads.size() ); 01018 01019 if( test_type ) 01020 iMeshP_getNumOfType( imesh, prtn, part, set, iBase_VERTEX, &count, &ierr ); 01021 else 01022 iMeshP_getNumOfTopo( imesh, prtn, part, set, iMesh_POINT, &count, &ierr );CHKERR; 01023 ASSERT( count == 0 ); 01024 01025 // compare local contents (using non-root set) 01026 01027 junk1 = 0; 01028 num_ent = 0; 01029 ptr = 0; 01030 iMeshP_getEntities( imesh, prtn, part, set, test_type ? iBase_FACE : iBase_ALL_TYPES, 01031 test_type ? iMesh_ALL_TOPOLOGIES : iMesh_QUADRILATERAL, &ptr, &junk1, &num_ent, &ierr );CHKERR; 01032 act_quads.resize( num_ent ); 01033 std::copy( ptr, ptr + num_ent, act_quads.begin() ); 01034 free( ptr ); 01035 std::sort( half_quads.begin(), half_quads.end() ); 01036 std::sort( act_quads.begin(), act_quads.end() ); 01037 ASSERT( act_quads == half_quads ); 01038 01039 return iBase_SUCCESS; 01040 } 01041 01042 /**\brief Test query by entity type 01043 * 01044 * Test: 01045 * - iMeshP_getNumOfTypeAll 01046 * - iMeshP_getNumOfType 01047 * - iMeshP_getEntities 01048 * - 01049 */ 01050 int test_get_by_type( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 01051 { 01052 int ierr; 01053 ierr = test_get_by_type_topo_all( imesh, prtn, true, map.num_parts() ); 01054 PCHECK; 01055 ierr = test_get_by_type_topo_local( imesh, prtn, true ); 01056 PCHECK; 01057 return 0; 01058 } 01059 01060 /**\brief Test query by entity topology 01061 * 01062 * Test: 01063 * - iMeshP_getNumOfTopoAll 01064 * - iMeshP_getNumOfTopo 01065 * - iMeshP_getEntities 01066 * - 01067 */ 01068 int test_get_by_topo( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 01069 { 01070 int ierr; 01071 ierr = test_get_by_type_topo_all( imesh, prtn, false, map.num_parts() ); 01072 PCHECK; 01073 ierr = test_get_by_type_topo_local( imesh, prtn, false ); 01074 PCHECK; 01075 return 0; 01076 } 01077 01078 /**\brief Test mapping from part id to part handle 01079 * 01080 * Test: 01081 * - iMeshP_getPartIdFromPartHandle 01082 * - iMeshP_getPartIdsFromPartHandlesArr 01083 * - iMeshP_getPartHandleFromPartId 01084 * - iMeshP_getPartHandlesFromPartsIdsArr 01085 */ 01086 int test_part_id_handle( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 01087 { 01088 // get local part ids 01089 int rank, ierr; 01090 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 01091 std::vector< iMeshP_Part > ids; 01092 map.part_id_from_rank( rank, ids ); 01093 01094 // check single-part functions and build list of part handles 01095 std::vector< iMeshP_PartHandle > handles( ids.size() ); 01096 size_t i; 01097 for( i = 0; i < ids.size(); ++i ) 01098 { 01099 iMeshP_getPartHandleFromPartId( imesh, prtn, ids[i], &handles[i], &ierr );CHKERR; 01100 iMeshP_Part id; 01101 iMeshP_getPartIdFromPartHandle( imesh, prtn, handles[i], &id, &ierr );CHKERR; 01102 if( id != ids[i] ) break; 01103 } 01104 ASSERT( i == ids.size() ); 01105 01106 // test iMeshP_getPartIdsFromPartHandlesArr 01107 std::vector< iMeshP_Part > ids2( ids.size() ); 01108 int junk1 = ids.size(), junk2 = 0; 01109 iMeshP_Part* ptr = &ids2[0]; 01110 iMeshP_getPartIdsFromPartHandlesArr( imesh, prtn, &handles[0], handles.size(), &ptr, &junk1, &junk2, &ierr ); 01111 PCHECK; 01112 ASSERT( ptr == &ids2[0] ); 01113 ASSERT( junk2 == (int)ids2.size() ); 01114 ASSERT( ids == ids2 ); 01115 01116 // test iMeshP_getPartHandlesFromPartsIdsArr 01117 std::vector< iMeshP_PartHandle > handles2( handles.size() ); 01118 junk1 = handles.size(); 01119 junk2 = 0; 01120 iMeshP_PartHandle* ptr2 = &handles2[0]; 01121 iMeshP_getPartHandlesFromPartsIdsArr( imesh, prtn, &ids[0], ids.size(), &ptr2, &junk1, &junk2, &ierr ); 01122 PCHECK; 01123 ASSERT( ptr2 == &handles2[0] ); 01124 ASSERT( junk2 == (int)handles2.size() ); 01125 ASSERT( handles == handles2 ); 01126 01127 return 0; 01128 } 01129 01130 /**\brief Test get part rank 01131 * 01132 * Tests: 01133 * - iMeshP_getRankOfPart 01134 * - iMeshP_getRankOfPartArr 01135 */ 01136 int test_part_rank( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 01137 { 01138 int ierr = 0, rank; 01139 std::vector< iMeshP_Part > invalid, failed; 01140 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 01141 01142 // test iMeshP_getRankOfPart 01143 for( size_t i = 0; i < map.get_parts().size(); ++i ) 01144 { 01145 int pr; 01146 iMeshP_getRankOfPart( imesh, prtn, map.get_parts()[i], &pr, &ierr ); 01147 if( iBase_SUCCESS != ierr ) 01148 failed.push_back( map.get_parts()[i] ); 01149 else if( pr != map.get_ranks()[i] ) 01150 invalid.push_back( map.get_parts()[i] ); 01151 } 01152 if( !failed.empty() ) 01153 { 01154 std::cerr << "Processor " << rank << ": iMeshP_getRankOfPart failed for " << failed.size() << " parts." 01155 << std::endl; 01156 ierr = iBase_FAILURE; 01157 } 01158 if( !invalid.empty() ) 01159 { 01160 std::cerr << "Processor " << rank << ": iMeshP_getRankOfPart was incorrect for " << invalid.size() << " parts." 01161 << std::endl; 01162 ierr = iBase_FAILURE; 01163 } 01164 PCHECK; 01165 01166 // test iMeshP_getRankOfPartArr 01167 std::vector< int > ranks( map.get_parts().size() ); 01168 int junk1 = ranks.size(), junk2, *ptr = &ranks[0]; 01169 iMeshP_getRankOfPartArr( imesh, prtn, &map.get_parts()[0], map.get_parts().size(), &ptr, &junk1, &junk2, &ierr ); 01170 PCHECK; 01171 assert( ptr == &ranks[0] ); 01172 assert( junk1 == (int)ranks.size() ); 01173 ASSERT( junk2 == (int)ranks.size() ); 01174 for( size_t i = 0; i < map.get_parts().size(); ++i ) 01175 { 01176 if( ranks[i] != map.get_ranks()[i] ) invalid.push_back( map.get_parts()[i] ); 01177 } 01178 if( !invalid.empty() ) 01179 { 01180 std::cerr << "Processor " << rank << ": iMeshP_getRankOfPartArr was incorrect for " << invalid.size() 01181 << " parts." << std::endl; 01182 ierr = iBase_FAILURE; 01183 } 01184 PCHECK; 01185 01186 return 0; 01187 } 01188 01189 // see create_mesh(..) 01190 static void get_part_neighbors( int logical_part_id, int num_parts, int neighbors[5], int& num_neighbors ) 01191 { 01192 num_neighbors = 0; 01193 if( logical_part_id + 1 < num_parts ) neighbors[num_neighbors++] = logical_part_id + 1; 01194 if( logical_part_id + 2 < num_parts ) neighbors[num_neighbors++] = logical_part_id + 2; 01195 if( logical_part_id % 2 ) 01196 { 01197 neighbors[num_neighbors++] = logical_part_id - 1; 01198 if( logical_part_id > 2 ) 01199 { 01200 neighbors[num_neighbors++] = logical_part_id - 3; 01201 neighbors[num_neighbors++] = logical_part_id - 2; 01202 } 01203 } 01204 else 01205 { 01206 if( logical_part_id + 3 < num_parts ) neighbors[num_neighbors++] = logical_part_id + 3; 01207 if( logical_part_id > 1 ) 01208 { 01209 neighbors[num_neighbors++] = logical_part_id - 1; 01210 neighbors[num_neighbors++] = logical_part_id - 2; 01211 } 01212 } 01213 } 01214 01215 /**\brief Test querying of part neighbors 01216 * 01217 * Test: 01218 * - iMeshP_getNumPartNbors 01219 * - iMeshP_getNumPartNborsArr 01220 * - iMeshP_getPartNbors 01221 * - iMeshP_getPartNborsArr 01222 */ 01223 int test_get_neighbors( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 01224 { 01225 int ierr, rank; 01226 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 01227 01228 std::vector< iMeshP_Part > local_parts; 01229 map.part_id_from_rank( rank, local_parts ); 01230 01231 // get handles for local parts 01232 std::vector< iMeshP_PartHandle > handles( local_parts.size() ); 01233 iMeshP_PartHandle* ptr = &handles[0]; 01234 int junk1 = handles.size(), junk2 = 0; 01235 iMeshP_getPartHandlesFromPartsIdsArr( imesh, prtn, &local_parts[0], local_parts.size(), &ptr, &junk1, &junk2, 01236 &ierr ); 01237 PCHECK; 01238 assert( ptr == &handles[0] ); 01239 assert( junk2 == (int)handles.size() ); 01240 01241 // get logical ids for local parts 01242 std::vector< int > logical_ids; 01243 map.local_id_from_rank( rank, logical_ids ); 01244 01245 // get neighbors for each local part 01246 std::vector< std::vector< iMeshP_Part > > neighbors( logical_ids.size() ); 01247 for( size_t i = 0; i < logical_ids.size(); ++i ) 01248 { 01249 int logical_neighbors[5], num_neighbors; 01250 get_part_neighbors( logical_ids[i], map.num_parts(), logical_neighbors, num_neighbors ); 01251 neighbors[i].resize( num_neighbors ); 01252 for( int j = 0; j < num_neighbors; ++j ) 01253 neighbors[i][j] = map.part_id_from_local_id( logical_neighbors[j] ); 01254 std::sort( neighbors[i].begin(), neighbors[i].end() ); 01255 } 01256 01257 // test iMeshP_getNumPartNbors 01258 std::vector< iMeshP_Part > invalid, failed; 01259 for( size_t i = 0; i < local_parts.size(); ++i ) 01260 { 01261 int count; 01262 iMeshP_getNumPartNbors( imesh, prtn, handles[i], iBase_VERTEX, &count, &ierr ); 01263 if( ierr ) 01264 failed.push_back( local_parts[i] ); 01265 else if( count != (int)neighbors[i].size() ) 01266 invalid.push_back( local_parts[i] ); 01267 } 01268 if( !failed.empty() ) 01269 { 01270 std::cerr << "Processor " << rank << ": iMeshP_getNumPartNbors failed for " << failed.size() << " parts." 01271 << std::endl; 01272 ierr = iBase_FAILURE; 01273 PCHECK; 01274 } 01275 if( !invalid.empty() ) 01276 { 01277 std::cerr << "Processor " << rank << ": iMeshP_getNumPartNbors was incorrect for " << invalid.size() 01278 << " parts." << std::endl; 01279 ierr = iBase_FAILURE; 01280 PCHECK; 01281 } 01282 01283 // test iMeshP_getPartNbors 01284 ierr = 0; 01285 for( size_t i = 0; i < local_parts.size(); ++i ) 01286 { 01287 int count, junk = 0, another_count; 01288 iMeshP_Part* list = 0; 01289 iMeshP_getPartNbors( imesh, prtn, handles[i], iBase_VERTEX, &another_count, &list, &junk, &count, &ierr ); 01290 assert( count == another_count ); 01291 if( ierr ) 01292 failed.push_back( local_parts[i] ); 01293 else 01294 { 01295 std::sort( list, list + count ); 01296 std::vector< iMeshP_Part > cpy( list, list + count ); 01297 if( cpy != neighbors[i] ) invalid.push_back( local_parts[i] ); 01298 free( list ); 01299 } 01300 } 01301 if( !failed.empty() ) 01302 { 01303 std::cerr << "Processor " << rank << ": iMeshP_getPartNbors failed for " << failed.size() << " parts." 01304 << std::endl; 01305 ierr = iBase_FAILURE; 01306 } 01307 if( !invalid.empty() ) 01308 { 01309 std::cerr << "Processor " << rank << ": iMeshP_getPartNbors was incorrect for " << invalid.size() << " parts." 01310 << std::endl; 01311 ierr = iBase_FAILURE; 01312 } 01313 PCHECK; 01314 01315 // test iMeshP_getNumPartNborsArr 01316 std::vector< int > count_vect( handles.size() ); 01317 int* count_arr = &count_vect[0]; 01318 junk1 = handles.size(); 01319 iMeshP_getNumPartNborsArr( imesh, prtn, &handles[0], handles.size(), iBase_VERTEX, &count_arr, &junk1, &junk2, 01320 &ierr ); 01321 PCHECK; 01322 assert( count_arr == &count_vect[0] ); 01323 assert( junk2 == (int)handles.size() ); 01324 for( size_t i = 0; i < local_parts.size(); ++i ) 01325 { 01326 if( count_arr[i] != (int)neighbors[i].size() ) invalid.push_back( local_parts[i] ); 01327 } 01328 if( !invalid.empty() ) 01329 { 01330 std::cerr << "Processor " << rank << ": iMeshP_getNumPartNborsArr was incorrect for " << invalid.size() 01331 << " parts." << std::endl; 01332 ierr = iBase_FAILURE; 01333 } 01334 PCHECK; 01335 01336 // test iMeshP_getPartNborsArr 01337 iMeshP_Part* nbor_arr = 0; 01338 junk1 = handles.size(), junk2 = 0; 01339 int junk3 = 0, nbor_size; 01340 iMeshP_getPartNborsArr( imesh, prtn, &handles[0], handles.size(), iBase_VERTEX, &count_arr, &junk1, &junk2, 01341 &nbor_arr, &junk3, &nbor_size, &ierr ); 01342 PCHECK; 01343 assert( count_arr == &count_vect[0] ); 01344 assert( junk2 == (int)handles.size() ); 01345 std::vector< iMeshP_Part > all_nbors( nbor_arr, nbor_arr + nbor_size ); 01346 free( nbor_arr ); 01347 std::vector< iMeshP_Part >::iterator j = all_nbors.begin(); 01348 bool bad_length = false; 01349 for( size_t i = 0; i < local_parts.size(); ++i ) 01350 { 01351 if( all_nbors.end() - j > count_arr[i] ) 01352 { 01353 bad_length = true; 01354 break; 01355 } 01356 if( count_arr[i] != (int)neighbors[i].size() ) 01357 { 01358 invalid.push_back( local_parts[i] ); 01359 } 01360 else 01361 { 01362 std::vector< iMeshP_Part >::iterator e = j + count_arr[i]; 01363 std::sort( j, e ); 01364 if( !std::equal( j, e, neighbors[i].begin() ) ) invalid.push_back( local_parts[i] ); 01365 } 01366 } 01367 if( bad_length ) 01368 { 01369 std::cerr << "Processor " << rank << ": iMeshP_getPartNborsArr had inconsistent result array lengths." 01370 << std::endl; 01371 ierr = iBase_FAILURE; 01372 } 01373 if( !invalid.empty() ) 01374 { 01375 std::cerr << "Processor " << rank << ": iMeshP_getPartNborsArr was incorrect for " << invalid.size() 01376 << " parts." << std::endl; 01377 ierr = iBase_FAILURE; 01378 } 01379 PCHECK; 01380 01381 return 0; 01382 } 01383 01384 // Determine the expected vertices on the interface between two parts. 01385 // Returns no vertices for non-adjacient parts and fails if both parts 01386 // are the same. 01387 // See create_mesh(..) for the assumed mesh. 01388 static int interface_verts( iMesh_Instance imesh, 01389 iMeshP_PartitionHandle prtn, 01390 iMeshP_PartHandle local_part, 01391 iMeshP_Part other_part, 01392 const PartMap& map, 01393 std::vector< iBase_EntityHandle >& vtx_handles ) 01394 { 01395 int ierr, rank; 01396 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 01397 01398 iMeshP_Part local_id; 01399 iMeshP_getPartIdFromPartHandle( imesh, prtn, local_part, &local_id, &ierr );CHKERR; 01400 01401 const int local_logical = map.local_id_from_part_id( local_id ); 01402 const int other_logical = map.local_id_from_part_id( other_part ); 01403 01404 // get grid of local vertices 01405 01406 iBase_EntityHandle verts[3][3]; 01407 const double xbase = ( local_id / 2 ) * 2; 01408 const double ybase = ( local_id % 2 ) * 2; 01409 01410 // get quads in partition 01411 iBase_EntityHandle quads[4], *ptr = quads; 01412 int junk1 = 4, junk2; 01413 iMesh_getEntities( imesh, local_part, iBase_FACE, iMesh_QUADRILATERAL, &ptr, &junk1, &junk2, &ierr );CHKERR; 01414 assert( ptr == quads ); 01415 assert( junk1 == 4 ); 01416 assert( junk2 == 4 ); 01417 01418 // get vertices in quads 01419 iBase_EntityHandle conn[16]; 01420 int offsets[5], *off_ptr = offsets, junk3 = 5, junk4; 01421 ptr = conn; 01422 junk1 = 16; 01423 iMesh_getEntArrAdj( imesh, quads, 4, iBase_VERTEX, &ptr, &junk1, &junk2, &off_ptr, &junk3, &junk4, &ierr );CHKERR; 01424 assert( ptr == conn ); 01425 assert( junk1 == 16 ); 01426 assert( junk2 == 16 ); 01427 assert( off_ptr == offsets ); 01428 assert( junk3 == 5 ); 01429 assert( junk4 == 5 ); 01430 01431 // make unique vertex list 01432 std::sort( conn, conn + 16 ); 01433 const int num_vtx = std::unique( conn, conn + 16 ) - conn; 01434 assert( 9 == num_vtx ); 01435 01436 // get vertex coords 01437 std::vector< double > coords( 27 ); 01438 ierr = get_coords( imesh, conn, 9, &coords[0] );CHKERR; 01439 01440 // use vertex coords to determine logical position 01441 for( int i = 0; i < num_vtx; ++i ) 01442 { 01443 int x = (int)round( coords[3 * i] - xbase ); 01444 int y = (int)round( coords[3 * i + 1] - ybase ); 01445 if( x < 0 || x > 2 || y < 0 || y > 2 ) 01446 { 01447 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01448 << " Invalid vertex coordinate: (" << coords[3 * i] << ", " << coords[3 * i + 1] << ", " 01449 << coords[3 * i + 2] << ")" << std::endl 01450 << " For logical partition " << local_id << std::endl; 01451 return iBase_FAILURE; 01452 } 01453 verts[x][y] = conn[i]; 01454 } 01455 01456 if( local_logical % 2 ) 01457 { 01458 switch( other_logical - local_logical ) 01459 { 01460 case 0: 01461 return iBase_FAILURE; 01462 case 1: // upper right 01463 vtx_handles.resize( 1 ); 01464 vtx_handles[0] = verts[2][0]; 01465 break; 01466 case 2: // right 01467 vtx_handles.resize( 3 ); 01468 std::copy( verts[2], verts[2] + 3, vtx_handles.begin() ); 01469 break; 01470 case -1: // above 01471 vtx_handles.resize( 3 ); 01472 vtx_handles[0] = verts[0][0]; 01473 vtx_handles[1] = verts[1][0]; 01474 vtx_handles[2] = verts[2][0]; 01475 break; 01476 case -2: // left 01477 vtx_handles.resize( 3 ); 01478 std::copy( verts[0], verts[0] + 3, vtx_handles.begin() ); 01479 break; 01480 case -3: // upper left 01481 vtx_handles.resize( 1 ); 01482 vtx_handles[0] = verts[0][0]; 01483 break; 01484 default: 01485 vtx_handles.clear(); 01486 break; 01487 } 01488 } 01489 else 01490 { 01491 switch( other_logical - local_logical ) 01492 { 01493 case 0: 01494 return iBase_FAILURE; 01495 case 1: // below 01496 vtx_handles.resize( 3 ); 01497 vtx_handles[0] = verts[0][2]; 01498 vtx_handles[1] = verts[1][2]; 01499 vtx_handles[2] = verts[2][2]; 01500 break; 01501 case 2: // right 01502 vtx_handles.resize( 3 ); 01503 std::copy( verts[2], verts[2] + 3, vtx_handles.begin() ); 01504 break; 01505 case 3: // lower right 01506 vtx_handles.resize( 1 ); 01507 vtx_handles[0] = verts[2][2]; 01508 break; 01509 case -1: // lower left 01510 vtx_handles.resize( 1 ); 01511 vtx_handles[0] = verts[0][2]; 01512 break; 01513 case -2: // left 01514 vtx_handles.resize( 3 ); 01515 std::copy( verts[0], verts[0] + 3, vtx_handles.begin() ); 01516 break; 01517 default: 01518 vtx_handles.clear(); 01519 break; 01520 } 01521 } 01522 01523 return iBase_SUCCESS; 01524 } 01525 01526 /**\brief Test querying of part boundary entities 01527 * 01528 * Test: 01529 * - iMeshP_getNumPartBdryEnts 01530 * - iMeshP_getPartBdryEnts 01531 */ 01532 int test_get_part_boundary( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 01533 { 01534 int ierr, rank; 01535 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 01536 01537 // get local part handles and part ids, and global part id list 01538 std::vector< iMeshP_PartHandle > local_handles; 01539 std::vector< iMeshP_Part > local_ids; 01540 std::vector< iMeshP_Part > all_parts = map.get_parts(); 01541 std::map< iMeshP_PartHandle, std::vector< iBase_EntityHandle > > part_bdry; 01542 ierr = get_local_parts( imesh, prtn, local_handles, &local_ids );CHKERR; 01543 01544 // for each combination of local part with any other part, 01545 // check for valid function values. 01546 std::vector< std::pair< iMeshP_Part, iMeshP_Part > > num_failed, num_error, list_failed, list_error, error; 01547 for( size_t i = 0; i < local_handles.size(); ++i ) 01548 { 01549 iMeshP_PartHandle local_handle = local_handles[i]; 01550 iMeshP_Part local_id = local_ids[i]; 01551 for( std::vector< iMeshP_Part >::iterator j = all_parts.begin(); j != all_parts.end(); ++j ) 01552 { 01553 iMeshP_Part other_id = *j; 01554 if( other_id == local_id ) continue; 01555 01556 std::pair< iMeshP_Part, iMeshP_Part > part_pair; 01557 part_pair.first = local_id; 01558 part_pair.second = other_id; 01559 01560 // get expected values 01561 std::vector< iBase_EntityHandle > shared_verts; 01562 ierr = interface_verts( imesh, prtn, local_handle, other_id, map, shared_verts ); 01563 if( ierr != iBase_SUCCESS ) 01564 { 01565 error.push_back( part_pair ); 01566 continue; 01567 } 01568 std::sort( shared_verts.begin(), shared_verts.end() ); 01569 01570 // test iMeshP_getNumPartBdryEnts 01571 int count; 01572 iMeshP_getNumPartBdryEnts( imesh, prtn, local_handle, iBase_VERTEX, iMesh_POINT, other_id, &count, &ierr ); 01573 if( iBase_SUCCESS != ierr ) 01574 num_error.push_back( part_pair ); 01575 else if( count != (int)shared_verts.size() ) 01576 num_failed.push_back( part_pair ); 01577 01578 // test iMeshP_getPartBdryEnts 01579 iBase_EntityHandle* ptr = 0; 01580 int junk = 0; 01581 iMeshP_getPartBdryEnts( imesh, prtn, local_handle, iBase_VERTEX, iMesh_POINT, other_id, &ptr, &junk, &count, 01582 &ierr ); 01583 if( iBase_SUCCESS != ierr ) 01584 list_error.push_back( part_pair ); 01585 else 01586 { 01587 std::copy( ptr, ptr + count, std::back_inserter( part_bdry[local_handles[i]] ) ); 01588 std::sort( ptr, ptr + count ); 01589 if( (int)shared_verts.size() != count || !std::equal( shared_verts.begin(), shared_verts.end(), ptr ) ) 01590 list_failed.push_back( part_pair ); 01591 free( ptr ); 01592 } 01593 } 01594 } 01595 01596 if( !error.empty() ) 01597 { 01598 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01599 << " Internal error for " << error.size() << " part pairs." << std::endl; 01600 ierr = iBase_FAILURE; 01601 } 01602 if( !num_error.empty() ) 01603 { 01604 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01605 << " iMeshP_getNumPartBdryEnts return error for " << num_error.size() << " part pairs." << std::endl; 01606 ierr = iBase_FAILURE; 01607 } 01608 if( !list_error.empty() ) 01609 { 01610 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01611 << " iMeshP_getPartBdryEnts return error for " << list_error.size() << " part pairs." << std::endl; 01612 ierr = iBase_FAILURE; 01613 } 01614 if( !num_failed.empty() ) 01615 { 01616 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01617 << " iMeshP_getNumPartBdryEnts return incorrect results for " << num_failed.size() << " part pairs." 01618 << std::endl; 01619 ierr = iBase_FAILURE; 01620 } 01621 if( !list_failed.empty() ) 01622 { 01623 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01624 << " iMeshP_getPartBdryEnts return incorrect results for " << list_failed.size() << " part pairs." 01625 << std::endl; 01626 ierr = iBase_FAILURE; 01627 } 01628 01629 if( iBase_SUCCESS != ierr ) return ierr; 01630 01631 // test with iMeshP_ALL_PARTS 01632 for( size_t i = 0; i < local_handles.size(); ++i ) 01633 { 01634 std::vector< iBase_EntityHandle >& exp_bdry = part_bdry[local_handles[i]]; 01635 std::sort( exp_bdry.begin(), exp_bdry.end() ); 01636 exp_bdry.erase( std::unique( exp_bdry.begin(), exp_bdry.end() ), exp_bdry.end() ); 01637 std::pair< iMeshP_Part, iMeshP_Part > part_pair; 01638 part_pair.first = local_ids[i]; 01639 part_pair.second = iMeshP_ALL_PARTS; 01640 01641 int num = 0; 01642 iMeshP_getNumPartBdryEnts( imesh, prtn, local_handles[i], iBase_VERTEX, iMesh_POINT, iMeshP_ALL_PARTS, &num, 01643 &ierr ); 01644 if( ierr ) 01645 num_error.push_back( part_pair ); 01646 else if( num != (int)exp_bdry.size() ) 01647 num_failed.push_back( part_pair ); 01648 01649 iBase_EntityHandle* bdry = 0; 01650 int junk = num = 0; 01651 iMeshP_getPartBdryEnts( imesh, prtn, local_handles[i], iBase_VERTEX, iMesh_POINT, iMeshP_ALL_PARTS, &bdry, 01652 &junk, &num, &ierr ); 01653 if( ierr ) 01654 list_error.push_back( part_pair ); 01655 else 01656 { 01657 std::sort( bdry, bdry + num ); 01658 if( num != (int)exp_bdry.size() || !std::equal( bdry, bdry + num, exp_bdry.begin() ) ) 01659 list_failed.push_back( part_pair ); 01660 free( bdry ); 01661 } 01662 } 01663 if( !num_error.empty() ) 01664 { 01665 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01666 << " iMeshP_getNumPartBdryEnts return error for " << num_error.size() << " part pairs." << std::endl; 01667 ierr = iBase_FAILURE; 01668 } 01669 if( !list_error.empty() ) 01670 { 01671 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01672 << " iMeshP_getPartBdryEnts return error for " << list_error.size() << " part pairs." << std::endl; 01673 ierr = iBase_FAILURE; 01674 } 01675 if( !num_failed.empty() ) 01676 { 01677 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01678 << " iMeshP_getNumPartBdryEnts return incorrect results for " << num_failed.size() << " part pairs." 01679 << std::endl; 01680 ierr = iBase_FAILURE; 01681 } 01682 if( !list_failed.empty() ) 01683 { 01684 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01685 << " iMeshP_getPartBdryEnts return incorrect results for " << list_failed.size() << " part pairs." 01686 << std::endl; 01687 ierr = iBase_FAILURE; 01688 } 01689 01690 return ierr; 01691 } 01692 01693 /**\brief Test querying of part boundary entities 01694 * 01695 * Test: 01696 * - iMeshP_initPartBdryEntIter 01697 * - iMeshP_initPartBdryEntArrIter 01698 */ 01699 int test_part_boundary_iter( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 01700 { 01701 int ierr, rank, has_data; 01702 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 01703 01704 // get local part handles and part ids, and global part id list 01705 std::vector< iMeshP_PartHandle > local_handles; 01706 std::vector< iMeshP_Part > local_ids; 01707 std::vector< iMeshP_Part > all_parts = map.get_parts(); 01708 ierr = get_local_parts( imesh, prtn, local_handles, &local_ids );CHKERR; 01709 01710 std::vector< std::pair< iMeshP_Part, iMeshP_Part > > single_failed, single_error, single_step_error, array_failed, 01711 array_error, array_step_error; 01712 for( size_t i = 0; i < local_handles.size(); ++i ) 01713 { 01714 iMeshP_PartHandle local_handle = local_handles[i]; 01715 iMeshP_Part local_id = local_ids[i]; 01716 for( std::vector< iMeshP_Part >::iterator j = all_parts.begin(); j != all_parts.end(); ++j ) 01717 { 01718 iMeshP_Part other_id = *j; 01719 if( other_id == local_id ) continue; 01720 01721 std::pair< iMeshP_Part, iMeshP_Part > part_pair; 01722 part_pair.first = local_id; 01723 part_pair.second = other_id; 01724 01725 // get expected values 01726 std::vector< iBase_EntityHandle > shared_verts; 01727 ierr = interface_verts( imesh, prtn, local_handle, other_id, map, shared_verts ); 01728 if( ierr != iBase_SUCCESS || 0 == shared_verts.size() ) continue; 01729 std::sort( shared_verts.begin(), shared_verts.end() ); 01730 01731 // test single entity iterator 01732 iBase_EntityIterator siter; 01733 iMeshP_initPartBdryEntIter( imesh, prtn, local_handle, iBase_VERTEX, iMesh_POINT, other_id, &siter, &ierr ); 01734 if( ierr != iBase_SUCCESS ) 01735 { 01736 single_error.push_back( part_pair ); 01737 } 01738 else 01739 { 01740 std::vector< iBase_EntityHandle > results; 01741 for( ;; ) 01742 { 01743 iBase_EntityHandle handle; 01744 iMesh_getNextEntIter( imesh, siter, &handle, &has_data, &ierr ); 01745 if( ierr != iBase_SUCCESS ) 01746 { 01747 single_step_error.push_back( part_pair ); 01748 break; 01749 } 01750 if( !has_data ) break; 01751 results.push_back( handle ); 01752 } 01753 01754 std::sort( results.begin(), results.end() ); 01755 if( results.size() != shared_verts.size() || 01756 !std::equal( results.begin(), results.end(), shared_verts.begin() ) ) 01757 single_failed.push_back( part_pair ); 01758 } 01759 iMesh_endEntIter( imesh, siter, &ierr ); 01760 01761 // test array iterator 01762 iBase_EntityArrIterator aiter; 01763 iMeshP_initPartBdryEntArrIter( imesh, prtn, local_handle, iBase_VERTEX, iMesh_POINT, shared_verts.size(), 01764 other_id, &aiter, &ierr ); 01765 if( ierr != iBase_SUCCESS ) 01766 { 01767 array_error.push_back( part_pair ); 01768 iMesh_endEntArrIter( imesh, aiter, &ierr ); 01769 continue; 01770 } 01771 iBase_EntityHandle results[5], *ptr = results; 01772 int junk = 5, count; 01773 iMesh_getNextEntArrIter( imesh, aiter, &ptr, &junk, &count, &has_data, &ierr ); 01774 if( ierr != iBase_SUCCESS || !has_data ) 01775 { 01776 array_step_error.push_back( part_pair ); 01777 iMesh_endEntArrIter( imesh, aiter, &ierr ); 01778 continue; 01779 } 01780 iMesh_endEntArrIter( imesh, aiter, &ierr ); 01781 assert( count <= 5 ); 01782 assert( ptr == results ); 01783 std::sort( ptr, ptr + count ); 01784 if( count != (int)shared_verts.size() || !std::equal( shared_verts.begin(), shared_verts.end(), results ) ) 01785 array_failed.push_back( part_pair ); 01786 } 01787 } 01788 01789 if( !single_error.empty() ) 01790 { 01791 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01792 << " iMeshP_initPartBdryEntIter return error for " << single_error.size() << " part pairs." 01793 << std::endl; 01794 ierr = iBase_FAILURE; 01795 } 01796 if( !single_step_error.empty() ) 01797 { 01798 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01799 << " iMesh_getNextEntIter return error for " << single_step_error.size() << " part pairs." 01800 << std::endl; 01801 ierr = iBase_FAILURE; 01802 } 01803 if( !single_failed.empty() ) 01804 { 01805 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01806 << " iMeshP_initPartBdryEntIter iterator iterated over invalid entities for " << single_failed.size() 01807 << " part pairs." << std::endl; 01808 ierr = iBase_FAILURE; 01809 } 01810 01811 if( !array_error.empty() ) 01812 { 01813 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01814 << " iMeshP_initPartBdryEntArrIter return error for " << array_error.size() << " part pairs." 01815 << std::endl; 01816 ierr = iBase_FAILURE; 01817 } 01818 if( !array_step_error.empty() ) 01819 { 01820 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01821 << " iMesh_getNextEntArrIter return error for " << array_step_error.size() << " part pairs." 01822 << std::endl; 01823 ierr = iBase_FAILURE; 01824 } 01825 if( !array_failed.empty() ) 01826 { 01827 std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl 01828 << " iMeshP_initPartBdryEntArrIter iterator iterated over invalid entities for " 01829 << array_failed.size() << " part pairs." << std::endl; 01830 ierr = iBase_FAILURE; 01831 } 01832 01833 return ierr; 01834 } 01835 01836 /**\brief Test adjacent entity query 01837 * 01838 * Test: 01839 * - iMeshP_getAdjEntities 01840 */ 01841 int test_get_adjacencies( iMesh_Instance /* imesh */, iMeshP_PartitionHandle /* prtn */, const PartMap& ) 01842 { 01843 return iBase_SUCCESS; 01844 } 01845 01846 /**\brief Test entity iterators 01847 * 01848 * Test: 01849 * - iMeshP_initEntIter 01850 * - iMeshP_initEntArrIter 01851 */ 01852 int test_entity_iterator( iMesh_Instance /*imesh */, iMeshP_PartitionHandle /*prtn*/, const PartMap& ) 01853 { 01854 return iBase_SUCCESS; 01855 } 01856 01857 /**\brief Test entity owner queries 01858 * 01859 * Test: 01860 * - iMeshP_getEntOwnerPart 01861 * - iMeshP_getEntOwnerPartArr 01862 * - iMeshP_isEntOwner 01863 * - iMeshP_isEntOwnerArr 01864 */ 01865 int test_entity_owner( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& /* map */ ) 01866 { 01867 int ierr, rank, size; 01868 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 01869 MPI_Comm_size( MPI_COMM_WORLD, &size ); 01870 01871 // get local part handles and part ids 01872 std::vector< iMeshP_PartHandle > local_handles; 01873 std::vector< iMeshP_Part > local_ids; 01874 ierr = get_local_parts( imesh, prtn, local_handles, &local_ids ); 01875 PCHECK; 01876 01877 // test iMeshP_getEntOwnerPart for quads in each part 01878 std::vector< iBase_EntityHandle > all_quads; 01879 std::vector< iMeshP_Part > quad_owners; 01880 int invalid_count = 0; 01881 for( size_t i = 0; i < local_handles.size(); ++i ) 01882 { 01883 std::vector< iBase_EntityHandle > quads; 01884 ierr = get_entities( imesh, local_handles[0], iBase_FACE, iMesh_QUADRILATERAL, quads ); 01885 if( ierr ) break; 01886 01887 for( size_t j = 0; j < quads.size(); ++j ) 01888 { 01889 all_quads.push_back( quads[j] ); 01890 quad_owners.push_back( local_ids[i] ); 01891 iMeshP_Part owner; 01892 iMeshP_getEntOwnerPart( imesh, prtn, quads[j], &owner, &ierr ); 01893 if( iBase_SUCCESS != ierr ) break; 01894 01895 if( owner != local_ids[i] ) ++invalid_count; 01896 } 01897 if( iBase_SUCCESS != ierr ) break; 01898 } 01899 PCHECK; 01900 ASSERT( 0 == invalid_count ); 01901 01902 // test iMeshP_getEntOwnerPartArr for quads in each part 01903 invalid_count = 0; 01904 for( size_t i = 0; i < local_handles.size(); ++i ) 01905 { 01906 std::vector< iBase_EntityHandle > quads; 01907 ierr = get_entities( imesh, local_handles[0], iBase_FACE, iMesh_QUADRILATERAL, quads ); 01908 if( ierr ) break; 01909 01910 std::vector< iMeshP_Part > owners( quads.size() ), expected( quads.size(), local_ids[i] ); 01911 int junk = owners.size(), count; 01912 iMeshP_Part* ptr = &owners[0]; 01913 iMeshP_getEntOwnerPartArr( imesh, prtn, &quads[0], quads.size(), &ptr, &junk, &count, &ierr ); 01914 if( ierr ) break; 01915 assert( ptr == &owners[0] ); 01916 assert( junk == (int)owners.size() ); 01917 assert( count == (int)quads.size() ); 01918 if( owners != expected ) ++invalid_count; 01919 } 01920 PCHECK; 01921 ASSERT( 0 == invalid_count ); 01922 01923 // get all vertices 01924 iBase_EntityHandle* vtx_arr = 0; 01925 int junk1 = 0, num_vtx; 01926 int *junk2 = 0, junk3 = 0, junk4; 01927 iMesh_getEntArrAdj( imesh, &all_quads[0], all_quads.size(), iBase_VERTEX, &vtx_arr, &junk1, &num_vtx, &junk2, 01928 &junk3, &junk4, &ierr ); 01929 PCHECK; 01930 free( junk2 ); 01931 std::sort( vtx_arr, vtx_arr + num_vtx ); 01932 num_vtx = std::unique( vtx_arr, vtx_arr + num_vtx ) - vtx_arr; 01933 std::vector< iBase_EntityHandle > all_verts( vtx_arr, vtx_arr + num_vtx ); 01934 free( vtx_arr ); 01935 01936 // check consistency between iMeshP_getEntOwnerPart and iMeshP_getEntOwnerPartArr 01937 // for all vertices 01938 std::vector< iMeshP_Part > vert_owners( all_verts.size() ); 01939 junk1 = vert_owners.size(); 01940 iMeshP_Part* junk5 = &vert_owners[0]; 01941 iMeshP_getEntOwnerPartArr( imesh, prtn, &all_verts[0], all_verts.size(), &junk5, &junk1, &junk3, &ierr ); 01942 PCHECK; 01943 assert( junk5 == &vert_owners[0] ); 01944 assert( junk1 == (int)vert_owners.size() ); 01945 assert( junk3 == (int)all_verts.size() ); 01946 01947 invalid_count = 0; 01948 for( size_t i = 0; i < all_verts.size(); ++i ) 01949 { 01950 iMeshP_Part owner; 01951 iMeshP_getEntOwnerPart( imesh, prtn, all_verts[i], &owner, &ierr ); 01952 if( iBase_SUCCESS != ierr || owner != vert_owners[i] ) ++invalid_count; 01953 } 01954 ASSERT( 0 == invalid_count ); 01955 01956 // get lists for all entities 01957 std::vector< iBase_EntityHandle > all_entities( all_verts ); 01958 std::copy( all_quads.begin(), all_quads.end(), std::back_inserter( all_entities ) ); 01959 std::vector< iMeshP_Part > all_owners( vert_owners ); 01960 std::copy( quad_owners.begin(), quad_owners.end(), std::back_inserter( all_owners ) ); 01961 01962 // check consistency of iMeshP_isEntOwner for all entities 01963 invalid_count = 0; 01964 ierr = iBase_SUCCESS; 01965 for( size_t i = 0; i < local_handles.size(); ++i ) 01966 { 01967 for( size_t j = 0; ierr == iBase_SUCCESS && j < all_entities.size(); ++j ) 01968 { 01969 int is_owner; 01970 iMeshP_isEntOwner( imesh, prtn, local_handles[i], all_entities[j], &is_owner, &ierr ); 01971 if( ierr != iBase_SUCCESS ) break; 01972 if( !is_owner == ( local_ids[i] == all_owners[j] ) ) ++invalid_count; 01973 } 01974 } 01975 PCHECK; 01976 ASSERT( 0 == invalid_count ); 01977 01978 // check consistency of iMeshP_isEntOwnerArr for all entities 01979 for( size_t i = 0; i < local_handles.size(); ++i ) 01980 { 01981 std::vector< int > is_owner_list( all_entities.size() ); 01982 junk1 = is_owner_list.size(); 01983 int* junk6 = &is_owner_list[0]; 01984 iMeshP_isEntOwnerArr( imesh, prtn, local_handles[i], &all_entities[0], all_entities.size(), &junk6, &junk1, 01985 &junk3, &ierr ); 01986 if( iBase_SUCCESS != ierr ) break; 01987 assert( junk6 == &is_owner_list[0] ); 01988 assert( junk1 == (int)is_owner_list.size() ); 01989 assert( junk3 == (int)all_entities.size() ); 01990 invalid_count = 0; 01991 for( size_t j = 0; j < all_entities.size(); ++j ) 01992 { 01993 if( !( is_owner_list[j] ) == ( local_ids[0] == all_owners[j] ) ) ++invalid_count; 01994 } 01995 } 01996 PCHECK; 01997 ASSERT( 0 == invalid_count ); 01998 01999 // check globally consistent owners for all vertices 02000 02001 // first communicate total number of vertex entries to be sent to root proc 02002 int local_count = all_verts.size(), global_count = 0; 02003 ierr = MPI_Reduce( &local_count, &global_count, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD );CHKERR; 02004 02005 // for each vertex, store { (x << 2) | y, owning part id } 02006 std::vector< int > vtxdata( 2 * all_verts.size() ); 02007 std::vector< double > coords( 3 * all_verts.size() ); 02008 ierr = get_coords( imesh, &all_verts[0], all_verts.size(), &coords[0] );CHKERR; 02009 for( size_t i = 0; i < all_verts.size(); ++i ) 02010 { 02011 int x = (int)round( coords[3 * i] ); 02012 int y = (int)round( coords[3 * i + 1] ); 02013 vtxdata[2 * i] = ( x << 3 ) | y; 02014 vtxdata[2 * i + 1] = vert_owners[i]; 02015 } 02016 02017 // collect all data on root procesor 02018 std::vector< int > all_data( 2 * global_count ); 02019 std::vector< int > displ( size ), counts( size ); 02020 for( int i = 0; i < size; i++ ) 02021 { 02022 counts[i] = vtxdata.size(); 02023 displ[i] = i * vtxdata.size(); 02024 } 02025 02026 // we could have used a simple gather, because all sequences are the same 02027 ierr = MPI_Gatherv( &vtxdata[0], vtxdata.size(), MPI_INT, &all_data[0], &counts[0], &displ[0], MPI_INT, 0, 02028 MPI_COMM_WORLD );CHKERR; 02029 02030 if( rank == 0 ) 02031 { 02032 // map from vertex tag to indices into data 02033 std::multimap< int, int > data_map; 02034 for( int i = 0; i < global_count; ++i ) 02035 { 02036 std::pair< int, int > p; 02037 p.first = all_data[2 * i]; 02038 p.second = i; 02039 data_map.insert( p ); 02040 } 02041 02042 // check consistent data for each vtx 02043 std::multimap< int, int >::const_iterator a, b; 02044 for( a = data_map.begin(); a != data_map.end(); a = b ) 02045 { 02046 for( b = a; b != data_map.end() && a->first == b->first; ++b ) 02047 { 02048 int idx1 = a->second; 02049 int idx2 = b->second; 02050 if( all_data[2 * idx1 + 1] == all_data[2 * idx2 + 1] ) continue; 02051 02052 ierr = iBase_FAILURE; 02053 02054 int proc1 = std::lower_bound( displ.begin(), displ.end(), 2 * idx1 ) - displ.begin(); 02055 if( displ[proc1] != 2 * idx1 ) ++proc1; 02056 int proc2 = std::lower_bound( displ.begin(), displ.end(), 2 * idx2 ) - displ.begin(); 02057 if( displ[proc2] != 2 * idx2 ) ++proc2; 02058 02059 std::cerr << "Error at " __FILE__ ":" << __LINE__ << " : " << std::endl 02060 << " For vertex at (" << ( a->first >> 2 ) << ", " << ( a->first & 3 ) << ") :" << std::endl 02061 << " Processor " << proc1 << " has " << all_data[2 * idx1 + 1] << " as the owning part" 02062 << std::endl 02063 << " Processor " << proc2 << " has " << all_data[2 * idx2 + 1] << " as the owning part" 02064 << std::endl; 02065 } 02066 } 02067 } 02068 02069 return ierr; 02070 } 02071 02072 static int get_part_boundary_verts( iMesh_Instance imesh, 02073 iMeshP_PartitionHandle prtn, 02074 const PartMap& map, 02075 iMeshP_PartHandle part, 02076 std::vector< iBase_EntityHandle >& boundary ) 02077 { 02078 int ierr, logical_id; 02079 ierr = map.part_from_coords( imesh, part, logical_id );CHKERR; 02080 02081 int neighbors[5], num_neighbors; 02082 get_part_neighbors( logical_id, map.get_parts().size(), neighbors, num_neighbors ); 02083 02084 for( int j = 0; j < num_neighbors; ++j ) 02085 { 02086 std::vector< iBase_EntityHandle > iface; 02087 ierr = interface_verts( imesh, prtn, part, neighbors[j], map, iface );CHKERR; 02088 std::copy( iface.begin(), iface.end(), std::back_inserter( boundary ) ); 02089 } 02090 02091 std::sort( boundary.begin(), boundary.end() ); 02092 boundary.erase( std::unique( boundary.begin(), boundary.end() ), boundary.end() ); 02093 return iBase_SUCCESS; 02094 } 02095 02096 /**\brief Test entity status 02097 * 02098 * Test: 02099 * - iMeshP_getEntStatus 02100 * - iMeshP_getEntStatusArr 02101 */ 02102 int test_entity_status( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 02103 { 02104 int ierr, rank, size; 02105 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 02106 MPI_Comm_size( MPI_COMM_WORLD, &size ); 02107 02108 // get local part handles 02109 std::vector< iMeshP_PartHandle > parts; 02110 ierr = get_local_parts( imesh, prtn, parts ); 02111 PCHECK; 02112 02113 // for each part 02114 int num_quad_ent_incorrect = 0, num_quad_ent_error = 0; 02115 int num_quad_arr_incorrect = 0, num_quad_arr_error = 0; 02116 int num_vert_ent_incorrect = 0, num_vert_ent_error = 0; 02117 int num_vert_arr_incorrect = 0, num_vert_arr_error = 0; 02118 for( size_t i = 0; i < parts.size(); ++i ) 02119 { 02120 const iMeshP_PartHandle part = parts[i]; 02121 02122 // get quads and vertices 02123 std::vector< iBase_EntityHandle > quads, verts; 02124 ierr = get_part_quads_and_verts( imesh, part, quads, verts ); 02125 if( ierr ) break; 02126 02127 // check quad status (no ghosting yet) 02128 for( size_t j = 0; j < quads.size(); ++j ) 02129 { 02130 int status; 02131 iMeshP_getEntStatus( imesh, prtn, part, quads[j], &status, &ierr ); 02132 if( ierr != iBase_SUCCESS ) 02133 { 02134 ++num_quad_ent_error; 02135 ierr = iBase_SUCCESS; 02136 continue; 02137 } 02138 02139 if( status != iMeshP_INTERNAL ) ++num_quad_ent_incorrect; 02140 } 02141 02142 // check quad status using iMeshP_getEntStatusArr 02143 std::vector< int > stat_list( quads.size() ); 02144 int* junk1 = &stat_list[0]; 02145 int junk2 = stat_list.size(), count; 02146 iMeshP_getEntStatusArr( imesh, prtn, part, &quads[0], quads.size(), &junk1, &junk2, &count, &ierr ); 02147 if( ierr != iBase_SUCCESS ) 02148 { 02149 ++num_quad_arr_error; 02150 ierr = iBase_SUCCESS; 02151 continue; 02152 } 02153 assert( junk1 == &stat_list[0] ); 02154 assert( junk2 == (int)stat_list.size() ); 02155 assert( count == (int)quads.size() ); 02156 for( size_t j = 0; j < quads.size(); ++j ) 02157 if( stat_list[j] != iMeshP_INTERNAL ) ++num_quad_arr_incorrect; 02158 02159 // figure out which vertices are on the boundary 02160 std::vector< iBase_EntityHandle > boundary; 02161 ierr = get_part_boundary_verts( imesh, prtn, map, part, boundary ); 02162 if( ierr ) break; 02163 std::sort( boundary.begin(), boundary.end() ); 02164 02165 // check vertex status (no ghosting yet) 02166 for( size_t j = 0; j < verts.size(); ++j ) 02167 { 02168 int status; 02169 iMeshP_getEntStatus( imesh, prtn, part, verts[j], &status, &ierr ); 02170 if( ierr != iBase_SUCCESS ) 02171 { 02172 ++num_vert_ent_error; 02173 ierr = iBase_SUCCESS; 02174 continue; 02175 } 02176 bool on_boundary = std::binary_search( boundary.begin(), boundary.end(), verts[j] ); 02177 if( status != ( on_boundary ? iMeshP_BOUNDARY : iMeshP_INTERNAL ) ) ++num_vert_ent_incorrect; 02178 } 02179 02180 // check vert status using iMeshP_getEntStatusArr 02181 stat_list.resize( verts.size() ); 02182 junk1 = &stat_list[0]; 02183 junk2 = stat_list.size(); 02184 iMeshP_getEntStatusArr( imesh, prtn, part, &verts[0], verts.size(), &junk1, &junk2, &count, &ierr ); 02185 if( ierr != iBase_SUCCESS ) 02186 { 02187 ++num_vert_arr_error; 02188 ierr = iBase_SUCCESS; 02189 continue; 02190 } 02191 assert( junk1 == &stat_list[0] ); 02192 assert( junk2 == (int)stat_list.size() ); 02193 assert( count == (int)verts.size() ); 02194 for( size_t j = 0; j < verts.size(); ++j ) 02195 { 02196 bool on_boundary = std::binary_search( boundary.begin(), boundary.end(), verts[j] ); 02197 if( stat_list[j] != ( on_boundary ? iMeshP_BOUNDARY : iMeshP_INTERNAL ) ) ++num_vert_arr_incorrect; 02198 } 02199 } 02200 PCHECK; // check if loop interrupted by any internal errors 02201 02202 ASSERT( 0 == num_quad_ent_error ); 02203 ASSERT( 0 == num_quad_arr_error ); 02204 ASSERT( 0 == num_vert_ent_error ); 02205 ASSERT( 0 == num_vert_arr_error ); 02206 ASSERT( 0 == num_quad_ent_incorrect ); 02207 ASSERT( 0 == num_quad_arr_incorrect ); 02208 ASSERT( 0 == num_vert_ent_incorrect ); 02209 ASSERT( 0 == num_vert_arr_incorrect ); 02210 02211 return iBase_SUCCESS; 02212 } 02213 02214 /**\brief Test information about entity copies for interface entities 02215 * 02216 * Test: 02217 * - iMeshP_getNumCopies 02218 * - iMeshP_getCopyParts 02219 */ 02220 int test_entity_copy_parts( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 02221 { 02222 int ierr, rank, size; 02223 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 02224 MPI_Comm_size( MPI_COMM_WORLD, &size ); 02225 02226 // get local part handles 02227 std::vector< iMeshP_PartHandle > parts; 02228 ierr = get_local_parts( imesh, prtn, parts ); 02229 PCHECK; 02230 ASSERT( !parts.empty() ); 02231 02232 // select a singe part to test 02233 const iMeshP_PartHandle part = parts[0]; 02234 int logical_id; 02235 ierr = map.part_from_coords( imesh, part, logical_id );CHKERR; 02236 const iMeshP_Part part_id = map.part_id_from_local_id( logical_id ); 02237 02238 // get vertices in part 02239 std::vector< iBase_EntityHandle > quads, verts; 02240 ierr = get_part_quads_and_verts( imesh, part, quads, verts ); 02241 PCHECK; 02242 02243 // get neighbors 02244 int neighbors[5], num_neighbors; 02245 get_part_neighbors( logical_id, map.get_parts().size(), neighbors, num_neighbors ); 02246 02247 // build map of sharing data for each vertex 02248 std::map< iBase_EntityHandle, std::vector< iMeshP_Part > > vert_sharing; 02249 for( int j = 0; j < num_neighbors; ++j ) 02250 { 02251 std::vector< iBase_EntityHandle > iface; 02252 ierr = interface_verts( imesh, prtn, part, neighbors[j], map, iface );CHKERR; 02253 for( size_t k = 0; k < iface.size(); ++k ) 02254 vert_sharing[iface[k]].push_back( map.part_id_from_local_id( neighbors[j] ) ); 02255 } 02256 02257 // test getNumCopies for each vertex 02258 std::map< iBase_EntityHandle, std::vector< iMeshP_Part > >::iterator i; 02259 int num_failed = 0, num_incorrect = 0; 02260 for( i = vert_sharing.begin(); i != vert_sharing.end(); ++i ) 02261 { 02262 int count; 02263 iBase_EntityHandle vtx = i->first; 02264 iMeshP_getNumCopies( imesh, prtn, vtx, &count, &ierr ); 02265 if( ierr ) 02266 ++num_failed; 02267 else if( (unsigned)count != i->second.size() + 1 ) // add one for the part we queried from 02268 ++num_incorrect; 02269 } 02270 ASSERT( 0 == num_failed ); 02271 ASSERT( 0 == num_incorrect ); 02272 02273 // get getCopyParts for each vertex 02274 num_failed = num_incorrect = 0; 02275 for( i = vert_sharing.begin(); i != vert_sharing.end(); ++i ) 02276 { 02277 iMeshP_Part* list = 0; 02278 int junk = 0, count; 02279 iMeshP_getCopyParts( imesh, prtn, i->first, &list, &junk, &count, &ierr ); 02280 if( iBase_SUCCESS != ierr ) 02281 { 02282 ++num_failed; 02283 continue; 02284 } 02285 if( (unsigned)count != i->second.size() + 1 ) 02286 { // add one for the part we queried from 02287 ++num_incorrect; 02288 free( list ); 02289 continue; 02290 } 02291 02292 std::vector< iMeshP_Part > expected( i->second ); 02293 expected.push_back( part_id ); 02294 std::sort( list, list + count ); 02295 std::sort( expected.begin(), expected.end() ); 02296 bool eq = std::equal( list, list + count, expected.begin() ); 02297 free( list ); 02298 if( !eq ) ++num_incorrect; 02299 } 02300 ASSERT( 0 == num_failed ); 02301 ASSERT( 0 == num_incorrect ); 02302 02303 return iBase_SUCCESS; 02304 } 02305 02306 // store remote handle data for a vertex 02307 struct VtxCopyData 02308 { 02309 std::vector< iMeshP_Part > parts; 02310 std::vector< iBase_EntityHandle > handles; 02311 }; 02312 02313 /**\brief Test information about entity copies for interface entities 02314 * 02315 * Test: 02316 * - iMeshP_getCopies 02317 * - iMeshP_getCopyOnPart 02318 * - iMeshP_getOwnerCopy 02319 */ 02320 int test_entity_copies( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& /* map */ ) 02321 { 02322 int ierr, rank, size; 02323 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 02324 MPI_Comm_size( MPI_COMM_WORLD, &size ); 02325 02326 // generate a unique ID for each vertex using the coordinates. 02327 // see create_mesh(..): each vertex has integer coordinates (x,y,0) 02328 // with x in [0,inf] and y in [0,4] 02329 // then to an Allgatherv to exchange handles for each processor 02330 02331 // cast everything to iBase_EntityHandle so we can pack it all in one communication 02332 MPI_Datatype tmp_type; 02333 if( sizeof( iBase_EntityHandle ) == sizeof( unsigned ) ) 02334 tmp_type = MPI_UNSIGNED; 02335 else if( sizeof( iBase_EntityHandle ) == sizeof( unsigned long ) ) 02336 tmp_type = MPI_UNSIGNED_LONG; 02337 else if( sizeof( iBase_EntityHandle ) == sizeof( unsigned long long ) ) 02338 tmp_type = MPI_UNSIGNED_LONG_LONG; 02339 else 02340 return iBase_FAILURE; 02341 const MPI_Datatype type = tmp_type; // make it const 02342 02343 // get local part handles 02344 std::vector< iMeshP_PartHandle > parts; 02345 ierr = get_local_parts( imesh, prtn, parts ); 02346 PCHECK; 02347 std::vector< iMeshP_Part > part_ids( parts.size() ); 02348 iMeshP_Part* junk1 = &part_ids[0]; 02349 int junk2 = part_ids.size(), junk3; 02350 iMeshP_getPartIdsFromPartHandlesArr( imesh, prtn, &parts[0], parts.size(), &junk1, &junk2, &junk3, &ierr ); 02351 PCHECK; 02352 assert( junk1 == &part_ids[0] ); 02353 assert( junk2 == (int)part_ids.size() ); 02354 assert( junk3 == (int)parts.size() ); 02355 02356 // build list of {vtx_id, part_id, handle} tuples to send 02357 // also build list of local vertex handles 02358 std::vector< iBase_EntityHandle > local_data, local_vertices; 02359 for( size_t i = 0; i < parts.size(); ++i ) 02360 { 02361 // get vertices 02362 std::vector< iBase_EntityHandle > quads, verts; 02363 ierr = get_part_quads_and_verts( imesh, parts[i], quads, verts ); 02364 if( ierr ) break; 02365 02366 // add all vertices to local_data 02367 for( size_t j = 0; j < verts.size(); ++j ) 02368 { 02369 int tag = 0; 02370 ierr = vertex_tag( imesh, verts[j], tag ); 02371 if( ierr ) break; 02372 long tmp_h = tag; 02373 local_data.push_back( (iBase_EntityHandle)tmp_h ); 02374 tmp_h = part_ids[i]; 02375 local_data.push_back( (iBase_EntityHandle)tmp_h ); 02376 local_data.push_back( verts[j] ); 02377 } 02378 if( ierr ) break; 02379 02380 std::copy( verts.begin(), verts.end(), std::back_inserter( local_vertices ) ); 02381 } 02382 02383 // build list of local vertices 02384 std::sort( local_vertices.begin(), local_vertices.end() ); 02385 local_vertices.erase( std::unique( local_vertices.begin(), local_vertices.end() ), local_vertices.end() ); 02386 std::vector< int > local_vtx_tags( local_vertices.size() );CHKERR; 02387 for( size_t i = 0; i < local_vertices.size(); ++i ) 02388 { 02389 ierr = vertex_tag( imesh, local_vertices[i], local_vtx_tags[i] ); 02390 if( ierr ) break; 02391 } 02392 CHKERR; 02393 02394 // communicate data 02395 std::vector< int > gcounts( size ), gdisp( size ); 02396 int local_data_size = local_data.size(); 02397 ierr = MPI_Allgather( &local_data_size, 1, MPI_INT, &gcounts[0], 1, MPI_INT, MPI_COMM_WORLD );CHKERR; 02398 gdisp[0] = 0; 02399 for( int i = 1; i < size; ++i ) 02400 gdisp[i] = gdisp[i - 1] + gcounts[i - 1]; 02401 std::vector< iBase_EntityHandle > global_data( gdisp[size - 1] + gcounts[size - 1] ); 02402 ierr = MPI_Allgatherv( &local_data[0], local_data_size, type, &global_data[0], &gcounts[0], &gdisp[0], type, 02403 MPI_COMM_WORLD );CHKERR; 02404 02405 // arrange global data in a more useful way 02406 std::map< int, VtxCopyData > vtx_sharing; 02407 assert( global_data.size() % 3 == 0 ); 02408 for( size_t i = 0; i < global_data.size(); i += 3 ) 02409 { 02410 int tag = (int)(size_t)global_data[i]; 02411 iMeshP_Part part = (iMeshP_Part)(size_t)global_data[i + 1]; 02412 iBase_EntityHandle handle = global_data[i + 2]; 02413 vtx_sharing[tag].parts.push_back( part ); 02414 vtx_sharing[tag].handles.push_back( handle ); 02415 } 02416 02417 // test iMeshP_getCopies for each local vertex 02418 int num_error = 0, num_incorrect = 0, junk4; 02419 for( size_t i = 0; i < local_vertices.size(); ++i ) 02420 { 02421 int num_copies = -1; 02422 // iMeshP_Part* part_ids = 0; 02423 iMeshP_Part* ptr_part_ids = 0; // Use ptr_part_ids to avoid shadowing std::vector<iMeshP_Part> part_ids 02424 iBase_EntityHandle* copies = 0; 02425 junk2 = junk3 = junk4 = 0; 02426 iMeshP_getCopies( imesh, prtn, local_vertices[i], &ptr_part_ids, &junk2, &num_copies, &copies, &junk3, &junk4, 02427 &ierr ); 02428 if( iBase_SUCCESS != ierr ) 02429 { 02430 ++num_error; 02431 continue; 02432 } 02433 assert( junk4 == num_copies ); 02434 02435 VtxCopyData& expected = vtx_sharing[local_vtx_tags[i]]; 02436 if( num_copies != (int)expected.parts.size() ) 02437 ++num_incorrect; 02438 else 02439 for( size_t j = 0; j < expected.parts.size(); ++j ) 02440 { 02441 int idx = std::find( ptr_part_ids, ptr_part_ids + num_copies, expected.parts[j] ) - ptr_part_ids; 02442 if( idx == num_copies || copies[idx] != expected.handles[j] ) 02443 { 02444 ++num_incorrect; 02445 break; 02446 } 02447 } 02448 free( ptr_part_ids ); 02449 free( copies ); 02450 } 02451 ASSERT( 0 == num_error ); 02452 ASSERT( 0 == num_incorrect ); 02453 02454 // test iMeshP_getCopyOnPart for each local vertex 02455 num_error = num_incorrect = 0; 02456 for( size_t i = 0; i < local_vertices.size(); ++i ) 02457 { 02458 VtxCopyData& expected = vtx_sharing[local_vtx_tags[i]]; 02459 for( size_t j = 0; j < expected.parts.size(); ++j ) 02460 { 02461 iBase_EntityHandle copy; 02462 iMeshP_getCopyOnPart( imesh, prtn, local_vertices[i], expected.parts[j], ©, &ierr ); 02463 if( iBase_SUCCESS != ierr ) 02464 ++num_error; 02465 else if( expected.handles[j] != copy ) 02466 ++num_incorrect; 02467 } 02468 } 02469 ASSERT( 0 == num_error ); 02470 ASSERT( 0 == num_incorrect ); 02471 02472 // test iMeshP_getOwnerCopy for each local vertex 02473 num_error = num_incorrect = 0; 02474 for( size_t i = 0; i < local_vertices.size(); ++i ) 02475 { 02476 VtxCopyData& expected = vtx_sharing[local_vtx_tags[i]]; 02477 iMeshP_Part owner_id = 0; 02478 iMeshP_getEntOwnerPart( imesh, prtn, local_vertices[i], &owner_id, &ierr ); 02479 if( iBase_SUCCESS != ierr ) continue; // not testing getEntOwnerPart here 02480 02481 size_t idx = std::find( expected.parts.begin(), expected.parts.end(), owner_id ) - expected.parts.begin(); 02482 if( idx == expected.parts.size() ) continue; // not testing getEntOwnerPart here 02483 02484 iMeshP_Part owner_id_2 = 0; 02485 iBase_EntityHandle copy = 0; 02486 iMeshP_getOwnerCopy( imesh, prtn, local_vertices[i], &owner_id_2, ©, &ierr ); 02487 if( iBase_SUCCESS != ierr ) 02488 ++num_error; 02489 else if( owner_id_2 != owner_id && copy != expected.handles[idx] ) 02490 ++num_incorrect; 02491 } 02492 ASSERT( 0 == num_error ); 02493 ASSERT( 0 == num_incorrect ); 02494 02495 return iBase_SUCCESS; 02496 } 02497 02498 int get_num_adj_quads( iMesh_Instance imesh, iBase_EntityHandle vtx, int& num ) 02499 { 02500 iBase_EntityHandle* list = 0; 02501 int ierr, junk = 0; 02502 iMesh_getEntAdj( imesh, vtx, iBase_FACE, &list, &junk, &num, &ierr ); 02503 if( iBase_SUCCESS == ierr ) free( list ); 02504 return ierr; 02505 } 02506 02507 int get_adj( iMesh_Instance imesh, iBase_EntityHandle ent, int type, std::vector< iBase_EntityHandle >& adj ) 02508 { 02509 iBase_EntityHandle* list = 0; 02510 int ierr, num, junk = 0; 02511 iMesh_getEntAdj( imesh, ent, type, &list, &junk, &num, &ierr ); 02512 if( iBase_SUCCESS == ierr ) 02513 { 02514 std::copy( list, list + num, std::back_inserter( adj ) ); 02515 free( list ); 02516 } 02517 return ierr; 02518 } 02519 02520 // assume regular quad mesh 02521 int get_boundary_vertices( iMesh_Instance imesh, std::vector< iBase_EntityHandle >& bdry ) 02522 { 02523 int ierr, n; 02524 iBase_EntitySetHandle root; 02525 iMesh_getRootSet( imesh, &root, &ierr );CHKERR; 02526 std::vector< iBase_EntityHandle > all_verts; 02527 ierr = get_entities( imesh, root, iBase_VERTEX, iMesh_POINT, all_verts );CHKERR; 02528 bdry.clear(); 02529 for( size_t i = 0; i < all_verts.size(); ++i ) 02530 { 02531 ierr = get_num_adj_quads( imesh, all_verts[i], n );CHKERR; 02532 if( n != 4 ) bdry.push_back( all_verts[i] ); 02533 } 02534 return iBase_SUCCESS; 02535 } 02536 02537 int check_one_layer( iMesh_Instance imesh, 02538 iBase_EntityHandle vtx, 02539 const std::vector< iBase_EntityHandle >& sorted_vertices ) 02540 { 02541 int ierr; 02542 if( std::binary_search( sorted_vertices.begin(), sorted_vertices.end(), vtx ) ) return iBase_SUCCESS; 02543 std::vector< iBase_EntityHandle > quads, verts; 02544 ierr = get_adj( imesh, vtx, iBase_FACE, quads );CHKERR; 02545 for( size_t i = 0; i < quads.size(); ++i ) 02546 { 02547 verts.clear(); 02548 ierr = get_adj( imesh, quads[i], iBase_VERTEX, verts );CHKERR; 02549 for( size_t j = 0; j < verts.size(); ++j ) 02550 { 02551 if( std::binary_search( sorted_vertices.begin(), sorted_vertices.end(), verts[j] ) ) return iBase_SUCCESS; 02552 } 02553 } 02554 02555 return iBase_FAILURE; 02556 } 02557 02558 // get number of adjacent quads to each vertex, both on the local 02559 // processor and in the entire mesh 02560 int get_num_adj_all( iMesh_Instance imesh, 02561 const std::vector< iBase_EntityHandle >& verts, 02562 std::vector< int >& num_local_adj, 02563 std::vector< int >& num_all_adj ) 02564 { 02565 int ierr, size; 02566 MPI_Comm_size( MPI_COMM_WORLD, &size ); 02567 02568 std::vector< int > vtx_tags( verts.size() ); 02569 num_local_adj.resize( verts.size() ); 02570 for( size_t i = 0; i < verts.size(); ++i ) 02571 { 02572 ierr = get_num_adj_quads( imesh, verts[i], num_local_adj[i] );CHKERR; 02573 ierr = vertex_tag( imesh, verts[i], vtx_tags[i] );CHKERR; 02574 } 02575 02576 std::vector< int > counts( size ), displ( size ); 02577 int num_vtx = verts.size(); 02578 ierr = MPI_Allgather( &num_vtx, 1, MPI_INT, &counts[0], 1, MPI_INT, MPI_COMM_WORLD );CHKERR; 02579 displ[0] = 0; 02580 for( int i = 1; i < size; ++i ) 02581 displ[i] = displ[i - 1] + counts[i - 1]; 02582 int total = displ[size - 1] + counts[size - 1]; 02583 std::vector< int > all_tags( total ), all_adj_counts( total ); 02584 ierr = MPI_Allgatherv( &vtx_tags[0], vtx_tags.size(), MPI_INT, &all_tags[0], &counts[0], &displ[0], MPI_INT, 02585 MPI_COMM_WORLD );CHKERR; 02586 ierr = MPI_Allgatherv( &num_local_adj[0], num_local_adj.size(), MPI_INT, &all_adj_counts[0], &counts[0], &displ[0], 02587 MPI_INT, MPI_COMM_WORLD );CHKERR; 02588 02589 num_all_adj.clear(); 02590 num_all_adj.resize( total, 0 ); 02591 for( int i = 0; i < total; ++i ) 02592 { 02593 std::vector< int >::iterator it = std::find( vtx_tags.begin(), vtx_tags.end(), all_tags[i] ); 02594 if( it == vtx_tags.end() ) continue; 02595 int idx = it - vtx_tags.begin(); 02596 num_all_adj[idx] += all_adj_counts[i]; 02597 } 02598 02599 return iBase_SUCCESS; 02600 } 02601 02602 /**\brief Test creation of ghost entities 02603 * 02604 * Test: 02605 * - iMeshP_createGhostEntsAll 02606 */ 02607 int test_create_ghost_ents( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& /* map */ ) 02608 { 02609 int ierr; 02610 02611 // get boundary vertices 02612 std::vector< iBase_EntityHandle > bdry; 02613 ierr = get_boundary_vertices( imesh, bdry ); 02614 PCHECK; 02615 // get counts of adjacent entities 02616 std::vector< int > num_local_adj, num_global_adj; 02617 ierr = get_num_adj_all( imesh, bdry, num_local_adj, num_global_adj ); 02618 PCHECK; 02619 // create one layer of ghost entities 02620 iMeshP_createGhostEntsAll( imesh, prtn, iBase_FACE, iBase_VERTEX, 1, 0, &ierr ); 02621 PCHECK; 02622 // check that each vertex has the correct number of adjacent entities 02623 int num_incorrect = 0; 02624 for( size_t i = 0; i < bdry.size(); ++i ) 02625 { 02626 int n; 02627 ierr = get_num_adj_quads( imesh, bdry[i], n ); 02628 if( iBase_SUCCESS != ierr || num_global_adj[i] != n ) ++num_incorrect; 02629 } 02630 ASSERT( 0 == num_incorrect ); 02631 // get new the new boundary 02632 std::vector< iBase_EntityHandle > new_bdry; 02633 ierr = get_boundary_vertices( imesh, new_bdry ); 02634 PCHECK; 02635 // check that each vertex on the new boundary is separated by 02636 // at most one layer from the old boundary 02637 std::sort( bdry.begin(), bdry.end() ); 02638 num_incorrect = 0; 02639 for( size_t i = 0; i < new_bdry.size(); ++i ) 02640 { 02641 ierr = check_one_layer( imesh, new_bdry[i], bdry ); 02642 if( ierr ) ++num_incorrect; 02643 } 02644 ASSERT( 0 == num_incorrect ); 02645 // make another layer of ghost entiites 02646 bdry.swap( new_bdry ); 02647 new_bdry.clear(); 02648 ierr = get_num_adj_all( imesh, bdry, num_local_adj, num_global_adj ); 02649 PCHECK; 02650 iMeshP_createGhostEntsAll( imesh, prtn, iBase_FACE, iBase_VERTEX, 2, 0, &ierr ); 02651 PCHECK; 02652 // check that each vertex has the correct number of adjacent entities 02653 num_incorrect = 0; 02654 for( size_t i = 0; i < bdry.size(); ++i ) 02655 { 02656 int n; 02657 ierr = get_num_adj_quads( imesh, bdry[i], n ); 02658 if( iBase_SUCCESS != ierr || num_global_adj[i] != n ) ++num_incorrect; 02659 } 02660 // check that each vertex on the new boundary is separated by 02661 // at most one layer from the old boundary 02662 std::sort( bdry.begin(), bdry.end() ); 02663 num_incorrect = 0; 02664 for( size_t i = 0; i < new_bdry.size(); ++i ) 02665 { 02666 ierr = check_one_layer( imesh, new_bdry[i], bdry ); 02667 if( ierr ) ++num_incorrect; 02668 } 02669 ASSERT( 0 == num_incorrect ); 02670 02671 return iBase_SUCCESS; 02672 } 02673 02674 /**\brief Test exchange entities 02675 * 02676 * Test: 02677 * - iMeshP_exchEntArrToPartsAll 02678 */ 02679 int test_exchange_ents( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map ) 02680 { 02681 int ierr, rank, size; 02682 int num_err = 0; 02683 iMeshP_RequestHandle request; 02684 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 02685 MPI_Comm_size( MPI_COMM_WORLD, &size ); 02686 02687 std::vector< iBase_EntityHandle > all_elems; 02688 std::vector< iMeshP_Part > all_ids; 02689 std::vector< iBase_EntityHandle > quads; 02690 02691 // get local part handles and part ids 02692 std::vector< iMeshP_PartHandle > local_handles; 02693 std::vector< iMeshP_Part > local_ids; 02694 ierr = get_local_parts( imesh, prtn, local_handles, &local_ids ); 02695 PCHECK; 02696 02697 // get loacal quads before exchange 02698 quads.clear(); 02699 ierr = get_entities( imesh, local_handles[0], iBase_FACE, iMesh_QUADRILATERAL, quads );CHKERR; 02700 int n_quads = quads.size(); 02701 02702 // send all elements in local processor to all other processors 02703 for( size_t i = 0; i < map.get_parts().size(); ++i ) 02704 { 02705 if( map.get_parts()[i] == (unsigned int)rank ) continue; // skip own rank 02706 02707 for( int j = 0; j < n_quads; j++ ) 02708 { 02709 all_elems.push_back( quads[j] ); 02710 all_ids.push_back( map.get_parts()[i] ); 02711 } 02712 } 02713 02714 // exchange entities 02715 iMeshP_exchEntArrToPartsAll( imesh, prtn, &all_elems[0], all_elems.size(), &all_ids[0], 0, 0, &request, &ierr ); 02716 if( iBase_SUCCESS != ierr ) ++num_err; 02717 02718 // get local quads after exchange 02719 quads.clear(); 02720 ierr = get_entities( imesh, local_handles[0], iBase_FACE, iMesh_QUADRILATERAL, quads );CHKERR; 02721 02722 // # of elements should be # of quads * # of processors 02723 ASSERT( quads.size() == (unsigned int)n_quads * size ); 02724 02725 ASSERT( 0 == num_err ); 02726 02727 return iBase_SUCCESS; 02728 } 02729 02730 /**\brief Test commuinication of tag data 02731 * 02732 * Test: 02733 * - iMeshP_pushTags 02734 * - iMeshP_pushTagsEnt 02735 */ 02736 int test_push_tag_data_common( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, int num_ghost_layers ) 02737 { 02738 const char* src_name = "test_src"; 02739 const char* dst_name = "test_dst"; 02740 int ierr, rank; 02741 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 02742 02743 if( num_ghost_layers ) 02744 { 02745 iMeshP_createGhostEntsAll( imesh, prtn, iBase_FACE, iBase_VERTEX, num_ghost_layers, 0, &ierr ); 02746 PCHECK; 02747 } 02748 02749 iBase_TagHandle src_tag, dst_tag; 02750 iMesh_createTag( imesh, src_name, 1, iBase_INTEGER, &src_tag, &ierr, strlen( src_name ) );CHKERR; 02751 iMesh_createTag( imesh, dst_name, 1, iBase_INTEGER, &dst_tag, &ierr, strlen( dst_name ) );CHKERR; 02752 02753 iBase_EntitySetHandle root; 02754 iMesh_getRootSet( imesh, &root, &ierr );CHKERR; 02755 02756 std::vector< iBase_EntityHandle > verts; 02757 ierr = get_entities( imesh, root, iBase_VERTEX, iMesh_POINT, verts );CHKERR; 02758 02759 // test iMeshP_pushTags 02760 // each processor writes its rank on all vertices 02761 // after push, each vertex should be tagged with the rank of its owner 02762 02763 std::vector< int > tag_vals( verts.size(), rank ); 02764 iMesh_setIntArrData( imesh, &verts[0], verts.size(), src_tag, &tag_vals[0], tag_vals.size(), &ierr );CHKERR; 02765 02766 iMeshP_pushTags( imesh, prtn, src_tag, dst_tag, iBase_VERTEX, iMesh_POINT, &ierr ); 02767 PCHECK; 02768 02769 tag_vals.clear(); 02770 tag_vals.resize( verts.size(), -1 ); 02771 iBase_TagHandle id_tag; 02772 iMesh_getTagHandle( imesh, "GLOBAL_ID", &id_tag, &ierr, strlen( "GLOBAL_ID" ) ); 02773 std::vector< int > ids( verts.size() ); 02774 int *junk1 = &ids[0], junk2 = ids.size(), junk3; 02775 iMesh_getIntArrData( imesh, &verts[0], verts.size(), id_tag, &junk1, &junk2, &junk3, &ierr ); 02776 PCHECK; 02777 int errcount = 0; 02778 for( size_t i = 0; i < verts.size(); ++i ) 02779 { 02780 iMesh_getIntData( imesh, verts[i], dst_tag, &tag_vals[i], &ierr ); 02781 if( ierr != iBase_SUCCESS ) 02782 { 02783 std::cerr << "Rank " << rank << " : getIntData failed for vertex " << ids[i] << std::endl; 02784 std::cerr.flush(); 02785 ++errcount; 02786 } 02787 } 02788 ASSERT( 0 == errcount ); 02789 02790 // int *junk1 = &tag_vals[0], junk2 = tag_vals.size(), junk3; 02791 // iMesh_getIntArrData( imesh, &verts[0], verts.size(), dst_tag, &junk1, &junk2, &junk3, &ierr 02792 // ); PCHECK; assert( junk1 == &tag_vals[0] ); assert( junk2 == (int)tag_vals.size() ); assert( 02793 // junk3 == (int)verts.size() ); 02794 02795 std::vector< int > expected( verts.size() ); 02796 std::vector< iMeshP_Part > parts( verts.size() ); 02797 iMeshP_Part* junk4 = &parts[0]; 02798 junk2 = parts.size(); 02799 iMeshP_getEntOwnerPartArr( imesh, prtn, &verts[0], verts.size(), &junk4, &junk2, &junk3, &ierr ); 02800 PCHECK; 02801 assert( junk4 == &parts[0] ); 02802 assert( junk2 == (int)parts.size() ); 02803 assert( junk3 == (int)verts.size() ); 02804 junk1 = &expected[0]; 02805 junk2 = expected.size(); 02806 iMeshP_getRankOfPartArr( imesh, prtn, &parts[0], parts.size(), &junk1, &junk2, &junk3, &ierr ); 02807 PCHECK; 02808 assert( junk1 == &expected[0] ); 02809 assert( junk2 == (int)expected.size() ); 02810 assert( junk3 == (int)parts.size() ); 02811 02812 ASSERT( tag_vals == expected ); 02813 02814 // test iMeshP_pushTagsEnt 02815 // write -1 on all vertices 02816 // For each vertex owned by this processor and shared with more than 02817 // two others, write the rank of the owning processor. 02818 02819 tag_vals.clear(); 02820 tag_vals.resize( verts.size(), -1 ); 02821 iMesh_setIntArrData( imesh, &verts[0], verts.size(), src_tag, &tag_vals[0], tag_vals.size(), &ierr ); 02822 PCHECK; 02823 tag_vals.resize( verts.size(), -1 ); 02824 iMesh_setIntArrData( imesh, &verts[0], verts.size(), dst_tag, &tag_vals[0], tag_vals.size(), &ierr ); 02825 PCHECK; 02826 02827 std::vector< iBase_EntityHandle > some; 02828 for( size_t i = 0; i < verts.size(); ++i ) 02829 { 02830 int num; 02831 iMeshP_getNumCopies( imesh, prtn, verts[i], &num, &ierr ); 02832 if( iBase_SUCCESS != ierr ) break; 02833 if( num > 2 ) 02834 some.push_back( verts[i] ); 02835 else 02836 expected[i] = -1; 02837 } 02838 02839 tag_vals.clear(); 02840 tag_vals.resize( some.size(), rank ); 02841 iMesh_setIntArrData( imesh, &some[0], some.size(), src_tag, &tag_vals[0], tag_vals.size(), &ierr ); 02842 PCHECK; 02843 02844 iMeshP_pushTagsEnt( imesh, prtn, src_tag, dst_tag, &some[0], some.size(), &ierr ); 02845 PCHECK; 02846 02847 tag_vals.clear(); 02848 tag_vals.resize( verts.size(), -1 ); 02849 junk1 = &tag_vals[0]; 02850 junk2 = tag_vals.size(); 02851 iMesh_getIntArrData( imesh, &verts[0], verts.size(), dst_tag, &junk1, &junk2, &junk3, &ierr );CHKERR; 02852 assert( junk1 == &tag_vals[0] ); 02853 assert( junk2 == (int)tag_vals.size() ); 02854 assert( junk3 == (int)verts.size() ); 02855 02856 ASSERT( tag_vals == expected ); 02857 return iBase_SUCCESS; 02858 } 02859 02860 /**\brief Test commuinication of tag data 02861 * 02862 * Test: 02863 * - iMeshP_pushTags 02864 * - iMeshP_pushTagsEnt 02865 */ 02866 int test_push_tag_data_iface( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& ) 02867 { 02868 return test_push_tag_data_common( imesh, prtn, 0 ); 02869 } 02870 02871 /**\brief Test commuinication of tag data 02872 * 02873 * Test: 02874 * - iMeshP_pushTags 02875 * - iMeshP_pushTagsEnt 02876 */ 02877 int test_push_tag_data_ghost( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& ) 02878 { 02879 return test_push_tag_data_common( imesh, prtn, 1 ); 02880 } 02881 02882 /************************************************************************** 02883 PartMap class 02884 **************************************************************************/ 02885 02886 int PartMap::build_map( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, int num_expected_parts ) 02887 { 02888 int ierr, rank, size; 02889 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 02890 MPI_Comm_size( MPI_COMM_WORLD, &size ); 02891 02892 // get local parts 02893 std::vector< iMeshP_PartHandle > local_parts; 02894 std::vector< iMeshP_Part > imesh_ids; 02895 ierr = get_local_parts( imesh, prtn, local_parts, &imesh_ids );CHKERR; 02896 02897 // get logical ids for local parts 02898 std::vector< int > local_ids( local_parts.size() ); 02899 for( size_t i = 0; i < local_parts.size(); ++i ) 02900 { 02901 ierr = part_from_coords( imesh, local_parts[i], local_ids[i] );CHKERR; 02902 } 02903 02904 // get total number of parts 02905 int num_global = 0, num_local = local_parts.size(); 02906 ierr = MPI_Allreduce( &num_local, &num_global, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD );CHKERR; 02907 if( num_global != num_expected_parts ) 02908 { 02909 std::cerr << "Invalid/unexpected global part count at " __FILE__ ":" << __LINE__ << " (proc " << rank 02910 << "): " << std::endl 02911 << " Expected: " << num_expected_parts << std::endl 02912 << " Actual: " << num_global << std::endl; 02913 return 1; 02914 } 02915 02916 // get counts and displacements for Allgatherv calls 02917 std::vector< int > dspls( size ), counts( size ); 02918 ierr = MPI_Allgather( &num_local, 1, MPI_INT, &counts[0], 1, MPI_INT, MPI_COMM_WORLD );CHKERR; 02919 dspls[0] = 0; 02920 for( int i = 1; i < size; ++i ) 02921 dspls[i] = dspls[i - 1] + counts[i - 1]; 02922 02923 // gather iMeshP_Part list from each processor 02924 std::vector< unsigned > global_part_ids( num_expected_parts ); 02925 assert( sizeof( iMeshP_Part ) == sizeof( int ) ); 02926 ierr = MPI_Allgatherv( &imesh_ids[0], num_local, MPI_UNSIGNED, &global_part_ids[0], &counts[0], &dspls[0], 02927 MPI_UNSIGNED, MPI_COMM_WORLD );CHKERR; 02928 02929 // gather local ids from each processor 02930 std::vector< int > global_id_list( num_expected_parts ); 02931 ierr = MPI_Allgatherv( &local_ids[0], num_local, MPI_INT, &global_id_list[0], &counts[0], &dspls[0], MPI_INT, 02932 MPI_COMM_WORLD );CHKERR; 02933 02934 // build owner list 02935 std::vector< int > global_owners( num_expected_parts ); 02936 for( int i = 0; i < size; ++i ) 02937 for( int j = 0; j < counts[i]; ++j ) 02938 global_owners[dspls[i] + j] = i; 02939 02940 // populate member lists 02941 sortedPartList = global_part_ids; 02942 std::sort( sortedPartList.begin(), sortedPartList.end() ); 02943 partLocalIds.resize( num_expected_parts ); 02944 partRanks.resize( num_expected_parts ); 02945 for( int i = 0; i < num_expected_parts; ++i ) 02946 { 02947 int idx = std::lower_bound( sortedPartList.begin(), sortedPartList.end(), global_part_ids[i] ) - 02948 sortedPartList.begin(); 02949 partLocalIds[idx] = global_id_list[i]; 02950 partRanks[idx] = global_owners[i]; 02951 } 02952 02953 // do some consistency checking 02954 if( std::unique( sortedPartList.begin(), sortedPartList.end() ) != sortedPartList.end() ) 02955 { 02956 if( rank == 0 ) 02957 { 02958 std::cerr << "ERROR: Duplicate iMeshP_Part values detected at " __FILE__ ":" << __LINE__ << std::endl; 02959 } 02960 return 1; 02961 } 02962 02963 // build revesre local id map and check for duplicates 02964 localIdReverseMap.clear(); 02965 localIdReverseMap.resize( num_expected_parts, -1 ); 02966 for( int i = 0; i < num_expected_parts; ++i ) 02967 { 02968 int idx = partLocalIds[i]; 02969 if( localIdReverseMap[idx] != -1 ) 02970 { 02971 if( rank == 0 ) 02972 { 02973 std::cerr << "ERROR: Part mesh has been duplicated in multiple parts." << std::endl 02974 << " Detected at " __FILE__ ":" << __LINE__ << std::endl 02975 << " See PartMap::part_from_coords" << std::endl; 02976 } 02977 return 1; 02978 } 02979 if( idx >= num_expected_parts ) 02980 { 02981 if( rank == 0 ) 02982 { 02983 std::cerr << "ERROR: Part mesh invalid/incorrect mesh." << std::endl 02984 << " Detected at " __FILE__ ":" << __LINE__ << std::endl 02985 << " See PartMap::part_from_coords" << std::endl; 02986 } 02987 return 1; 02988 } 02989 02990 localIdReverseMap[idx] = i; 02991 } 02992 02993 return 0; 02994 } 02995 02996 void PartMap::part_id_from_rank( int rank, std::vector< iMeshP_Part >& parts ) const 02997 { 02998 for( size_t i = 0; i < sortedPartList.size(); ++i ) 02999 if( partRanks[i] == rank ) parts.push_back( sortedPartList[i] ); 03000 } 03001 03002 void PartMap::local_id_from_rank( int rank, std::vector< int >& ids ) const 03003 { 03004 for( size_t i = 0; i < sortedPartList.size(); ++i ) 03005 if( partRanks[i] == rank ) ids.push_back( partLocalIds[i] ); 03006 } 03007 03008 int PartMap::part_from_coords( iMesh_Instance imesh, iMeshP_PartHandle part, int& id ) 03009 { 03010 int ierr, rank; 03011 MPI_Comm_rank( MPI_COMM_WORLD, &rank ); 03012 03013 // get elements 03014 const int num_elem = 4; 03015 iBase_EntityHandle array[num_elem]; 03016 iBase_EntityHandle* ptr = array; 03017 int junk1 = num_elem, n = -1; 03018 iMesh_getEntities( imesh, part, iBase_FACE, iMesh_QUADRILATERAL, &ptr, &junk1, &n, &ierr );CHKERR; 03019 assert( ptr == array ); 03020 assert( junk1 == num_elem ); 03021 if( n != num_elem ) 03022 { 03023 std::cerr << "Internal error at " __FILE__ ":" << __LINE__ << " (proc " << rank 03024 << "): Expected all parts to have " << num_elem << " elements. Found one with " << n << std::endl; 03025 return 1; 03026 } 03027 03028 // get vertices 03029 iBase_EntityHandle adj_array[4 * num_elem]; 03030 int junk2, junk3, offset_array[5]; 03031 ptr = adj_array; 03032 junk1 = sizeof( adj_array ) / sizeof( adj_array[0] ); 03033 junk2 = sizeof( offset_array ) / sizeof( offset_array[0] ); 03034 int* ptr2 = offset_array; 03035 iMesh_getEntArrAdj( imesh, array, num_elem, iBase_VERTEX, &ptr, &junk1, &n, &ptr2, &junk2, &junk3, &ierr );CHKERR; 03036 assert( ptr == adj_array ); 03037 assert( ptr2 == offset_array ); 03038 assert( junk1 == sizeof( adj_array ) / sizeof( adj_array[0] ) ); 03039 assert( junk2 == sizeof( offset_array ) / sizeof( offset_array[0] ) ); 03040 assert( n == 4 * num_elem ); 03041 assert( offset_array[0] == 0 ); 03042 for( int i = 1; i < junk3; ++i ) 03043 assert( offset_array[i] - offset_array[i - 1] == 4 ); 03044 03045 // find center vertex 03046 iBase_EntityHandle vtx; 03047 bool all_match; 03048 for( int i = 0; i < 4; ++i ) 03049 { 03050 vtx = adj_array[i]; 03051 all_match = true; 03052 for( int j = 1; j < 4; ++j ) 03053 { 03054 iBase_EntityHandle* mvtx = adj_array + 4 * j; 03055 int k; 03056 for( k = 0; k < 4; ++k ) 03057 if( mvtx[k] == vtx ) break; 03058 if( k == 4 ) all_match = false; 03059 } 03060 if( all_match ) break; 03061 } 03062 assert( all_match ); 03063 03064 // get center vertex coordinates 03065 double x, y, z; 03066 iMesh_getVtxCoord( imesh, vtx, &x, &y, &z, &ierr );CHKERR; 03067 assert( 0.0 == z ); 03068 const int xi = ( (int)round( x ) - 1 ) / 2; 03069 const int yi = ( (int)round( y ) - 1 ) / 2; 03070 assert( xi >= 0 ); 03071 assert( yi >= 0 ); 03072 assert( fabs( x - 2 * xi - 1 ) < 1e-12 ); 03073 assert( fabs( y - 2 * yi - 1 ) < 1e-12 ); 03074 03075 id = 2 * xi + yi; 03076 return 0; 03077 }