MOAB: Mesh Oriented datABase
(version 5.4.1)
|
00001 #include <iostream> 00002 #include <moab/Core.hpp> 00003 #include "moab/MergeMesh.hpp" 00004 #include "moab/ProgOptions.hpp" 00005 00006 #ifdef MOAB_HAVE_MPI 00007 #include "moab_mpi.h" 00008 #include "moab/ParallelComm.hpp" 00009 #include "moab/ParallelMergeMesh.hpp" 00010 #endif 00011 00012 using namespace moab; 00013 00014 const char BRIEF_DESC[] = "Merges mesh files or entities in a mesh file. Use available options as desired."; 00015 std::ostringstream LONG_DESC; 00016 00017 int main( int argc, char* argv[] ) 00018 { 00019 #ifdef MOAB_HAVE_MPI 00020 int myID = 0, numprocs = 1; 00021 MPI_Init( &argc, &argv ); 00022 MPI_Comm_rank( MPI_COMM_WORLD, &myID ); 00023 MPI_Comm_size( MPI_COMM_WORLD, &numprocs ); 00024 #endif 00025 bool fsimple = true; // default 00026 bool ffile = false; // parmerge 00027 bool fall = false; // merge all 00028 std::string mtag = ""; // tag based merge 00029 std::string input_file, output_file; 00030 double merge_tol = 1.0e-4; 00031 00032 LONG_DESC << "mbmerge tool has the ability to merge nodes in a mesh. For skin-based merge with " 00033 "multiple" 00034 "files parallel options is also supported." 00035 << std::endl 00036 << "If no method is specified, the default is simple merge. Simple merge case gets " 00037 "all the skins available" 00038 << " in the mesh file and merges the nodes to obtain a conformal mesh. Options to " 00039 "merge all duplicate nodes" 00040 << " and merging based on a specific tag on the nodes are also supported." << std::endl; 00041 00042 ProgOptions opts( LONG_DESC.str(), BRIEF_DESC ); 00043 00044 opts.addOpt< void >( "all,a", "merge all including interior.", &fall ); 00045 opts.addOpt< std::string >( "mergetag name,t", "merge based on nodes that have a specific tag name assigned", 00046 &mtag ); 00047 opts.addOpt< double >( "mergetolerance,e", "merge tolerance, default is 1e-4", &merge_tol ); 00048 opts.addOpt< void >( "simple,s", "simple merge, merge based on skins provided as in the input mesh (Default)", 00049 &fsimple ); 00050 opts.addRequiredArg< std::string >( "input_file", "Input file to be merged", &input_file ); 00051 opts.addRequiredArg< std::string >( "output_file", "Output mesh file name with extension", &output_file ); 00052 #ifdef MOAB_HAVE_MPI 00053 opts.addOpt< void >( "file,f", 00054 "files based merge using skin for individual meshes. Input file is a file containing names" 00055 "of mesh files to be merged (each on a different line). Only works with parallel build", 00056 &ffile ); 00057 #endif 00058 opts.parseCommandLine( argc, argv ); 00059 00060 moab::Core* mb = new moab::Core(); 00061 moab::ErrorCode rval; 00062 00063 if( mtag != "" ) 00064 { 00065 rval = mb->load_mesh( input_file.c_str() ); 00066 if( rval != moab::MB_SUCCESS ) 00067 { 00068 std::cerr << "Error Opening Mesh File " << input_file << std::endl; 00069 return 1; 00070 } 00071 else 00072 { 00073 std::cout << "Read input mesh file: " << input_file << std::endl; 00074 } 00075 int dim = 0; 00076 moab::Range verts; 00077 mb->get_entities_by_dimension( 0, dim, verts ); 00078 if( rval != moab::MB_SUCCESS ) 00079 { 00080 std::cerr << "failed to get entities by dimension" << std::endl; 00081 return 1; 00082 } 00083 Tag tag_for_merge; 00084 rval = mb->tag_get_handle( mtag.c_str(), tag_for_merge ); 00085 if( rval != moab::MB_SUCCESS ) 00086 { 00087 std::cerr << "unable to get tag: " << mtag << " specified" << std::endl; 00088 return 1; 00089 } 00090 MergeMesh mm( mb ); 00091 rval = mm.merge_using_integer_tag( verts, tag_for_merge ); 00092 if( rval != moab::MB_SUCCESS ) 00093 { 00094 std::cerr << "error in routine merge using integer tag" << std::endl; 00095 return 1; 00096 } 00097 rval = mb->write_file( output_file.c_str() ); 00098 if( rval != moab::MB_SUCCESS ) 00099 { 00100 std::cerr << "Error Writing Mesh File " << output_file << std::endl; 00101 return 1; 00102 } 00103 else 00104 { 00105 std::cout << "Wrote output mesh file: " << output_file << std::endl; 00106 } 00107 } 00108 00109 else if( fall == true ) 00110 { 00111 00112 rval = mb->load_mesh( input_file.c_str() ); 00113 if( rval != moab::MB_SUCCESS ) 00114 { 00115 std::cerr << "Error Opening Mesh File " << input_file << std::endl; 00116 return 1; 00117 } 00118 else 00119 { 00120 std::cout << "Read input mesh file: " << input_file << std::endl; 00121 } 00122 MergeMesh mm( mb ); 00123 rval = mm.merge_all( 0, merge_tol ); // root set 00124 if( rval != moab::MB_SUCCESS ) 00125 { 00126 std::cerr << "error in merge_all routine" << std::endl; 00127 return 1; 00128 } 00129 rval = mb->write_file( output_file.c_str() ); 00130 if( rval != moab::MB_SUCCESS ) 00131 { 00132 std::cerr << "Error Writing Mesh File " << output_file << std::endl; 00133 return 1; 00134 } 00135 else 00136 { 00137 std::cout << "Wrote output mesh file: " << output_file << std::endl; 00138 } 00139 } 00140 else if( ffile == true ) 00141 { 00142 /* 00143 Parmerge 00144 Takes multiple mesh files and merges them into a single output file. 00145 This is a driver for ParallelMergeMesh 00146 Does not currently work if #procs > #meshfiles 00147 00148 <input_file> is text file containing each mesh file on a line 00149 i.e. 00150 00151 /my/path/file1 00152 /my/path/file2 00153 ... 00154 /my/path/fileN 00155 00156 <output_file> file is a single file where the entire mesh is written to 00157 It must be of type ".h5m" 00158 00159 <tolerance> is the merging tolerance 00160 00161 Typical usage of: 00162 mpirun -n <#procs> parmerge <input_file> <output_file> <tolerance> 00163 */ 00164 // Read in files from input files 00165 // Round robin distribution of reading meshes 00166 #ifdef MOAB_HAVE_MPI 00167 std::ifstream file( input_file.c_str() ); 00168 if( file.is_open() ) 00169 { 00170 std::string line; 00171 int count = 0; 00172 // Read each line 00173 while( file.good() ) 00174 { 00175 getline( file, line ); 00176 if( myID == count && line != "" ) 00177 { 00178 // Read in the file 00179 rval = mb->load_mesh( line.c_str() ); 00180 if( rval != moab::MB_SUCCESS ) 00181 { 00182 std::cerr << "Error Opening Mesh File " << line << std::endl; 00183 MPI_Abort( MPI_COMM_WORLD, 1 ); 00184 file.close(); 00185 return 1; 00186 } 00187 else 00188 { 00189 std::cout << "Read input mesh file: " << line << std::endl; 00190 } 00191 } 00192 count = ( count + 1 ) % numprocs; 00193 } 00194 file.close(); 00195 } 00196 else 00197 { 00198 std::cerr << "Error Opening Input File " << input_file << std::endl; 00199 MPI_Abort( MPI_COMM_WORLD, 1 ); 00200 return 1; 00201 } 00202 00203 // Get a pcomm object 00204 moab::ParallelComm* pc = new moab::ParallelComm( mb, MPI_COMM_WORLD ); 00205 00206 // Call the resolve parallel function 00207 moab::ParallelMergeMesh pm( pc, merge_tol ); 00208 rval = pm.merge(); 00209 if( rval != moab::MB_SUCCESS ) 00210 { 00211 std::cerr << "Merge Failed" << std::endl; 00212 MPI_Abort( MPI_COMM_WORLD, 1 ); 00213 return 1; 00214 } 00215 00216 // Write out the file 00217 rval = mb->write_file( output_file.c_str(), 0, "PARALLEL=WRITE_PART" ); 00218 if( rval != moab::MB_SUCCESS ) 00219 { 00220 std::cerr << "Writing output file failed. Code:"; 00221 // Temporary File error info. 00222 std::cerr << mb->get_error_string( rval ) << std::endl; 00223 std::string foo = ""; 00224 mb->get_last_error( foo ); 00225 std::cerr << "File Error: " << foo << std::endl; 00226 return 1; 00227 } 00228 else if( myID == 0 ) 00229 { 00230 std::cout << "Wrote output mesh file: " << output_file << std::endl; 00231 } 00232 00233 // The barrier may be necessary to stop items from being deleted when needed 00234 // But probably not necessary 00235 MPI_Barrier( MPI_COMM_WORLD ); 00236 00237 delete pc; 00238 #endif 00239 } 00240 else if( fsimple == true ) 00241 { 00242 rval = mb->load_mesh( input_file.c_str() ); 00243 if( rval != moab::MB_SUCCESS ) 00244 { 00245 std::cerr << "Error Opening Mesh File " << input_file << std::endl; 00246 return 1; 00247 } 00248 else 00249 { 00250 std::cout << "Read input mesh file: " << input_file << std::endl; 00251 } 00252 int dim = 3; 00253 moab::Range ents; 00254 mb->get_entities_by_dimension( 0, dim, ents ); 00255 if( rval != moab::MB_SUCCESS ) 00256 { 00257 std::cerr << "error getting entities by dimension" << std::endl; 00258 return 1; 00259 } 00260 MergeMesh mm( mb ); 00261 rval = mm.merge_entities( ents, merge_tol ); 00262 if( rval != moab::MB_SUCCESS ) 00263 { 00264 std::cerr << "error in merge entities routine" << std::endl; 00265 return 1; 00266 } 00267 rval = mb->write_file( output_file.c_str() ); 00268 if( rval != moab::MB_SUCCESS ) 00269 { 00270 std::cerr << " Writing Mesh File " << output_file << std::endl; 00271 return 1; 00272 } 00273 else 00274 { 00275 std::cout << "Wrote output mesh file: " << output_file << std::endl; 00276 } 00277 } 00278 else 00279 { 00280 std::cerr << " Unhandled option " << std::endl; 00281 return 1; 00282 } 00283 00284 delete mb; 00285 #ifdef MOAB_HAVE_MPI 00286 MPI_Finalize(); 00287 #endif 00288 00289 return 0; 00290 }