Loading [MathJax]/extensions/tex2jax.js
Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
merge.cpp File Reference
#include <iostream>
#include <moab/Core.hpp>
#include "moab/MergeMesh.hpp"
#include "moab/ProgOptions.hpp"
+ Include dependency graph for merge.cpp:

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

const char BRIEF_DESC [] = "Merges mesh files or entities in a mesh file. Use available options as desired."
 
std::ostringstream LONG_DESC
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 17 of file merge.cpp.

18 { 19 #ifdef MOAB_HAVE_MPI 20  int myID = 0, numprocs = 1; 21  MPI_Init( &argc, &argv ); 22  MPI_Comm_rank( MPI_COMM_WORLD, &myID ); 23  MPI_Comm_size( MPI_COMM_WORLD, &numprocs ); 24 #endif 25  bool fsimple = true; // default 26  bool ffile = false; // parmerge 27  bool fall = false; // merge all 28  std::string mtag = ""; // tag based merge 29  std::string input_file, output_file; 30  double merge_tol = 1.0e-4; 31  32  LONG_DESC << "mbmerge tool has the ability to merge nodes in a mesh. For skin-based merge with " 33  "multiple" 34  "files parallel options is also supported." 35  << std::endl 36  << "If no method is specified, the default is simple merge. Simple merge case gets " 37  "all the skins available" 38  << " in the mesh file and merges the nodes to obtain a conformal mesh. Options to " 39  "merge all duplicate nodes" 40  << " and merging based on a specific tag on the nodes are also supported." << std::endl; 41  42  ProgOptions opts( LONG_DESC.str(), BRIEF_DESC ); 43  44  opts.addOpt< void >( "all,a", "merge all including interior.", &fall ); 45  opts.addOpt< std::string >( "mergetag name,t", "merge based on nodes that have a specific tag name assigned", 46  &mtag ); 47  opts.addOpt< double >( "mergetolerance,e", "merge tolerance, default is 1e-4", &merge_tol ); 48  opts.addOpt< void >( "simple,s", "simple merge, merge based on skins provided as in the input mesh (Default)", 49  &fsimple ); 50  opts.addRequiredArg< std::string >( "input_file", "Input file to be merged", &input_file ); 51  opts.addRequiredArg< std::string >( "output_file", "Output mesh file name with extension", &output_file ); 52 #ifdef MOAB_HAVE_MPI 53  opts.addOpt< void >( "file,f", 54  "files based merge using skin for individual meshes. Input file is a file containing names" 55  "of mesh files to be merged (each on a different line). Only works with parallel build", 56  &ffile ); 57 #endif 58  opts.parseCommandLine( argc, argv ); 59  60  moab::Core* mb = new moab::Core(); 61  moab::ErrorCode rval; 62  63  if( mtag != "" ) 64  { 65  rval = mb->load_mesh( input_file.c_str() ); 66  if( rval != moab::MB_SUCCESS ) 67  { 68  std::cerr << "Error Opening Mesh File " << input_file << std::endl; 69  return 1; 70  } 71  else 72  { 73  std::cout << "Read input mesh file: " << input_file << std::endl; 74  } 75  int dim = 0; 76  moab::Range verts; 77  mb->get_entities_by_dimension( 0, dim, verts ); 78  if( rval != moab::MB_SUCCESS ) 79  { 80  std::cerr << "failed to get entities by dimension" << std::endl; 81  return 1; 82  } 83  Tag tag_for_merge; 84  rval = mb->tag_get_handle( mtag.c_str(), tag_for_merge ); 85  if( rval != moab::MB_SUCCESS ) 86  { 87  std::cerr << "unable to get tag: " << mtag << " specified" << std::endl; 88  return 1; 89  } 90  MergeMesh mm( mb ); 91  rval = mm.merge_using_integer_tag( verts, tag_for_merge ); 92  if( rval != moab::MB_SUCCESS ) 93  { 94  std::cerr << "error in routine merge using integer tag" << std::endl; 95  return 1; 96  } 97  rval = mb->write_file( output_file.c_str() ); 98  if( rval != moab::MB_SUCCESS ) 99  { 100  std::cerr << "Error Writing Mesh File " << output_file << std::endl; 101  return 1; 102  } 103  else 104  { 105  std::cout << "Wrote output mesh file: " << output_file << std::endl; 106  } 107  } 108  109  else if( fall == true ) 110  { 111  112  rval = mb->load_mesh( input_file.c_str() ); 113  if( rval != moab::MB_SUCCESS ) 114  { 115  std::cerr << "Error Opening Mesh File " << input_file << std::endl; 116  return 1; 117  } 118  else 119  { 120  std::cout << "Read input mesh file: " << input_file << std::endl; 121  } 122  MergeMesh mm( mb ); 123  rval = mm.merge_all( 0, merge_tol ); // root set 124  if( rval != moab::MB_SUCCESS ) 125  { 126  std::cerr << "error in merge_all routine" << std::endl; 127  return 1; 128  } 129  rval = mb->write_file( output_file.c_str() ); 130  if( rval != moab::MB_SUCCESS ) 131  { 132  std::cerr << "Error Writing Mesh File " << output_file << std::endl; 133  return 1; 134  } 135  else 136  { 137  std::cout << "Wrote output mesh file: " << output_file << std::endl; 138  } 139  } 140  else if( ffile == true ) 141  { 142  /* 143  Parmerge 144  Takes multiple mesh files and merges them into a single output file. 145  This is a driver for ParallelMergeMesh 146  Does not currently work if #procs > #meshfiles 147  148  <input_file> is text file containing each mesh file on a line 149  i.e. 150  151  /my/path/file1 152  /my/path/file2 153  ... 154  /my/path/fileN 155  156  <output_file> file is a single file where the entire mesh is written to 157  It must be of type ".h5m" 158  159  <tolerance> is the merging tolerance 160  161  Typical usage of: 162  mpirun -n <#procs> parmerge <input_file> <output_file> <tolerance> 163  */ 164  // Read in files from input files 165  // Round robin distribution of reading meshes 166 #ifdef MOAB_HAVE_MPI 167  std::ifstream file( input_file.c_str() ); 168  if( file.is_open() ) 169  { 170  std::string line; 171  int count = 0; 172  // Read each line 173  while( file.good() ) 174  { 175  getline( file, line ); 176  if( myID == count && line != "" ) 177  { 178  // Read in the file 179  rval = mb->load_mesh( line.c_str() ); 180  if( rval != moab::MB_SUCCESS ) 181  { 182  std::cerr << "Error Opening Mesh File " << line << std::endl; 183  MPI_Abort( MPI_COMM_WORLD, 1 ); 184  file.close(); 185  return 1; 186  } 187  else 188  { 189  std::cout << "Read input mesh file: " << line << std::endl; 190  } 191  } 192  count = ( count + 1 ) % numprocs; 193  } 194  file.close(); 195  } 196  else 197  { 198  std::cerr << "Error Opening Input File " << input_file << std::endl; 199  MPI_Abort( MPI_COMM_WORLD, 1 ); 200  return 1; 201  } 202  203  // Get a pcomm object 204  moab::ParallelComm* pc = new moab::ParallelComm( mb, MPI_COMM_WORLD ); 205  206  // Call the resolve parallel function 207  moab::ParallelMergeMesh pm( pc, merge_tol ); 208  rval = pm.merge(); 209  if( rval != moab::MB_SUCCESS ) 210  { 211  std::cerr << "Merge Failed" << std::endl; 212  MPI_Abort( MPI_COMM_WORLD, 1 ); 213  return 1; 214  } 215  216  // Write out the file 217  rval = mb->write_file( output_file.c_str(), 0, "PARALLEL=WRITE_PART" ); 218  if( rval != moab::MB_SUCCESS ) 219  { 220  std::cerr << "Writing output file failed. Code:"; 221  // Temporary File error info. 222  std::cerr << mb->get_error_string( rval ) << std::endl; 223  std::string foo = ""; 224  mb->get_last_error( foo ); 225  std::cerr << "File Error: " << foo << std::endl; 226  return 1; 227  } 228  else if( myID == 0 ) 229  { 230  std::cout << "Wrote output mesh file: " << output_file << std::endl; 231  } 232  233  // The barrier may be necessary to stop items from being deleted when needed 234  // But probably not necessary 235  MPI_Barrier( MPI_COMM_WORLD ); 236  237  delete pc; 238 #endif 239  } 240  else if( fsimple == true ) 241  { 242  rval = mb->load_mesh( input_file.c_str() ); 243  if( rval != moab::MB_SUCCESS ) 244  { 245  std::cerr << "Error Opening Mesh File " << input_file << std::endl; 246  return 1; 247  } 248  else 249  { 250  std::cout << "Read input mesh file: " << input_file << std::endl; 251  } 252  int dim = 3; 253  moab::Range ents; 254  mb->get_entities_by_dimension( 0, dim, ents ); 255  if( rval != moab::MB_SUCCESS ) 256  { 257  std::cerr << "error getting entities by dimension" << std::endl; 258  return 1; 259  } 260  MergeMesh mm( mb ); 261  rval = mm.merge_entities( ents, merge_tol ); 262  if( rval != moab::MB_SUCCESS ) 263  { 264  std::cerr << "error in merge entities routine" << std::endl; 265  return 1; 266  } 267  rval = mb->write_file( output_file.c_str() ); 268  if( rval != moab::MB_SUCCESS ) 269  { 270  std::cerr << " Writing Mesh File " << output_file << std::endl; 271  return 1; 272  } 273  else 274  { 275  std::cout << "Wrote output mesh file: " << output_file << std::endl; 276  } 277  } 278  else 279  { 280  std::cerr << " Unhandled option " << std::endl; 281  return 1; 282  } 283  284  delete mb; 285 #ifdef MOAB_HAVE_MPI 286  MPI_Finalize(); 287 #endif 288  289  return 0; 290 }

References ProgOptions::addOpt(), ProgOptions::addRequiredArg(), BRIEF_DESC, dim, ErrorCode, moab::Core::get_entities_by_dimension(), moab::Core::get_error_string(), moab::Core::get_last_error(), input_file, moab::Core::load_mesh(), LONG_DESC, mb, MB_SUCCESS, moab::ParallelMergeMesh::merge(), moab::MergeMesh::merge_all(), moab::MergeMesh::merge_entities(), moab::MergeMesh::merge_using_integer_tag(), ProgOptions::parseCommandLine(), moab::Core::tag_get_handle(), and moab::Core::write_file().

Variable Documentation

◆ BRIEF_DESC

const char BRIEF_DESC[] = "Merges mesh files or entities in a mesh file. Use available options as desired."

Definition at line 14 of file merge.cpp.

Referenced by main().

◆ LONG_DESC

std::ostringstream LONG_DESC

Definition at line 15 of file merge.cpp.

Referenced by main().