MOAB: Mesh Oriented datABase  (version 5.5.0)
MOAB_iMeshP_unit_tests.cpp
Go to the documentation of this file.
1 #include "iMeshP.h"
2 #include "moab_mpi.h"
3 #include <iostream>
4 #include <algorithm>
5 #include <vector>
6 #include <sstream>
7 #include <cassert>
8 #include <cmath>
9 #include <map>
10 #include <cstring>
11 #include <cstdio> // remove()
12 
13 #if !defined( _MSC_VER ) && !defined( __MINGW32__ )
14 #include <unistd.h>
15 #endif
16 
17 #define STRINGIFY_( X ) #X
18 #define STRINGIFY( X ) STRINGIFY_( X )
19 const char* const FILENAME = "iMeshP_test_file";
20 
21 /**************************************************************************
22  Error Checking
23  **************************************************************************/
24 
25 #define CHKERR \
26  do \
27  { \
28  if( ierr ) \
29  { \
30  std::cerr << "Error code " << ierr << " at " << __FILE__ << ":" << __LINE__ << std::endl; \
31  return ierr; \
32  } \
33  } while( false )
34 
35 #define PCHECK \
36  do \
37  { \
38  ierr = is_any_proc_error( ierr ); \
39  CHKERR; \
40  } while( false )
41 
42 // Use my_rank instead of rank to avoid shadowed declaration
43 #define ASSERT( A ) \
44  do \
45  { \
46  if( is_any_proc_error( !( A ) ) ) \
47  { \
48  int my_rank = 0; \
49  MPI_Comm_rank( MPI_COMM_WORLD, &my_rank ); \
50  if( 0 == my_rank ) \
51  std::cerr << "Failed assertion: " #A << std::endl << " at " __FILE__ ":" << __LINE__ << std::endl; \
52  return 1; \
53  } \
54  } while( false )
55 
56 // Test if is_my_error is non-zero on any processor in MPI_COMM_WORLD
57 int is_any_proc_error( int is_my_error )
58 {
59  int result;
60  int err = MPI_Allreduce( &is_my_error, &result, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD );
61  return err || result;
62 }
63 
64 /**************************************************************************
65  Test Declarations
66  **************************************************************************/
67 
68 class PartMap;
69 
70 /**\brief Consistency check for parallel load
71  *
72  * All other tests depend on this one.
73  */
74 int test_load( iMesh_Instance, iMeshP_PartitionHandle prtn, PartMap& map, int comm_size );
75 
76 /**\brief Test partition query methods
77  *
78  * Test:
79  * - iMeshP_getPartitionComm
80  * - iMeshP_getNumPartitions
81  * - iMeshP_getPartitions
82  */
84 
85 /**\brief Test part quyery methods
86  *
87  * Test:
88  * - iMeshP_getNumGlobalParts
89  * - iMeshP_getNumLocalParts
90  * - iMeshP_getLocalParts
91  */
93 
94 /**\brief Test query by entity type
95  *
96  * Test:
97  * - iMeshP_getNumOfTypeAll
98  * - iMeshP_getNumOfType
99  * - iMeshP_getEntities
100  * -
101  */
103 
104 /**\brief Test query by entity topology
105  *
106  * Test:
107  * - iMeshP_getNumOfTopoAll
108  * - iMeshP_getNumOfTopo
109  * - iMeshP_getEntities
110  * -
111  */
113 
114 /**\brief Test mapping from part id to part handle
115  *
116  * Test:
117  * - iMeshP_getPartIdFromPartHandle
118  * - iMeshP_getPartIdsFromPartHandlesArr
119  * - iMeshP_getPartHandleFromPartId
120  * - iMeshP_getPartHandlesFromPartsIdsArr
121  * - iMeshP_getRankOfPart
122  * - iMeshP_getRankOfPartArr
123  */
125 
126 /**\brief Test get part rank
127  *
128  * Tests:
129  * - iMeshP_getRankOfPart
130  * - iMeshP_getRankOfPartArr
131  */
133 
134 /**\brief Test querying of part neighbors
135  *
136  * Test:
137  * - iMeshP_getNumPartNbors
138  * - iMeshP_getNumPartNborsArr
139  * - iMeshP_getPartNbors
140  * - iMeshP_getPartNborsArr
141  */
143 
144 /**\brief Test querying of part boundary entities
145  *
146  * Test:
147  * - iMeshP_getNumPartBdryEnts
148  * - iMeshP_getPartBdryEnts
149  */
151 
152 /**\brief Test querying of part boundary entities
153  *
154  * Test:
155  * - iMeshP_initPartBdryEntIter
156  * - iMeshP_initPartBdryEntArrIter
157  */
159 
160 /**\brief Test adjacent entity query
161  *
162  * Test:
163  * - iMeshP_getAdjEntities
164  */
166 
167 /**\brief Test entity iterators
168  *
169  * Test:
170  * - iMeshP_initEntIter
171  * - iMeshP_initEntArrIter
172  */
174 
175 /**\brief Test entity owner queries
176  *
177  * Test:
178  * - iMeshP_getEntOwnerPart
179  * - iMeshP_getEntOwnerPartArr
180  * - iMeshP_isEntOwner
181  * - iMeshP_isEntOwnerArr
182  */
184 
185 /**\brief Test entity status
186  *
187  * Test:
188  * - iMeshP_getEntStatus
189  * - iMeshP_getEntStatusArr
190  */
192 
193 /**\brief Test information about entity copies for interface entities
194  *
195  * Test:
196  * - iMeshP_getNumCopies
197  * - iMeshP_getCopyParts
198  */
200 
201 /**\brief Test information about entity copies for interface entities
202  *
203  * Test:
204  * - iMeshP_getCopies
205  * - iMeshP_getCopyOnPart
206  * - iMeshP_getOwnerCopy
207  */
209 
210 /**\brief Test creation of ghost entities
211  *
212  * Test:
213  * - iMeshP_createGhostEntsAll
214  */
216 
217 /**\brief Test commuinication of tag data
218  *
219  * Test:
220  * - iMeshP_pushTags
221  * - iMeshP_pushTagsEnt
222  */
225 
226 int test_exchange_ents( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, const PartMap& map );
227 
228 /**************************************************************************
229  Helper Funcions
230  **************************************************************************/
231 
232 class PartMap
233 {
234  public:
235  int num_parts() const
236  {
237  return sortedPartList.size();
238  }
239 
240  iMeshP_Part part_id_from_local_id( int local_id ) const
241  {
242  return sortedPartList[idx_from_local_id( local_id )];
243  }
244 
246  {
247  return partLocalIds[idx_from_part_id( part )];
248  }
249 
250  int rank_from_part_id( iMeshP_Part part ) const
251  {
252  return partRanks[idx_from_part_id( part )];
253  }
254 
255  int rank_from_local_id( int id ) const
256  {
257  return partRanks[idx_from_local_id( id )];
258  }
259 
260  int count_from_rank( int rank ) const
261  {
262  return std::count( partRanks.begin(), partRanks.end(), rank );
263  }
264 
265  void part_id_from_rank( int rank, std::vector< iMeshP_Part >& parts ) const;
266 
267  void local_id_from_rank( int rank, std::vector< int >& ids ) const;
268 
269  const std::vector< iMeshP_Part >& get_parts() const
270  {
271  return sortedPartList;
272  }
273 
274  const std::vector< int >& get_ranks() const
275  {
276  return partRanks;
277  }
278 
279  int build_map( iMesh_Instance imesh, iMeshP_PartitionHandle partition, int num_expected_parts );
280 
281  static int part_from_coords( iMesh_Instance imesh, iMeshP_PartHandle part, int& id_out );
282 
283  private:
284  inline int idx_from_part_id( iMeshP_Part id ) const
285  {
286  return std::lower_bound( sortedPartList.begin(), sortedPartList.end(), id ) - sortedPartList.begin();
287  }
288  inline int idx_from_local_id( int id ) const
289  {
290  return localIdReverseMap[id];
291  }
292 
293  std::vector< iMeshP_Part > sortedPartList;
294  std::vector< int > partRanks;
295  std::vector< int > partLocalIds;
296  std::vector< int > localIdReverseMap;
297 };
298 
299 /**\brief Create mesh for use in parallel tests */
300 int create_mesh( const char* filename, int num_parts );
301 
302 int create_mesh_in_memory( int rank, int size, iMesh_Instance imesh, iMeshP_PartitionHandle& prtn, PartMap& map );
303 
304 /**\brief get unique identifier for each vertex */
305 int vertex_tag( iMesh_Instance imesh, iBase_EntityHandle vertex, int& tag );
306 
309  std::vector< iMeshP_PartHandle >& handles,
310  std::vector< iMeshP_Part >* ids = 0 )
311 {
312  iMeshP_PartHandle* arr = 0;
313  int ierr, alloc = 0, size;
314  iMeshP_getLocalParts( instance, prtn, &arr, &alloc, &size, &ierr );CHKERR;
315  handles.resize( size );
316  std::copy( arr, arr + size, handles.begin() );
317  free( arr );
318  if( !ids ) return iBase_SUCCESS;
319 
320  ids->resize( size );
321  alloc = size;
322  iMeshP_Part* ptr = &( *ids )[0];
323  iMeshP_getPartIdsFromPartHandlesArr( instance, prtn, &handles[0], handles.size(), &ptr, &alloc, &size, &ierr );CHKERR;
324  assert( size == (int)ids->size() );
325  assert( ptr == &( *ids )[0] );
326  return iBase_SUCCESS;
327 }
328 
329 static int get_entities( iMesh_Instance imesh,
331  iBase_EntityType type,
333  std::vector< iBase_EntityHandle >& entities )
334 {
335  iBase_EntityHandle* array = 0;
336  int junk = 0, size = 0, err;
337  iMesh_getEntities( imesh, set, type, topo, &array, &junk, &size, &err );
338  if( !err )
339  {
340  entities.clear();
341  entities.resize( size );
342  std::copy( array, array + size, entities.begin() );
343  free( array );
344  }
345  return err;
346 }
347 
349  iMeshP_PartHandle part,
350  std::vector< iBase_EntityHandle >& elems,
351  std::vector< iBase_EntityHandle >& verts )
352 {
353  int ierr = get_entities( imesh, part, iBase_FACE, iMesh_QUADRILATERAL, elems );CHKERR;
354 
355  verts.resize( 4 * elems.size() );
356  std::vector< int > junk( elems.size() + 1 );
357  int junk1 = verts.size(), count, junk2 = junk.size(), junk3;
358  iBase_EntityHandle* junk4 = &verts[0];
359  int* junk5 = &junk[0];
360  iMesh_getEntArrAdj( imesh, &elems[0], elems.size(), iBase_VERTEX, &junk4, &junk1, &count, &junk5, &junk2, &junk3,
361  &ierr );CHKERR;
362  assert( junk1 == (int)verts.size() );
363  assert( count == (int)( 4 * elems.size() ) );
364  assert( junk2 == (int)junk.size() );
365  assert( junk4 == &verts[0] );
366  assert( junk5 == &junk[0] );
367  std::sort( verts.begin(), verts.end() );
368  verts.erase( std::unique( verts.begin(), verts.end() ), verts.end() );
369  return iBase_SUCCESS;
370 }
371 
372 static int get_coords( iMesh_Instance imesh, const iBase_EntityHandle* verts, int num_verts, double* coords )
373 {
374  double* junk1 = coords;
375  int junk2 = 3 * num_verts;
376  int junk3;
377  int ierr;
378  iMesh_getVtxArrCoords( imesh, verts, num_verts, iBase_INTERLEAVED, &junk1, &junk2, &junk3, &ierr );
379  if( iBase_SUCCESS != ierr ) return ierr;
380  assert( junk1 == coords );
381  assert( junk2 == 3 * num_verts );
382  assert( junk3 == 3 * num_verts );
383  return iBase_SUCCESS;
384 }
385 
386 /**************************************************************************
387  Main Method
388  **************************************************************************/
389 
390 #define RUN_TEST( A ) run_test( &( A ), #A )
391 
392 int run_test( int ( *func )( iMesh_Instance, iMeshP_PartitionHandle, const PartMap& ), const char* func_name )
393 {
394  int rank, size, ierr;
395  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
396  MPI_Comm_size( MPI_COMM_WORLD, &size );
397  iMesh_Instance imesh;
398  iMesh_newMesh( 0, &imesh, &ierr, 0 );
399  PCHECK;
400 
402  iMeshP_createPartitionAll( imesh, MPI_COMM_WORLD, &prtn, &ierr );
403  PCHECK;
404 
405  PartMap map;
406 
407 #ifdef MOAB_HAVE_HDF5
408  if( rank == 0 )
409  {
411  }
412  MPI_Bcast( &ierr, 1, MPI_INT, 0, MPI_COMM_WORLD );
413  if( ierr )
414  {
415  if( rank == 0 )
416  {
417  std::cerr << "Failed to create input test file on root processor. Aborting." << std::endl;
418  }
419  abort();
420  }
421 
422  ierr = test_load( imesh, prtn, map, size );
423  if( ierr )
424  {
425  if( rank == 0 )
426  {
427  std::cerr << "Failed to load input mesh." << std::endl
428  << "Cannot run further tests." << std::endl
429  << "ABORTING" << std::endl;
430  }
431  abort();
432  }
433 #else
434  // so we have MPI and no HDF5; in order to run the test we need to create the
435  // model in memory, and then call sync to resolve shared ents, as if it was read
436  ierr = create_mesh_in_memory( rank, size, imesh, prtn, map );
437  MPI_Bcast( &ierr, 1, MPI_INT, 0, MPI_COMM_WORLD );
438  if( ierr )
439  {
440  if( rank == 0 )
441  {
442  std::cerr << "Failed to create mesh. Aborting." << std::endl;
443  }
444  abort();
445  }
446 
447 #endif
448  int result = ( *func )( imesh, prtn, map );
449  int is_err = is_any_proc_error( result );
450  if( rank == 0 )
451  {
452  if( is_err )
453  std::cout << func_name << " : FAILED!!" << std::endl;
454  else
455  std::cout << func_name << " : success" << std::endl;
456  }
457 
458  iMesh_dtor( imesh, &ierr );CHKERR;
459  return is_err;
460 }
461 
462 int main( int argc, char* argv[] )
463 {
464  MPI_Init( &argc, &argv );
465  int size, rank;
466  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
467  MPI_Comm_size( MPI_COMM_WORLD, &size );
468 
469  if( argc > 2 && !strcmp( argv[1], "-p" ) )
470  {
471 #if !defined( _MSC_VER ) && !defined( __MINGW32__ )
472  std::cout << "Processor " << rank << " of " << size << " with PID " << getpid() << std::endl;
473  std::cout.flush();
474 #endif
475  // loop forever on requested processor, giving the user time
476  // to attach a debugger. Once the debugger in attached, user
477  // can change 'pause'. E.g. on gdb do "set var pause = 0"
478  if( atoi( argv[2] ) == rank )
479  {
480  volatile int pause = 1;
481  while( pause )
482  ;
483  }
484  MPI_Barrier( MPI_COMM_WORLD );
485  }
486 
487  int num_errors = 0;
488  num_errors += RUN_TEST( test_get_partitions );
489  num_errors += RUN_TEST( test_get_parts );
490  num_errors += RUN_TEST( test_get_by_type );
491  num_errors += RUN_TEST( test_get_by_topo );
492  num_errors += RUN_TEST( test_part_id_handle );
493  num_errors += RUN_TEST( test_part_rank );
494  num_errors += RUN_TEST( test_get_neighbors );
495  num_errors += RUN_TEST( test_get_part_boundary );
496  num_errors += RUN_TEST( test_part_boundary_iter );
497  num_errors += RUN_TEST( test_get_adjacencies );
498  num_errors += RUN_TEST( test_entity_iterator );
499  num_errors += RUN_TEST( test_entity_owner );
500  num_errors += RUN_TEST( test_entity_status );
501  num_errors += RUN_TEST( test_entity_copy_parts );
502  num_errors += RUN_TEST( test_entity_copies );
503  num_errors += RUN_TEST( test_push_tag_data_iface );
504  num_errors += RUN_TEST( test_push_tag_data_ghost );
505  num_errors += RUN_TEST( test_create_ghost_ents );
506  num_errors += RUN_TEST( test_exchange_ents );
507 
508  // wait until all procs are done before writing summary data
509  std::cout.flush();
510  MPI_Barrier( MPI_COMM_WORLD );
511 
512 #ifdef MOAB_HAVE_HDF5
513  // clean up output file
514  if( rank == 0 ) remove( FILENAME );
515 #endif
516 
517  if( rank == 0 )
518  {
519  if( !num_errors )
520  std::cout << "All tests passed" << std::endl;
521  else
522  std::cout << num_errors << " TESTS FAILED!" << std::endl;
523  }
524 
525  MPI_Finalize();
526  return num_errors;
527 }
528 
529 // Create a mesh
530 //
531 //
532 // Groups of four quads will be arranged into parts as follows:
533 // +------+------+------+------+------+-----
534 // | | |
535 // | | |
536 // + Part 0 + Part 2 + Part 4
537 // | | |
538 // | | |
539 // +------+------+------+------+------+-----
540 // | | |
541 // | | |
542 // + Part 1 + Part 3 + Part 5
543 // | | |
544 // | | |
545 // +------+------+------+------+------+-----
546 //
547 // Vertices will be enumerated as follows:
548 // 1------6-----11-----16-----21-----26-----
549 // | | |
550 // | | |
551 // 2 7 12 17 22 27
552 // | | |
553 // | | |
554 // 3------8-----13-----18-----23-----28-----
555 // | | |
556 // | | |
557 // 4 9 14 19 24 29
558 // | | |
559 // | | |
560 // 5-----10-----15-----20-----25-----30-----
561 //
562 // Element IDs will be [4*rank+1,4*rank+5]
563 template < int size >
564 struct EHARR
565 {
568  {
569  return h[i];
570  }
571  operator iBase_EntityHandle*()
572  {
573  return h;
574  }
575 };
576 int create_mesh( const char* filename, int num_parts )
577 {
578  const char* tagname = "GLOBAL_ID";
579  int ierr;
580 
581  iMesh_Instance imesh;
582  iMesh_newMesh( 0, &imesh, &ierr, 0 );CHKERR;
583 
584  const int num_full_cols = 2 * ( num_parts / 2 );
585  const int need_half_cols = num_parts % 2;
586  const int num_cols = num_full_cols + 2 * need_half_cols;
587  const int num_vtx = 5 + 5 * num_cols - 4 * ( num_parts % 2 );
588  std::vector< EHARR< 5 > > vertices( num_cols + 1 );
589  std::vector< EHARR< 4 > > elements( num_cols );
590  std::vector< int > vertex_ids( num_vtx );
591  std::vector< iBase_EntityHandle > vertex_list( num_vtx );
592  for( int i = 0; i < num_vtx; ++i )
593  vertex_ids[i] = i + 1;
594 
595  // create vertices
596  int vl_pos = 0;
597  for( int i = 0; i <= num_cols; ++i )
598  {
599  double coords[15] = { static_cast< double >( i ), 0, 0, static_cast< double >( i ), 1, 0,
600  static_cast< double >( i ), 2, 0, static_cast< double >( i ), 3, 0,
601  static_cast< double >( i ), 4, 0 };
602  iBase_EntityHandle* ptr = vertices[i];
603  const int n = ( num_full_cols == num_cols || i <= num_full_cols ) ? 5 : 3;
604  int junk1 = n, junk2 = n;
605  iMesh_createVtxArr( imesh, n, iBase_INTERLEAVED, coords, 3 * n, &ptr, &junk1, &junk2, &ierr );CHKERR;
606  assert( ptr == vertices[i] );
607  assert( junk1 == n );
608  assert( junk2 == n );
609  for( int j = 0; j < n; ++j )
610  vertex_list[vl_pos++] = vertices[i][j];
611  }
612 
613  // create elements
614  for( int i = 0; i < num_cols; ++i )
615  {
616  iBase_EntityHandle conn[16];
617  for( int j = 0; j < 4; ++j )
618  {
619  conn[4 * j] = vertices[i][j];
620  conn[4 * j + 1] = vertices[i][j + 1];
621  conn[4 * j + 2] = vertices[i + 1][j + 1];
622  conn[4 * j + 3] = vertices[i + 1][j];
623  }
624  iBase_EntityHandle* ptr = elements[i];
625  const int n = ( i < num_full_cols ) ? 4 : 2;
626  int junk1 = n, junk2 = n, junk3 = n, junk4 = n;
627  int stat[4];
628  int* ptr2 = stat;
629  iMesh_createEntArr( imesh, iMesh_QUADRILATERAL, conn, 4 * n, &ptr, &junk1, &junk2, &ptr2, &junk3, &junk4,
630  &ierr );CHKERR;
631  assert( ptr == elements[i] );
632  assert( junk1 == n );
633  assert( junk2 == n );
634  assert( ptr2 == stat );
635  assert( junk3 == n );
636  assert( junk4 == n );
637  }
638 
639  // create partition
640  iMeshP_PartitionHandle partition;
641  iMeshP_createPartitionAll( imesh, MPI_COMM_SELF, &partition, &ierr );CHKERR;
642  for( int i = 0; i < num_parts; ++i )
643  {
644  iMeshP_PartHandle part;
645  iMeshP_createPart( imesh, partition, &part, &ierr );CHKERR;
646  iBase_EntityHandle quads[] = { elements[2 * ( i / 2 )][2 * ( i % 2 )],
647  elements[2 * ( i / 2 ) + 1][2 * ( i % 2 )],
648  elements[2 * ( i / 2 )][2 * ( i % 2 ) + 1],
649  elements[2 * ( i / 2 ) + 1][2 * ( i % 2 ) + 1] };
650  iMesh_addEntArrToSet( imesh, quads, 4, part, &ierr );CHKERR;
651  }
652 
653  // assign global ids to vertices
655  iMesh_getTagHandle( imesh, tagname, &id_tag, &ierr, strlen( tagname ) );
656  if( iBase_SUCCESS == ierr )
657  {
658  int tag_size, tag_type;
659  iMesh_getTagSizeValues( imesh, id_tag, &tag_size, &ierr );CHKERR;
660  if( tag_size != 1 ) return iBase_TAG_ALREADY_EXISTS;
661  iMesh_getTagType( imesh, id_tag, &tag_type, &ierr );CHKERR;
662  if( tag_type != iBase_INTEGER ) return iBase_TAG_ALREADY_EXISTS;
663  }
664  else
665  {
666  iMesh_createTag( imesh, tagname, 1, iBase_INTEGER, &id_tag, &ierr, strlen( tagname ) );CHKERR;
667  }
668  iMesh_setIntArrData( imesh, &vertex_list[0], num_vtx, id_tag, &vertex_ids[0], num_vtx, &ierr );CHKERR;
669 
670  // write file
672  iMesh_getRootSet( imesh, &root_set, &ierr );
673  iMeshP_saveAll( imesh, partition, root_set, filename, 0, &ierr, strlen( filename ), 0 );CHKERR;
674 
675  iMesh_dtor( imesh, &ierr );CHKERR;
676 
677  return 0;
678 }
680  int num_parts,
681  iMesh_Instance imesh,
682  iMeshP_PartitionHandle& partition,
683  PartMap& map )
684 {
685  const char* tagname = "GLOBAL_ID";
686  int ierr;
687 
688  const int num_cols = 2;
689  const int num_vtx = 9;
690  // we are on the top or botton row
691  int bottom = rank % 2; // 0 2
692  // 1 3
693  std::vector< EHARR< 3 > > vertices( 3 );
694  std::vector< EHARR< 2 > > elements( 2 ); // 4 elements per process
695  std::vector< int > vertex_ids( num_vtx );
696  std::vector< iBase_EntityHandle > vertex_list( num_vtx );
697  int start = 1 + 2 * bottom + 10 * ( rank / 2 );
698 
699  for( int i = 0; i < 3; ++i )
700  for( int j = 0; j < 3; ++j )
701  {
702  vertex_ids[i + 3 * j] = start + i + 5 * j;
703  }
704 
705  // create vertices
706  int vl_pos = 0;
707  int startI = 2 * ( rank / 2 ); // so it will be 0, 0, 2, 2, ...)
708  for( int i = 0; i <= 2; ++i )
709  {
710  double coords[9] = {
711  static_cast< double >( i + startI ), 2. * bottom, 0,
712  static_cast< double >( i + startI ), 1 + 2. * bottom, 0,
713  static_cast< double >( i + startI ), 2 + 2. * bottom, 0,
714  };
715  iBase_EntityHandle* ptr = vertices[i];
716  const int n = 3;
717  int junk1 = n, junk2 = n;
718  iMesh_createVtxArr( imesh, n, iBase_INTERLEAVED, coords, 3 * n, &ptr, &junk1, &junk2, &ierr );CHKERR;
719  assert( ptr == vertices[i] );
720  assert( junk1 == n );
721  assert( junk2 == n );
722  for( int j = 0; j < n; ++j )
723  vertex_list[vl_pos++] = vertices[i][j];
724  }
725 
726  // create elements
727  for( int i = 0; i < num_cols; ++i )
728  {
729  iBase_EntityHandle conn[8];
730  for( int j = 0; j < 2; ++j )
731  {
732  conn[4 * j] = vertices[i][j];
733  conn[4 * j + 1] = vertices[i][j + 1];
734  conn[4 * j + 2] = vertices[i + 1][j + 1];
735  conn[4 * j + 3] = vertices[i + 1][j];
736  }
737  iBase_EntityHandle* ptr = elements[i];
738  const int n = 2;
739  int junk1 = n, junk2 = n, junk3 = n, junk4 = n;
740  int stat[4];
741  int* ptr2 = stat;
742  iMesh_createEntArr( imesh, iMesh_QUADRILATERAL, conn, 4 * n, &ptr, &junk1, &junk2, &ptr2, &junk3, &junk4,
743  &ierr );CHKERR;
744  assert( ptr == elements[i] );
745  assert( junk1 == n );
746  assert( junk2 == n );
747  assert( ptr2 == stat );
748  assert( junk3 == n );
749  assert( junk4 == n );
750  }
751 
752  // create partition
753  iMeshP_createPartitionAll( imesh, MPI_COMM_WORLD, &partition, &ierr );CHKERR;
754 
755  iMeshP_PartHandle part;
756  iMeshP_createPart( imesh, partition, &part, &ierr );CHKERR;
757  iBase_EntityHandle quads[] = { elements[0][0], elements[0][1], elements[1][0], elements[1][1] };
758  iMesh_addEntArrToSet( imesh, quads, 4, part, &ierr );CHKERR;
759 
760  // assign global ids to vertices
762  iMesh_getTagHandle( imesh, tagname, &id_tag, &ierr, strlen( tagname ) );
763  if( iBase_SUCCESS == ierr )
764  {
765  int tag_size, tag_type;
766  iMesh_getTagSizeValues( imesh, id_tag, &tag_size, &ierr );CHKERR;
767  if( tag_size != 1 ) return iBase_TAG_ALREADY_EXISTS;
768  iMesh_getTagType( imesh, id_tag, &tag_type, &ierr );CHKERR;
769  if( tag_type != iBase_INTEGER ) return iBase_TAG_ALREADY_EXISTS;
770  }
771  else
772  {
773  iMesh_createTag( imesh, tagname, 1, iBase_INTEGER, &id_tag, &ierr, strlen( tagname ) );CHKERR;
774  }
775  iMesh_setIntArrData( imesh, &vertex_list[0], num_vtx, id_tag, &vertex_ids[0], num_vtx, &ierr );CHKERR;
776 
777  // some mesh sync
778  iMeshP_syncPartitionAll( imesh, partition, &ierr );CHKERR;
779  iMeshP_syncMeshAll( imesh, partition, &ierr );CHKERR;
780 
781  ierr = map.build_map( imesh, partition, num_parts );CHKERR;
782  return 0;
783 }
784 // generate unique for each vertex from coordinates.
785 // Assume integer coordinate values with x in [0,inf] and y in [0,4]
786 // as generated by create_mean(..).
788 {
789  int ierr;
790  double x, y, z;
791  iMesh_getVtxCoord( imesh, vertex, &x, &y, &z, &ierr );CHKERR;
792 
793  int xc = (int)round( x );
794  int yc = (int)round( y );
795  tag = 5 * xc + yc + 1;
796  return ierr;
797 }
798 
799 /**************************************************************************
800  Test Implementations
801  **************************************************************************/
802 
803 int test_load( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, PartMap& map, int proc_size )
804 {
805  int ierr;
806 
808  iMesh_getRootSet( imesh, &root_set, &ierr );
809  const char* opt = "moab:PARTITION=PARALLEL_PARTITION";
810  iMeshP_loadAll( imesh, prtn, root_set, FILENAME, opt, &ierr, strlen( FILENAME ), strlen( opt ) );
811  PCHECK;
812 
813  ierr = map.build_map( imesh, prtn, proc_size );CHKERR;
814  return iBase_SUCCESS;
815 }
816 
817 /**\brief Test partition query methods
818  *
819  * Test:
820  * - iMeshP_getPartitionComm
821  * - iMeshP_getNumPartitions
822  * - iMeshP_getPartitions
823  */
825 {
826  int ierr;
827 
828  // test iMeshP_getPartitionCom
829  MPI_Comm comm = MPI_COMM_SELF;
830  iMeshP_getPartitionComm( imesh, prtn, &comm, &ierr );
831  PCHECK;
832  ASSERT( comm == MPI_COMM_WORLD );
833 
834  // test iMeshP_getPartitions
835  iMeshP_PartitionHandle* array = 0;
836  int alloc = 0, size = -1;
837  iMeshP_getPartitions( imesh, &array, &alloc, &size, &ierr );
838  PCHECK;
839  ASSERT( array != 0 );
840  ASSERT( alloc == size );
841  ASSERT( size > 0 );
842  int idx = std::find( array, array + size, prtn ) - array;
843  free( array );
844  ASSERT( idx < size );
845 
846  // test iMesP_getNumPartitions
847  int size2 = -1;
848  iMeshP_getNumPartitions( imesh, &size2, &ierr );
849  PCHECK;
850  ASSERT( size2 == size );
851  return 0;
852 }
853 
854 /**\brief Test part quyery methods
855  *
856  * Test:
857  * - iMeshP_getNumGlobalParts
858  * - iMeshP_getNumLocalParts
859  * - iMeshP_getLocalParts
860  */
862 {
863  int size, rank, ierr;
864  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
865  MPI_Comm_size( MPI_COMM_WORLD, &size );
866 
867  int num_part_g;
868  iMeshP_getNumGlobalParts( imesh, prtn, &num_part_g, &ierr );
869  PCHECK;
870  ASSERT( num_part_g == map.num_parts() );
871 
872  int num_part_l;
873  iMeshP_getNumLocalParts( imesh, prtn, &num_part_l, &ierr );
874  PCHECK;
875  ASSERT( num_part_l == map.count_from_rank( rank ) );
876 
877  std::vector< iMeshP_PartHandle > parts( num_part_l );
878  iMeshP_PartHandle* ptr = &parts[0];
879  int junk1 = num_part_l, count = -1;
880  iMeshP_getLocalParts( imesh, prtn, &ptr, &junk1, &count, &ierr );
881  PCHECK;
882  assert( ptr == &parts[0] );
883  assert( junk1 == num_part_l );
884  ASSERT( count == num_part_l );
885 
886  return iBase_SUCCESS;
887 }
888 
889 static int test_get_by_type_topo_all( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, bool test_type, int num_parts )
890 {
891  // calculate number of quads and vertices in entire mesh
892  // from number of parts (see create_mesh(..) function.)
893  const int expected_global_quad_count = 4 * num_parts;
894  const int num_col = 2 * ( num_parts / 2 + num_parts % 2 );
895  const int expected_global_vtx_count = num_parts == 1 ? 9 : num_parts % 2 ? 1 + 5 * num_col : 5 + 5 * num_col;
896 
897  // test getNumOf*All for root set
898  int ierr, count;
900  iMesh_getRootSet( imesh, &root, &ierr );
901  if( test_type )
902  iMeshP_getNumOfTypeAll( imesh, prtn, root, iBase_VERTEX, &count, &ierr );
903  else
904  iMeshP_getNumOfTopoAll( imesh, prtn, root, iMesh_POINT, &count, &ierr );
905  PCHECK;
906  ASSERT( count == expected_global_vtx_count );
907  if( test_type )
908  iMeshP_getNumOfTypeAll( imesh, prtn, root, iBase_FACE, &count, &ierr );
909  else
910  iMeshP_getNumOfTopoAll( imesh, prtn, root, iMesh_QUADRILATERAL, &count, &ierr );
911  PCHECK;
912  ASSERT( count == expected_global_quad_count );
913 
914  // create an entity set containing half of the quads
915  std::vector< iBase_EntityHandle > all_quads, half_quads;
916  ierr = get_entities( imesh, root, iBase_FACE, iMesh_QUADRILATERAL, all_quads );
917  assert( 0 == all_quads.size() % 2 );
918  half_quads.resize( all_quads.size() / 2 );
919  for( size_t i = 0; i < all_quads.size() / 2; ++i )
920  half_quads[i] = all_quads[2 * i];
922  iMesh_createEntSet( imesh, 1, &set, &ierr );CHKERR;
923  iMesh_addEntArrToSet( imesh, &half_quads[0], half_quads.size(), set, &ierr );CHKERR;
924 
925  // test getNumOf*All with defined set
926  if( test_type )
927  iMeshP_getNumOfTypeAll( imesh, prtn, set, iBase_VERTEX, &count, &ierr );
928  else
929  iMeshP_getNumOfTopoAll( imesh, prtn, set, iMesh_POINT, &count, &ierr );
930  PCHECK;
931  ASSERT( count == 0 );
932  if( test_type )
933  iMeshP_getNumOfTypeAll( imesh, prtn, set, iBase_FACE, &count, &ierr );
934  else
935  iMeshP_getNumOfTopoAll( imesh, prtn, set, iMesh_QUADRILATERAL, &count, &ierr );
936  PCHECK;
937  ASSERT( count == expected_global_quad_count / 2 );
938 
939  return 0;
940 }
941 
942 static int test_get_by_type_topo_local( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, bool test_type )
943 {
944  int ierr;
946  iMesh_getRootSet( imesh, &root, &ierr );
947 
948  // select a single part
949  std::vector< iMeshP_PartHandle > parts;
950  ierr = get_local_parts( imesh, prtn, parts );CHKERR;
951  iMeshP_PartHandle part = parts.front();
952 
953  // get the entities contained in the part
954  std::vector< iBase_EntityHandle > part_quads, part_all;
955  ierr = get_entities( imesh, part, iBase_FACE, iMesh_QUADRILATERAL, part_quads );CHKERR;
956  ierr = get_entities( imesh, part, iBase_ALL_TYPES, iMesh_ALL_TOPOLOGIES, part_all );CHKERR;
957 
958  // compare local counts (using root set)
959 
960  int count;
961  if( test_type )
962  iMeshP_getNumOfType( imesh, prtn, part, root, iBase_FACE, &count, &ierr );
963  else
964  iMeshP_getNumOfTopo( imesh, prtn, part, root, iMesh_QUADRILATERAL, &count, &ierr );CHKERR;
965  ASSERT( count == (int)part_quads.size() );
966 
967  if( test_type )
968  iMeshP_getNumOfType( imesh, prtn, part, root, iBase_ALL_TYPES, &count, &ierr );
969  else
970  iMeshP_getNumOfTopo( imesh, prtn, part, root, iMesh_ALL_TOPOLOGIES, &count, &ierr );CHKERR;
971  ASSERT( count == (int)part_all.size() );
972 
973  // compare local contents (using root set)
974 
975  iBase_EntityHandle* ptr = 0;
976  int num_ent, junk1 = 0;
977  iMeshP_getEntities( imesh, prtn, part, root, test_type ? iBase_FACE : iBase_ALL_TYPES,
978  test_type ? iMesh_ALL_TOPOLOGIES : iMesh_QUADRILATERAL, &ptr, &junk1, &num_ent, &ierr );CHKERR;
979  std::vector< iBase_EntityHandle > act_quads( ptr, ptr + num_ent );
980  free( ptr );
981  junk1 = num_ent = 0;
982  ptr = 0;
983  iMeshP_getEntities( imesh, prtn, part, root, iBase_ALL_TYPES, iMesh_ALL_TOPOLOGIES, &ptr, &junk1, &num_ent, &ierr );CHKERR;
984  std::vector< iBase_EntityHandle > act_all( ptr, ptr + num_ent );
985  free( ptr );
986  std::sort( part_quads.begin(), part_quads.end() );
987  std::sort( part_all.begin(), part_all.end() );
988  std::sort( act_quads.begin(), act_quads.end() );
989  std::sort( act_all.begin(), act_all.end() );
990  ASSERT( part_quads == act_quads );
991  ASSERT( part_all == act_all );
992 
993  // create an entity set containing half of the quads from the part
994  std::vector< iBase_EntityHandle > half_quads( part_quads.size() / 2 );
995  for( size_t i = 0; i < half_quads.size(); ++i )
996  half_quads[i] = part_quads[2 * i];
998  iMesh_createEntSet( imesh, 1, &set, &ierr );CHKERR;
999  iMesh_addEntArrToSet( imesh, &half_quads[0], half_quads.size(), set, &ierr );CHKERR;
1000 
1001  // check if there exists any quads not in the part that we
1002  // can add to the set
1003  std::vector< iBase_EntityHandle > all_quads, other_quads;
1004  ierr = get_entities( imesh, root, iBase_FACE, iMesh_QUADRILATERAL, all_quads );CHKERR;
1005  std::sort( all_quads.begin(), all_quads.end() );
1006  std::sort( part_quads.begin(), part_quads.end() );
1007  std::set_difference( all_quads.begin(), all_quads.end(), part_quads.begin(), part_quads.end(),
1008  std::back_inserter( other_quads ) );
1009  iMesh_addEntArrToSet( imesh, &other_quads[0], other_quads.size(), set, &ierr );CHKERR;
1010 
1011  // compare local counts (using non-root set)
1012 
1013  if( test_type )
1014  iMeshP_getNumOfType( imesh, prtn, part, set, iBase_FACE, &count, &ierr );
1015  else
1016  iMeshP_getNumOfTopo( imesh, prtn, part, set, iMesh_QUADRILATERAL, &count, &ierr );CHKERR;
1017  ASSERT( count == (int)half_quads.size() );
1018 
1019  if( test_type )
1020  iMeshP_getNumOfType( imesh, prtn, part, set, iBase_VERTEX, &count, &ierr );
1021  else
1022  iMeshP_getNumOfTopo( imesh, prtn, part, set, iMesh_POINT, &count, &ierr );CHKERR;
1023  ASSERT( count == 0 );
1024 
1025  // compare local contents (using non-root set)
1026 
1027  junk1 = 0;
1028  num_ent = 0;
1029  ptr = 0;
1030  iMeshP_getEntities( imesh, prtn, part, set, test_type ? iBase_FACE : iBase_ALL_TYPES,
1031  test_type ? iMesh_ALL_TOPOLOGIES : iMesh_QUADRILATERAL, &ptr, &junk1, &num_ent, &ierr );CHKERR;
1032  act_quads.resize( num_ent );
1033  std::copy( ptr, ptr + num_ent, act_quads.begin() );
1034  free( ptr );
1035  std::sort( half_quads.begin(), half_quads.end() );
1036  std::sort( act_quads.begin(), act_quads.end() );
1037  ASSERT( act_quads == half_quads );
1038 
1039  return iBase_SUCCESS;
1040 }
1041 
1042 /**\brief Test query by entity type
1043  *
1044  * Test:
1045  * - iMeshP_getNumOfTypeAll
1046  * - iMeshP_getNumOfType
1047  * - iMeshP_getEntities
1048  * -
1049  */
1051 {
1052  int ierr;
1053  ierr = test_get_by_type_topo_all( imesh, prtn, true, map.num_parts() );
1054  PCHECK;
1055  ierr = test_get_by_type_topo_local( imesh, prtn, true );
1056  PCHECK;
1057  return 0;
1058 }
1059 
1060 /**\brief Test query by entity topology
1061  *
1062  * Test:
1063  * - iMeshP_getNumOfTopoAll
1064  * - iMeshP_getNumOfTopo
1065  * - iMeshP_getEntities
1066  * -
1067  */
1069 {
1070  int ierr;
1071  ierr = test_get_by_type_topo_all( imesh, prtn, false, map.num_parts() );
1072  PCHECK;
1073  ierr = test_get_by_type_topo_local( imesh, prtn, false );
1074  PCHECK;
1075  return 0;
1076 }
1077 
1078 /**\brief Test mapping from part id to part handle
1079  *
1080  * Test:
1081  * - iMeshP_getPartIdFromPartHandle
1082  * - iMeshP_getPartIdsFromPartHandlesArr
1083  * - iMeshP_getPartHandleFromPartId
1084  * - iMeshP_getPartHandlesFromPartsIdsArr
1085  */
1087 {
1088  // get local part ids
1089  int rank, ierr;
1090  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
1091  std::vector< iMeshP_Part > ids;
1092  map.part_id_from_rank( rank, ids );
1093 
1094  // check single-part functions and build list of part handles
1095  std::vector< iMeshP_PartHandle > handles( ids.size() );
1096  size_t i;
1097  for( i = 0; i < ids.size(); ++i )
1098  {
1099  iMeshP_getPartHandleFromPartId( imesh, prtn, ids[i], &handles[i], &ierr );CHKERR;
1100  iMeshP_Part id;
1101  iMeshP_getPartIdFromPartHandle( imesh, prtn, handles[i], &id, &ierr );CHKERR;
1102  if( id != ids[i] ) break;
1103  }
1104  ASSERT( i == ids.size() );
1105 
1106  // test iMeshP_getPartIdsFromPartHandlesArr
1107  std::vector< iMeshP_Part > ids2( ids.size() );
1108  int junk1 = ids.size(), junk2 = 0;
1109  iMeshP_Part* ptr = &ids2[0];
1110  iMeshP_getPartIdsFromPartHandlesArr( imesh, prtn, &handles[0], handles.size(), &ptr, &junk1, &junk2, &ierr );
1111  PCHECK;
1112  ASSERT( ptr == &ids2[0] );
1113  ASSERT( junk2 == (int)ids2.size() );
1114  ASSERT( ids == ids2 );
1115 
1116  // test iMeshP_getPartHandlesFromPartsIdsArr
1117  std::vector< iMeshP_PartHandle > handles2( handles.size() );
1118  junk1 = handles.size();
1119  junk2 = 0;
1120  iMeshP_PartHandle* ptr2 = &handles2[0];
1121  iMeshP_getPartHandlesFromPartsIdsArr( imesh, prtn, &ids[0], ids.size(), &ptr2, &junk1, &junk2, &ierr );
1122  PCHECK;
1123  ASSERT( ptr2 == &handles2[0] );
1124  ASSERT( junk2 == (int)handles2.size() );
1125  ASSERT( handles == handles2 );
1126 
1127  return 0;
1128 }
1129 
1130 /**\brief Test get part rank
1131  *
1132  * Tests:
1133  * - iMeshP_getRankOfPart
1134  * - iMeshP_getRankOfPartArr
1135  */
1137 {
1138  int ierr = 0, rank;
1139  std::vector< iMeshP_Part > invalid, failed;
1140  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
1141 
1142  // test iMeshP_getRankOfPart
1143  for( size_t i = 0; i < map.get_parts().size(); ++i )
1144  {
1145  int pr;
1146  iMeshP_getRankOfPart( imesh, prtn, map.get_parts()[i], &pr, &ierr );
1147  if( iBase_SUCCESS != ierr )
1148  failed.push_back( map.get_parts()[i] );
1149  else if( pr != map.get_ranks()[i] )
1150  invalid.push_back( map.get_parts()[i] );
1151  }
1152  if( !failed.empty() )
1153  {
1154  std::cerr << "Processor " << rank << ": iMeshP_getRankOfPart failed for " << failed.size() << " parts."
1155  << std::endl;
1156  ierr = iBase_FAILURE;
1157  }
1158  if( !invalid.empty() )
1159  {
1160  std::cerr << "Processor " << rank << ": iMeshP_getRankOfPart was incorrect for " << invalid.size() << " parts."
1161  << std::endl;
1162  ierr = iBase_FAILURE;
1163  }
1164  PCHECK;
1165 
1166  // test iMeshP_getRankOfPartArr
1167  std::vector< int > ranks( map.get_parts().size() );
1168  int junk1 = ranks.size(), junk2, *ptr = &ranks[0];
1169  iMeshP_getRankOfPartArr( imesh, prtn, &map.get_parts()[0], map.get_parts().size(), &ptr, &junk1, &junk2, &ierr );
1170  PCHECK;
1171  assert( ptr == &ranks[0] );
1172  assert( junk1 == (int)ranks.size() );
1173  ASSERT( junk2 == (int)ranks.size() );
1174  for( size_t i = 0; i < map.get_parts().size(); ++i )
1175  {
1176  if( ranks[i] != map.get_ranks()[i] ) invalid.push_back( map.get_parts()[i] );
1177  }
1178  if( !invalid.empty() )
1179  {
1180  std::cerr << "Processor " << rank << ": iMeshP_getRankOfPartArr was incorrect for " << invalid.size()
1181  << " parts." << std::endl;
1182  ierr = iBase_FAILURE;
1183  }
1184  PCHECK;
1185 
1186  return 0;
1187 }
1188 
1189 // see create_mesh(..)
1190 static void get_part_neighbors( int logical_part_id, int num_parts, int neighbors[5], int& num_neighbors )
1191 {
1192  num_neighbors = 0;
1193  if( logical_part_id + 1 < num_parts ) neighbors[num_neighbors++] = logical_part_id + 1;
1194  if( logical_part_id + 2 < num_parts ) neighbors[num_neighbors++] = logical_part_id + 2;
1195  if( logical_part_id % 2 )
1196  {
1197  neighbors[num_neighbors++] = logical_part_id - 1;
1198  if( logical_part_id > 2 )
1199  {
1200  neighbors[num_neighbors++] = logical_part_id - 3;
1201  neighbors[num_neighbors++] = logical_part_id - 2;
1202  }
1203  }
1204  else
1205  {
1206  if( logical_part_id + 3 < num_parts ) neighbors[num_neighbors++] = logical_part_id + 3;
1207  if( logical_part_id > 1 )
1208  {
1209  neighbors[num_neighbors++] = logical_part_id - 1;
1210  neighbors[num_neighbors++] = logical_part_id - 2;
1211  }
1212  }
1213 }
1214 
1215 /**\brief Test querying of part neighbors
1216  *
1217  * Test:
1218  * - iMeshP_getNumPartNbors
1219  * - iMeshP_getNumPartNborsArr
1220  * - iMeshP_getPartNbors
1221  * - iMeshP_getPartNborsArr
1222  */
1224 {
1225  int ierr, rank;
1226  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
1227 
1228  std::vector< iMeshP_Part > local_parts;
1229  map.part_id_from_rank( rank, local_parts );
1230 
1231  // get handles for local parts
1232  std::vector< iMeshP_PartHandle > handles( local_parts.size() );
1233  iMeshP_PartHandle* ptr = &handles[0];
1234  int junk1 = handles.size(), junk2 = 0;
1235  iMeshP_getPartHandlesFromPartsIdsArr( imesh, prtn, &local_parts[0], local_parts.size(), &ptr, &junk1, &junk2,
1236  &ierr );
1237  PCHECK;
1238  assert( ptr == &handles[0] );
1239  assert( junk2 == (int)handles.size() );
1240 
1241  // get logical ids for local parts
1242  std::vector< int > logical_ids;
1243  map.local_id_from_rank( rank, logical_ids );
1244 
1245  // get neighbors for each local part
1246  std::vector< std::vector< iMeshP_Part > > neighbors( logical_ids.size() );
1247  for( size_t i = 0; i < logical_ids.size(); ++i )
1248  {
1249  int logical_neighbors[5], num_neighbors;
1250  get_part_neighbors( logical_ids[i], map.num_parts(), logical_neighbors, num_neighbors );
1251  neighbors[i].resize( num_neighbors );
1252  for( int j = 0; j < num_neighbors; ++j )
1253  neighbors[i][j] = map.part_id_from_local_id( logical_neighbors[j] );
1254  std::sort( neighbors[i].begin(), neighbors[i].end() );
1255  }
1256 
1257  // test iMeshP_getNumPartNbors
1258  std::vector< iMeshP_Part > invalid, failed;
1259  for( size_t i = 0; i < local_parts.size(); ++i )
1260  {
1261  int count;
1262  iMeshP_getNumPartNbors( imesh, prtn, handles[i], iBase_VERTEX, &count, &ierr );
1263  if( ierr )
1264  failed.push_back( local_parts[i] );
1265  else if( count != (int)neighbors[i].size() )
1266  invalid.push_back( local_parts[i] );
1267  }
1268  if( !failed.empty() )
1269  {
1270  std::cerr << "Processor " << rank << ": iMeshP_getNumPartNbors failed for " << failed.size() << " parts."
1271  << std::endl;
1272  ierr = iBase_FAILURE;
1273  PCHECK;
1274  }
1275  if( !invalid.empty() )
1276  {
1277  std::cerr << "Processor " << rank << ": iMeshP_getNumPartNbors was incorrect for " << invalid.size()
1278  << " parts." << std::endl;
1279  ierr = iBase_FAILURE;
1280  PCHECK;
1281  }
1282 
1283  // test iMeshP_getPartNbors
1284  ierr = 0;
1285  for( size_t i = 0; i < local_parts.size(); ++i )
1286  {
1287  int count, junk = 0, another_count;
1288  iMeshP_Part* list = 0;
1289  iMeshP_getPartNbors( imesh, prtn, handles[i], iBase_VERTEX, &another_count, &list, &junk, &count, &ierr );
1290  assert( count == another_count );
1291  if( ierr )
1292  failed.push_back( local_parts[i] );
1293  else
1294  {
1295  std::sort( list, list + count );
1296  std::vector< iMeshP_Part > cpy( list, list + count );
1297  if( cpy != neighbors[i] ) invalid.push_back( local_parts[i] );
1298  free( list );
1299  }
1300  }
1301  if( !failed.empty() )
1302  {
1303  std::cerr << "Processor " << rank << ": iMeshP_getPartNbors failed for " << failed.size() << " parts."
1304  << std::endl;
1305  ierr = iBase_FAILURE;
1306  }
1307  if( !invalid.empty() )
1308  {
1309  std::cerr << "Processor " << rank << ": iMeshP_getPartNbors was incorrect for " << invalid.size() << " parts."
1310  << std::endl;
1311  ierr = iBase_FAILURE;
1312  }
1313  PCHECK;
1314 
1315  // test iMeshP_getNumPartNborsArr
1316  std::vector< int > count_vect( handles.size() );
1317  int* count_arr = &count_vect[0];
1318  junk1 = handles.size();
1319  iMeshP_getNumPartNborsArr( imesh, prtn, &handles[0], handles.size(), iBase_VERTEX, &count_arr, &junk1, &junk2,
1320  &ierr );
1321  PCHECK;
1322  assert( count_arr == &count_vect[0] );
1323  assert( junk2 == (int)handles.size() );
1324  for( size_t i = 0; i < local_parts.size(); ++i )
1325  {
1326  if( count_arr[i] != (int)neighbors[i].size() ) invalid.push_back( local_parts[i] );
1327  }
1328  if( !invalid.empty() )
1329  {
1330  std::cerr << "Processor " << rank << ": iMeshP_getNumPartNborsArr was incorrect for " << invalid.size()
1331  << " parts." << std::endl;
1332  ierr = iBase_FAILURE;
1333  }
1334  PCHECK;
1335 
1336  // test iMeshP_getPartNborsArr
1337  iMeshP_Part* nbor_arr = 0;
1338  junk1 = handles.size(), junk2 = 0;
1339  int junk3 = 0, nbor_size;
1340  iMeshP_getPartNborsArr( imesh, prtn, &handles[0], handles.size(), iBase_VERTEX, &count_arr, &junk1, &junk2,
1341  &nbor_arr, &junk3, &nbor_size, &ierr );
1342  PCHECK;
1343  assert( count_arr == &count_vect[0] );
1344  assert( junk2 == (int)handles.size() );
1345  std::vector< iMeshP_Part > all_nbors( nbor_arr, nbor_arr + nbor_size );
1346  free( nbor_arr );
1347  std::vector< iMeshP_Part >::iterator j = all_nbors.begin();
1348  bool bad_length = false;
1349  for( size_t i = 0; i < local_parts.size(); ++i )
1350  {
1351  if( all_nbors.end() - j > count_arr[i] )
1352  {
1353  bad_length = true;
1354  break;
1355  }
1356  if( count_arr[i] != (int)neighbors[i].size() )
1357  {
1358  invalid.push_back( local_parts[i] );
1359  }
1360  else
1361  {
1362  std::vector< iMeshP_Part >::iterator e = j + count_arr[i];
1363  std::sort( j, e );
1364  if( !std::equal( j, e, neighbors[i].begin() ) ) invalid.push_back( local_parts[i] );
1365  }
1366  }
1367  if( bad_length )
1368  {
1369  std::cerr << "Processor " << rank << ": iMeshP_getPartNborsArr had inconsistent result array lengths."
1370  << std::endl;
1371  ierr = iBase_FAILURE;
1372  }
1373  if( !invalid.empty() )
1374  {
1375  std::cerr << "Processor " << rank << ": iMeshP_getPartNborsArr was incorrect for " << invalid.size()
1376  << " parts." << std::endl;
1377  ierr = iBase_FAILURE;
1378  }
1379  PCHECK;
1380 
1381  return 0;
1382 }
1383 
1384 // Determine the expected vertices on the interface between two parts.
1385 // Returns no vertices for non-adjacient parts and fails if both parts
1386 // are the same.
1387 // See create_mesh(..) for the assumed mesh.
1390  iMeshP_PartHandle local_part,
1391  iMeshP_Part other_part,
1392  const PartMap& map,
1393  std::vector< iBase_EntityHandle >& vtx_handles )
1394 {
1395  int ierr, rank;
1396  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
1397 
1398  iMeshP_Part local_id;
1399  iMeshP_getPartIdFromPartHandle( imesh, prtn, local_part, &local_id, &ierr );CHKERR;
1400 
1401  const int local_logical = map.local_id_from_part_id( local_id );
1402  const int other_logical = map.local_id_from_part_id( other_part );
1403 
1404  // get grid of local vertices
1405 
1406  iBase_EntityHandle verts[3][3];
1407  const double xbase = ( local_id / 2 ) * 2;
1408  const double ybase = ( local_id % 2 ) * 2;
1409 
1410  // get quads in partition
1411  iBase_EntityHandle quads[4], *ptr = quads;
1412  int junk1 = 4, junk2;
1413  iMesh_getEntities( imesh, local_part, iBase_FACE, iMesh_QUADRILATERAL, &ptr, &junk1, &junk2, &ierr );CHKERR;
1414  assert( ptr == quads );
1415  assert( junk1 == 4 );
1416  assert( junk2 == 4 );
1417 
1418  // get vertices in quads
1419  iBase_EntityHandle conn[16];
1420  int offsets[5], *off_ptr = offsets, junk3 = 5, junk4;
1421  ptr = conn;
1422  junk1 = 16;
1423  iMesh_getEntArrAdj( imesh, quads, 4, iBase_VERTEX, &ptr, &junk1, &junk2, &off_ptr, &junk3, &junk4, &ierr );CHKERR;
1424  assert( ptr == conn );
1425  assert( junk1 == 16 );
1426  assert( junk2 == 16 );
1427  assert( off_ptr == offsets );
1428  assert( junk3 == 5 );
1429  assert( junk4 == 5 );
1430 
1431  // make unique vertex list
1432  std::sort( conn, conn + 16 );
1433  const int num_vtx = std::unique( conn, conn + 16 ) - conn;
1434  assert( 9 == num_vtx );
1435 
1436  // get vertex coords
1437  std::vector< double > coords( 27 );
1438  ierr = get_coords( imesh, conn, 9, &coords[0] );CHKERR;
1439 
1440  // use vertex coords to determine logical position
1441  for( int i = 0; i < num_vtx; ++i )
1442  {
1443  int x = (int)round( coords[3 * i] - xbase );
1444  int y = (int)round( coords[3 * i + 1] - ybase );
1445  if( x < 0 || x > 2 || y < 0 || y > 2 )
1446  {
1447  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1448  << " Invalid vertex coordinate: (" << coords[3 * i] << ", " << coords[3 * i + 1] << ", "
1449  << coords[3 * i + 2] << ")" << std::endl
1450  << " For logical partition " << local_id << std::endl;
1451  return iBase_FAILURE;
1452  }
1453  verts[x][y] = conn[i];
1454  }
1455 
1456  if( local_logical % 2 )
1457  {
1458  switch( other_logical - local_logical )
1459  {
1460  case 0:
1461  return iBase_FAILURE;
1462  case 1: // upper right
1463  vtx_handles.resize( 1 );
1464  vtx_handles[0] = verts[2][0];
1465  break;
1466  case 2: // right
1467  vtx_handles.resize( 3 );
1468  std::copy( verts[2], verts[2] + 3, vtx_handles.begin() );
1469  break;
1470  case -1: // above
1471  vtx_handles.resize( 3 );
1472  vtx_handles[0] = verts[0][0];
1473  vtx_handles[1] = verts[1][0];
1474  vtx_handles[2] = verts[2][0];
1475  break;
1476  case -2: // left
1477  vtx_handles.resize( 3 );
1478  std::copy( verts[0], verts[0] + 3, vtx_handles.begin() );
1479  break;
1480  case -3: // upper left
1481  vtx_handles.resize( 1 );
1482  vtx_handles[0] = verts[0][0];
1483  break;
1484  default:
1485  vtx_handles.clear();
1486  break;
1487  }
1488  }
1489  else
1490  {
1491  switch( other_logical - local_logical )
1492  {
1493  case 0:
1494  return iBase_FAILURE;
1495  case 1: // below
1496  vtx_handles.resize( 3 );
1497  vtx_handles[0] = verts[0][2];
1498  vtx_handles[1] = verts[1][2];
1499  vtx_handles[2] = verts[2][2];
1500  break;
1501  case 2: // right
1502  vtx_handles.resize( 3 );
1503  std::copy( verts[2], verts[2] + 3, vtx_handles.begin() );
1504  break;
1505  case 3: // lower right
1506  vtx_handles.resize( 1 );
1507  vtx_handles[0] = verts[2][2];
1508  break;
1509  case -1: // lower left
1510  vtx_handles.resize( 1 );
1511  vtx_handles[0] = verts[0][2];
1512  break;
1513  case -2: // left
1514  vtx_handles.resize( 3 );
1515  std::copy( verts[0], verts[0] + 3, vtx_handles.begin() );
1516  break;
1517  default:
1518  vtx_handles.clear();
1519  break;
1520  }
1521  }
1522 
1523  return iBase_SUCCESS;
1524 }
1525 
1526 /**\brief Test querying of part boundary entities
1527  *
1528  * Test:
1529  * - iMeshP_getNumPartBdryEnts
1530  * - iMeshP_getPartBdryEnts
1531  */
1533 {
1534  int ierr, rank;
1535  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
1536 
1537  // get local part handles and part ids, and global part id list
1538  std::vector< iMeshP_PartHandle > local_handles;
1539  std::vector< iMeshP_Part > local_ids;
1540  std::vector< iMeshP_Part > all_parts = map.get_parts();
1541  std::map< iMeshP_PartHandle, std::vector< iBase_EntityHandle > > part_bdry;
1542  ierr = get_local_parts( imesh, prtn, local_handles, &local_ids );CHKERR;
1543 
1544  // for each combination of local part with any other part,
1545  // check for valid function values.
1546  std::vector< std::pair< iMeshP_Part, iMeshP_Part > > num_failed, num_error, list_failed, list_error, error;
1547  for( size_t i = 0; i < local_handles.size(); ++i )
1548  {
1549  iMeshP_PartHandle local_handle = local_handles[i];
1550  iMeshP_Part local_id = local_ids[i];
1551  for( std::vector< iMeshP_Part >::iterator j = all_parts.begin(); j != all_parts.end(); ++j )
1552  {
1553  iMeshP_Part other_id = *j;
1554  if( other_id == local_id ) continue;
1555 
1556  std::pair< iMeshP_Part, iMeshP_Part > part_pair;
1557  part_pair.first = local_id;
1558  part_pair.second = other_id;
1559 
1560  // get expected values
1561  std::vector< iBase_EntityHandle > shared_verts;
1562  ierr = interface_verts( imesh, prtn, local_handle, other_id, map, shared_verts );
1563  if( ierr != iBase_SUCCESS )
1564  {
1565  error.push_back( part_pair );
1566  continue;
1567  }
1568  std::sort( shared_verts.begin(), shared_verts.end() );
1569 
1570  // test iMeshP_getNumPartBdryEnts
1571  int count;
1572  iMeshP_getNumPartBdryEnts( imesh, prtn, local_handle, iBase_VERTEX, iMesh_POINT, other_id, &count, &ierr );
1573  if( iBase_SUCCESS != ierr )
1574  num_error.push_back( part_pair );
1575  else if( count != (int)shared_verts.size() )
1576  num_failed.push_back( part_pair );
1577 
1578  // test iMeshP_getPartBdryEnts
1579  iBase_EntityHandle* ptr = 0;
1580  int junk = 0;
1581  iMeshP_getPartBdryEnts( imesh, prtn, local_handle, iBase_VERTEX, iMesh_POINT, other_id, &ptr, &junk, &count,
1582  &ierr );
1583  if( iBase_SUCCESS != ierr )
1584  list_error.push_back( part_pair );
1585  else
1586  {
1587  std::copy( ptr, ptr + count, std::back_inserter( part_bdry[local_handles[i]] ) );
1588  std::sort( ptr, ptr + count );
1589  if( (int)shared_verts.size() != count || !std::equal( shared_verts.begin(), shared_verts.end(), ptr ) )
1590  list_failed.push_back( part_pair );
1591  free( ptr );
1592  }
1593  }
1594  }
1595 
1596  if( !error.empty() )
1597  {
1598  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1599  << " Internal error for " << error.size() << " part pairs." << std::endl;
1600  ierr = iBase_FAILURE;
1601  }
1602  if( !num_error.empty() )
1603  {
1604  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1605  << " iMeshP_getNumPartBdryEnts return error for " << num_error.size() << " part pairs." << std::endl;
1606  ierr = iBase_FAILURE;
1607  }
1608  if( !list_error.empty() )
1609  {
1610  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1611  << " iMeshP_getPartBdryEnts return error for " << list_error.size() << " part pairs." << std::endl;
1612  ierr = iBase_FAILURE;
1613  }
1614  if( !num_failed.empty() )
1615  {
1616  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1617  << " iMeshP_getNumPartBdryEnts return incorrect results for " << num_failed.size() << " part pairs."
1618  << std::endl;
1619  ierr = iBase_FAILURE;
1620  }
1621  if( !list_failed.empty() )
1622  {
1623  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1624  << " iMeshP_getPartBdryEnts return incorrect results for " << list_failed.size() << " part pairs."
1625  << std::endl;
1626  ierr = iBase_FAILURE;
1627  }
1628 
1629  if( iBase_SUCCESS != ierr ) return ierr;
1630 
1631  // test with iMeshP_ALL_PARTS
1632  for( size_t i = 0; i < local_handles.size(); ++i )
1633  {
1634  std::vector< iBase_EntityHandle >& exp_bdry = part_bdry[local_handles[i]];
1635  std::sort( exp_bdry.begin(), exp_bdry.end() );
1636  exp_bdry.erase( std::unique( exp_bdry.begin(), exp_bdry.end() ), exp_bdry.end() );
1637  std::pair< iMeshP_Part, iMeshP_Part > part_pair;
1638  part_pair.first = local_ids[i];
1639  part_pair.second = iMeshP_ALL_PARTS;
1640 
1641  int num = 0;
1642  iMeshP_getNumPartBdryEnts( imesh, prtn, local_handles[i], iBase_VERTEX, iMesh_POINT, iMeshP_ALL_PARTS, &num,
1643  &ierr );
1644  if( ierr )
1645  num_error.push_back( part_pair );
1646  else if( num != (int)exp_bdry.size() )
1647  num_failed.push_back( part_pair );
1648 
1649  iBase_EntityHandle* bdry = 0;
1650  int junk = num = 0;
1651  iMeshP_getPartBdryEnts( imesh, prtn, local_handles[i], iBase_VERTEX, iMesh_POINT, iMeshP_ALL_PARTS, &bdry,
1652  &junk, &num, &ierr );
1653  if( ierr )
1654  list_error.push_back( part_pair );
1655  else
1656  {
1657  std::sort( bdry, bdry + num );
1658  if( num != (int)exp_bdry.size() || !std::equal( bdry, bdry + num, exp_bdry.begin() ) )
1659  list_failed.push_back( part_pair );
1660  free( bdry );
1661  }
1662  }
1663  if( !num_error.empty() )
1664  {
1665  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1666  << " iMeshP_getNumPartBdryEnts return error for " << num_error.size() << " part pairs." << std::endl;
1667  ierr = iBase_FAILURE;
1668  }
1669  if( !list_error.empty() )
1670  {
1671  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1672  << " iMeshP_getPartBdryEnts return error for " << list_error.size() << " part pairs." << std::endl;
1673  ierr = iBase_FAILURE;
1674  }
1675  if( !num_failed.empty() )
1676  {
1677  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1678  << " iMeshP_getNumPartBdryEnts return incorrect results for " << num_failed.size() << " part pairs."
1679  << std::endl;
1680  ierr = iBase_FAILURE;
1681  }
1682  if( !list_failed.empty() )
1683  {
1684  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1685  << " iMeshP_getPartBdryEnts return incorrect results for " << list_failed.size() << " part pairs."
1686  << std::endl;
1687  ierr = iBase_FAILURE;
1688  }
1689 
1690  return ierr;
1691 }
1692 
1693 /**\brief Test querying of part boundary entities
1694  *
1695  * Test:
1696  * - iMeshP_initPartBdryEntIter
1697  * - iMeshP_initPartBdryEntArrIter
1698  */
1700 {
1701  int ierr, rank, has_data;
1702  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
1703 
1704  // get local part handles and part ids, and global part id list
1705  std::vector< iMeshP_PartHandle > local_handles;
1706  std::vector< iMeshP_Part > local_ids;
1707  std::vector< iMeshP_Part > all_parts = map.get_parts();
1708  ierr = get_local_parts( imesh, prtn, local_handles, &local_ids );CHKERR;
1709 
1710  std::vector< std::pair< iMeshP_Part, iMeshP_Part > > single_failed, single_error, single_step_error, array_failed,
1711  array_error, array_step_error;
1712  for( size_t i = 0; i < local_handles.size(); ++i )
1713  {
1714  iMeshP_PartHandle local_handle = local_handles[i];
1715  iMeshP_Part local_id = local_ids[i];
1716  for( std::vector< iMeshP_Part >::iterator j = all_parts.begin(); j != all_parts.end(); ++j )
1717  {
1718  iMeshP_Part other_id = *j;
1719  if( other_id == local_id ) continue;
1720 
1721  std::pair< iMeshP_Part, iMeshP_Part > part_pair;
1722  part_pair.first = local_id;
1723  part_pair.second = other_id;
1724 
1725  // get expected values
1726  std::vector< iBase_EntityHandle > shared_verts;
1727  ierr = interface_verts( imesh, prtn, local_handle, other_id, map, shared_verts );
1728  if( ierr != iBase_SUCCESS || 0 == shared_verts.size() ) continue;
1729  std::sort( shared_verts.begin(), shared_verts.end() );
1730 
1731  // test single entity iterator
1732  iBase_EntityIterator siter;
1733  iMeshP_initPartBdryEntIter( imesh, prtn, local_handle, iBase_VERTEX, iMesh_POINT, other_id, &siter, &ierr );
1734  if( ierr != iBase_SUCCESS )
1735  {
1736  single_error.push_back( part_pair );
1737  }
1738  else
1739  {
1740  std::vector< iBase_EntityHandle > results;
1741  for( ;; )
1742  {
1743  iBase_EntityHandle handle;
1744  iMesh_getNextEntIter( imesh, siter, &handle, &has_data, &ierr );
1745  if( ierr != iBase_SUCCESS )
1746  {
1747  single_step_error.push_back( part_pair );
1748  break;
1749  }
1750  if( !has_data ) break;
1751  results.push_back( handle );
1752  }
1753 
1754  std::sort( results.begin(), results.end() );
1755  if( results.size() != shared_verts.size() ||
1756  !std::equal( results.begin(), results.end(), shared_verts.begin() ) )
1757  single_failed.push_back( part_pair );
1758  }
1759  iMesh_endEntIter( imesh, siter, &ierr );
1760 
1761  // test array iterator
1763  iMeshP_initPartBdryEntArrIter( imesh, prtn, local_handle, iBase_VERTEX, iMesh_POINT, shared_verts.size(),
1764  other_id, &aiter, &ierr );
1765  if( ierr != iBase_SUCCESS )
1766  {
1767  array_error.push_back( part_pair );
1768  iMesh_endEntArrIter( imesh, aiter, &ierr );
1769  continue;
1770  }
1771  iBase_EntityHandle results[5], *ptr = results;
1772  int junk = 5, count;
1773  iMesh_getNextEntArrIter( imesh, aiter, &ptr, &junk, &count, &has_data, &ierr );
1774  if( ierr != iBase_SUCCESS || !has_data )
1775  {
1776  array_step_error.push_back( part_pair );
1777  iMesh_endEntArrIter( imesh, aiter, &ierr );
1778  continue;
1779  }
1780  iMesh_endEntArrIter( imesh, aiter, &ierr );
1781  assert( count <= 5 );
1782  assert( ptr == results );
1783  std::sort( ptr, ptr + count );
1784  if( count != (int)shared_verts.size() || !std::equal( shared_verts.begin(), shared_verts.end(), results ) )
1785  array_failed.push_back( part_pair );
1786  }
1787  }
1788 
1789  if( !single_error.empty() )
1790  {
1791  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1792  << " iMeshP_initPartBdryEntIter return error for " << single_error.size() << " part pairs."
1793  << std::endl;
1794  ierr = iBase_FAILURE;
1795  }
1796  if( !single_step_error.empty() )
1797  {
1798  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1799  << " iMesh_getNextEntIter return error for " << single_step_error.size() << " part pairs."
1800  << std::endl;
1801  ierr = iBase_FAILURE;
1802  }
1803  if( !single_failed.empty() )
1804  {
1805  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1806  << " iMeshP_initPartBdryEntIter iterator iterated over invalid entities for " << single_failed.size()
1807  << " part pairs." << std::endl;
1808  ierr = iBase_FAILURE;
1809  }
1810 
1811  if( !array_error.empty() )
1812  {
1813  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1814  << " iMeshP_initPartBdryEntArrIter return error for " << array_error.size() << " part pairs."
1815  << std::endl;
1816  ierr = iBase_FAILURE;
1817  }
1818  if( !array_step_error.empty() )
1819  {
1820  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1821  << " iMesh_getNextEntArrIter return error for " << array_step_error.size() << " part pairs."
1822  << std::endl;
1823  ierr = iBase_FAILURE;
1824  }
1825  if( !array_failed.empty() )
1826  {
1827  std::cerr << "Processor " << rank << ": Error at " __FILE__ ":" << __LINE__ << std::endl
1828  << " iMeshP_initPartBdryEntArrIter iterator iterated over invalid entities for "
1829  << array_failed.size() << " part pairs." << std::endl;
1830  ierr = iBase_FAILURE;
1831  }
1832 
1833  return ierr;
1834 }
1835 
1836 /**\brief Test adjacent entity query
1837  *
1838  * Test:
1839  * - iMeshP_getAdjEntities
1840  */
1842 {
1843  return iBase_SUCCESS;
1844 }
1845 
1846 /**\brief Test entity iterators
1847  *
1848  * Test:
1849  * - iMeshP_initEntIter
1850  * - iMeshP_initEntArrIter
1851  */
1853 {
1854  return iBase_SUCCESS;
1855 }
1856 
1857 /**\brief Test entity owner queries
1858  *
1859  * Test:
1860  * - iMeshP_getEntOwnerPart
1861  * - iMeshP_getEntOwnerPartArr
1862  * - iMeshP_isEntOwner
1863  * - iMeshP_isEntOwnerArr
1864  */
1866 {
1867  int ierr, rank, size;
1868  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
1869  MPI_Comm_size( MPI_COMM_WORLD, &size );
1870 
1871  // get local part handles and part ids
1872  std::vector< iMeshP_PartHandle > local_handles;
1873  std::vector< iMeshP_Part > local_ids;
1874  ierr = get_local_parts( imesh, prtn, local_handles, &local_ids );
1875  PCHECK;
1876 
1877  // test iMeshP_getEntOwnerPart for quads in each part
1878  std::vector< iBase_EntityHandle > all_quads;
1879  std::vector< iMeshP_Part > quad_owners;
1880  int invalid_count = 0;
1881  for( size_t i = 0; i < local_handles.size(); ++i )
1882  {
1883  std::vector< iBase_EntityHandle > quads;
1884  ierr = get_entities( imesh, local_handles[0], iBase_FACE, iMesh_QUADRILATERAL, quads );
1885  if( ierr ) break;
1886 
1887  for( size_t j = 0; j < quads.size(); ++j )
1888  {
1889  all_quads.push_back( quads[j] );
1890  quad_owners.push_back( local_ids[i] );
1891  iMeshP_Part owner;
1892  iMeshP_getEntOwnerPart( imesh, prtn, quads[j], &owner, &ierr );
1893  if( iBase_SUCCESS != ierr ) break;
1894 
1895  if( owner != local_ids[i] ) ++invalid_count;
1896  }
1897  if( iBase_SUCCESS != ierr ) break;
1898  }
1899  PCHECK;
1900  ASSERT( 0 == invalid_count );
1901 
1902  // test iMeshP_getEntOwnerPartArr for quads in each part
1903  invalid_count = 0;
1904  for( size_t i = 0; i < local_handles.size(); ++i )
1905  {
1906  std::vector< iBase_EntityHandle > quads;
1907  ierr = get_entities( imesh, local_handles[0], iBase_FACE, iMesh_QUADRILATERAL, quads );
1908  if( ierr ) break;
1909 
1910  std::vector< iMeshP_Part > owners( quads.size() ), expected( quads.size(), local_ids[i] );
1911  int junk = owners.size(), count;
1912  iMeshP_Part* ptr = &owners[0];
1913  iMeshP_getEntOwnerPartArr( imesh, prtn, &quads[0], quads.size(), &ptr, &junk, &count, &ierr );
1914  if( ierr ) break;
1915  assert( ptr == &owners[0] );
1916  assert( junk == (int)owners.size() );
1917  assert( count == (int)quads.size() );
1918  if( owners != expected ) ++invalid_count;
1919  }
1920  PCHECK;
1921  ASSERT( 0 == invalid_count );
1922 
1923  // get all vertices
1924  iBase_EntityHandle* vtx_arr = 0;
1925  int junk1 = 0, num_vtx;
1926  int *junk2 = 0, junk3 = 0, junk4;
1927  iMesh_getEntArrAdj( imesh, &all_quads[0], all_quads.size(), iBase_VERTEX, &vtx_arr, &junk1, &num_vtx, &junk2,
1928  &junk3, &junk4, &ierr );
1929  PCHECK;
1930  free( junk2 );
1931  std::sort( vtx_arr, vtx_arr + num_vtx );
1932  num_vtx = std::unique( vtx_arr, vtx_arr + num_vtx ) - vtx_arr;
1933  std::vector< iBase_EntityHandle > all_verts( vtx_arr, vtx_arr + num_vtx );
1934  free( vtx_arr );
1935 
1936  // check consistency between iMeshP_getEntOwnerPart and iMeshP_getEntOwnerPartArr
1937  // for all vertices
1938  std::vector< iMeshP_Part > vert_owners( all_verts.size() );
1939  junk1 = vert_owners.size();
1940  iMeshP_Part* junk5 = &vert_owners[0];
1941  iMeshP_getEntOwnerPartArr( imesh, prtn, &all_verts[0], all_verts.size(), &junk5, &junk1, &junk3, &ierr );
1942  PCHECK;
1943  assert( junk5 == &vert_owners[0] );
1944  assert( junk1 == (int)vert_owners.size() );
1945  assert( junk3 == (int)all_verts.size() );
1946 
1947  invalid_count = 0;
1948  for( size_t i = 0; i < all_verts.size(); ++i )
1949  {
1950  iMeshP_Part owner;
1951  iMeshP_getEntOwnerPart( imesh, prtn, all_verts[i], &owner, &ierr );
1952  if( iBase_SUCCESS != ierr || owner != vert_owners[i] ) ++invalid_count;
1953  }
1954  ASSERT( 0 == invalid_count );
1955 
1956  // get lists for all entities
1957  std::vector< iBase_EntityHandle > all_entities( all_verts );
1958  std::copy( all_quads.begin(), all_quads.end(), std::back_inserter( all_entities ) );
1959  std::vector< iMeshP_Part > all_owners( vert_owners );
1960  std::copy( quad_owners.begin(), quad_owners.end(), std::back_inserter( all_owners ) );
1961 
1962  // check consistency of iMeshP_isEntOwner for all entities
1963  invalid_count = 0;
1964  ierr = iBase_SUCCESS;
1965  for( size_t i = 0; i < local_handles.size(); ++i )
1966  {
1967  for( size_t j = 0; ierr == iBase_SUCCESS && j < all_entities.size(); ++j )
1968  {
1969  int is_owner;
1970  iMeshP_isEntOwner( imesh, prtn, local_handles[i], all_entities[j], &is_owner, &ierr );
1971  if( ierr != iBase_SUCCESS ) break;
1972  if( !is_owner == ( local_ids[i] == all_owners[j] ) ) ++invalid_count;
1973  }
1974  }
1975  PCHECK;
1976  ASSERT( 0 == invalid_count );
1977 
1978  // check consistency of iMeshP_isEntOwnerArr for all entities
1979  for( size_t i = 0; i < local_handles.size(); ++i )
1980  {
1981  std::vector< int > is_owner_list( all_entities.size() );
1982  junk1 = is_owner_list.size();
1983  int* junk6 = &is_owner_list[0];
1984  iMeshP_isEntOwnerArr( imesh, prtn, local_handles[i], &all_entities[0], all_entities.size(), &junk6, &junk1,
1985  &junk3, &ierr );
1986  if( iBase_SUCCESS != ierr ) break;
1987  assert( junk6 == &is_owner_list[0] );
1988  assert( junk1 == (int)is_owner_list.size() );
1989  assert( junk3 == (int)all_entities.size() );
1990  invalid_count = 0;
1991  for( size_t j = 0; j < all_entities.size(); ++j )
1992  {
1993  if( !( is_owner_list[j] ) == ( local_ids[0] == all_owners[j] ) ) ++invalid_count;
1994  }
1995  }
1996  PCHECK;
1997  ASSERT( 0 == invalid_count );
1998 
1999  // check globally consistent owners for all vertices
2000 
2001  // first communicate total number of vertex entries to be sent to root proc
2002  int local_count = all_verts.size(), global_count = 0;
2003  ierr = MPI_Reduce( &local_count, &global_count, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD );CHKERR;
2004 
2005  // for each vertex, store { (x << 2) | y, owning part id }
2006  std::vector< int > vtxdata( 2 * all_verts.size() );
2007  std::vector< double > coords( 3 * all_verts.size() );
2008  ierr = get_coords( imesh, &all_verts[0], all_verts.size(), &coords[0] );CHKERR;
2009  for( size_t i = 0; i < all_verts.size(); ++i )
2010  {
2011  int x = (int)round( coords[3 * i] );
2012  int y = (int)round( coords[3 * i + 1] );
2013  vtxdata[2 * i] = ( x << 3 ) | y;
2014  vtxdata[2 * i + 1] = vert_owners[i];
2015  }
2016 
2017  // collect all data on root procesor
2018  std::vector< int > all_data( 2 * global_count );
2019  std::vector< int > displ( size ), counts( size );
2020  for( int i = 0; i < size; i++ )
2021  {
2022  counts[i] = vtxdata.size();
2023  displ[i] = i * vtxdata.size();
2024  }
2025 
2026  // we could have used a simple gather, because all sequences are the same
2027  ierr = MPI_Gatherv( &vtxdata[0], vtxdata.size(), MPI_INT, &all_data[0], &counts[0], &displ[0], MPI_INT, 0,
2029 
2030  if( rank == 0 )
2031  {
2032  // map from vertex tag to indices into data
2033  std::multimap< int, int > data_map;
2034  for( int i = 0; i < global_count; ++i )
2035  {
2036  std::pair< int, int > p;
2037  p.first = all_data[2 * i];
2038  p.second = i;
2039  data_map.insert( p );
2040  }
2041 
2042  // check consistent data for each vtx
2043  std::multimap< int, int >::const_iterator a, b;
2044  for( a = data_map.begin(); a != data_map.end(); a = b )
2045  {
2046  for( b = a; b != data_map.end() && a->first == b->first; ++b )
2047  {
2048  int idx1 = a->second;
2049  int idx2 = b->second;
2050  if( all_data[2 * idx1 + 1] == all_data[2 * idx2 + 1] ) continue;
2051 
2052  ierr = iBase_FAILURE;
2053 
2054  int proc1 = std::lower_bound( displ.begin(), displ.end(), 2 * idx1 ) - displ.begin();
2055  if( displ[proc1] != 2 * idx1 ) ++proc1;
2056  int proc2 = std::lower_bound( displ.begin(), displ.end(), 2 * idx2 ) - displ.begin();
2057  if( displ[proc2] != 2 * idx2 ) ++proc2;
2058 
2059  std::cerr << "Error at " __FILE__ ":" << __LINE__ << " : " << std::endl
2060  << " For vertex at (" << ( a->first >> 2 ) << ", " << ( a->first & 3 ) << ") :" << std::endl
2061  << " Processor " << proc1 << " has " << all_data[2 * idx1 + 1] << " as the owning part"
2062  << std::endl
2063  << " Processor " << proc2 << " has " << all_data[2 * idx2 + 1] << " as the owning part"
2064  << std::endl;
2065  }
2066  }
2067  }
2068 
2069  return ierr;
2070 }
2071 
2074  const PartMap& map,
2075  iMeshP_PartHandle part,
2076  std::vector< iBase_EntityHandle >& boundary )
2077 {
2078  int ierr, logical_id;
2079  ierr = map.part_from_coords( imesh, part, logical_id );CHKERR;
2080 
2081  int neighbors[5], num_neighbors;
2082  get_part_neighbors( logical_id, map.get_parts().size(), neighbors, num_neighbors );
2083 
2084  for( int j = 0; j < num_neighbors; ++j )
2085  {
2086  std::vector< iBase_EntityHandle > iface;
2087  ierr = interface_verts( imesh, prtn, part, neighbors[j], map, iface );CHKERR;
2088  std::copy( iface.begin(), iface.end(), std::back_inserter( boundary ) );
2089  }
2090 
2091  std::sort( boundary.begin(), boundary.end() );
2092  boundary.erase( std::unique( boundary.begin(), boundary.end() ), boundary.end() );
2093  return iBase_SUCCESS;
2094 }
2095 
2096 /**\brief Test entity status
2097  *
2098  * Test:
2099  * - iMeshP_getEntStatus
2100  * - iMeshP_getEntStatusArr
2101  */
2103 {
2104  int ierr, rank, size;
2105  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
2106  MPI_Comm_size( MPI_COMM_WORLD, &size );
2107 
2108  // get local part handles
2109  std::vector< iMeshP_PartHandle > parts;
2110  ierr = get_local_parts( imesh, prtn, parts );
2111  PCHECK;
2112 
2113  // for each part
2114  int num_quad_ent_incorrect = 0, num_quad_ent_error = 0;
2115  int num_quad_arr_incorrect = 0, num_quad_arr_error = 0;
2116  int num_vert_ent_incorrect = 0, num_vert_ent_error = 0;
2117  int num_vert_arr_incorrect = 0, num_vert_arr_error = 0;
2118  for( size_t i = 0; i < parts.size(); ++i )
2119  {
2120  const iMeshP_PartHandle part = parts[i];
2121 
2122  // get quads and vertices
2123  std::vector< iBase_EntityHandle > quads, verts;
2124  ierr = get_part_quads_and_verts( imesh, part, quads, verts );
2125  if( ierr ) break;
2126 
2127  // check quad status (no ghosting yet)
2128  for( size_t j = 0; j < quads.size(); ++j )
2129  {
2130  int status;
2131  iMeshP_getEntStatus( imesh, prtn, part, quads[j], &status, &ierr );
2132  if( ierr != iBase_SUCCESS )
2133  {
2134  ++num_quad_ent_error;
2135  ierr = iBase_SUCCESS;
2136  continue;
2137  }
2138 
2139  if( status != iMeshP_INTERNAL ) ++num_quad_ent_incorrect;
2140  }
2141 
2142  // check quad status using iMeshP_getEntStatusArr
2143  std::vector< int > stat_list( quads.size() );
2144  int* junk1 = &stat_list[0];
2145  int junk2 = stat_list.size(), count;
2146  iMeshP_getEntStatusArr( imesh, prtn, part, &quads[0], quads.size(), &junk1, &junk2, &count, &ierr );
2147  if( ierr != iBase_SUCCESS )
2148  {
2149  ++num_quad_arr_error;
2150  ierr = iBase_SUCCESS;
2151  continue;
2152  }
2153  assert( junk1 == &stat_list[0] );
2154  assert( junk2 == (int)stat_list.size() );
2155  assert( count == (int)quads.size() );
2156  for( size_t j = 0; j < quads.size(); ++j )
2157  if( stat_list[j] != iMeshP_INTERNAL ) ++num_quad_arr_incorrect;
2158 
2159  // figure out which vertices are on the boundary
2160  std::vector< iBase_EntityHandle > boundary;
2161  ierr = get_part_boundary_verts( imesh, prtn, map, part, boundary );
2162  if( ierr ) break;
2163  std::sort( boundary.begin(), boundary.end() );
2164 
2165  // check vertex status (no ghosting yet)
2166  for( size_t j = 0; j < verts.size(); ++j )
2167  {
2168  int status;
2169  iMeshP_getEntStatus( imesh, prtn, part, verts[j], &status, &ierr );
2170  if( ierr != iBase_SUCCESS )
2171  {
2172  ++num_vert_ent_error;
2173  ierr = iBase_SUCCESS;
2174  continue;
2175  }
2176  bool on_boundary = std::binary_search( boundary.begin(), boundary.end(), verts[j] );
2177  if( status != ( on_boundary ? iMeshP_BOUNDARY : iMeshP_INTERNAL ) ) ++num_vert_ent_incorrect;
2178  }
2179 
2180  // check vert status using iMeshP_getEntStatusArr
2181  stat_list.resize( verts.size() );
2182  junk1 = &stat_list[0];
2183  junk2 = stat_list.size();
2184  iMeshP_getEntStatusArr( imesh, prtn, part, &verts[0], verts.size(), &junk1, &junk2, &count, &ierr );
2185  if( ierr != iBase_SUCCESS )
2186  {
2187  ++num_vert_arr_error;
2188  ierr = iBase_SUCCESS;
2189  continue;
2190  }
2191  assert( junk1 == &stat_list[0] );
2192  assert( junk2 == (int)stat_list.size() );
2193  assert( count == (int)verts.size() );
2194  for( size_t j = 0; j < verts.size(); ++j )
2195  {
2196  bool on_boundary = std::binary_search( boundary.begin(), boundary.end(), verts[j] );
2197  if( stat_list[j] != ( on_boundary ? iMeshP_BOUNDARY : iMeshP_INTERNAL ) ) ++num_vert_arr_incorrect;
2198  }
2199  }
2200  PCHECK; // check if loop interrupted by any internal errors
2201 
2202  ASSERT( 0 == num_quad_ent_error );
2203  ASSERT( 0 == num_quad_arr_error );
2204  ASSERT( 0 == num_vert_ent_error );
2205  ASSERT( 0 == num_vert_arr_error );
2206  ASSERT( 0 == num_quad_ent_incorrect );
2207  ASSERT( 0 == num_quad_arr_incorrect );
2208  ASSERT( 0 == num_vert_ent_incorrect );
2209  ASSERT( 0 == num_vert_arr_incorrect );
2210 
2211  return iBase_SUCCESS;
2212 }
2213 
2214 /**\brief Test information about entity copies for interface entities
2215  *
2216  * Test:
2217  * - iMeshP_getNumCopies
2218  * - iMeshP_getCopyParts
2219  */
2221 {
2222  int ierr, rank, size;
2223  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
2224  MPI_Comm_size( MPI_COMM_WORLD, &size );
2225 
2226  // get local part handles
2227  std::vector< iMeshP_PartHandle > parts;
2228  ierr = get_local_parts( imesh, prtn, parts );
2229  PCHECK;
2230  ASSERT( !parts.empty() );
2231 
2232  // select a singe part to test
2233  const iMeshP_PartHandle part = parts[0];
2234  int logical_id;
2235  ierr = map.part_from_coords( imesh, part, logical_id );CHKERR;
2236  const iMeshP_Part part_id = map.part_id_from_local_id( logical_id );
2237 
2238  // get vertices in part
2239  std::vector< iBase_EntityHandle > quads, verts;
2240  ierr = get_part_quads_and_verts( imesh, part, quads, verts );
2241  PCHECK;
2242 
2243  // get neighbors
2244  int neighbors[5], num_neighbors;
2245  get_part_neighbors( logical_id, map.get_parts().size(), neighbors, num_neighbors );
2246 
2247  // build map of sharing data for each vertex
2248  std::map< iBase_EntityHandle, std::vector< iMeshP_Part > > vert_sharing;
2249  for( int j = 0; j < num_neighbors; ++j )
2250  {
2251  std::vector< iBase_EntityHandle > iface;
2252  ierr = interface_verts( imesh, prtn, part, neighbors[j], map, iface );CHKERR;
2253  for( size_t k = 0; k < iface.size(); ++k )
2254  vert_sharing[iface[k]].push_back( map.part_id_from_local_id( neighbors[j] ) );
2255  }
2256 
2257  // test getNumCopies for each vertex
2258  std::map< iBase_EntityHandle, std::vector< iMeshP_Part > >::iterator i;
2259  int num_failed = 0, num_incorrect = 0;
2260  for( i = vert_sharing.begin(); i != vert_sharing.end(); ++i )
2261  {
2262  int count;
2263  iBase_EntityHandle vtx = i->first;
2264  iMeshP_getNumCopies( imesh, prtn, vtx, &count, &ierr );
2265  if( ierr )
2266  ++num_failed;
2267  else if( (unsigned)count != i->second.size() + 1 ) // add one for the part we queried from
2268  ++num_incorrect;
2269  }
2270  ASSERT( 0 == num_failed );
2271  ASSERT( 0 == num_incorrect );
2272 
2273  // get getCopyParts for each vertex
2274  num_failed = num_incorrect = 0;
2275  for( i = vert_sharing.begin(); i != vert_sharing.end(); ++i )
2276  {
2277  iMeshP_Part* list = 0;
2278  int junk = 0, count;
2279  iMeshP_getCopyParts( imesh, prtn, i->first, &list, &junk, &count, &ierr );
2280  if( iBase_SUCCESS != ierr )
2281  {
2282  ++num_failed;
2283  continue;
2284  }
2285  if( (unsigned)count != i->second.size() + 1 )
2286  { // add one for the part we queried from
2287  ++num_incorrect;
2288  free( list );
2289  continue;
2290  }
2291 
2292  std::vector< iMeshP_Part > expected( i->second );
2293  expected.push_back( part_id );
2294  std::sort( list, list + count );
2295  std::sort( expected.begin(), expected.end() );
2296  bool eq = std::equal( list, list + count, expected.begin() );
2297  free( list );
2298  if( !eq ) ++num_incorrect;
2299  }
2300  ASSERT( 0 == num_failed );
2301  ASSERT( 0 == num_incorrect );
2302 
2303  return iBase_SUCCESS;
2304 }
2305 
2306 // store remote handle data for a vertex
2308 {
2309  std::vector< iMeshP_Part > parts;
2310  std::vector< iBase_EntityHandle > handles;
2311 };
2312 
2313 /**\brief Test information about entity copies for interface entities
2314  *
2315  * Test:
2316  * - iMeshP_getCopies
2317  * - iMeshP_getCopyOnPart
2318  * - iMeshP_getOwnerCopy
2319  */
2321 {
2322  int ierr, rank, size;
2323  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
2324  MPI_Comm_size( MPI_COMM_WORLD, &size );
2325 
2326  // generate a unique ID for each vertex using the coordinates.
2327  // see create_mesh(..): each vertex has integer coordinates (x,y,0)
2328  // with x in [0,inf] and y in [0,4]
2329  // then to an Allgatherv to exchange handles for each processor
2330 
2331  // cast everything to iBase_EntityHandle so we can pack it all in one communication
2332  MPI_Datatype tmp_type;
2333  if( sizeof( iBase_EntityHandle ) == sizeof( unsigned ) )
2334  tmp_type = MPI_UNSIGNED;
2335  else if( sizeof( iBase_EntityHandle ) == sizeof( unsigned long ) )
2336  tmp_type = MPI_UNSIGNED_LONG;
2337  else if( sizeof( iBase_EntityHandle ) == sizeof( unsigned long long ) )
2338  tmp_type = MPI_UNSIGNED_LONG_LONG;
2339  else
2340  return iBase_FAILURE;
2341  const MPI_Datatype type = tmp_type; // make it const
2342 
2343  // get local part handles
2344  std::vector< iMeshP_PartHandle > parts;
2345  ierr = get_local_parts( imesh, prtn, parts );
2346  PCHECK;
2347  std::vector< iMeshP_Part > part_ids( parts.size() );
2348  iMeshP_Part* junk1 = &part_ids[0];
2349  int junk2 = part_ids.size(), junk3;
2350  iMeshP_getPartIdsFromPartHandlesArr( imesh, prtn, &parts[0], parts.size(), &junk1, &junk2, &junk3, &ierr );
2351  PCHECK;
2352  assert( junk1 == &part_ids[0] );
2353  assert( junk2 == (int)part_ids.size() );
2354  assert( junk3 == (int)parts.size() );
2355 
2356  // build list of {vtx_id, part_id, handle} tuples to send
2357  // also build list of local vertex handles
2358  std::vector< iBase_EntityHandle > local_data, local_vertices;
2359  for( size_t i = 0; i < parts.size(); ++i )
2360  {
2361  // get vertices
2362  std::vector< iBase_EntityHandle > quads, verts;
2363  ierr = get_part_quads_and_verts( imesh, parts[i], quads, verts );
2364  if( ierr ) break;
2365 
2366  // add all vertices to local_data
2367  for( size_t j = 0; j < verts.size(); ++j )
2368  {
2369  int tag = 0;
2370  ierr = vertex_tag( imesh, verts[j], tag );
2371  if( ierr ) break;
2372  long tmp_h = tag;
2373  local_data.push_back( (iBase_EntityHandle)tmp_h );
2374  tmp_h = part_ids[i];
2375  local_data.push_back( (iBase_EntityHandle)tmp_h );
2376  local_data.push_back( verts[j] );
2377  }
2378  if( ierr ) break;
2379 
2380  std::copy( verts.begin(), verts.end(), std::back_inserter( local_vertices ) );
2381  }
2382 
2383  // build list of local vertices
2384  std::sort( local_vertices.begin(), local_vertices.end() );
2385  local_vertices.erase( std::unique( local_vertices.begin(), local_vertices.end() ), local_vertices.end() );
2386  std::vector< int > local_vtx_tags( local_vertices.size() );CHKERR;
2387  for( size_t i = 0; i < local_vertices.size(); ++i )
2388  {
2389  ierr = vertex_tag( imesh, local_vertices[i], local_vtx_tags[i] );
2390  if( ierr ) break;
2391  }
2392  CHKERR;
2393 
2394  // communicate data
2395  std::vector< int > gcounts( size ), gdisp( size );
2396  int local_data_size = local_data.size();
2397  ierr = MPI_Allgather( &local_data_size, 1, MPI_INT, &gcounts[0], 1, MPI_INT, MPI_COMM_WORLD );CHKERR;
2398  gdisp[0] = 0;
2399  for( int i = 1; i < size; ++i )
2400  gdisp[i] = gdisp[i - 1] + gcounts[i - 1];
2401  std::vector< iBase_EntityHandle > global_data( gdisp[size - 1] + gcounts[size - 1] );
2402  ierr = MPI_Allgatherv( &local_data[0], local_data_size, type, &global_data[0], &gcounts[0], &gdisp[0], type,
2404 
2405  // arrange global data in a more useful way
2406  std::map< int, VtxCopyData > vtx_sharing;
2407  assert( global_data.size() % 3 == 0 );
2408  for( size_t i = 0; i < global_data.size(); i += 3 )
2409  {
2410  int tag = (int)(size_t)global_data[i];
2411  iMeshP_Part part = (iMeshP_Part)(size_t)global_data[i + 1];
2412  iBase_EntityHandle handle = global_data[i + 2];
2413  vtx_sharing[tag].parts.push_back( part );
2414  vtx_sharing[tag].handles.push_back( handle );
2415  }
2416 
2417  // test iMeshP_getCopies for each local vertex
2418  int num_error = 0, num_incorrect = 0, junk4;
2419  for( size_t i = 0; i < local_vertices.size(); ++i )
2420  {
2421  int num_copies = -1;
2422  // iMeshP_Part* part_ids = 0;
2423  iMeshP_Part* ptr_part_ids = 0; // Use ptr_part_ids to avoid shadowing std::vector<iMeshP_Part> part_ids
2424  iBase_EntityHandle* copies = 0;
2425  junk2 = junk3 = junk4 = 0;
2426  iMeshP_getCopies( imesh, prtn, local_vertices[i], &ptr_part_ids, &junk2, &num_copies, &copies, &junk3, &junk4,
2427  &ierr );
2428  if( iBase_SUCCESS != ierr )
2429  {
2430  ++num_error;
2431  continue;
2432  }
2433  assert( junk4 == num_copies );
2434 
2435  VtxCopyData& expected = vtx_sharing[local_vtx_tags[i]];
2436  if( num_copies != (int)expected.parts.size() )
2437  ++num_incorrect;
2438  else
2439  for( size_t j = 0; j < expected.parts.size(); ++j )
2440  {
2441  int idx = std::find( ptr_part_ids, ptr_part_ids + num_copies, expected.parts[j] ) - ptr_part_ids;
2442  if( idx == num_copies || copies[idx] != expected.handles[j] )
2443  {
2444  ++num_incorrect;
2445  break;
2446  }
2447  }
2448  free( ptr_part_ids );
2449  free( copies );
2450  }
2451  ASSERT( 0 == num_error );
2452  ASSERT( 0 == num_incorrect );
2453 
2454  // test iMeshP_getCopyOnPart for each local vertex
2455  num_error = num_incorrect = 0;
2456  for( size_t i = 0; i < local_vertices.size(); ++i )
2457  {
2458  VtxCopyData& expected = vtx_sharing[local_vtx_tags[i]];
2459  for( size_t j = 0; j < expected.parts.size(); ++j )
2460  {
2461  iBase_EntityHandle copy;
2462  iMeshP_getCopyOnPart( imesh, prtn, local_vertices[i], expected.parts[j], &copy, &ierr );
2463  if( iBase_SUCCESS != ierr )
2464  ++num_error;
2465  else if( expected.handles[j] != copy )
2466  ++num_incorrect;
2467  }
2468  }
2469  ASSERT( 0 == num_error );
2470  ASSERT( 0 == num_incorrect );
2471 
2472  // test iMeshP_getOwnerCopy for each local vertex
2473  num_error = num_incorrect = 0;
2474  for( size_t i = 0; i < local_vertices.size(); ++i )
2475  {
2476  VtxCopyData& expected = vtx_sharing[local_vtx_tags[i]];
2477  iMeshP_Part owner_id = 0;
2478  iMeshP_getEntOwnerPart( imesh, prtn, local_vertices[i], &owner_id, &ierr );
2479  if( iBase_SUCCESS != ierr ) continue; // not testing getEntOwnerPart here
2480 
2481  size_t idx = std::find( expected.parts.begin(), expected.parts.end(), owner_id ) - expected.parts.begin();
2482  if( idx == expected.parts.size() ) continue; // not testing getEntOwnerPart here
2483 
2484  iMeshP_Part owner_id_2 = 0;
2485  iBase_EntityHandle copy = 0;
2486  iMeshP_getOwnerCopy( imesh, prtn, local_vertices[i], &owner_id_2, &copy, &ierr );
2487  if( iBase_SUCCESS != ierr )
2488  ++num_error;
2489  else if( owner_id_2 != owner_id && copy != expected.handles[idx] )
2490  ++num_incorrect;
2491  }
2492  ASSERT( 0 == num_error );
2493  ASSERT( 0 == num_incorrect );
2494 
2495  return iBase_SUCCESS;
2496 }
2497 
2499 {
2500  iBase_EntityHandle* list = 0;
2501  int ierr, junk = 0;
2502  iMesh_getEntAdj( imesh, vtx, iBase_FACE, &list, &junk, &num, &ierr );
2503  if( iBase_SUCCESS == ierr ) free( list );
2504  return ierr;
2505 }
2506 
2507 int get_adj( iMesh_Instance imesh, iBase_EntityHandle ent, int type, std::vector< iBase_EntityHandle >& adj )
2508 {
2509  iBase_EntityHandle* list = 0;
2510  int ierr, num, junk = 0;
2511  iMesh_getEntAdj( imesh, ent, type, &list, &junk, &num, &ierr );
2512  if( iBase_SUCCESS == ierr )
2513  {
2514  std::copy( list, list + num, std::back_inserter( adj ) );
2515  free( list );
2516  }
2517  return ierr;
2518 }
2519 
2520 // assume regular quad mesh
2521 int get_boundary_vertices( iMesh_Instance imesh, std::vector< iBase_EntityHandle >& bdry )
2522 {
2523  int ierr, n;
2524  iBase_EntitySetHandle root;
2525  iMesh_getRootSet( imesh, &root, &ierr );CHKERR;
2526  std::vector< iBase_EntityHandle > all_verts;
2527  ierr = get_entities( imesh, root, iBase_VERTEX, iMesh_POINT, all_verts );CHKERR;
2528  bdry.clear();
2529  for( size_t i = 0; i < all_verts.size(); ++i )
2530  {
2531  ierr = get_num_adj_quads( imesh, all_verts[i], n );CHKERR;
2532  if( n != 4 ) bdry.push_back( all_verts[i] );
2533  }
2534  return iBase_SUCCESS;
2535 }
2536 
2539  const std::vector< iBase_EntityHandle >& sorted_vertices )
2540 {
2541  int ierr;
2542  if( std::binary_search( sorted_vertices.begin(), sorted_vertices.end(), vtx ) ) return iBase_SUCCESS;
2543  std::vector< iBase_EntityHandle > quads, verts;
2544  ierr = get_adj( imesh, vtx, iBase_FACE, quads );CHKERR;
2545  for( size_t i = 0; i < quads.size(); ++i )
2546  {
2547  verts.clear();
2548  ierr = get_adj( imesh, quads[i], iBase_VERTEX, verts );CHKERR;
2549  for( size_t j = 0; j < verts.size(); ++j )
2550  {
2551  if( std::binary_search( sorted_vertices.begin(), sorted_vertices.end(), verts[j] ) ) return iBase_SUCCESS;
2552  }
2553  }
2554 
2555  return iBase_FAILURE;
2556 }
2557 
2558 // get number of adjacent quads to each vertex, both on the local
2559 // processor and in the entire mesh
2561  const std::vector< iBase_EntityHandle >& verts,
2562  std::vector< int >& num_local_adj,
2563  std::vector< int >& num_all_adj )
2564 {
2565  int ierr, size;
2566  MPI_Comm_size( MPI_COMM_WORLD, &size );
2567 
2568  std::vector< int > vtx_tags( verts.size() );
2569  num_local_adj.resize( verts.size() );
2570  for( size_t i = 0; i < verts.size(); ++i )
2571  {
2572  ierr = get_num_adj_quads( imesh, verts[i], num_local_adj[i] );CHKERR;
2573  ierr = vertex_tag( imesh, verts[i], vtx_tags[i] );CHKERR;
2574  }
2575 
2576  std::vector< int > counts( size ), displ( size );
2577  int num_vtx = verts.size();
2578  ierr = MPI_Allgather( &num_vtx, 1, MPI_INT, &counts[0], 1, MPI_INT, MPI_COMM_WORLD );CHKERR;
2579  displ[0] = 0;
2580  for( int i = 1; i < size; ++i )
2581  displ[i] = displ[i - 1] + counts[i - 1];
2582  int total = displ[size - 1] + counts[size - 1];
2583  std::vector< int > all_tags( total ), all_adj_counts( total );
2584  ierr = MPI_Allgatherv( &vtx_tags[0], vtx_tags.size(), MPI_INT, &all_tags[0], &counts[0], &displ[0], MPI_INT,
2586  ierr = MPI_Allgatherv( &num_local_adj[0], num_local_adj.size(), MPI_INT, &all_adj_counts[0], &counts[0], &displ[0],
2587  MPI_INT, MPI_COMM_WORLD );CHKERR;
2588 
2589  num_all_adj.clear();
2590  num_all_adj.resize( total, 0 );
2591  for( int i = 0; i < total; ++i )
2592  {
2593  std::vector< int >::iterator it = std::find( vtx_tags.begin(), vtx_tags.end(), all_tags[i] );
2594  if( it == vtx_tags.end() ) continue;
2595  int idx = it - vtx_tags.begin();
2596  num_all_adj[idx] += all_adj_counts[i];
2597  }
2598 
2599  return iBase_SUCCESS;
2600 }
2601 
2602 /**\brief Test creation of ghost entities
2603  *
2604  * Test:
2605  * - iMeshP_createGhostEntsAll
2606  */
2608 {
2609  int ierr;
2610 
2611  // get boundary vertices
2612  std::vector< iBase_EntityHandle > bdry;
2613  ierr = get_boundary_vertices( imesh, bdry );
2614  PCHECK;
2615  // get counts of adjacent entities
2616  std::vector< int > num_local_adj, num_global_adj;
2617  ierr = get_num_adj_all( imesh, bdry, num_local_adj, num_global_adj );
2618  PCHECK;
2619  // create one layer of ghost entities
2620  iMeshP_createGhostEntsAll( imesh, prtn, iBase_FACE, iBase_VERTEX, 1, 0, &ierr );
2621  PCHECK;
2622  // check that each vertex has the correct number of adjacent entities
2623  int num_incorrect = 0;
2624  for( size_t i = 0; i < bdry.size(); ++i )
2625  {
2626  int n;
2627  ierr = get_num_adj_quads( imesh, bdry[i], n );
2628  if( iBase_SUCCESS != ierr || num_global_adj[i] != n ) ++num_incorrect;
2629  }
2630  ASSERT( 0 == num_incorrect );
2631  // get new the new boundary
2632  std::vector< iBase_EntityHandle > new_bdry;
2633  ierr = get_boundary_vertices( imesh, new_bdry );
2634  PCHECK;
2635  // check that each vertex on the new boundary is separated by
2636  // at most one layer from the old boundary
2637  std::sort( bdry.begin(), bdry.end() );
2638  num_incorrect = 0;
2639  for( size_t i = 0; i < new_bdry.size(); ++i )
2640  {
2641  ierr = check_one_layer( imesh, new_bdry[i], bdry );
2642  if( ierr ) ++num_incorrect;
2643  }
2644  ASSERT( 0 == num_incorrect );
2645  // make another layer of ghost entiites
2646  bdry.swap( new_bdry );
2647  new_bdry.clear();
2648  ierr = get_num_adj_all( imesh, bdry, num_local_adj, num_global_adj );
2649  PCHECK;
2650  iMeshP_createGhostEntsAll( imesh, prtn, iBase_FACE, iBase_VERTEX, 2, 0, &ierr );
2651  PCHECK;
2652  // check that each vertex has the correct number of adjacent entities
2653  num_incorrect = 0;
2654  for( size_t i = 0; i < bdry.size(); ++i )
2655  {
2656  int n;
2657  ierr = get_num_adj_quads( imesh, bdry[i], n );
2658  if( iBase_SUCCESS != ierr || num_global_adj[i] != n ) ++num_incorrect;
2659  }
2660  // check that each vertex on the new boundary is separated by
2661  // at most one layer from the old boundary
2662  std::sort( bdry.begin(), bdry.end() );
2663  num_incorrect = 0;
2664  for( size_t i = 0; i < new_bdry.size(); ++i )
2665  {
2666  ierr = check_one_layer( imesh, new_bdry[i], bdry );
2667  if( ierr ) ++num_incorrect;
2668  }
2669  ASSERT( 0 == num_incorrect );
2670 
2671  return iBase_SUCCESS;
2672 }
2673 
2674 /**\brief Test exchange entities
2675  *
2676  * Test:
2677  * - iMeshP_exchEntArrToPartsAll
2678  */
2680 {
2681  int ierr, rank, size;
2682  int num_err = 0;
2683  iMeshP_RequestHandle request;
2684  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
2685  MPI_Comm_size( MPI_COMM_WORLD, &size );
2686 
2687  std::vector< iBase_EntityHandle > all_elems;
2688  std::vector< iMeshP_Part > all_ids;
2689  std::vector< iBase_EntityHandle > quads;
2690 
2691  // get local part handles and part ids
2692  std::vector< iMeshP_PartHandle > local_handles;
2693  std::vector< iMeshP_Part > local_ids;
2694  ierr = get_local_parts( imesh, prtn, local_handles, &local_ids );
2695  PCHECK;
2696 
2697  // get loacal quads before exchange
2698  quads.clear();
2699  ierr = get_entities( imesh, local_handles[0], iBase_FACE, iMesh_QUADRILATERAL, quads );CHKERR;
2700  int n_quads = quads.size();
2701 
2702  // send all elements in local processor to all other processors
2703  for( size_t i = 0; i < map.get_parts().size(); ++i )
2704  {
2705  if( map.get_parts()[i] == (unsigned int)rank ) continue; // skip own rank
2706 
2707  for( int j = 0; j < n_quads; j++ )
2708  {
2709  all_elems.push_back( quads[j] );
2710  all_ids.push_back( map.get_parts()[i] );
2711  }
2712  }
2713 
2714  // exchange entities
2715  iMeshP_exchEntArrToPartsAll( imesh, prtn, &all_elems[0], all_elems.size(), &all_ids[0], 0, 0, &request, &ierr );
2716  if( iBase_SUCCESS != ierr ) ++num_err;
2717 
2718  // get local quads after exchange
2719  quads.clear();
2720  ierr = get_entities( imesh, local_handles[0], iBase_FACE, iMesh_QUADRILATERAL, quads );CHKERR;
2721 
2722  // # of elements should be # of quads * # of processors
2723  ASSERT( quads.size() == (unsigned int)n_quads * size );
2724 
2725  ASSERT( 0 == num_err );
2726 
2727  return iBase_SUCCESS;
2728 }
2729 
2730 /**\brief Test commuinication of tag data
2731  *
2732  * Test:
2733  * - iMeshP_pushTags
2734  * - iMeshP_pushTagsEnt
2735  */
2736 int test_push_tag_data_common( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, int num_ghost_layers )
2737 {
2738  const char* src_name = "test_src";
2739  const char* dst_name = "test_dst";
2740  int ierr, rank;
2741  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
2742 
2743  if( num_ghost_layers )
2744  {
2745  iMeshP_createGhostEntsAll( imesh, prtn, iBase_FACE, iBase_VERTEX, num_ghost_layers, 0, &ierr );
2746  PCHECK;
2747  }
2748 
2749  iBase_TagHandle src_tag, dst_tag;
2750  iMesh_createTag( imesh, src_name, 1, iBase_INTEGER, &src_tag, &ierr, strlen( src_name ) );CHKERR;
2751  iMesh_createTag( imesh, dst_name, 1, iBase_INTEGER, &dst_tag, &ierr, strlen( dst_name ) );CHKERR;
2752 
2753  iBase_EntitySetHandle root;
2754  iMesh_getRootSet( imesh, &root, &ierr );CHKERR;
2755 
2756  std::vector< iBase_EntityHandle > verts;
2757  ierr = get_entities( imesh, root, iBase_VERTEX, iMesh_POINT, verts );CHKERR;
2758 
2759  // test iMeshP_pushTags
2760  // each processor writes its rank on all vertices
2761  // after push, each vertex should be tagged with the rank of its owner
2762 
2763  std::vector< int > tag_vals( verts.size(), rank );
2764  iMesh_setIntArrData( imesh, &verts[0], verts.size(), src_tag, &tag_vals[0], tag_vals.size(), &ierr );CHKERR;
2765 
2766  iMeshP_pushTags( imesh, prtn, src_tag, dst_tag, iBase_VERTEX, iMesh_POINT, &ierr );
2767  PCHECK;
2768 
2769  tag_vals.clear();
2770  tag_vals.resize( verts.size(), -1 );
2772  iMesh_getTagHandle( imesh, "GLOBAL_ID", &id_tag, &ierr, strlen( "GLOBAL_ID" ) );
2773  std::vector< int > ids( verts.size() );
2774  int *junk1 = &ids[0], junk2 = ids.size(), junk3;
2775  iMesh_getIntArrData( imesh, &verts[0], verts.size(), id_tag, &junk1, &junk2, &junk3, &ierr );
2776  PCHECK;
2777  int errcount = 0;
2778  for( size_t i = 0; i < verts.size(); ++i )
2779  {
2780  iMesh_getIntData( imesh, verts[i], dst_tag, &tag_vals[i], &ierr );
2781  if( ierr != iBase_SUCCESS )
2782  {
2783  std::cerr << "Rank " << rank << " : getIntData failed for vertex " << ids[i] << std::endl;
2784  std::cerr.flush();
2785  ++errcount;
2786  }
2787  }
2788  ASSERT( 0 == errcount );
2789 
2790  // int *junk1 = &tag_vals[0], junk2 = tag_vals.size(), junk3;
2791  // iMesh_getIntArrData( imesh, &verts[0], verts.size(), dst_tag, &junk1, &junk2, &junk3, &ierr
2792  // ); PCHECK; assert( junk1 == &tag_vals[0] ); assert( junk2 == (int)tag_vals.size() ); assert(
2793  // junk3 == (int)verts.size() );
2794 
2795  std::vector< int > expected( verts.size() );
2796  std::vector< iMeshP_Part > parts( verts.size() );
2797  iMeshP_Part* junk4 = &parts[0];
2798  junk2 = parts.size();
2799  iMeshP_getEntOwnerPartArr( imesh, prtn, &verts[0], verts.size(), &junk4, &junk2, &junk3, &ierr );
2800  PCHECK;
2801  assert( junk4 == &parts[0] );
2802  assert( junk2 == (int)parts.size() );
2803  assert( junk3 == (int)verts.size() );
2804  junk1 = &expected[0];
2805  junk2 = expected.size();
2806  iMeshP_getRankOfPartArr( imesh, prtn, &parts[0], parts.size(), &junk1, &junk2, &junk3, &ierr );
2807  PCHECK;
2808  assert( junk1 == &expected[0] );
2809  assert( junk2 == (int)expected.size() );
2810  assert( junk3 == (int)parts.size() );
2811 
2812  ASSERT( tag_vals == expected );
2813 
2814  // test iMeshP_pushTagsEnt
2815  // write -1 on all vertices
2816  // For each vertex owned by this processor and shared with more than
2817  // two others, write the rank of the owning processor.
2818 
2819  tag_vals.clear();
2820  tag_vals.resize( verts.size(), -1 );
2821  iMesh_setIntArrData( imesh, &verts[0], verts.size(), src_tag, &tag_vals[0], tag_vals.size(), &ierr );
2822  PCHECK;
2823  tag_vals.resize( verts.size(), -1 );
2824  iMesh_setIntArrData( imesh, &verts[0], verts.size(), dst_tag, &tag_vals[0], tag_vals.size(), &ierr );
2825  PCHECK;
2826 
2827  std::vector< iBase_EntityHandle > some;
2828  for( size_t i = 0; i < verts.size(); ++i )
2829  {
2830  int num;
2831  iMeshP_getNumCopies( imesh, prtn, verts[i], &num, &ierr );
2832  if( iBase_SUCCESS != ierr ) break;
2833  if( num > 2 )
2834  some.push_back( verts[i] );
2835  else
2836  expected[i] = -1;
2837  }
2838 
2839  tag_vals.clear();
2840  tag_vals.resize( some.size(), rank );
2841  iMesh_setIntArrData( imesh, &some[0], some.size(), src_tag, &tag_vals[0], tag_vals.size(), &ierr );
2842  PCHECK;
2843 
2844  iMeshP_pushTagsEnt( imesh, prtn, src_tag, dst_tag, &some[0], some.size(), &ierr );
2845  PCHECK;
2846 
2847  tag_vals.clear();
2848  tag_vals.resize( verts.size(), -1 );
2849  junk1 = &tag_vals[0];
2850  junk2 = tag_vals.size();
2851  iMesh_getIntArrData( imesh, &verts[0], verts.size(), dst_tag, &junk1, &junk2, &junk3, &ierr );CHKERR;
2852  assert( junk1 == &tag_vals[0] );
2853  assert( junk2 == (int)tag_vals.size() );
2854  assert( junk3 == (int)verts.size() );
2855 
2856  ASSERT( tag_vals == expected );
2857  return iBase_SUCCESS;
2858 }
2859 
2860 /**\brief Test commuinication of tag data
2861  *
2862  * Test:
2863  * - iMeshP_pushTags
2864  * - iMeshP_pushTagsEnt
2865  */
2867 {
2868  return test_push_tag_data_common( imesh, prtn, 0 );
2869 }
2870 
2871 /**\brief Test commuinication of tag data
2872  *
2873  * Test:
2874  * - iMeshP_pushTags
2875  * - iMeshP_pushTagsEnt
2876  */
2878 {
2879  return test_push_tag_data_common( imesh, prtn, 1 );
2880 }
2881 
2882 /**************************************************************************
2883  PartMap class
2884  **************************************************************************/
2885 
2886 int PartMap::build_map( iMesh_Instance imesh, iMeshP_PartitionHandle prtn, int num_expected_parts )
2887 {
2888  int ierr, rank, size;
2889  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
2890  MPI_Comm_size( MPI_COMM_WORLD, &size );
2891 
2892  // get local parts
2893  std::vector< iMeshP_PartHandle > local_parts;
2894  std::vector< iMeshP_Part > imesh_ids;
2895  ierr = get_local_parts( imesh, prtn, local_parts, &imesh_ids );CHKERR;
2896 
2897  // get logical ids for local parts
2898  std::vector< int > local_ids( local_parts.size() );
2899  for( size_t i = 0; i < local_parts.size(); ++i )
2900  {
2901  ierr = part_from_coords( imesh, local_parts[i], local_ids[i] );CHKERR;
2902  }
2903 
2904  // get total number of parts
2905  int num_global = 0, num_local = local_parts.size();
2906  ierr = MPI_Allreduce( &num_local, &num_global, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD );CHKERR;
2907  if( num_global != num_expected_parts )
2908  {
2909  std::cerr << "Invalid/unexpected global part count at " __FILE__ ":" << __LINE__ << " (proc " << rank
2910  << "): " << std::endl
2911  << " Expected: " << num_expected_parts << std::endl
2912  << " Actual: " << num_global << std::endl;
2913  return 1;
2914  }
2915 
2916  // get counts and displacements for Allgatherv calls
2917  std::vector< int > dspls( size ), counts( size );
2918  ierr = MPI_Allgather( &num_local, 1, MPI_INT, &counts[0], 1, MPI_INT, MPI_COMM_WORLD );CHKERR;
2919  dspls[0] = 0;
2920  for( int i = 1; i < size; ++i )
2921  dspls[i] = dspls[i - 1] + counts[i - 1];
2922 
2923  // gather iMeshP_Part list from each processor
2924  std::vector< unsigned > global_part_ids( num_expected_parts );
2925  assert( sizeof( iMeshP_Part ) == sizeof( int ) );
2926  ierr = MPI_Allgatherv( &imesh_ids[0], num_local, MPI_UNSIGNED, &global_part_ids[0], &counts[0], &dspls[0],
2927  MPI_UNSIGNED, MPI_COMM_WORLD );CHKERR;
2928 
2929  // gather local ids from each processor
2930  std::vector< int > global_id_list( num_expected_parts );
2931  ierr = MPI_Allgatherv( &local_ids[0], num_local, MPI_INT, &global_id_list[0], &counts[0], &dspls[0], MPI_INT,
2933 
2934  // build owner list
2935  std::vector< int > global_owners( num_expected_parts );
2936  for( int i = 0; i < size; ++i )
2937  for( int j = 0; j < counts[i]; ++j )
2938  global_owners[dspls[i] + j] = i;
2939 
2940  // populate member lists
2941  sortedPartList = global_part_ids;
2942  std::sort( sortedPartList.begin(), sortedPartList.end() );
2943  partLocalIds.resize( num_expected_parts );
2944  partRanks.resize( num_expected_parts );
2945  for( int i = 0; i < num_expected_parts; ++i )
2946  {
2947  int idx = std::lower_bound( sortedPartList.begin(), sortedPartList.end(), global_part_ids[i] ) -
2948  sortedPartList.begin();
2949  partLocalIds[idx] = global_id_list[i];
2950  partRanks[idx] = global_owners[i];
2951  }
2952 
2953  // do some consistency checking
2954  if( std::unique( sortedPartList.begin(), sortedPartList.end() ) != sortedPartList.end() )
2955  {
2956  if( rank == 0 )
2957  {
2958  std::cerr << "ERROR: Duplicate iMeshP_Part values detected at " __FILE__ ":" << __LINE__ << std::endl;
2959  }
2960  return 1;
2961  }
2962 
2963  // build revesre local id map and check for duplicates
2964  localIdReverseMap.clear();
2965  localIdReverseMap.resize( num_expected_parts, -1 );
2966  for( int i = 0; i < num_expected_parts; ++i )
2967  {
2968  int idx = partLocalIds[i];
2969  if( localIdReverseMap[idx] != -1 )
2970  {
2971  if( rank == 0 )
2972  {
2973  std::cerr << "ERROR: Part mesh has been duplicated in multiple parts." << std::endl
2974  << " Detected at " __FILE__ ":" << __LINE__ << std::endl
2975  << " See PartMap::part_from_coords" << std::endl;
2976  }
2977  return 1;
2978  }
2979  if( idx >= num_expected_parts )
2980  {
2981  if( rank == 0 )
2982  {
2983  std::cerr << "ERROR: Part mesh invalid/incorrect mesh." << std::endl
2984  << " Detected at " __FILE__ ":" << __LINE__ << std::endl
2985  << " See PartMap::part_from_coords" << std::endl;
2986  }
2987  return 1;
2988  }
2989 
2990  localIdReverseMap[idx] = i;
2991  }
2992 
2993  return 0;
2994 }
2995 
2996 void PartMap::part_id_from_rank( int rank, std::vector< iMeshP_Part >& parts ) const
2997 {
2998  for( size_t i = 0; i < sortedPartList.size(); ++i )
2999  if( partRanks[i] == rank ) parts.push_back( sortedPartList[i] );
3000 }
3001 
3002 void PartMap::local_id_from_rank( int rank, std::vector< int >& ids ) const
3003 {
3004  for( size_t i = 0; i < sortedPartList.size(); ++i )
3005  if( partRanks[i] == rank ) ids.push_back( partLocalIds[i] );
3006 }
3007 
3009 {
3010  int ierr, rank;
3011  MPI_Comm_rank( MPI_COMM_WORLD, &rank );
3012 
3013  // get elements
3014  const int num_elem = 4;
3015  iBase_EntityHandle array[num_elem];
3016  iBase_EntityHandle* ptr = array;
3017  int junk1 = num_elem, n = -1;
3018  iMesh_getEntities( imesh, part, iBase_FACE, iMesh_QUADRILATERAL, &ptr, &junk1, &n, &ierr );CHKERR;
3019  assert( ptr == array );
3020  assert( junk1 == num_elem );
3021  if( n != num_elem )
3022  {
3023  std::cerr << "Internal error at " __FILE__ ":" << __LINE__ << " (proc " << rank
3024  << "): Expected all parts to have " << num_elem << " elements. Found one with " << n << std::endl;
3025  return 1;
3026  }
3027 
3028  // get vertices
3029  iBase_EntityHandle adj_array[4 * num_elem];
3030  int junk2, junk3, offset_array[5];
3031  ptr = adj_array;
3032  junk1 = sizeof( adj_array ) / sizeof( adj_array[0] );
3033  junk2 = sizeof( offset_array ) / sizeof( offset_array[0] );
3034  int* ptr2 = offset_array;
3035  iMesh_getEntArrAdj( imesh, array, num_elem, iBase_VERTEX, &ptr, &junk1, &n, &ptr2, &junk2, &junk3, &ierr );CHKERR;
3036  assert( ptr == adj_array );
3037  assert( ptr2 == offset_array );
3038  assert( junk1 == sizeof( adj_array ) / sizeof( adj_array[0] ) );
3039  assert( junk2 == sizeof( offset_array ) / sizeof( offset_array[0] ) );
3040  assert( n == 4 * num_elem );
3041  assert( offset_array[0] == 0 );
3042  for( int i = 1; i < junk3; ++i )
3043  assert( offset_array[i] - offset_array[i - 1] == 4 );
3044 
3045  // find center vertex
3047  bool all_match;
3048  for( int i = 0; i < 4; ++i )
3049  {
3050  vtx = adj_array[i];
3051  all_match = true;
3052  for( int j = 1; j < 4; ++j )
3053  {
3054  iBase_EntityHandle* mvtx = adj_array + 4 * j;
3055  int k;
3056  for( k = 0; k < 4; ++k )
3057  if( mvtx[k] == vtx ) break;
3058  if( k == 4 ) all_match = false;
3059  }
3060  if( all_match ) break;
3061  }
3062  assert( all_match );
3063 
3064  // get center vertex coordinates
3065  double x, y, z;
3066  iMesh_getVtxCoord( imesh, vtx, &x, &y, &z, &ierr );CHKERR;
3067  assert( 0.0 == z );
3068  const int xi = ( (int)round( x ) - 1 ) / 2;
3069  const int yi = ( (int)round( y ) - 1 ) / 2;
3070  assert( xi >= 0 );
3071  assert( yi >= 0 );
3072  assert( fabs( x - 2 * xi - 1 ) < 1e-12 );
3073  assert( fabs( y - 2 * yi - 1 ) < 1e-12 );
3074 
3075  id = 2 * xi + yi;
3076  return 0;
3077 }