Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
convert.cpp
Go to the documentation of this file.
1 /**
2  * MOAB, a Mesh-Oriented datABase, is a software component for creating,
3  * storing and accessing finite element mesh data.
4  *
5  * Copyright 2004 Sandia Corporation. Under the terms of Contract
6  * DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government
7  * retains certain rights in this software.
8  *
9  * This library is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public
11  * License as published by the Free Software Foundation; either
12  * version 2.1 of the License, or (at your option) any later version.
13  *
14  */
15 
16 // If Microsoft compiler, then WIN32
17 #ifndef WIN32
18 #define WIN32
19 #endif
20 
21 #include "moab/Core.hpp"
22 #include "moab/Range.hpp"
23 #include "MBTagConventions.hpp"
24 #include "moab/ReaderWriterSet.hpp"
25 #include "moab/ReorderTool.hpp"
26 #include <iostream>
27 #include <fstream>
28 #include <sstream>
29 #include <iomanip>
30 #include <set>
31 #include <list>
32 #include <cstdlib>
33 #include <algorithm>
34 #ifndef WIN32
35 #include <sys/times.h>
36 #include <limits.h>
37 #include <unistd.h>
38 #endif
39 #include <ctime>
40 #ifdef MOAB_HAVE_MPI
41 #include "moab/ParallelComm.hpp"
42 #endif
43 
44 #ifdef MOAB_HAVE_TEMPESTREMAP
46 
47 constexpr bool offlineGenerator = true;
48 #endif
49 
50 #include <cstdio>
51 
52 /* Exit values */
53 #define USAGE_ERROR 1
54 #define READ_ERROR 2
55 #define WRITE_ERROR 3
56 #define OTHER_ERROR 4
57 #define ENT_NOT_FOUND 5
58 
59 using namespace moab;
60 
61 static void print_usage( const char* name, std::ostream& stream )
62 {
63  stream << "Usage: " << name
64  << " [-a <sat_file>|-A] [-t] [subset options] [-f format] <input_file> [<input_file2> "
65  "...] <output_file>"
66  << std::endl
67  << "\t-f <format> - Specify output file format" << std::endl
68  << "\t-a <acis_file> - ACIS SAT file dumped by .cub reader (same as \"-o "
69  "SAT_FILE=acis_file\""
70  << std::endl
71  << "\t-A - .cub file reader should not dump a SAT file (depricated default)" << std::endl
72  << "\t-o option - Specify write option." << std::endl
73  << "\t-O option - Specify read option." << std::endl
74  << "\t-t - Time read and write of files." << std::endl
75  << "\t-g - Enable verbose/debug output." << std::endl
76  << "\t-h - Print this help text and exit." << std::endl
77  << "\t-l - List available file formats and exit." << std::endl
78  << "\t-I <dim> - Generate internal entities of specified dimension." << std::endl
79 #ifdef MOAB_HAVE_MPI
80  << "\t-P - Append processor ID to output file name" << std::endl
81  << "\t-p - Replace '%' with processor ID in input and output file name" << std::endl
82  << "\t-M[0|1|2] - Read/write in parallel, optionally also doing "
83  "resolve_shared_ents (1) and exchange_ghosts (2)"
84  << std::endl
85  << "\t-z <file> - Read metis partition information corresponding to an MPAS grid "
86  "file and create h5m partition file"
87  << std::endl
88 #endif
89 
90 #ifdef MOAB_HAVE_TEMPESTREMAP
91  << "\t-B - Use TempestRemap exodus file reader and convert to MOAB format" << std::endl
92  << "\t-b - Convert MOAB mesh to TempestRemap exodus file writer" << std::endl
93  << "\t-S - Scale climate mesh to unit sphere" << std::endl
94  << "\t-i - Name of the global DoF tag to use with mbtempest" << std::endl
95  << "\t-r - Order of field DoF (discretization) data; FV=1,SE=[1,N]" << std::endl
96 #endif
97  << "\t-- - treat all subsequent options as file names" << std::endl
98  << "\t (allows file names beginning with '-')" << std::endl
99  << " subset options: " << std::endl
100  << "\tEach of the following options should be followed by " << std::endl
101  << "\ta list of ids. IDs may be separated with commas. " << std::endl
102  << "\tRanges of IDs may be specified with a '-' between " << std::endl
103  << "\ttwo values. The list may not contain spaces." << std::endl
104  << "\t-v - Volume" << std::endl
105  << "\t-s - Surface" << std::endl
106  << "\t-c - Curve" << std::endl
107  << "\t-V - Vertex" << std::endl
108  << "\t-m - Material set (block)" << std::endl
109  << "\t-d - Dirichlet set (nodeset)" << std::endl
110  << "\t-n - Neumann set (sideset)" << std::endl
111  << "\t-D - Parallel partitioning set (PARALLEL_PARTITION)" << std::endl
112  << "\tThe presence of one or more of the following flags limits " << std::endl
113  << "\tthe exported mesh to only elements of the corresponding " << std::endl
114  << "\tdimension. Vertices are always exported." << std::endl
115  << "\t-1 - Edges " << std::endl
116  << "\t-2 - Tri, Quad, Polygon " << std::endl
117  << "\t-3 - Tet, Hex, Prism, etc. " << std::endl;
118 }
119 
120 static void print_help( const char* name )
121 {
122  std::cout << " This program can be used to convert between mesh file\n"
123  " formats, extract a subset of a mesh file to a separate\n"
124  " file, or both. The type of file to write is determined\n"
125  " from the file extension (e.g. \".vtk\") portion of the\n"
126  " output file name.\n"
127  " \n"
128  " While MOAB strives to export and import all data from\n"
129  " each supported file format, most file formats do\n"
130  " not support MOAB's entire data model. Thus MOAB cannot\n"
131  " guarantee lossless conversion for any file formats\n"
132  " other than the native HDF5 representation.\n"
133  "\n";
134 
135  print_usage( name, std::cout );
136  exit( 0 );
137 }
138 
139 static void usage_error( const char* name )
140 {
141  print_usage( name, std::cerr );
142 #ifdef MOAB_HAVE_MPI
143  MPI_Finalize();
144 #endif
145  exit( USAGE_ERROR );
146 }
147 
148 static void list_formats( Interface* );
149 static bool parse_id_list( const char* string, std::set< int >& );
150 static void print_id_list( const char*, std::ostream& stream, const std::set< int >& list );
151 static void reset_times();
152 static void write_times( std::ostream& stream );
153 static void remove_entities_from_sets( Interface* gMB, Range& dead_entities, Range& empty_sets );
154 static void remove_from_vector( std::vector< EntityHandle >& vect, const Range& ents_to_remove );
155 static bool make_opts_string( std::vector< std::string > options, std::string& result );
156 static std::string percent_subst( const std::string& s, int val );
157 
158 static int process_partition_file( Interface* gMB, std::string& metis_partition_file );
159 
160 int main( int argc, char* argv[] )
161 {
162  int proc_id = 0;
163 #ifdef MOAB_HAVE_MPI
164  MPI_Init( &argc, &argv );
165  MPI_Comm_rank( MPI_COMM_WORLD, &proc_id );
166 #endif
167 
168 #ifdef MOAB_HAVE_TEMPESTREMAP
169  bool tempestin = false, tempestout = false;
170 #endif
171 
172  Core core;
173  Interface* gMB = &core;
174  ErrorCode result;
175  Range range;
176 
177 #if( defined( MOAB_HAVE_MPI ) && defined( MOAB_HAVE_TEMPESTREMAP ) )
178  moab::ParallelComm* pcomm = new moab::ParallelComm( gMB, MPI_COMM_WORLD, 0 );
179 #endif
180 
181  bool append_rank = false;
182  bool percent_rank_subst = false;
183  bool file_written = false;
184  int i, dim;
185  std::list< std::string >::iterator j;
186  bool dims[4] = { false, false, false, false };
187  const char* format = NULL; // output file format
188  std::list< std::string > in; // input file name list
189  std::string out; // output file name
190  bool verbose = false;
191  std::set< int > geom[4], mesh[4]; // user-specified IDs
192  std::vector< EntityHandle > set_list; // list of user-specified sets to write
193  std::vector< std::string > write_opts, read_opts;
194  std::string metis_partition_file;
195 #ifdef MOAB_HAVE_TEMPESTREMAP
196  std::string globalid_tag_name;
197  int spectral_order = 1;
198  bool unitscaling = false;
199 #endif
200 
201  const char* const mesh_tag_names[] = { DIRICHLET_SET_TAG_NAME, NEUMANN_SET_TAG_NAME, MATERIAL_SET_TAG_NAME,
203  const char* const geom_names[] = { "VERTEX", "CURVE", "SURFACE", "VOLUME" };
204 
205  // scan arguments
206  bool do_flag = true;
207  bool print_times = false;
208  bool generate[] = { false, false, false };
209  bool pval;
210  bool parallel = false, resolve_shared = false, exchange_ghosts = false;
211  for( i = 1; i < argc; i++ )
212  {
213  if( !argv[i][0] ) usage_error( argv[0] );
214 
215  if( do_flag && argv[i][0] == '-' )
216  {
217  if( !argv[i][1] || ( argv[i][1] != 'M' && argv[i][2] ) ) usage_error( argv[0] );
218 
219  switch( argv[i][1] )
220  {
221  // do flag arguments:
222  case '-':
223  do_flag = false;
224  break;
225  case 'g':
226  verbose = true;
227  break;
228  case 't':
229  print_times = true;
230  break;
231  case 'A':
232  break;
233  case 'h':
234  case 'H':
235  print_help( argv[0] );
236  break;
237  case 'l':
238  list_formats( gMB );
239  break;
240 #ifdef MOAB_HAVE_MPI
241  case 'P':
242  append_rank = true;
243  break;
244  case 'p':
245  percent_rank_subst = true;
246  break;
247  case 'M':
248  parallel = true;
249  if( argv[i][2] == '1' || argv[i][2] == '2' ) resolve_shared = true;
250  if( argv[i][2] == '2' ) exchange_ghosts = true;
251  break;
252 #endif
253 #ifdef MOAB_HAVE_TEMPESTREMAP
254  case 'B':
255  tempestin = true;
256  break;
257  case 'b':
258  tempestout = true;
259  break;
260  case 'S':
261  unitscaling = true;
262  break;
263 #endif
264  case '1':
265  case '2':
266  case '3':
267  dims[argv[i][1] - '0'] = true;
268  break;
269  // do options that require additional args:
270  default:
271  ++i;
272  if( i == argc || argv[i][0] == '-' )
273  {
274  std::cerr << "Expected argument following " << argv[i - 1] << std::endl;
275  usage_error( argv[0] );
276  }
277  if( argv[i - 1][1] == 'I' )
278  {
279  dim = atoi( argv[i] );
280  if( dim < 1 || dim > 2 )
281  {
282  std::cerr << "Invalid dimension value following -I" << std::endl;
283  usage_error( argv[0] );
284  }
285  generate[dim] = true;
286  continue;
287  }
288  pval = false;
289  switch( argv[i - 1][1] )
290  {
291  case 'a':
292  read_opts.push_back( std::string( "SAT_FILE=" ) + argv[i] );
293  pval = true;
294  break;
295  case 'f':
296  format = argv[i];
297  pval = true;
298  break;
299  case 'o':
300  write_opts.push_back( argv[i] );
301  pval = true;
302  break;
303  case 'O':
304  read_opts.push_back( argv[i] );
305  pval = true;
306  break;
307 #ifdef MOAB_HAVE_TEMPESTREMAP
308  case 'i':
309  globalid_tag_name = std::string( argv[i] );
310  pval = true;
311  break;
312  case 'r':
313  spectral_order = atoi( argv[i] );
314  pval = true;
315  break;
316 #endif
317  case 'v':
318  pval = parse_id_list( argv[i], geom[3] );
319  break;
320  case 's':
321  pval = parse_id_list( argv[i], geom[2] );
322  break;
323  case 'c':
324  pval = parse_id_list( argv[i], geom[1] );
325  break;
326  case 'V':
327  pval = parse_id_list( argv[i], geom[0] );
328  break;
329  case 'D':
330  pval = parse_id_list( argv[i], mesh[3] );
331  break;
332  case 'm':
333  pval = parse_id_list( argv[i], mesh[2] );
334  break;
335  case 'n':
336  pval = parse_id_list( argv[i], mesh[1] );
337  break;
338  case 'd':
339  pval = parse_id_list( argv[i], mesh[0] );
340  break;
341  case 'z':
342  metis_partition_file = argv[i];
343  pval = true;
344  break;
345  default:
346  std::cerr << "Invalid option: " << argv[i] << std::endl;
347  }
348 
349  if( !pval )
350  {
351  std::cerr << "Invalid flag or flag value: " << argv[i - 1] << " " << argv[i] << std::endl;
352  usage_error( argv[0] );
353  }
354  }
355  }
356  // do file names
357  else
358  {
359  in.push_back( argv[i] );
360  }
361  }
362  if( in.size() < 2 )
363  {
364  std::cerr << "No output file name specified." << std::endl;
365  usage_error( argv[0] );
366  }
367  // output file name is the last one specified
368  out = in.back();
369  in.pop_back();
370 
371  if( append_rank )
372  {
373  std::ostringstream mod;
374  mod << out << "." << proc_id;
375  out = mod.str();
376  }
377 
378  if( percent_rank_subst )
379  {
380  for( j = in.begin(); j != in.end(); ++j )
381  *j = percent_subst( *j, proc_id );
382  out = percent_subst( out, proc_id );
383  }
384 
385  // construct options string from individual options
386  std::string read_options, write_options;
387  if( parallel )
388  {
389  read_opts.push_back( "PARALLEL=READ_PART" );
390  read_opts.push_back( "PARTITION=PARALLEL_PARTITION" );
391  if( !append_rank && !percent_rank_subst ) write_opts.push_back( "PARALLEL=WRITE_PART" );
392  }
393  if( resolve_shared ) read_opts.push_back( "PARALLEL_RESOLVE_SHARED_ENTS" );
394  if( exchange_ghosts ) read_opts.push_back( "PARALLEL_GHOSTS=3.0.1" );
395 
396  if( !make_opts_string( read_opts, read_options ) || !make_opts_string( write_opts, write_options ) )
397  {
398 #ifdef MOAB_HAVE_MPI
399  MPI_Finalize();
400 #endif
401  return USAGE_ERROR;
402  }
403 
404  if( !metis_partition_file.empty() )
405  {
406  if( ( in.size() != 1 ) || ( proc_id != 0 ) )
407  {
408  std::cerr << " mpas partition allows only one input file, in serial conversion\n";
409 #ifdef MOAB_HAVE_MPI
410  MPI_Finalize();
411 #endif
412  return USAGE_ERROR;
413  }
414  }
415 
416  Tag id_tag = gMB->globalId_tag();
417 
418  // Read the input file.
419 #ifdef MOAB_HAVE_TEMPESTREMAP
420  if( tempestin && in.size() > 1 )
421  {
422  std::cerr << " we can read only one tempest files at a time\n";
423 #ifdef MOAB_HAVE_MPI
424  MPI_Finalize();
425 #endif
426  return USAGE_ERROR;
427  }
428  TempestRemapper* remapper = NULL;
429  if( tempestin or tempestout )
430  {
431 #ifdef MOAB_HAVE_MPI
432  remapper = new moab::TempestRemapper( gMB, pcomm, offlineGenerator );
433 #else
434  remapper = new moab::TempestRemapper( gMB, offlineGenerator );
435 #endif
436  }
437 
438  bool use_overlap_context = false;
439  Tag srcParentTag, tgtParentTag;
440 
441 #endif
442  for( j = in.begin(); j != in.end(); ++j )
443  {
444  std::string inFileName = *j;
445 
446  reset_times();
447 
448 #ifdef MOAB_HAVE_TEMPESTREMAP
449  if( remapper )
450  {
451  remapper->meshValidate = false;
452  // remapper->constructEdgeMap = true;
453  remapper->initialize();
454 
455  if( tempestin )
456  {
458 
459  // convert
460  result = remapper->LoadMesh( moab::Remapper::SourceMesh, inFileName, moab::TempestRemapper::DEFAULT );MB_CHK_ERR( result );
461 
462  Mesh* tempestMesh = remapper->GetMesh( moab::Remapper::SourceMesh );
463  tempestMesh->RemoveZeroEdges();
464  tempestMesh->RemoveCoincidentNodes();
465 
466  // Load the meshes and validate
467  result = remapper->ConvertTempestMesh( moab::Remapper::SourceMesh );
468 
469  if( unitscaling )
470  {
471  result = moab::IntxUtils::ScaleToRadius( gMB, srcmesh, 1.0 );MB_CHK_ERR( result );
472  }
473 
474  // Check if we are converting a RLL grid
475  NcFile ncInput( inFileName.c_str(), NcFile::ReadOnly );
476 
477  NcError error_temp( NcError::silent_nonfatal );
478  // get the attribute
479  NcAtt* attRectilinear = ncInput.get_att( "rectilinear" );
480 
481  // If rectilinear attribute present, mark it
482  std::vector< int > vecDimSizes( 3, 0 );
483  Tag rectilinearTag;
484  // Tag data contains: guessed mesh type, mesh size1, mesh size 2
485  // Example: CS(0)/ICO(1)/ICOD(2), num_elements, num_nodes
486  // : RLL(3), num_lat, num_lon
487  result = gMB->tag_get_handle( "ClimateMetadata", 3, MB_TYPE_INTEGER, rectilinearTag,
488  MB_TAG_SPARSE | MB_TAG_CREAT, vecDimSizes.data() );MB_CHK_SET_ERR( result, "can't create rectilinear sizes tag" );
489 
490  if( attRectilinear != nullptr )
491  {
492  // Obtain rectilinear attributes (dimension sizes)
493  NcAtt* attRectilinearDim0Size = ncInput.get_att( "rectilinear_dim0_size" );
494  NcAtt* attRectilinearDim1Size = ncInput.get_att( "rectilinear_dim1_size" );
495 
496  if( attRectilinearDim0Size == nullptr )
497  {
498  _EXCEPTIONT( "Missing attribute \"rectilinear_dim0_size\"" );
499  }
500  if( attRectilinearDim1Size == nullptr )
501  {
502  _EXCEPTIONT( "Missing attribute \"rectilinear_dim1_size\"" );
503  }
504 
505  int nDim0Size = attRectilinearDim0Size->as_int( 0 );
506  int nDim1Size = attRectilinearDim1Size->as_int( 0 );
507 
508  // Obtain rectilinear attributes (dimension names)
509  NcAtt* attRectilinearDim0Name = ncInput.get_att( "rectilinear_dim0_name" );
510  NcAtt* attRectilinearDim1Name = ncInput.get_att( "rectilinear_dim1_name" );
511 
512  if( attRectilinearDim0Name == nullptr )
513  {
514  _EXCEPTIONT( "Missing attribute \"rectilinear_dim0_name\"" );
515  }
516  if( attRectilinearDim1Name == nullptr )
517  {
518  _EXCEPTIONT( "Missing attribute \"rectilinear_dim1_name\"" );
519  }
520 
521  std::string strDim0Name = attRectilinearDim0Name->as_string( 0 );
522  std::string strDim1Name = attRectilinearDim1Name->as_string( 0 );
523 
524  std::map< std::string, int > vecDimNameSizes;
525  // Push rectilinear attributes into array
526  vecDimNameSizes[strDim0Name] = nDim0Size;
527  vecDimNameSizes[strDim1Name] = nDim1Size;
528  vecDimSizes[0] = static_cast< int >( moab::TempestRemapper::RLL );
529  vecDimSizes[1] = vecDimNameSizes["lat"];
530  vecDimSizes[2] = vecDimNameSizes["lon"];
531 
532  printf( "Rectilinear RLL mesh size: (lat) %d X (lon) %d\n", vecDimSizes[1], vecDimSizes[2] );
533 
534  moab::EntityHandle mSet = 0;
535  // mSet = remapper->GetMeshSet( moab::Remapper::SourceMesh );
536  result = gMB->tag_set_data( rectilinearTag, &mSet, 1, vecDimSizes.data() );MB_CHK_ERR( result );
537  }
538  else
539  {
540  const Range& elems = remapper->GetMeshEntities( moab::Remapper::SourceMesh );
541  bool isQuads = elems.all_of_type( moab::MBQUAD );
542  bool isTris = elems.all_of_type( moab::MBTRI );
543  // vecDimSizes[0] = static_cast< int >( remapper->GetMeshType( moab::Remapper::SourceMesh );
544  vecDimSizes[0] = ( isQuads ? static_cast< int >( moab::TempestRemapper::CS )
545  : ( isTris ? static_cast< int >( moab::TempestRemapper::ICO )
546  : static_cast< int >( moab::TempestRemapper::ICOD ) ) );
547  vecDimSizes[1] = elems.size();
548  vecDimSizes[2] = remapper->GetMeshVertices( moab::Remapper::SourceMesh ).size();
549 
550  switch( vecDimSizes[0] )
551  {
552  case 0:
553  printf( "Cubed-Sphere mesh: %d (elems), %d (nodes)\n", vecDimSizes[1], vecDimSizes[2] );
554  break;
555  case 2:
556  printf( "Icosahedral (triangular) mesh: %d (elems), %d (nodes)\n", vecDimSizes[1],
557  vecDimSizes[2] );
558  break;
559  case 3:
560  default:
561  printf( "Polygonal mesh: %d (elems), %d (nodes)\n", vecDimSizes[1], vecDimSizes[2] );
562  break;
563  }
564 
565  moab::EntityHandle mSet = 0;
566  // mSet = remapper->GetMeshSet( moab::Remapper::SourceMesh );
567  result = gMB->tag_set_data( rectilinearTag, &mSet, 1, vecDimSizes.data() );MB_CHK_ERR( result );
568  }
569 
570  ncInput.close();
571 
572  const size_t nOverlapFaces = tempestMesh->faces.size();
573  if( tempestMesh->vecSourceFaceIx.size() == nOverlapFaces &&
574  tempestMesh->vecSourceFaceIx.size() == nOverlapFaces )
575  {
576  int defaultInt = -1;
577  use_overlap_context = true;
578  // Check if our MOAB mesh has RED and BLUE tags; this would indicate we are
579  // converting an overlap grid
580  result = gMB->tag_get_handle( "TargetParent", 1, MB_TYPE_INTEGER, tgtParentTag,
581  MB_TAG_DENSE | MB_TAG_CREAT, &defaultInt );MB_CHK_SET_ERR( result, "can't create target parent tag" );
582 
583  result = gMB->tag_get_handle( "SourceParent", 1, MB_TYPE_INTEGER, srcParentTag,
584  MB_TAG_DENSE | MB_TAG_CREAT, &defaultInt );MB_CHK_SET_ERR( result, "can't create source parent tag" );
585 
586  const Range& faces = remapper->GetMeshEntities( moab::Remapper::SourceMesh );
587 
588  std::vector< int > gids( faces.size() ), srcpar( faces.size() ), tgtpar( faces.size() );
589  result = gMB->tag_get_data( id_tag, faces, &gids[0] );MB_CHK_ERR( result );
590 
591  for( unsigned ii = 0; ii < faces.size(); ++ii )
592  {
593  srcpar[ii] = tempestMesh->vecSourceFaceIx[gids[ii] - 1];
594  tgtpar[ii] = tempestMesh->vecTargetFaceIx[gids[ii] - 1];
595  }
596 
597  result = gMB->tag_set_data( srcParentTag, faces, &srcpar[0] );MB_CHK_ERR( result );
598  result = gMB->tag_set_data( tgtParentTag, faces, &tgtpar[0] );MB_CHK_ERR( result );
599 
600  srcpar.clear();
601  tgtpar.clear();
602  gids.clear();
603  }
604  }
605  else if( tempestout )
606  {
609 
610  // load the mesh in MOAB format
611  std::vector< int > metadata( 2 );
612  result = remapper->LoadNativeMesh( *j, srcmesh, metadata );MB_CHK_ERR( result );
613 
614  if( unitscaling )
615  {
616  result = moab::IntxUtils::ScaleToRadius( gMB, srcmesh, 1.0 );MB_CHK_ERR( result );
617  }
618 
619  // Check if our MOAB mesh has RED and BLUE tags; this would indicate we are converting
620  // an overlap grid
621  ErrorCode rval1 = gMB->tag_get_handle( "SourceParent", srcParentTag );
622  ErrorCode rval2 = gMB->tag_get_handle( "TargetParent", tgtParentTag );
623  if( rval1 == MB_SUCCESS && rval2 == MB_SUCCESS )
624  {
625  use_overlap_context = true;
626  ovmesh = srcmesh;
627 
628  Tag countTag;
629  result = gMB->tag_get_handle( "Counting", countTag );MB_CHK_ERR( result );
630 
631  // Load the meshes and validate
632  Tag order;
633  ReorderTool reorder_tool( &core );
634  result = reorder_tool.handle_order_from_int_tag( srcParentTag, -1, order );MB_CHK_ERR( result );
635  result = reorder_tool.reorder_entities( order );MB_CHK_ERR( result );
636  result = gMB->tag_delete( order );MB_CHK_ERR( result );
637  result = remapper->ConvertMeshToTempest( moab::Remapper::OverlapMesh );MB_CHK_ERR( result );
638  }
639  else
640  {
641  if( metadata[0] == static_cast< int >( moab::TempestRemapper::RLL ) )
642  {
643  assert( metadata.size() );
644  std::cout << "Converting a RLL mesh with rectilinear dimension: " << metadata[0] << " X "
645  << metadata[1] << std::endl;
646  }
647 
648  // Convert the mesh and validate
649  result = remapper->ConvertMeshToTempest( moab::Remapper::SourceMesh );MB_CHK_ERR( result );
650  }
651  }
652  }
653  else
654  result = gMB->load_file( j->c_str(), 0, read_options.c_str() );
655 #else
656  result = gMB->load_file( j->c_str(), 0, read_options.c_str() );
657 #endif
658  if( MB_SUCCESS != result )
659  {
660  std::cerr << "Failed to load \"" << *j << "\"." << std::endl;
661  std::cerr << "Error code: " << gMB->get_error_string( result ) << " (" << result << ")" << std::endl;
662  std::string message;
663  if( MB_SUCCESS == gMB->get_last_error( message ) && !message.empty() )
664  std::cerr << "Error message: " << message << std::endl;
665 #ifdef MOAB_HAVE_MPI
666  MPI_Finalize();
667 #endif
668  return READ_ERROR;
669  }
670  if( !proc_id ) std::cerr << "Read \"" << *j << "\"" << std::endl;
671  if( print_times && !proc_id ) write_times( std::cout );
672  }
673 
674  // Determine if the user has specified any geometry sets to write
675  bool have_geom = false;
676  for( dim = 0; dim <= 3; ++dim )
677  {
678  if( !geom[dim].empty() ) have_geom = true;
679  if( verbose ) print_id_list( geom_names[dim], std::cout, geom[dim] );
680  }
681 
682  // True if the user has specified any sets to write
683  bool have_sets = have_geom;
684 
685  // Get geometry tags
686  Tag dim_tag;
687  if( have_geom )
688  {
689  if( id_tag == 0 )
690  {
691  std::cerr << "No ID tag defined." << std::endl;
692  have_geom = false;
693  }
694  result = gMB->tag_get_handle( GEOM_DIMENSION_TAG_NAME, 1, MB_TYPE_INTEGER, dim_tag );
695  if( MB_SUCCESS != result )
696  {
697  std::cerr << "No geometry tag defined." << std::endl;
698  have_geom = false;
699  }
700  }
701 
702  // Get geometry sets
703  if( have_geom )
704  {
705  int id_val;
706  Tag tags[] = { id_tag, dim_tag };
707  const void* vals[] = { &id_val, &dim };
708  for( dim = 0; dim <= 3; ++dim )
709  {
710  int init_count = set_list.size();
711  for( std::set< int >::iterator iter = geom[dim].begin(); iter != geom[dim].end(); ++iter )
712  {
713  id_val = *iter;
714  range.clear();
715  result = gMB->get_entities_by_type_and_tag( 0, MBENTITYSET, tags, vals, 2, range );
716  if( MB_SUCCESS != result || range.empty() )
717  {
718  range.clear();
719  std::cerr << geom_names[dim] << " " << id_val << " not found.\n";
720  }
721  std::copy( range.begin(), range.end(), std::back_inserter( set_list ) );
722  }
723 
724  if( verbose )
725  std::cout << "Found " << ( set_list.size() - init_count ) << ' ' << geom_names[dim] << " sets"
726  << std::endl;
727  }
728  }
729 
730  // Get mesh groupings
731  for( i = 0; i < 4; ++i )
732  {
733  if( verbose ) print_id_list( mesh_tag_names[i], std::cout, mesh[i] );
734 
735  if( mesh[i].empty() ) continue;
736  have_sets = true;
737 
738  // Get tag
739  Tag tag;
740  result = gMB->tag_get_handle( mesh_tag_names[i], 1, MB_TYPE_INTEGER, tag );
741  if( MB_SUCCESS != result )
742  {
743  std::cerr << "Tag not found: " << mesh_tag_names[i] << std::endl;
744  continue;
745  }
746 
747  // get entity sets
748  int init_count = set_list.size();
749  for( std::set< int >::iterator iter = mesh[i].begin(); iter != mesh[i].end(); ++iter )
750  {
751  range.clear();
752  const void* vals[] = { &*iter };
753  result = gMB->get_entities_by_type_and_tag( 0, MBENTITYSET, &tag, vals, 1, range );
754  if( MB_SUCCESS != result || range.empty() )
755  {
756  range.clear();
757  std::cerr << mesh_tag_names[i] << " " << *iter << " not found.\n";
758  }
759  std::copy( range.begin(), range.end(), std::back_inserter( set_list ) );
760  }
761 
762  if( verbose )
763  std::cout << "Found " << ( set_list.size() - init_count ) << ' ' << mesh_tag_names[i] << " sets"
764  << std::endl;
765  }
766 
767  // Check if output is limited to certain dimensions of elements
768  bool bydim = false;
769  for( dim = 1; dim < 4; ++dim )
770  if( dims[dim] ) bydim = true;
771 
772  // Check conflicting input
773  if( bydim )
774  {
775  if( generate[1] && !dims[1] )
776  {
777  std::cerr << "Warning: Request to generate 1D internal entities but not export them." << std::endl;
778  generate[1] = false;
779  }
780  if( generate[2] && !dims[2] )
781  {
782  std::cerr << "Warning: Request to generate 2D internal entities but not export them." << std::endl;
783  generate[2] = false;
784  }
785  }
786 
787  // Generate any internal entities
788  if( generate[1] || generate[2] )
789  {
790  EntityHandle all_mesh = 0;
791  const EntityHandle* sets = &all_mesh;
792  int num_sets = 1;
793  if( have_sets )
794  {
795  num_sets = set_list.size();
796  sets = &set_list[0];
797  }
798  for( i = 0; i < num_sets; ++i )
799  {
800  Range dim3, dim2, adj;
801  gMB->get_entities_by_dimension( sets[i], 3, dim3, true );
802  if( generate[1] )
803  {
804  gMB->get_entities_by_dimension( sets[i], 2, dim2, true );
805  gMB->get_adjacencies( dim3, 1, true, adj, Interface::UNION );
806  gMB->get_adjacencies( dim2, 1, true, adj, Interface::UNION );
807  }
808  if( generate[2] )
809  {
810  gMB->get_adjacencies( dim3, 2, true, adj, Interface::UNION );
811  }
812  if( sets[i] ) gMB->add_entities( sets[i], adj );
813  }
814  }
815 
816  // Delete any entities not of the dimensions to be exported
817  if( bydim )
818  {
819  // Get list of dead elements
820  Range dead_entities, tmp_range;
821  for( dim = 1; dim <= 3; ++dim )
822  {
823  if( dims[dim] ) continue;
824  gMB->get_entities_by_dimension( 0, dim, tmp_range );
825  dead_entities.merge( tmp_range );
826  }
827  // Remove dead entities from all sets, and add all
828  // empty sets to list of dead entities.
829  Range empty_sets;
830  remove_entities_from_sets( gMB, dead_entities, empty_sets );
831  while( !empty_sets.empty() )
832  {
833  if( !set_list.empty() ) remove_from_vector( set_list, empty_sets );
834  dead_entities.merge( empty_sets );
835  tmp_range.clear();
836  remove_entities_from_sets( gMB, empty_sets, tmp_range );
837  empty_sets = subtract( tmp_range, dead_entities );
838  }
839  // Destroy dead entities
840  gMB->delete_entities( dead_entities );
841  }
842 
843  // If user specified sets to write, but none were found, exit.
844  if( have_sets && set_list.empty() )
845  {
846  std::cerr << "Nothing to write." << std::endl;
847 #ifdef MOAB_HAVE_MPI
848  MPI_Finalize();
849 #endif
850  return ENT_NOT_FOUND;
851  }
852 
853  // interpret the mpas partition file created by gpmetis
854  if( !metis_partition_file.empty() )
855  {
856  int err = process_partition_file( gMB, metis_partition_file );
857  if( err )
858  {
859  std::cerr << "Failed to process partition file \"" << metis_partition_file << "\"." << std::endl;
860 #ifdef MOAB_HAVE_MPI
861  MPI_Finalize();
862 #endif
863  return WRITE_ERROR;
864  }
865  }
866  if( verbose )
867  {
868  if( have_sets )
869  std::cout << "Found " << set_list.size() << " specified sets to write (total)." << std::endl;
870  else
871  std::cout << "No sets specifed. Writing entire mesh." << std::endl;
872  }
873 
874  // Write the output file
875  reset_times();
876 #ifdef MOAB_HAVE_TEMPESTREMAP
877  if( remapper )
878  {
879  Range faces;
880  Mesh* tempestMesh =
881  remapper->GetMesh( ( use_overlap_context ? moab::Remapper::OverlapMesh : moab::Remapper::SourceMesh ) );
882  moab::EntityHandle& srcmesh =
883  remapper->GetMeshSet( ( use_overlap_context ? moab::Remapper::OverlapMesh : moab::Remapper::SourceMesh ) );
884  result = gMB->get_entities_by_dimension( srcmesh, 2, faces );MB_CHK_ERR( result );
885  int ntot_elements = 0, nelements = faces.size();
886 #ifdef MOAB_HAVE_MPI
887  int ierr = MPI_Allreduce( &nelements, &ntot_elements, 1, MPI_INT, MPI_SUM, pcomm->comm() );
888  if( ierr != 0 ) MB_CHK_SET_ERR( MB_FAILURE, "MPI_Allreduce failed to get total source elements" );
889 #else
890  ntot_elements = nelements;
891 #endif
892 
893  Tag gidTag = gMB->globalId_tag();
894  std::vector< int > gids( faces.size() );
895  result = gMB->tag_get_data( gidTag, faces, &gids[0] );MB_CHK_ERR( result );
896 
897  if( faces.size() > 1 && gids[0] == gids[1] && !use_overlap_context )
898  {
899 #ifdef MOAB_HAVE_MPI
900  result = pcomm->assign_global_ids( srcmesh, 2, 1, false );MB_CHK_ERR( result );
901 #else
902  result = remapper->assign_vertex_element_IDs( gidTag, srcmesh, 2, 1 );MB_CHK_ERR( result );
903  result = remapper->assign_vertex_element_IDs( gidTag, srcmesh, 0, 1 );MB_CHK_ERR( result );
904 #endif
905  }
906 
907  // VSM: If user requested explicitly for some metadata, we need to generate the DoF ID tag
908  // and set the appropriate numbering based on specified discretization order
909  // Useful only for SE meshes with GLL DoFs
910  if( spectral_order > 1 && globalid_tag_name.size() > 1 )
911  {
912  result = remapper->GenerateMeshMetadata( *tempestMesh, ntot_elements, faces, NULL, globalid_tag_name,
913  spectral_order );MB_CHK_ERR( result );
914  }
915 
916  if( tempestout )
917  {
918  // Check if our MOAB mesh has RED and BLUE tags; this would indicate we are converting an
919  // overlap grid
920  if( use_overlap_context && false )
921  {
922  const int nOverlapFaces = faces.size();
923  // Overlap mesh: resize the source and target connection arrays
924  tempestMesh->vecSourceFaceIx.resize( nOverlapFaces ); // 0-based indices corresponding to source mesh
925  tempestMesh->vecTargetFaceIx.resize( nOverlapFaces ); // 0-based indices corresponding to target mesh
926  result = gMB->tag_get_data( srcParentTag, faces, &tempestMesh->vecSourceFaceIx[0] );MB_CHK_ERR( result );
927  result = gMB->tag_get_data( tgtParentTag, faces, &tempestMesh->vecTargetFaceIx[0] );MB_CHK_ERR( result );
928  }
929  // Write out the mesh using TempestRemap
930  tempestMesh->Write( out, NcFile::Netcdf4 );
931  file_written = true;
932  }
933  delete remapper; // cleanup
934  }
935 #endif
936 
937  if( !file_written )
938  {
939  if( have_sets )
940  result = gMB->write_file( out.c_str(), format, write_options.c_str(), &set_list[0], set_list.size() );
941  else
942  result = gMB->write_file( out.c_str(), format, write_options.c_str() );
943  if( MB_SUCCESS != result )
944  {
945  std::cerr << "Failed to write \"" << out << "\"." << std::endl;
946  std::cerr << "Error code: " << gMB->get_error_string( result ) << " (" << result << ")" << std::endl;
947  std::string message;
948  if( MB_SUCCESS == gMB->get_last_error( message ) && !message.empty() )
949  std::cerr << "Error message: " << message << std::endl;
950 #ifdef MOAB_HAVE_MPI
951  MPI_Finalize();
952 #endif
953  return WRITE_ERROR;
954  }
955  }
956 
957  if( !proc_id ) std::cerr << "Wrote \"" << out << "\"" << std::endl;
958  if( print_times && !proc_id ) write_times( std::cout );
959 
960 #ifdef MOAB_HAVE_MPI
961  MPI_Finalize();
962 #endif
963  return 0;
964 }
965 
966 bool parse_id_list( const char* string, std::set< int >& results )
967 {
968  bool okay = true;
969  char* mystr = strdup( string );
970  for( const char* ptr = strtok( mystr, "," ); ptr; ptr = strtok( 0, "," ) )
971  {
972  char* endptr;
973  long val = strtol( ptr, &endptr, 0 );
974  if( endptr == ptr || val <= 0 )
975  {
976  std::cerr << "Not a valid id: " << ptr << std::endl;
977  okay = false;
978  break;
979  }
980 
981  long val2 = val;
982  if( *endptr == '-' )
983  {
984  const char* sptr = endptr + 1;
985  val2 = strtol( sptr, &endptr, 0 );
986  if( endptr == sptr || val2 <= 0 )
987  {
988  std::cerr << "Not a valid id: " << sptr << std::endl;
989  okay = false;
990  break;
991  }
992  if( val2 < val )
993  {
994  std::cerr << "Invalid id range: " << ptr << std::endl;
995  okay = false;
996  break;
997  }
998  }
999 
1000  if( *endptr )
1001  {
1002  std::cerr << "Unexpected character: " << *endptr << std::endl;
1003  okay = false;
1004  break;
1005  }
1006 
1007  for( ; val <= val2; ++val )
1008  if( !results.insert( (int)val ).second ) std::cerr << "Warning: duplicate Id: " << val << std::endl;
1009  }
1010 
1011  free( mystr );
1012  return okay;
1013 }
1014 
1015 void print_id_list( const char* head, std::ostream& stream, const std::set< int >& list )
1016 {
1017  stream << head << ": ";
1018 
1019  if( list.empty() )
1020  {
1021  stream << "(none)" << std::endl;
1022  return;
1023  }
1024 
1025  int start, prev;
1026  std::set< int >::const_iterator iter = list.begin();
1027  start = prev = *( iter++ );
1028  for( ;; )
1029  {
1030  if( iter == list.end() || *iter != 1 + prev )
1031  {
1032  stream << start;
1033  if( prev != start ) stream << '-' << prev;
1034  if( iter == list.end() ) break;
1035  stream << ", ";
1036  start = *iter;
1037  }
1038  prev = *( iter++ );
1039  }
1040 
1041  stream << std::endl;
1042 }
1043 
1044 static void print_time( int clk_per_sec, const char* prefix, clock_t ticks, std::ostream& stream )
1045 {
1046  ticks *= clk_per_sec / 100;
1047  clock_t centi = ticks % 100;
1048  clock_t seconds = ticks / 100;
1049  stream << prefix;
1050  if( seconds < 120 )
1051  {
1052  stream << ( ticks / 100 ) << "." << centi << "s" << std::endl;
1053  }
1054  else
1055  {
1056  clock_t minutes = ( seconds / 60 ) % 60;
1057  clock_t hours = ( seconds / 3600 );
1058  seconds %= 60;
1059  if( hours ) stream << hours << "h";
1060  if( minutes ) stream << minutes << "m";
1061  if( seconds || centi ) stream << seconds << "." << centi << "s";
1062  stream << " (" << ( ticks / 100 ) << "." << centi << "s)" << std::endl;
1063  }
1064 }
1065 
1067 
1068 #ifdef WIN32
1069 
1071 {
1072  abs_time = clock();
1073 }
1074 
1075 void write_times( std::ostream& stream )
1076 {
1077  clock_t abs_tm = clock();
1078  print_time( CLOCKS_PER_SEC, " ", abs_tm - abs_time, stream );
1079  abs_time = abs_tm;
1080 }
1081 
1082 #else
1083 
1084 void reset_times()
1085 {
1086  tms timebuf;
1087  abs_time = times( &timebuf );
1088  usr_time = timebuf.tms_utime;
1089  sys_time = timebuf.tms_stime;
1090 }
1091 
1092 void write_times( std::ostream& stream )
1093 {
1094  clock_t usr_tm, sys_tm, abs_tm;
1095  tms timebuf;
1096  abs_tm = times( &timebuf );
1097  usr_tm = timebuf.tms_utime;
1098  sys_tm = timebuf.tms_stime;
1099  print_time( sysconf( _SC_CLK_TCK ), " real: ", abs_tm - abs_time, stream );
1100  print_time( sysconf( _SC_CLK_TCK ), " user: ", usr_tm - usr_time, stream );
1101  print_time( sysconf( _SC_CLK_TCK ), " system: ", sys_tm - sys_time, stream );
1102  abs_time = abs_tm;
1103  usr_time = usr_tm;
1104  sys_time = sys_tm;
1105 }
1106 
1107 #endif
1108 
1109 bool make_opts_string( std::vector< std::string > options, std::string& opts )
1110 {
1111  opts.clear();
1112  if( options.empty() ) return true;
1113 
1114  // choose a separator character
1115  std::vector< std::string >::const_iterator i;
1116  char separator = '\0';
1117  const char* alt_separators = ";+,:\t\n";
1118  for( const char* sep_ptr = alt_separators; *sep_ptr; ++sep_ptr )
1119  {
1120  bool seen = false;
1121  for( i = options.begin(); i != options.end(); ++i )
1122  if( i->find( *sep_ptr, 0 ) != std::string::npos )
1123  {
1124  seen = true;
1125  break;
1126  }
1127  if( !seen )
1128  {
1129  separator = *sep_ptr;
1130  break;
1131  }
1132  }
1133  if( !separator )
1134  {
1135  std::cerr << "Error: cannot find separator character for options string" << std::endl;
1136  return false;
1137  }
1138  if( separator != ';' )
1139  {
1140  opts = ";";
1141  opts += separator;
1142  }
1143 
1144  // concatenate options
1145  i = options.begin();
1146  opts += *i;
1147  for( ++i; i != options.end(); ++i )
1148  {
1149  opts += separator;
1150  opts += *i;
1151  }
1152 
1153  return true;
1154 }
1155 
1157 {
1158  const char iface_name[] = "ReaderWriterSet";
1159  ErrorCode err;
1160  ReaderWriterSet* set = 0;
1162  std::ostream& str = std::cout;
1163 
1164  // get ReaderWriterSet
1165  err = gMB->query_interface( set );
1166  if( err != MB_SUCCESS || !set )
1167  {
1168  std::cerr << "Internal error: Interface \"" << iface_name << "\" not available.\n";
1169  exit( OTHER_ERROR );
1170  }
1171 
1172  // get field width for format description
1173  size_t w = 0;
1174  for( i = set->begin(); i != set->end(); ++i )
1175  if( i->description().length() > w ) w = i->description().length();
1176 
1177  // write table header
1178  str << "Format " << std::setw( w ) << std::left << "Description"
1179  << " Read Write File Name Suffixes\n"
1180  << "------ " << std::setw( w ) << std::setfill( '-' ) << "" << std::setfill( ' ' )
1181  << " ---- ----- ------------------\n";
1182 
1183  // write table data
1184  for( i = set->begin(); i != set->end(); ++i )
1185  {
1186  std::vector< std::string > ext;
1187  i->get_extensions( ext );
1188  str << std::setw( 6 ) << i->name() << " " << std::setw( w ) << std::left << i->description() << " "
1189  << ( i->have_reader() ? " yes" : " no" ) << " " << ( i->have_writer() ? " yes" : " no" ) << " ";
1190  for( std::vector< std::string >::iterator j = ext.begin(); j != ext.end(); ++j )
1191  str << " " << *j;
1192  str << std::endl;
1193  }
1194  str << std::endl;
1195 
1196  gMB->release_interface( set );
1197  exit( 0 );
1198 }
1199 
1200 void remove_entities_from_sets( Interface* gMB, Range& dead_entities, Range& empty_sets )
1201 {
1202  empty_sets.clear();
1203  Range sets;
1204  gMB->get_entities_by_type( 0, MBENTITYSET, sets );
1205  for( Range::iterator i = sets.begin(); i != sets.end(); ++i )
1206  {
1207  Range set_contents;
1208  gMB->get_entities_by_handle( *i, set_contents, false );
1209  set_contents = intersect( set_contents, dead_entities );
1210  gMB->remove_entities( *i, set_contents );
1211  set_contents.clear();
1212  gMB->get_entities_by_handle( *i, set_contents, false );
1213  if( set_contents.empty() ) empty_sets.insert( *i );
1214  }
1215 }
1216 
1217 void remove_from_vector( std::vector< EntityHandle >& vect, const Range& ents_to_remove )
1218 {
1220  std::vector< EntityHandle >::iterator j;
1221  for( i = ents_to_remove.begin(); i != ents_to_remove.end(); ++i )
1222  {
1223  j = std::find( vect.begin(), vect.end(), *i );
1224  if( j != vect.end() ) vect.erase( j );
1225  }
1226 }
1227 
1228 std::string percent_subst( const std::string& s, int val )
1229 {
1230  if( s.empty() ) return s;
1231 
1232  size_t j = s.find( '%' );
1233  if( j == std::string::npos ) return s;
1234 
1235  std::ostringstream st;
1236  st << s.substr( 0, j );
1237  st << val;
1238 
1239  size_t i;
1240  while( ( i = s.find( '%', j + 1 ) ) != std::string::npos )
1241  {
1242  st << s.substr( j, i - j );
1243  st << val;
1244  j = i;
1245  }
1246  st << s.substr( j + 1 );
1247  return st.str();
1248 }
1249 
1250 int process_partition_file( Interface* mb, std::string& metis_partition_file )
1251 {
1252  // how many faces in the file ? how do we make sure it is an mpas file?
1253  // mpas atmosphere files can be downloaded from here
1254  // https://mpas-dev.github.io/atmosphere/atmosphere_meshes.html
1255  Range faces;
1256  ErrorCode rval = mb->get_entities_by_dimension( 0, 2, faces );MB_CHK_ERR( rval );
1257  std::cout << " MPAS model has " << faces.size() << " polygons\n";
1258 
1259  // read the partition file
1260  std::ifstream partfile;
1261  partfile.open( metis_partition_file.c_str() );
1262  std::string line;
1263  std::vector< int > parts;
1264  parts.resize( faces.size(), -1 );
1265  int i = 0;
1266  if( partfile.is_open() )
1267  {
1268  while( getline( partfile, line ) )
1269  {
1270  // cout << line << '\n';
1271  parts[i++] = atoi( line.c_str() );
1272  if( i > (int)faces.size() )
1273  {
1274  std::cerr << " too many lines in partition file \n. bail out \n";
1275  return 1;
1276  }
1277  }
1278  partfile.close();
1279  }
1280  std::vector< int >::iterator pmax = max_element( parts.begin(), parts.end() );
1281  std::vector< int >::iterator pmin = min_element( parts.begin(), parts.end() );
1282  if( *pmin <= -1 )
1283  {
1284  std::cerr << " partition file is incomplete, *pmin is -1 !! \n";
1285  return 1;
1286  }
1287  std::cout << " partitions range: " << *pmin << " " << *pmax << "\n";
1288  Tag part_set_tag;
1289  int dum_id = -1;
1290  rval = mb->tag_get_handle( "PARALLEL_PARTITION", 1, MB_TYPE_INTEGER, part_set_tag, MB_TAG_SPARSE | MB_TAG_CREAT,
1291  &dum_id );MB_CHK_ERR( rval );
1292 
1293  // get any sets already with this tag, and clear them
1294  // remove the parallel partition sets if they exist
1295  Range tagged_sets;
1296  rval = mb->get_entities_by_type_and_tag( 0, MBENTITYSET, &part_set_tag, NULL, 1, tagged_sets, Interface::UNION );MB_CHK_ERR( rval );
1297  if( !tagged_sets.empty() )
1298  {
1299  rval = mb->clear_meshset( tagged_sets );MB_CHK_ERR( rval );
1300  rval = mb->tag_delete_data( part_set_tag, tagged_sets );MB_CHK_ERR( rval );
1301  }
1302  Tag gid;
1303  rval = mb->tag_get_handle( "GLOBAL_ID", gid );MB_CHK_ERR( rval );
1304  int num_sets = *pmax + 1;
1305  if( *pmin != 0 )
1306  {
1307  std::cout << " problem reading parts; min is not 0 \n";
1308  return 1;
1309  }
1310  for( i = 0; i < num_sets; i++ )
1311  {
1312  EntityHandle new_set;
1313  rval = mb->create_meshset( MESHSET_SET, new_set );MB_CHK_ERR( rval );
1314  tagged_sets.insert( new_set );
1315  }
1316  int* dum_ids = new int[num_sets];
1317  for( i = 0; i < num_sets; i++ )
1318  dum_ids[i] = i;
1319 
1320  rval = mb->tag_set_data( part_set_tag, tagged_sets, dum_ids );MB_CHK_ERR( rval );
1321  delete[] dum_ids;
1322 
1323  std::vector< int > gids;
1324  int num_faces = (int)faces.size();
1325  gids.resize( num_faces );
1326  rval = mb->tag_get_data( gid, faces, &gids[0] );MB_CHK_ERR( rval );
1327 
1328  for( int j = 0; j < num_faces; j++ )
1329  {
1330  int eid = gids[j];
1331  EntityHandle eh = faces[j];
1332  int partition = parts[eid - 1];
1333  if( partition < 0 || partition >= num_sets )
1334  {
1335  std::cout << " wrong partition number \n";
1336  return 1;
1337  }
1338  rval = mb->add_entities( tagged_sets[partition], &eh, 1 );MB_CHK_ERR( rval );
1339  }
1340  return 0;
1341 }