MOAB: Mesh Oriented datABase  (version 5.5.0)
parmerge.cpp File Reference
#include "moab/ParallelMergeMesh.hpp"
#include "moab/Core.hpp"
#include "moab/Range.hpp"
#include "moab_mpi.h"
#include <iostream>
#include <string>
#include <cstdlib>
#include "moab/ParallelComm.hpp"
#include "MBParallelConventions.h"
#include <fstream>
#include <sstream>
#include "moab/Skinner.hpp"
+ Include dependency graph for parmerge.cpp:

Go to the source code of this file.

Macros

#define PrintInfo   false
 

Functions

void print_output (moab::ParallelComm *pc, moab::Core *mb, int numprocs, int myID, bool perform)
 
int main (int argc, char *argv[])
 

Macro Definition Documentation

◆ PrintInfo

#define PrintInfo   false

Definition at line 21 of file parmerge.cpp.

Function Documentation

◆ main()

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

Definition at line 48 of file parmerge.cpp.

49 {
50  // Check argument count
51  if( argc != 4 )
52  {
53  std::cerr << "Usage: " << argv[0] << " <inputfile> <outputfile> <tolerance>" << std::endl;
54  return 1;
55  }
56  // Check the output file extension
57  std::string outfile( argv[2] );
58  if( outfile.compare( outfile.size() - 4, 4, ".h5m" ) != 0 )
59  {
60  std::cerr << "Invalid Parallel Output File" << std::endl;
61  return 1;
62  }
63 
64  // Read in tolerance
65  double epsilon;
66  if( !( std::istringstream( argv[3] ) >> epsilon ) )
67  {
68  std::cerr << "Unable to parse tolerance" << std::endl;
69  return 1;
70  }
71 
72  // Initialize MPI
73  int numprocs, myID;
74  MPI_Init( &argc, &argv );
75  MPI_Comm_size( MPI_COMM_WORLD, &numprocs );
76  MPI_Comm_rank( MPI_COMM_WORLD, &myID );
77 
78  // Read in files from input files
79  // Round robin distribution of reading meshes
80  moab::Core* mb = new moab::Core();
81  moab::ErrorCode rval;
82  std::ifstream file( argv[1] );
83  if( file.is_open() )
84  {
85  std::string line;
86  int count = 0;
87  // Read each line
88  while( file.good() )
89  {
90  getline( file, line );
91  if( myID == count && line != "" )
92  {
93  // Read in the file
94  rval = mb->load_mesh( line.c_str() );
95  if( rval != moab::MB_SUCCESS )
96  {
97  std::cerr << "Error Opening Mesh File " << line << std::endl;
98  MPI_Abort( MPI_COMM_WORLD, 1 );
99  file.close();
100  return 1;
101  }
102  }
103  count = ( count + 1 ) % numprocs;
104  }
105  file.close();
106  }
107  else
108  {
109  std::cerr << "Error Opening Input File " << argv[1] << std::endl;
110  MPI_Abort( MPI_COMM_WORLD, 1 );
111  return 1;
112  }
113 
114  // Get a pcomm object
116 
117  // Call the resolve parallel function
118  moab::ParallelMergeMesh pm( pc, epsilon );
119  rval = pm.merge();
120  if( rval != moab::MB_SUCCESS )
121  {
122  std::cerr << "Merge Failed" << std::endl;
123  MPI_Abort( MPI_COMM_WORLD, 1 );
124  return 1;
125  }
126 
127  print_output( pc, mb, myID, numprocs, PrintInfo );
128 
129  // Write out the file
130  rval = mb->write_file( outfile.c_str(), 0, "PARALLEL=WRITE_PART" );
131  if( rval != moab::MB_SUCCESS )
132  {
133  std::cerr << "Writing output file failed. Code:";
134  // Temporary File error info.
135  std::cerr << mb->get_error_string( rval ) << std::endl;
136  std::string foo = "";
137  mb->get_last_error( foo );
138  std::cerr << "File Error: " << foo << std::endl;
139  return 1;
140  }
141 
142  // The barrier may be necessary to stop items from being deleted when needed
143  // But probably not necessary
144  MPI_Barrier( MPI_COMM_WORLD );
145 
146  delete pc;
147  delete mb;
148  MPI_Finalize();
149 
150  return 0;
151 }

References ErrorCode, mb, MB_SUCCESS, moab::ParallelMergeMesh::merge(), MPI_COMM_WORLD, outfile, print_output(), and PrintInfo.

◆ print_output()

void print_output ( moab::ParallelComm pc,
moab::Core mb,
int  numprocs,
int  myID,
bool  perform 
)

Definition at line 155 of file parmerge.cpp.

156 {
157  moab::Range ents, skin;
158  int o_ct = 0, no_ct = 0, tmp = 0, o_tot = 0, no_tot = 0;
159  if( perform )
160  {
161  if( myID == 0 ) std::cout << "------------------------------------------" << std::endl;
162  // Check the count of total vertices
163  mb->get_entities_by_dimension( 0, 0, ents );
164  for( moab::Range::iterator rit = ents.begin(); rit != ents.end(); ++rit )
165  {
166  pc->get_owner( *rit, tmp );
167  if( tmp == myID )
168  {
169  o_ct++;
170  }
171  }
172  MPI_Reduce( &o_ct, &o_tot, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD );
173  if( myID == 0 )
174  {
175  std::cout << "There are " << o_tot << " vertices." << std::endl;
176  std::cout << "------------------------------------------" << std::endl;
177  }
178  // Check the count of owned and not owned skin faces.
179  // owned-not owned == total skin faces
180  moab::Skinner skinner( mb );
181  o_ct = 0;
182  no_ct = 0;
183  o_tot = 0;
184  no_tot = 0;
185  skin.clear();
186  ents.clear();
187  mb->get_entities_by_dimension( 0, 3, ents );
188  skinner.find_skin( 0, ents, 2, skin );
189  for( moab::Range::iterator s_rit = skin.begin(); s_rit != skin.end(); ++s_rit )
190  {
191  pc->get_owner( *s_rit, tmp );
192  if( tmp == myID )
193  {
194  o_ct++;
195  }
196  else
197  {
198  no_ct++;
199  }
200  }
201  MPI_Reduce( &o_ct, &o_tot, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD );
202  MPI_Reduce( &no_ct, &no_tot, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD );
203  if( myID == 0 )
204  {
205  std::cout << "There are " << o_tot << " owned skin faces." << std::endl;
206  std::cout << "There are " << no_tot << " not owned skin faces." << std::endl;
207  std::cout << "The difference (Global Skin Faces) is " << ( o_tot - no_tot ) << "." << std::endl;
208  std::cout << "------------------------------------------" << std::endl;
209  }
210  }
211 }

References moab::Range::begin(), moab::Range::clear(), moab::Range::end(), moab::Skinner::find_skin(), moab::ParallelComm::get_owner(), mb, MPI_COMM_WORLD, and skin().

Referenced by main().