Mesh Oriented datABase  (version 5.6.0)
An array-based unstructured mesh library
ReadCCMIO.cpp
Go to the documentation of this file.
1 #include <cstdlib> // For exit()
2 #include <vector>
3 #include <map>
4 #include <iostream>
5 #include <string>
6 #include <algorithm>
7 
8 #include "moab/CN.hpp"
9 #include "moab/Range.hpp"
10 #include "moab/Interface.hpp"
11 #include "MBTagConventions.hpp"
12 #include "Internals.hpp"
13 #include "moab/ReadUtilIface.hpp"
14 #include "moab/FileOptions.hpp"
15 #include "ReadCCMIO.hpp"
16 #include "moab/MeshTopoUtil.hpp"
17 
18 #include "ccmio.h"
19 
20 /*
21  * CCMIO file structure
22  *
23  * Root
24  * State(kCCMIOState)
25  * Processor*
26  * Vertices
27  * ->ReadVerticesx, ReadMap
28  * Topology
29  * Boundary faces*(kCCMIOBoundaryFaces)
30  * ->ReadFaces, ReadFaceCells, ReadMap
31  * Internal faces(kCCMIOInternalFaces)
32  * Cells (kCCMIOCells)
33  * ->ReadCells (mapID), ReadMap, ReadCells (cellTypes)
34  * Solution
35  * Phase
36  * Field
37  * FieldData
38  * Problem(kCCMIOProblemDescription)
39  * CellType* (kCCMIOCellType)
40  * Index (GetEntityIndex), MaterialId(ReadOpti), MaterialType(ReadOptstr),
41  * PorosityId(ReadOpti), SpinId(ReadOpti), GroupId(ReadOpti)
42  *
43  * MaterialType (CCMIOReadOptstr in readexample)
44  * constants (see readexample)
45  * lagrangian data (CCMIOReadLagrangianData)
46  * vertices label (CCMIOEntityDescription)
47  * restart info: char solver[], iteratoins, time, char timeUnits[], angle
48  * (CCMIOReadRestartInfo, kCCMIORestartData), reference data?
49  * phase:
50  * field: char name[], dims, CCMIODataType datatype, char units[]
51  * dims = kCCMIOScalar (CCMIOReadFieldDataf),
52  * kCCMIOVector (CCMIOReadMultiDimensionalFieldData),
53  * kCCMIOTensor
54  * MonitoringSets: num, name (CellSet, VertexSet, BoundarySet, BlockSet, SplineSet, CoupleSet)
55  * CCMIOGetProstarSet, CCMIOReadOpt1i,
56  */
57 
59 {
68  kCellType
69 };
70 
71 namespace moab
72 {
73 
74 static int const kNValues = 10; // Number of values of each element to print
75 static char const kDefaultState[] = "default";
76 static char const kUnitsName[] = "Units";
77 static int const kVertOffset = 2;
78 static int const kCellInc = 4;
79 
80 #define CHK_SET_CCMERR( ccm_err_code, ccm_err_msg ) \
81  { \
82  if( kCCMIONoErr != ( ccm_err_code ) && kCCMIONoFileErr != ( ccm_err_code ) && \
83  kCCMIONoNodeErr != ( ccm_err_code ) ) \
84  MB_SET_ERR( MB_FAILURE, ccm_err_msg ); \
85  }
86 
88 {
89  return new ReadCCMIO( iface );
90 }
91 
93  : mMaterialIdTag( 0 ), mMaterialTypeTag( 0 ), mRadiationTag( 0 ), mPorosityIdTag( 0 ), mSpinIdTag( 0 ),
94  mGroupIdTag( 0 ), mColorIdxTag( 0 ), mProcessorIdTag( 0 ), mLightMaterialTag( 0 ), mFreeSurfaceMaterialTag( 0 ),
95  mThicknessTag( 0 ), mProstarRegionNumberTag( 0 ), mBoundaryTypeTag( 0 ), mCreatingProgramTag( 0 ), mbImpl( impl ),
96  hasSolution( false )
97 {
98  assert( impl != NULL );
99 
101 
102  // Initialize in case tag_get_handle fails below
103  mMaterialSetTag = 0;
104  mDirichletSetTag = 0;
105  mNeumannSetTag = 0;
106  mHasMidNodesTag = 0;
107  mGlobalIdTag = 0;
108  mNameTag = 0;
109 
110  //! Get and cache predefined tag handles
111  const int negone = -1;
113  MB_TAG_CREAT | MB_TAG_SPARSE, &negone );
114  MB_CHK_SET_ERR_RET( result, "Failed to get MATERIAL_SET tag" );
115 
117  MB_TAG_CREAT | MB_TAG_SPARSE, &negone );
118  MB_CHK_SET_ERR_RET( result, "Failed to get DIRICHLET_SET tag" );
119 
121  MB_TAG_CREAT | MB_TAG_SPARSE, &negone );
122  MB_CHK_SET_ERR_RET( result, "Failed to get NEUMANN_SET tag" );
123 
124  const int negonearr[] = { -1, -1, -1, -1 };
126  MB_TAG_CREAT | MB_TAG_SPARSE, negonearr );
127  MB_CHK_SET_ERR_RET( result, "Failed to get HAS_MID_NODES tag" );
128 
129  mGlobalIdTag = impl->globalId_tag();
130 
131  result =
133  MB_CHK_SET_ERR_RET( result, "Failed to get NAME tag" );
134 }
135 
137 {
139 }
140 
141 ErrorCode ReadCCMIO::load_file( const char* file_name,
142  const EntityHandle* file_set,
143  const FileOptions& /* opts */,
144  const ReaderIface::SubsetList* subset_list,
145  const Tag* /* file_id_tag */ )
146 {
147  CCMIOID rootID, problemID, stateID, processorID, verticesID, topologyID, solutionID;
148  CCMIOError error = kCCMIONoErr;
149 
150  if( subset_list )
151  {
152  MB_SET_ERR( MB_UNSUPPORTED_OPERATION, "Reading subset of files not supported for CCMOI data" );
153  }
154 
155  CCMIOOpenFile( &error, file_name, kCCMIORead, &rootID );
156  CHK_SET_CCMERR( error, "Problem opening file" );
157 
158  // Get the file state
159  MB_CHK_SET_ERR( get_state( rootID, problemID, stateID ), "Failed to get state" );
160 
161  // Get processors
162  std::vector< CCMIOSize_t > procs;
163  bool has_solution = false;
164  MB_CHK_SET_ERR( get_processors( stateID, processorID, verticesID, topologyID, solutionID, procs, has_solution ),
165  "Failed to get processors" );
166 
167  std::vector< CCMIOSize_t >::iterator vit;
168  Range new_ents, *new_ents_ptr = NULL;
169  if( file_set ) new_ents_ptr = &new_ents;
170 
171  for( vit = procs.begin(); vit != procs.end(); ++vit )
172  {
173  MB_CHK_SET_ERR( read_processor( stateID, problemID, processorID, verticesID, topologyID, *vit, new_ents_ptr ),
174  "Failed to read processors" );
175  }
176 
177  // Load some meta-data
178  MB_CHK_SET_ERR( load_metadata( rootID, problemID, stateID, processorID, file_set ),
179  "Failed to load some meta-data" );
180 
181  // Now, put all this into the file set, if there is one
182  if( file_set )
183  {
184  MB_CHK_SET_ERR( mbImpl->add_entities( *file_set, new_ents ), "Failed to add new entities to file set" );
185  }
186 
187  return rval;
188 }
189 
190 ErrorCode ReadCCMIO::get_state( CCMIOID rootID, CCMIOID& problemID, CCMIOID& stateID )
191 {
192  CCMIOError error = kCCMIONoErr;
193 
194  // First try default
195  CCMIOGetState( &error, rootID, "default", &problemID, &stateID );
196  if( kCCMIONoErr != error )
197  {
198  CCMIOSize_t i = CCMIOSIZEC( 0 );
199  CCMIOError tmp_error = kCCMIONoErr;
200  CCMIONextEntity( &tmp_error, rootID, kCCMIOState, &i, &stateID );
201  if( kCCMIONoErr == tmp_error ) CCMIONextEntity( &error, rootID, kCCMIOProblemDescription, &i, &problemID );
202  }
203  CHK_SET_CCMERR( error, "Couldn't find state" );
204 
205  return MB_SUCCESS;
206 }
207 
209  CCMIOID problemID,
210  CCMIOID /* stateID */,
211  CCMIOID processorID,
212  const EntityHandle* file_set )
213 {
214  // Read the simulation title.
215  CCMIOError error = kCCMIONoErr;
216  ErrorCode rval = MB_SUCCESS;
217  CCMIONode rootNode;
218  if( kCCMIONoErr == CCMIOGetEntityNode( &error, rootID, &rootNode ) )
219  {
220  char* name = NULL;
221  CCMIOGetTitle( &error, rootNode, &name );
222 
223  if( NULL != name && strlen( name ) != 0 )
224  {
225  // Make a tag for it and tag the read set
226  Tag simname;
227  MB_CHK_SET_ERR( mbImpl->tag_get_handle( "Title", strlen( name ), MB_TYPE_OPAQUE, simname,
229  "Simulation name tag not found or created" );
230  EntityHandle set = file_set ? *file_set : 0;
231  MB_CHK_SET_ERR( mbImpl->tag_set_data( simname, &set, 1, name ), "Problem setting simulation name tag" );
232  }
233  if( name ) free( name );
234  }
235 
236  // Creating program
237  EntityHandle dumh = ( file_set ? *file_set : 0 );
238  MB_CHK_SET_ERR( get_str_option( "CreatingProgram", dumh, mCreatingProgramTag, processorID ),
239  "Trouble getting CreatingProgram tag" );
240 
241  MB_CHK_SET_ERR( load_matset_data( problemID ), "Failure loading matset data" );
242 
243  MB_CHK_SET_ERR( load_neuset_data( problemID ), "Failure loading neuset data" );
244 
245  return rval;
246 }
247 
249 {
250  // Make sure there are matsets
251  if( newMatsets.empty() ) return MB_SUCCESS;
252 
253  // ... walk through each cell type
254  CCMIOSize_t i = CCMIOSIZEC( 0 );
255  CCMIOID next;
256  CCMIOError error = kCCMIONoErr;
257 
258  while( CCMIONextEntity( NULL, problemID, kCCMIOCellType, &i, &next ) == kCCMIONoErr )
259  {
260  // Get index, corresponding set, and label with material set tag
261  int mindex;
262  CCMIOGetEntityIndex( &error, next, &mindex );
263  std::map< int, EntityHandle >::iterator mit = newMatsets.find( mindex );
264  if( mit == newMatsets.end() )
265  // No actual faces for this matset; continue to next
266  continue;
267 
268  EntityHandle dum_ent = mit->second;
269  MB_CHK_SET_ERR( mbImpl->tag_set_data( mMaterialSetTag, &dum_ent, 1, &mindex ),
270  "Trouble setting material set tag" );
271 
272  // Set name
273  CCMIOSize_t len;
274  CCMIOEntityLabel( &error, next, &len, NULL );
275  std::vector< char > opt_string2( GETINT32( len ) + 1, '\0' );
276  CCMIOEntityLabel( &error, next, NULL, &opt_string2[0] );
277  if( opt_string2.size() >= NAME_TAG_SIZE )
278  opt_string2[NAME_TAG_SIZE - 1] = '\0';
279  else
280  ( opt_string2.resize( NAME_TAG_SIZE, '\0' ) );
281  MB_CHK_SET_ERR( mbImpl->tag_set_data( mNameTag, &dum_ent, 1, &opt_string2[0] ),
282  "Trouble setting name tag for material set" );
283 
284  // Material id
285  MB_CHK_SET_ERR( get_int_option( "MaterialId", dum_ent, mMaterialIdTag, next ),
286  "Trouble getting MaterialId tag" );
287 
288  MB_CHK_SET_ERR( get_str_option( "MaterialType", dum_ent, mMaterialTypeTag, next ),
289  "Trouble getting MaterialType tag" );
290 
291  MB_CHK_SET_ERR( get_int_option( "Radiation", dum_ent, mRadiationTag, next ),
292  "Trouble getting Radiation option" );
293 
294  MB_CHK_SET_ERR( get_int_option( "PorosityId", dum_ent, mPorosityIdTag, next ),
295  "Trouble getting PorosityId option" );
296 
297  MB_CHK_SET_ERR( get_int_option( "SpinId", dum_ent, mSpinIdTag, next ), "Trouble getting SpinId option" );
298 
299  MB_CHK_SET_ERR( get_int_option( "GroupId", dum_ent, mGroupIdTag, next ), "Trouble getting GroupId option" );
300 
301  MB_CHK_SET_ERR( get_int_option( "ColorIdx", dum_ent, mColorIdxTag, next ), "Trouble getting ColorIdx option" );
302 
303  MB_CHK_SET_ERR( get_int_option( "ProcessorId", dum_ent, mProcessorIdTag, next ),
304  "Trouble getting ProcessorId option" );
305 
306  MB_CHK_SET_ERR( get_int_option( "LightMaterial", dum_ent, mLightMaterialTag, next ),
307  "Trouble getting LightMaterial option" );
308 
309  MB_CHK_SET_ERR( get_int_option( "FreeSurfaceMaterial", dum_ent, mFreeSurfaceMaterialTag, next ),
310  "Trouble getting FreeSurfaceMaterial option" );
311 
312  MB_CHK_SET_ERR( get_dbl_option( "Thickness", dum_ent, mThicknessTag, next ),
313  "Trouble getting Thickness option" );
314  }
315 
316  return MB_SUCCESS;
317 }
318 
319 ErrorCode ReadCCMIO::get_int_option( const char* opt_str, EntityHandle seth, Tag& tag, CCMIOID node )
320 {
321  int idum;
322  ErrorCode rval;
323  if( kCCMIONoErr == CCMIOReadOpti( NULL, node, opt_str, &idum ) )
324  {
325  if( !tag )
326  {
328  "Failed to get tag handle" );
329  }
330 
331  MB_CHK_SET_ERR( mbImpl->tag_set_data( tag, &seth, 1, &idum ), "Failed to set tag data" );
332  }
333 
334  return MB_SUCCESS;
335 }
336 
337 ErrorCode ReadCCMIO::get_dbl_option( const char* opt_str, EntityHandle seth, Tag& tag, CCMIOID node )
338 {
339  float fdum;
340  if( kCCMIONoErr == CCMIOReadOptf( NULL, node, opt_str, &fdum ) )
341  {
342  ErrorCode rval;
343  if( !tag )
344  {
346  "Failed to get tag handle" );
347  }
348 
349  double dum_dbl = fdum;
350  MB_CHK_SET_ERR( mbImpl->tag_set_data( tag, &seth, 1, &dum_dbl ), "Failed to set tag data" );
351  }
352 
353  return MB_SUCCESS;
354 }
355 
357  EntityHandle seth,
358  Tag& tag,
359  CCMIOID node,
360  const char* other_tag_name )
361 {
362  int len;
363  CCMIOError error = kCCMIONoErr;
364  std::vector< char > opt_string;
365  if( kCCMIONoErr != CCMIOReadOptstr( NULL, node, opt_str, &len, NULL ) ) return MB_SUCCESS;
366 
367  opt_string.resize( len );
368  CCMIOReadOptstr( &error, node, opt_str, &len, &opt_string[0] );
369  ErrorCode rval = MB_SUCCESS;
370  if( !tag )
371  {
372  MB_CHK_SET_ERR( mbImpl->tag_get_handle( other_tag_name ? other_tag_name : opt_str, NAME_TAG_SIZE,
374  "Failed to get tag handle" );
375  }
376 
377  if( opt_string.size() > NAME_TAG_SIZE )
378  opt_string[NAME_TAG_SIZE - 1] = '\0';
379  else
380  ( opt_string.resize( NAME_TAG_SIZE, '\0' ) );
381 
382  MB_CHK_SET_ERR( mbImpl->tag_set_data( tag, &seth, 1, &opt_string[0] ), "Failed to set tag data" );
383 
384  return MB_SUCCESS;
385 }
386 
388 {
389  CCMIOSize_t i = CCMIOSIZEC( 0 );
390  CCMIOID next;
391 
392  // Make sure there are matsets
393  if( newNeusets.empty() ) return MB_SUCCESS;
394 
395  while( CCMIONextEntity( NULL, problemID, kCCMIOBoundaryRegion, &i, &next ) == kCCMIONoErr )
396  {
397  // Get index, corresponding set, and label with neumann set tag
398  int mindex;
399  CCMIOError error = kCCMIONoErr;
400  CCMIOGetEntityIndex( &error, next, &mindex );
401  std::map< int, EntityHandle >::iterator mit = newNeusets.find( mindex );
402  if( mit == newNeusets.end() )
403  // No actual faces for this neuset; continue to next
404  continue;
405 
406  EntityHandle dum_ent = mit->second;
407  MB_CHK_SET_ERR( mbImpl->tag_set_data( mNeumannSetTag, &dum_ent, 1, &mindex ),
408  "Trouble setting neumann set tag" );
409 
410  // Set name
411  MB_CHK_SET_ERR( get_str_option( "BoundaryName", dum_ent, mNameTag, next, NAME_TAG_NAME ),
412  "Trouble creating BoundaryName tag" );
413 
414  // BoundaryType
415  MB_CHK_SET_ERR( get_str_option( "BoundaryType", dum_ent, mBoundaryTypeTag, next ),
416  "Trouble creating BoundaryType tag" );
417 
418  // ProstarRegionNumber
419  MB_CHK_SET_ERR( get_int_option( "ProstarRegionNumber", dum_ent, mProstarRegionNumberTag, next ),
420  "Trouble creating ProstarRegionNumber tag" );
421  }
422 
423  return MB_SUCCESS;
424 }
425 
426 ErrorCode ReadCCMIO::read_processor( CCMIOID /* stateID */,
427  CCMIOID problemID,
428  CCMIOID processorID,
429  CCMIOID verticesID,
430  CCMIOID topologyID,
431  CCMIOSize_t proc,
432  Range* new_ents )
433 {
434  ErrorCode rval;
435 
436  // vert_map fields: s: none, i: gid, ul: vert handle, r: none
437  // TupleList vert_map(0, 1, 1, 0, 0);
438  TupleList vert_map;
439  MB_CHK_SET_ERR( read_vertices( proc, processorID, verticesID, topologyID, new_ents, vert_map ),
440  "Failed to read vertices" );
441 
442  MB_CHK_SET_ERR( read_cells( proc, problemID, verticesID, topologyID, vert_map, new_ents ), "Failed to read cells" );
443 
444  return rval;
445 }
446 
447 ErrorCode ReadCCMIO::read_cells( CCMIOSize_t /* proc */,
448  CCMIOID problemID,
449  CCMIOID /* verticesID */,
450  CCMIOID topologyID,
451  TupleList& vert_map,
452  Range* new_ents )
453 {
454  // Read the faces.
455  // face_map fields: s:forward/reverse, i: cell id, ul: face handle, r: none
456  ErrorCode rval;
457 #ifdef TUPLE_LIST
458  TupleList face_map( 1, 1, 1, 0, 0 );
459 #else
460  TupleList face_map;
461  SenseList sense_map;
462 #endif
463  MB_CHK_SET_ERR( read_all_faces( topologyID, vert_map, face_map,
464 #ifndef TUPLE_LIST
465  sense_map,
466 #endif
467  new_ents ),
468  "Failed to read all cells" );
469 
470  // Read the cell topology types, if any exist in the file
471  std::map< int, int > cell_topo_types;
472  MB_CHK_SET_ERR( read_topology_types( topologyID, cell_topo_types ), "Problem reading cell topo types" );
473 
474  // Now construct the cells; sort the face map by cell ids first
475 #ifdef TUPLE_LIST
476  MB_CHK_SET_ERR( face_map.sort( 1 ), "Couldn't sort face map by cell id" );
477 #endif
478  std::vector< EntityHandle > new_cells;
479  MB_CHK_SET_ERR( construct_cells( face_map,
480 #ifndef TUPLE_LIST
481  sense_map,
482 #endif
483  vert_map, cell_topo_types, new_cells ),
484  "Failed to construct cells" );
485  if( new_ents )
486  {
487  Range::iterator rit = new_ents->end();
488  std::vector< EntityHandle >::reverse_iterator vit;
489  for( vit = new_cells.rbegin(); vit != new_cells.rend(); ++vit )
490  rit = new_ents->insert( rit, *vit );
491  }
492 
493  MB_CHK_SET_ERR( read_gids_and_types( problemID, topologyID, new_cells ), "Failed to read gids and types" );
494 
495  return MB_SUCCESS;
496 }
497 
498 ErrorCode ReadCCMIO::read_topology_types( CCMIOID& topologyID, std::map< int, int >& cell_topo_types )
499 {
500  CCMIOError error = kCCMIONoErr;
501  CCMIOID cellID, mapID;
502  CCMIOSize_t ncells;
503  CCMIOGetEntity( &error, topologyID, kCCMIOCells, 0, &cellID );
504  CCMIOEntitySize( &error, cellID, &ncells, NULL );
505  int num_cells = GETINT32( ncells );
506 
507  // First, do a dummy read to see if we even have topo types in this mesh
508  int dum_int;
509  CCMIOReadOpt1i( &error, cellID, "CellTopologyType", &dum_int, CCMIOINDEXC( kCCMIOStart ),
510  CCMIOINDEXC( kCCMIOStart ) + 1 );
511  if( kCCMIONoErr != error ) return MB_SUCCESS;
512 
513  // OK, we have topo types; first get the map node
514  std::vector< int > dum_ints( num_cells );
515  CCMIOReadCells( &error, cellID, &mapID, &dum_ints[0], CCMIOINDEXC( kCCMIOStart ), CCMIOINDEXC( kCCMIOStart ) + 1 );
516  CHK_SET_CCMERR( error, "Failed to get the map node" );
517 
518  // Now read the map
519  CCMIOReadMap( &error, mapID, &dum_ints[0], CCMIOINDEXC( kCCMIOStart ), CCMIOINDEXC( kCCMIOEnd ) );
520  CHK_SET_CCMERR( error, "Failed to get cell ids" );
521  int i;
522  for( i = 0; i < num_cells; i++ )
523  cell_topo_types[dum_ints[i]] = 0;
524 
525  // Now read the cell topo types for real, reusing cell_topo_types
526  std::vector< int > topo_types( num_cells );
527  CCMIOReadOpt1i( &error, cellID, "CellTopologyType", &topo_types[0], CCMIOINDEXC( kCCMIOStart ),
528  CCMIOINDEXC( kCCMIOEnd ) );
529  CHK_SET_CCMERR( error, "Failed to get cell topo types" );
530  for( i = 0; i < num_cells; i++ )
531  cell_topo_types[dum_ints[i]] = topo_types[i];
532 
533  return MB_SUCCESS;
534 }
535 
536 ErrorCode ReadCCMIO::read_gids_and_types( CCMIOID /* problemID */,
537  CCMIOID topologyID,
538  std::vector< EntityHandle >& cells )
539 {
540  // Get the cells entity and number of cells
541  CCMIOSize_t dum_cells;
542  int num_cells;
543  CCMIOError error = kCCMIONoErr;
544  CCMIOID cellsID, mapID;
545  CCMIOGetEntity( &error, topologyID, kCCMIOCells, 0, &cellsID );
546  CCMIOEntitySize( &error, cellsID, &dum_cells, NULL );
547  num_cells = GETINT32( dum_cells );
548 
549  // Check the number of cells against how many are in the cell array
550  if( num_cells != (int)cells.size() ) MB_SET_ERR( MB_FAILURE, "Number of cells doesn't agree" );
551 
552  // Read the gid map and set global ids
553  std::vector< int > cell_gids( num_cells );
554  CCMIOReadCells( &error, cellsID, &mapID, NULL, CCMIOINDEXC( kCCMIOStart ), CCMIOINDEXC( kCCMIOEnd ) );
555  CHK_SET_CCMERR( error, "Couldn't read cells" );
556  CCMIOReadMap( &error, mapID, &cell_gids[0], CCMIOINDEXC( kCCMIOStart ), CCMIOINDEXC( kCCMIOEnd ) );
557  CHK_SET_CCMERR( error, "Couldn't read cell id map" );
558 
559  MB_CHK_SET_ERR( mbImpl->tag_set_data( mGlobalIdTag, &cells[0], cells.size(), &cell_gids[0] ),
560  "Couldn't set gids tag" );
561 
562  // Now read cell material types; reuse cell_gids
563  CCMIOReadCells( &error, cellsID, NULL, &cell_gids[0], CCMIOINDEXC( kCCMIOStart ), CCMIOINDEXC( kCCMIOEnd ) );
564  CHK_SET_CCMERR( error, "Trouble reading cell types" );
565 
566  // Create the matsets
567  std::map< int, Range > matset_ents;
568  for( int i = 0; i < num_cells; i++ )
569  matset_ents[cell_gids[i]].insert( cells[i] );
570 
571  for( std::map< int, Range >::iterator mit = matset_ents.begin(); mit != matset_ents.end(); ++mit )
572  {
573  EntityHandle matset;
574  MB_CHK_SET_ERR( mbImpl->create_meshset( MESHSET_SET, matset ), "Couldn't create material set" );
575  newMatsets[mit->first] = matset;
576 
577  MB_CHK_SET_ERR( mbImpl->add_entities( matset, mit->second ), "Couldn't add entities to material set" );
578  }
579 
580  return MB_SUCCESS;
581 }
582 
584 #ifndef TUPLE_LIST
585  SenseList& sense_map,
586 #endif
587  TupleList& /* vert_map */,
588  std::map< int, int >& cell_topo_types,
589  std::vector< EntityHandle >& new_cells )
590 {
591  std::vector< EntityHandle > facehs;
592  std::vector< int > senses;
593  EntityHandle cell;
594  ErrorCode tmp_rval, rval = MB_SUCCESS;
595  EntityType this_type = MBMAXTYPE;
596  bool has_mid_nodes = false;
597 #ifdef TUPLE_LIST
598  unsigned int i = 0;
599  while( i < face_map.n )
600  {
601  // Pull out face handles bounding the same cell
602  facehs.clear();
603  int this_id = face_map.get_int( i );
604  unsigned int inext = i;
605  while( face_map.get_int( inext ) == this_id && inext <= face_map.n )
606  {
607  inext++;
608  EntityHandle face = face_map.get_ulong( inext );
609  facehs.push_back( face );
610  senses.push_back( face_map.get_short( inext ) );
611  }
612  this_type = MBMAXTYPE;
613  has_mid_nodes = false;
614 #else
615  std::map< int, std::vector< EntityHandle > >::iterator fmit;
616  std::map< int, std::vector< int > >::iterator smit;
617  std::map< int, int >::iterator typeit;
618  for( fmit = face_map.begin(), smit = sense_map.begin(); fmit != face_map.end(); ++fmit, ++smit )
619  {
620  // Pull out face handles bounding the same cell
621  facehs.clear();
622  int this_id = ( *fmit ).first;
623  facehs = ( *fmit ).second;
624  senses.clear();
625  senses = ( *smit ).second;
626  typeit = cell_topo_types.find( this_id );
627  if( typeit != cell_topo_types.end() )
628  {
629  rval = ccmio_to_moab_type( typeit->second, this_type, has_mid_nodes );
630  }
631  else
632  {
633  this_type = MBMAXTYPE;
634  has_mid_nodes = false;
635  }
636 #endif
637  tmp_rval = create_cell_from_faces( facehs, senses, this_type, has_mid_nodes, cell );
638  if( MB_SUCCESS != tmp_rval )
639  rval = tmp_rval;
640  else
641  {
642  new_cells.push_back( cell );
643  // Tag cell with global id
644  tmp_rval = mbImpl->tag_set_data( mGlobalIdTag, &cell, 1, &this_id );
645  if( MB_SUCCESS != tmp_rval ) rval = tmp_rval;
646  }
647  }
648 
649  return rval;
650 }
651 
652 ErrorCode ReadCCMIO::ccmio_to_moab_type( int ccm_type, EntityType& moab_type, bool& has_mid_nodes )
653 {
654  switch( ccm_type )
655  {
656  case 1:
657  moab_type = MBVERTEX;
658  break;
659  case 2:
660  case 28:
661  moab_type = MBEDGE;
662  break;
663  case 29:
664  moab_type = MBMAXTYPE;
665  break;
666  case 3:
667  case 4:
668  moab_type = MBQUAD;
669  break;
670  case 11:
671  case 21:
672  moab_type = MBHEX;
673  break;
674  case 12:
675  case 22:
676  moab_type = MBPRISM;
677  break;
678  case 13:
679  case 23:
680  moab_type = MBTET;
681  break;
682  case 14:
683  case 24:
684  moab_type = MBPYRAMID;
685  break;
686  case 255:
687  moab_type = MBPOLYHEDRON;
688  break;
689  default:
690  moab_type = MBMAXTYPE;
691  }
692 
693  switch( ccm_type )
694  {
695  case 28:
696  case 4:
697  case 21:
698  case 22:
699  case 23:
700  case 24:
701  has_mid_nodes = true;
702  break;
703  default:
704  break;
705  }
706 
707  return MB_SUCCESS;
708 }
709 
710 ErrorCode ReadCCMIO::create_cell_from_faces( std::vector< EntityHandle >& facehs,
711  std::vector< int >& senses,
712  EntityType this_type,
713  bool /* has_mid_nodes */,
714  EntityHandle& cell )
715 {
716  ErrorCode rval;
717 
718  // Test up front to see if they're one type
719  EntityType face_type = mbImpl->type_from_handle( facehs[0] );
720  bool same_type = true;
721  for( std::vector< EntityHandle >::iterator vit = facehs.begin(); vit != facehs.end(); ++vit )
722  {
723  if( face_type != mbImpl->type_from_handle( *vit ) )
724  {
725  same_type = false;
726  break;
727  }
728  }
729 
730  std::vector< EntityHandle > verts;
731  EntityType input_type = this_type;
732  std::vector< EntityHandle > storage;
733  MeshTopoUtil mtu( mbImpl );
734 
735  // Preset this to maxtype, so we get an affirmative choice in loop
736  this_type = MBMAXTYPE;
737 
738  if( ( MBTET == input_type || MBMAXTYPE == input_type ) && same_type && face_type == MBTRI && facehs.size() == 4 )
739  {
740  // Try to get proper connectivity for tet
741 
742  // Get connectivity of first face, and reverse it if sense is forward, since
743  // base face always points into entity
744  MB_CHK_SET_ERR( mbImpl->get_connectivity( &facehs[0], 1, verts ), "Couldn't get connectivity" );
745  if( senses[0] > 0 ) std::reverse( verts.begin(), verts.end() );
746 
747  // Get the 4th vertex through the next tri
748  const EntityHandle* conn;
749  int conn_size;
750  MB_CHK_SET_ERR( mbImpl->get_connectivity( facehs[1], conn, conn_size, true, &storage ),
751  "Couldn't get connectivity" );
752  int i = 0;
753  while( std::find( verts.begin(), verts.end(), conn[i] ) != verts.end() && i < conn_size )
754  i++;
755 
756  // If i is not at the end of the verts, found the apex; otherwise fall back to polyhedron
757  if( conn_size != i )
758  {
759  this_type = MBTET;
760  verts.push_back( conn[i] );
761  }
762  }
763  else if( ( MBHEX == input_type || MBMAXTYPE == input_type ) && same_type && MBQUAD == face_type &&
764  facehs.size() == 6 )
765  {
766  // Build hex from quads
767  // Algorithm:
768  // - verts = vertices from 1st quad
769  // - Find quad q1 sharing verts[0] and verts[1]
770  // - Find quad q2 sharing other 2 verts in q1
771  // - Find v1 = opposite vert from verts[1] in q1 , v2 = opposite from verts[0]
772  // - Get i = offset of v1 in verts2 of q2, rotate verts2 by i
773  // - If verts2[(i + 1) % 4] != v2, flip verts2 by switching verts2[1] and verts2[3]
774  // - append verts2 to verts
775 
776  // Get the other vertices for this hex; need to find the quad with no common vertices
777  Range tmp_faces, tmp_verts;
778  // Get connectivity of first face, and reverse it if sense is forward, since
779  // base face always points into entity
780  MB_CHK_SET_ERR( mbImpl->get_connectivity( &facehs[0], 1, verts ), "Couldn't get connectivity" );
781  if( senses[0] > 0 ) std::reverse( verts.begin(), verts.end() );
782 
783  // Get q1, which shares 2 vertices with q0
784  std::copy( facehs.begin(), facehs.end(), range_inserter( tmp_faces ) );
785  rval = mbImpl->get_adjacencies( &verts[0], 2, 2, false, tmp_faces );
786  if( MB_SUCCESS != rval || tmp_faces.size() != 2 ) MB_SET_ERR( MB_FAILURE, "Couldn't get adj face" );
787  tmp_faces.erase( facehs[0] );
788  EntityHandle q1 = *tmp_faces.begin();
789  // Get other 2 verts of q1
790  MB_CHK_SET_ERR( mbImpl->get_connectivity( &q1, 1, tmp_verts ), "Couldn't get adj verts" );
791  tmp_verts.erase( verts[0] );
792  tmp_verts.erase( verts[1] );
793  // Get q2
794  std::copy( facehs.begin(), facehs.end(), range_inserter( tmp_faces ) );
795  rval = mbImpl->get_adjacencies( tmp_verts, 2, false, tmp_faces );
796  if( MB_SUCCESS != rval || tmp_faces.size() != 2 ) MB_SET_ERR( MB_FAILURE, "Couldn't get adj face" );
797  tmp_faces.erase( q1 );
798  EntityHandle q2 = *tmp_faces.begin();
799  // Get verts in q2
800  MB_CHK_SET_ERR( mbImpl->get_connectivity( &q2, 1, storage ), "Couldn't get adj vertices" );
801 
802  // Get verts in q1 opposite from v[1] and v[0] in q0
803  EntityHandle v0 = 0, v1 = 0;
804  MB_CHK_SET_ERR( mtu.opposite_entity( q1, verts[1], v0 ), "Couldn't get the opposite side entity" );
805  MB_CHK_SET_ERR( mtu.opposite_entity( q1, verts[0], v1 ), "Couldn't get the opposite side entity" );
806  if( v0 && v1 )
807  {
808  // Offset of v0 in q2, then rotate and flip
809  unsigned int ioff = std::find( storage.begin(), storage.end(), v0 ) - storage.begin();
810  if( 4 == ioff ) MB_SET_ERR( MB_FAILURE, "Trouble finding offset" );
811 
812  if( storage[( ioff + 1 ) % 4] != v1 )
813  {
814  std::reverse( storage.begin(), storage.end() );
815  ioff = std::find( storage.begin(), storage.end(), v0 ) - storage.begin();
816  }
817  if( 0 != ioff ) std::rotate( storage.begin(), storage.begin() + ioff, storage.end() );
818 
819  // Copy into verts, and make hex
820  std::copy( storage.begin(), storage.end(), std::back_inserter( verts ) );
821  this_type = MBHEX;
822  }
823  }
824 
825  if( MBMAXTYPE == this_type && facehs.size() == 5 )
826  {
827  // Some preliminaries
828  std::vector< EntityHandle > tris, quads;
829  for( unsigned int i = 0; i < 5; i++ )
830  {
831  if( MBTRI == mbImpl->type_from_handle( facehs[i] ) )
832  tris.push_back( facehs[i] );
833  else if( MBQUAD == mbImpl->type_from_handle( facehs[i] ) )
834  quads.push_back( facehs[i] );
835  }
836 
837  // Check for prisms
838  if( 2 == tris.size() && 3 == quads.size() )
839  {
840  // OK, we have the right number of tris and quads; try to find the proper verts
841 
842  // Get connectivity of first tri, and reverse if necessary
843  int index = std::find( facehs.begin(), facehs.end(), tris[0] ) - facehs.begin();
844  MB_CHK_SET_ERR( mbImpl->get_connectivity( &tris[0], 1, verts ), "Couldn't get connectivity" );
845  if( senses[index] > 0 ) std::reverse( verts.begin(), verts.end() );
846 
847  // Now align vertices of other tri, through a quad, similar to how we did hexes
848  // Get q1, which shares 2 vertices with t0
849  Range tmp_faces, tmp_verts;
850  std::copy( facehs.begin(), facehs.end(), range_inserter( tmp_faces ) );
851  rval = mbImpl->get_adjacencies( &verts[0], 2, 2, false, tmp_faces );
852  if( MB_SUCCESS != rval || tmp_faces.size() != 2 ) MB_SET_ERR( MB_FAILURE, "Couldn't get adj face" );
853  tmp_faces.erase( tris[0] );
854  EntityHandle q1 = *tmp_faces.begin();
855  // Get verts in q1
856  MB_CHK_SET_ERR( mbImpl->get_connectivity( &q1, 1, storage ), "Couldn't get adj vertices" );
857 
858  // Get verts in q1 opposite from v[1] and v[0] in q0
859  EntityHandle v0 = 0, v1 = 0;
860  MB_CHK_SET_ERR( mtu.opposite_entity( q1, verts[1], v0 ), "Couldn't get the opposite side entity" );
861  MB_CHK_SET_ERR( mtu.opposite_entity( q1, verts[0], v1 ), "Couldn't get the opposite side entity" );
862  if( v0 && v1 )
863  {
864  // Offset of v0 in t2, then rotate and flip
865  storage.clear();
866  MB_CHK_SET_ERR( mbImpl->get_connectivity( &tris[1], 1, storage ), "Couldn't get connectivity" );
867 
868  index = std::find( facehs.begin(), facehs.end(), tris[1] ) - facehs.begin();
869  if( senses[index] < 0 ) std::reverse( storage.begin(), storage.end() );
870  unsigned int ioff = std::find( storage.begin(), storage.end(), v0 ) - storage.begin();
871  if( 3 == ioff ) MB_SET_ERR( MB_FAILURE, "Trouble finding offset" );
872  for( unsigned int i = 0; i < 3; i++ )
873  verts.push_back( storage[( ioff + i ) % 3] );
874 
875  this_type = MBPRISM;
876  }
877  }
878  else if( tris.size() == 4 && quads.size() == 1 )
879  {
880  // Check for pyramid
881  // Get connectivity of first tri, and reverse if necessary
882  int index = std::find( facehs.begin(), facehs.end(), quads[0] ) - facehs.begin();
883  MB_CHK_SET_ERR( mbImpl->get_connectivity( &quads[0], 1, verts ), "Couldn't get connectivity" );
884  if( senses[index] > 0 ) std::reverse( verts.begin(), verts.end() );
885 
886  // Get apex node
887  MB_CHK_SET_ERR( mbImpl->get_connectivity( &tris[0], 1, storage ), "Couldn't get connectivity" );
888  for( unsigned int i = 0; i < 3; i++ )
889  {
890  if( std::find( verts.begin(), verts.end(), storage[i] ) == verts.end() )
891  {
892  verts.push_back( storage[i] );
893  break;
894  }
895  }
896 
897  if( 5 == verts.size() ) this_type = MBPYRAMID;
898  }
899  else
900  {
901  // Dummy else clause to stop in debugger
902  this_type = MBMAXTYPE;
903  }
904  }
905 
906  if( MBMAXTYPE != input_type && input_type != this_type && this_type != MBMAXTYPE )
907  std::cerr << "Warning: types disagree (cell_topo_type = " << CN::EntityTypeName( input_type )
908  << ", faces indicate type " << CN::EntityTypeName( this_type ) << std::endl;
909 
910  if( MBMAXTYPE != input_type && this_type == MBMAXTYPE && input_type != MBPOLYHEDRON )
911  std::cerr << "Warning: couldn't find proper connectivity for specified topo_type = "
912  << CN::EntityTypeName( input_type ) << std::endl;
913 
914  // Now make the element; if we fell back to polyhedron, use faces, otherwise use verts
915  if( MBPOLYHEDRON == this_type || MBMAXTYPE == this_type )
916  {
917  MB_CHK_SET_ERR( mbImpl->create_element( MBPOLYHEDRON, &facehs[0], facehs.size(), cell ),
918  "create_element failed" );
919  }
920  else
921  {
922  MB_CHK_SET_ERR( mbImpl->create_element( this_type, &verts[0], verts.size(), cell ), "create_element failed" );
923  }
924 
925  return MB_SUCCESS;
926 }
927 
929  TupleList& vert_map,
930  TupleList& face_map,
931 #ifndef TUPLE_LIST
932  SenseList& sense_map,
933 #endif
934  Range* new_faces )
935 {
936  CCMIOSize_t index;
937  CCMIOID faceID;
938  ErrorCode rval;
939  CCMIOError error = kCCMIONoErr;
940 
941  // Get total # internal/bdy faces, size the face map accordingly
942 #ifdef TUPLE_LIST
943  index = CCMIOSIZEC( 0 );
944  int nbdy_faces = 0;
945  CCMIOSize_t nf;
946  error = kCCMIONoErr;
947  while( kCCMIONoErr == CCMIONextEntity( NULL, topologyID, kCCMIOBoundaryFaces, &index, &faceID ) )
948  {
949  CCMIOEntitySize( &error, faceID, &nf, NULL );
950  nbdy_faces += nf;
951  }
952  CCMIOGetEntity( &error, topologyID, kCCMIOInternalFaces, 0, &faceID );
953  CCMIOEntitySize( &error, faceID, &nf, NULL );
954 
955  int nint_faces = nf;
956  face_map.resize( 2 * nint_faces + nbdy_faces );
957 #endif
958 
959  // Get multiple blocks of bdy faces
960  index = CCMIOSIZEC( 0 );
961  while( kCCMIONoErr == CCMIONextEntity( NULL, topologyID, kCCMIOBoundaryFaces, &index, &faceID ) )
962  {
963  MB_CHK_SET_ERR( read_faces( faceID, kCCMIOBoundaryFaces, vert_map, face_map,
964 #ifndef TUPLE_LIST
965  sense_map,
966 #endif
967  new_faces ),
968  "Trouble reading boundary faces" );
969  }
970 
971  // Now get internal faces
972  CCMIOGetEntity( &error, topologyID, kCCMIOInternalFaces, 0, &faceID );
973  CHK_SET_CCMERR( error, "Couldn't get internal faces" );
974 
975  MB_CHK_SET_ERR( read_faces( faceID, kCCMIOInternalFaces, vert_map, face_map,
976 #ifndef TUPLE_LIST
977  sense_map,
978 #endif
979  new_faces ),
980  "Trouble reading internal faces" );
981 
982  return rval;
983 }
984 
986  CCMIOEntity bdy_or_int,
987  TupleList& vert_map,
988  TupleList& face_map,
989 #ifndef TUPLE_LIST
990  SenseList& sense_map,
991 #endif
992  Range* new_faces )
993 {
994  if( kCCMIOInternalFaces != bdy_or_int && kCCMIOBoundaryFaces != bdy_or_int )
995  MB_SET_ERR( MB_FAILURE, "Face type isn't boundary or internal" );
996 
997  CCMIOSize_t dum_faces;
998  CCMIOError error = kCCMIONoErr;
999  CCMIOEntitySize( &error, faceID, &dum_faces, NULL );
1000  int num_faces = GETINT32( dum_faces );
1001 
1002  // Get the size of the face connectivity array (not really a straight connect
1003  // array, has n, connect(n), ...)
1004  CCMIOSize_t farray_size = CCMIOSIZEC( 0 );
1005  CCMIOReadFaces( &error, faceID, bdy_or_int, NULL, &farray_size, NULL, CCMIOINDEXC( kCCMIOStart ),
1006  CCMIOINDEXC( kCCMIOEnd ) );
1007  CHK_SET_CCMERR( error, "Trouble reading face connectivity length" );
1008 
1009  // Allocate vectors for holding farray and cells for each face; use new for finer
1010  // control of de-allocation
1011  int num_sides = ( kCCMIOInternalFaces == bdy_or_int ? 2 : 1 );
1012  int* farray = new int[GETINT32( farray_size )];
1013 
1014  // Read farray and make the faces
1015  CCMIOID mapID;
1016  CCMIOReadFaces( &error, faceID, bdy_or_int, &mapID, NULL, farray, CCMIOINDEXC( kCCMIOStart ),
1017  CCMIOINDEXC( kCCMIOEnd ) );
1018  CHK_SET_CCMERR( error, "Trouble reading face connectivity" );
1019 
1020  std::vector< EntityHandle > face_handles;
1021  MB_CHK_SET_ERR( make_faces( farray, vert_map, face_handles, num_faces ), "Failed to make the faces" );
1022 
1023  // Read face cells and make tuples
1024  int* face_cells;
1025  if( num_sides * num_faces < farray_size )
1026  face_cells = new int[num_sides * num_faces];
1027  else
1028  face_cells = farray;
1029  CCMIOReadFaceCells( &error, faceID, bdy_or_int, face_cells, CCMIOINDEXC( kCCMIOStart ), CCMIOINDEXC( kCCMIOEnd ) );
1030  CHK_SET_CCMERR( error, "Trouble reading face cells" );
1031 
1032  int* tmp_ptr = face_cells;
1033  for( unsigned int i = 0; i < face_handles.size(); i++ )
1034  {
1035 #ifdef TUPLE_LIST
1036  short forward = 1, reverse = -1;
1037  face_map.push_back( &forward, tmp_ptr++, &face_handles[i], NULL );
1038  if( 2 == num_sides ) face_map.push_back( &reverse, tmp_ptr++, &face_handles[i], NULL );
1039 #else
1040  face_map[*tmp_ptr].push_back( face_handles[i] );
1041  sense_map[*tmp_ptr++].push_back( 1 );
1042  if( 2 == num_sides )
1043  {
1044  face_map[*tmp_ptr].push_back( face_handles[i] );
1045  sense_map[*tmp_ptr++].push_back( -1 );
1046  }
1047 #endif
1048  }
1049 
1050  // Now read & set face gids, reuse face_cells 'cuz we know it's big enough
1051  CCMIOReadMap( &error, mapID, face_cells, CCMIOINDEXC( kCCMIOStart ), CCMIOINDEXC( kCCMIOEnd ) );
1052  CHK_SET_CCMERR( error, "Trouble reading face gids" );
1053 
1054  MB_CHK_SET_ERR( mbImpl->tag_set_data( mGlobalIdTag, &face_handles[0], face_handles.size(), face_cells ),
1055  "Couldn't set face global ids" );
1056 
1057  // Make a neumann set for these faces if they're all in a boundary face set
1058  if( kCCMIOBoundaryFaces == bdy_or_int )
1059  {
1060  EntityHandle neuset;
1061  MB_CHK_SET_ERR( mbImpl->create_meshset( MESHSET_SET, neuset ), "Failed to create neumann set" );
1062 
1063  // Don't trust entity index passed in
1064  int index;
1065  CCMIOGetEntityIndex( &error, faceID, &index );
1066  newNeusets[index] = neuset;
1067 
1068  MB_CHK_SET_ERR( mbImpl->add_entities( neuset, &face_handles[0], face_handles.size() ),
1069  "Failed to add faces to neumann set" );
1070 
1071  // Now tag as neumann set; will add id later
1072  int dum_val = 0;
1073  MB_CHK_SET_ERR( mbImpl->tag_set_data( mNeumannSetTag, &neuset, 1, &dum_val ), "Failed to tag neumann set" );
1074  }
1075 
1076  if( new_faces )
1077  {
1078  std::sort( face_handles.begin(), face_handles.end() );
1079  std::copy( face_handles.rbegin(), face_handles.rend(), range_inserter( *new_faces ) );
1080  }
1081 
1082  return MB_SUCCESS;
1083 }
1084 
1086  TupleList& vert_map,
1087  std::vector< EntityHandle >& new_faces,
1088  int num_faces )
1089 {
1090  std::vector< EntityHandle > verts;
1091  ErrorCode tmp_rval = MB_SUCCESS, rval = MB_SUCCESS;
1092 
1093  for( int i = 0; i < num_faces; i++ )
1094  {
1095  int num_verts = *farray++;
1096  verts.resize( num_verts );
1097 
1098  // Fill in connectivity by looking up by gid in vert tuple_list
1099  for( int j = 0; j < num_verts; j++ )
1100  {
1101 #ifdef TUPLE_LIST
1102  int tindex = vert_map.find( 1, farray[j] );
1103  if( -1 == tindex )
1104  {
1105  tmp_rval = MB_FAILURE;
1106  break;
1107  }
1108  verts[j] = vert_map.get_ulong( tindex, 0 );
1109 #else
1110  verts[j] = ( vert_map[farray[j]] )[0];
1111 #endif
1112  }
1113  farray += num_verts;
1114 
1115  if( MB_SUCCESS == tmp_rval )
1116  {
1117  // Make face
1118  EntityType ftype = ( 3 == num_verts ? MBTRI : ( 4 == num_verts ? MBQUAD : MBPOLYGON ) );
1119  EntityHandle faceh;
1120  tmp_rval = mbImpl->create_element( ftype, &verts[0], num_verts, faceh );
1121  if( faceh ) new_faces.push_back( faceh );
1122  }
1123 
1124  if( MB_SUCCESS != tmp_rval ) rval = tmp_rval;
1125  }
1126 
1127  return rval;
1128 }
1129 
1130 ErrorCode ReadCCMIO::read_vertices( CCMIOSize_t /* proc */,
1131  CCMIOID /* processorID */,
1132  CCMIOID verticesID,
1133  CCMIOID /* topologyID */,
1134  Range* verts,
1135  TupleList& vert_map )
1136 {
1137  CCMIOError error = kCCMIONoErr;
1138 
1139  // Pre-read the number of vertices, so we can pre-allocate & read directly in
1140  CCMIOSize_t nverts = CCMIOSIZEC( 0 );
1141  CCMIOEntitySize( &error, verticesID, &nverts, NULL );
1142  CHK_SET_CCMERR( error, "Couldn't get number of vertices" );
1143 
1144  // Get # dimensions
1145  CCMIOSize_t dims;
1146  float scale;
1147  CCMIOReadVerticesf( &error, verticesID, &dims, NULL, NULL, NULL, CCMIOINDEXC( 0 ), CCMIOINDEXC( 1 ) );
1148  CHK_SET_CCMERR( error, "Couldn't get number of dimensions" );
1149 
1150  // Allocate vertex space
1151  EntityHandle node_handle = 0;
1152  std::vector< double* > arrays;
1153  readMeshIface->get_node_coords( 3, GETINT32( nverts ), MB_START_ID, node_handle, arrays );
1154 
1155  // Read vertex coords
1156  CCMIOID mapID;
1157  std::vector< double > tmp_coords( GETINT32( dims ) * GETINT32( nverts ) );
1158  CCMIOReadVerticesd( &error, verticesID, &dims, &scale, &mapID, &tmp_coords[0], CCMIOINDEXC( 0 ),
1159  CCMIOINDEXC( 0 + nverts ) );
1160  CHK_SET_CCMERR( error, "Trouble reading vertex coordinates" );
1161 
1162  // Copy interleaved coords into moab blocked coordinate space
1163  int i = 0, threei = 0;
1164  for( ; i < nverts; i++ )
1165  {
1166  arrays[0][i] = tmp_coords[threei++];
1167  arrays[1][i] = tmp_coords[threei++];
1168  if( 3 == GETINT32( dims ) )
1169  arrays[2][i] = tmp_coords[threei++];
1170  else
1171  arrays[2][i] = 0.0;
1172  }
1173 
1174  // Scale, if necessary
1175  if( 1.0 != scale )
1176  {
1177  for( i = 0; i < nverts; i++ )
1178  {
1179  arrays[0][i] *= scale;
1180  arrays[1][i] *= scale;
1181  if( 3 == GETINT32( dims ) ) arrays[2][i] *= scale;
1182  }
1183  }
1184 
1185  // Read gids for vertices
1186  std::vector< int > gids( GETINT32( nverts ) );
1187  CCMIOReadMap( &error, mapID, &gids[0], CCMIOINDEXC( kCCMIOStart ), CCMIOINDEXC( kCCMIOEnd ) );
1188  CHK_SET_CCMERR( error, "Trouble reading vertex global ids" );
1189 
1190  // Put new vertex handles into range, and set gids for them
1191  Range new_verts( node_handle, node_handle + nverts - 1 );
1192  MB_CHK_SET_ERR( mbImpl->tag_set_data( mGlobalIdTag, new_verts, &gids[0] ), "Couldn't set gids on vertices" );
1193 
1194  // Pack vert_map with global ids and handles for these vertices
1195 #ifdef TUPLE_LIST
1196  vert_map.resize( GETINT32( nverts ) );
1197  for( i = 0; i < GETINT32( nverts ); i++ )
1198  {
1199  vert_map.push_back( NULL, &gids[i], &node_handle, NULL );
1200 #else
1201  for( i = 0; i < GETINT32( nverts ); i++ )
1202  {
1203  ( vert_map[gids[i]] ).push_back( node_handle );
1204 #endif
1205  node_handle += 1;
1206  }
1207 
1208  if( verts ) verts->merge( new_verts );
1209 
1210  return MB_SUCCESS;
1211 }
1212 
1214  CCMIOID& processorID,
1215  CCMIOID& verticesID,
1216  CCMIOID& topologyID,
1217  CCMIOID& solutionID,
1218  std::vector< CCMIOSize_t >& procs,
1219  bool& /* has_solution */ )
1220 {
1221  CCMIOSize_t proc = CCMIOSIZEC( 0 );
1222  CCMIOError error = kCCMIONoErr;
1223 
1224  CCMIONextEntity( &error, stateID, kCCMIOProcessor, &proc, &processorID );
1225  CHK_SET_CCMERR( error, "CCMIONextEntity() failed" );
1226  if( CCMIOReadProcessor( NULL, processorID, &verticesID, &topologyID, NULL, &solutionID ) != kCCMIONoErr )
1227  {
1228  // Maybe no solution; try again
1229  CCMIOReadProcessor( &error, processorID, &verticesID, &topologyID, NULL, NULL );
1230  CHK_SET_CCMERR( error, "CCMIOReadProcessor() failed" );
1231  hasSolution = false;
1232  }
1233 
1234  procs.push_back( proc );
1235 
1236  return MB_SUCCESS;
1237 }
1238 
1239 ErrorCode ReadCCMIO::read_tag_values( const char* /* file_name */,
1240  const char* /* tag_name */,
1241  const FileOptions& /* opts */,
1242  std::vector< int >& /* tag_values_out */,
1243  const SubsetList* /* subset_list */ )
1244 {
1245  return MB_FAILURE;
1246 }
1247 
1248 } // namespace moab