Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
mbslavepart.cpp File Reference
#include <iostream>
#include <exception>
#include <cmath>
#include <vector>
#include <string>
#include "moab/ProgOptions.hpp"
#include "moab/Core.hpp"
#include "moab/AdaptiveKDTree.hpp"
#include "moab/BVHTree.hpp"
#include "moab/IntxMesh/IntxUtils.hpp"
+ Include dependency graph for mbslavepart.cpp:

Go to the source code of this file.

Functions

template<typename K , typename V >
static V get_map_value (const std::map< K, V > &m, const K &key, const V &defval)
 
int main (int argc, char *argv[])
 

Function Documentation

◆ get_map_value()

template<typename K , typename V >
static V get_map_value ( const std::map< K, V > &  m,
const K &  key,
const V &  defval 
)
static

Definition at line 29 of file mbslavepart.cpp.

30 {
31  typename std::map< K, V >::const_iterator it = m.find( key );
32  if( it == m.end() )
33  {
34  return defval;
35  }
36  else
37  {
38  return it->second;
39  }
40 }

Referenced by main().

◆ main()

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

Definition at line 42 of file mbslavepart.cpp.

43 {
44  int nprocs = 1, dimension = 3;
45 #ifdef MOAB_HAVE_MPI
46  int proc_id = 0;
47  MPI_Init( &argc, &argv );
48  MPI_Comm_rank( MPI_COMM_WORLD, &proc_id );
49  MPI_Comm_size( MPI_COMM_WORLD, &nprocs );
50 #endif
51 
52  int defaultpart = 0;
53  double tolerance = 1e-2, treetolerance = 1e-13, btolerance = 1e-4;
54  std::string masterfile, slavefile, outfile( "slavemesh.h5m" );
55  bool keepsparts = false;
56  bool use_spherical = false;
57  ProgOptions opts;
58 
59  opts.addOpt< std::string >( "master,m", "Master mesh filename", &masterfile );
60  opts.addOpt< std::string >( "slave,s", "Slave mesh filename", &slavefile );
61  opts.addOpt< std::string >( "output,o", "Output partitioned mesh filename", &outfile );
62  opts.addOpt< int >( "dim,d", "Dimension of entities to use for partitioning", &dimension );
63  opts.addOpt< int >( "defaultpart,p", "Default partition number if target element is not found on source grid",
64  &defaultpart );
65  opts.addOpt< double >( "eps,e", "Tolerance for the point search", &tolerance );
66  opts.addOpt< double >( "beps,b", "Tolerance for the bounding box search", &btolerance );
67  opts.addOpt< void >( "keep,K",
68  "Keep the existing partitions in the slave mesh (use PARALLEL_PARTITION_SLAVE instead)",
69  &keepsparts );
70  opts.addOpt< void >( "spherical", "Hint that the meshes are defined on a spherical surface (Climate problems)",
71  &use_spherical );
72  opts.parseCommandLine( argc, argv );
73 
74  if( masterfile.empty() || slavefile.empty() )
75  {
76  opts.printHelp();
77 #ifdef MOAB_HAVE_MPI
78  MPI_Finalize();
79 #endif
80  exit( 1 );
81  }
82 
84  Core* mbCore = new Core();
85 
86  // Set the read options for parallel file loading
87  const std::string partition_set_name = "PARALLEL_PARTITION";
88  const std::string global_id_name = "GLOBAL_ID";
89  std::vector< std::string > read_opts, write_opts;
90  std::string read_options( "" ), write_options( "" );
91 
92  if( nprocs > 1 )
93  {
94  read_options = "PARALLEL=READ_PART;PARTITION=" + partition_set_name + ";PARALLEL_RESOLVE_SHARED_ENTS";
95  write_options = "PARALLEL=WRITE_PART";
96  }
97 
98  EntityHandle masterfileset, slavefileset;
101 
102  // Load file
103  error = mbCore->load_file( masterfile.c_str(), &masterfileset, read_options.c_str() );MB_CHK_ERR( error );
104  error = mbCore->load_file( slavefile.c_str(), &slavefileset, read_options.c_str() );MB_CHK_ERR( error );
105  // if (error != MB_SUCCESS && size > 1)
106  // {
107  // std::string newread_options = "PARALLEL=BCAST_DELETE;PARALLEL_RESOLVE_SHARED_ENTS";
108  // error = mbCore->load_file(slavefile.c_str(), &slavefileset, newread_options.c_str());
109  // }
110  // else MB_CHK_ERR(error);
111 
112  Tag gidtag = 0, parttag = 0, sparttag = 0;
113  int dum_id = -1;
114  error = mbCore->tag_get_handle( partition_set_name.c_str(), 1, MB_TYPE_INTEGER, parttag,
115  MB_TAG_SPARSE | MB_TAG_CREAT, &dum_id );MB_CHK_ERR( error );
116  gidtag = mbCore->globalId_tag();
117  if( keepsparts )
118  {
119  error = mbCore->tag_get_handle( std::string( partition_set_name + "_SLAVE" ).c_str(), 1, MB_TYPE_INTEGER,
120  sparttag, MB_TAG_CREAT | MB_TAG_SPARSE, &dum_id );MB_CHK_ERR( error );
121  }
122 
123  Range melems, msets, selems, ssets;
124 
125  // Get the partition sets on the master mesh
126  std::map< int, int > mpartvals;
127  error = mbCore->get_entities_by_type_and_tag( masterfileset, MBENTITYSET, &parttag, NULL, 1, msets,
129  if( msets.size() == 0 )
130  {
131  std::cout << "No partition sets found in the master mesh. Quitting..." << std::endl;
132  exit( 1 );
133  }
134 
135  for( unsigned i = 0; i < msets.size(); ++i )
136  {
137  EntityHandle mset = msets[i];
138 
139  moab::Range msetelems;
140  error = mbCore->get_entities_by_dimension( mset, dimension, msetelems );MB_CHK_ERR( error );
141  melems.merge( msetelems );
142 
143  int partID;
144  error = mbCore->tag_get_data( parttag, &mset, 1, &partID );MB_CHK_ERR( error );
145 
146  // Get the global ID and use that as the indicator
147  std::vector< int > gidMelems( msetelems.size() );
148  error = mbCore->tag_get_data( gidtag, msetelems, gidMelems.data() );MB_CHK_ERR( error );
149 
150  for( unsigned j = 0; j < msetelems.size(); ++j )
151  mpartvals[gidMelems[j]] = partID;
152  // mpartvals[msetelems[j]]=partID;
153 #ifdef VERBOSE
154  std::cout << "Part " << partID << " has " << msetelems.size() << " elements." << std::endl;
155 #endif
156  }
157 
158  // Get information about the slave file set
159  error = mbCore->get_entities_by_type_and_tag( slavefileset, MBENTITYSET, &parttag, NULL, 1, ssets,
161  // TODO: expand and add other dimensional elements
162  error = mbCore->get_entities_by_dimension( slavefileset, dimension, selems );MB_CHK_ERR( error );
163 
164  std::cout << "Master (elements, parts) : (" << melems.size() << ", " << msets.size()
165  << "), Slave (elements, parts) : (" << selems.size() << ", " << ssets.size() << ")" << std::endl;
166 
167  // double master_radius = 1.0;
168  double slave_radius = 1.0;
169  std::vector< double > mastercoords;
170  Range masterverts, slaveverts;
171  {
172  error = mbCore->get_entities_by_dimension( masterfileset, 0, masterverts );MB_CHK_ERR( error );
173  error = mbCore->get_entities_by_dimension( slavefileset, 0, slaveverts );MB_CHK_ERR( error );
174  }
175  if( use_spherical )
176  {
177  double points[6];
178  EntityHandle mfrontback[2] = { masterverts[0], masterverts[masterverts.size() - 1] };
179  error = mbCore->get_coords( &mfrontback[0], 2, points );MB_CHK_ERR( error );
180  // master_radius = 0.5 * ( std::sqrt( points[0] * points[0] + points[1] * points[1] + points[2] * points[2] ) +
181  // std::sqrt( points[3] * points[3] + points[4] * points[4] + points[5] * points[5] ) );
182  EntityHandle sfrontback[2] = { slaveverts[0], slaveverts[slaveverts.size() - 1] };
183  error = mbCore->get_coords( &sfrontback[0], 2, points );MB_CHK_ERR( error );
184  slave_radius = 0.5 * ( std::sqrt( points[0] * points[0] + points[1] * points[1] + points[2] * points[2] ) +
185  std::sqrt( points[3] * points[3] + points[4] * points[4] + points[5] * points[5] ) );
186  // Let us rescale both master and slave meshes to a unit sphere
187  error = moab::IntxUtils::ScaleToRadius( mbCore, masterfileset, 1.0 );MB_CHK_ERR( error );
188  error = moab::IntxUtils::ScaleToRadius( mbCore, slavefileset, 1.0 );MB_CHK_ERR( error );
189  }
190 
191  try
192  {
193  std::map< int, moab::Range > spartvals;
194  int npoints_notfound = 0;
195  {
196  FileOptions fopts( ( use_spherical ? "SPHERICAL" : "" ) );
197  moab::EntityHandle tree_root = masterfileset;
198 
199  moab::AdaptiveKDTree tree( mbCore );
200  error = tree.build_tree( melems, &tree_root, &fopts );
201 
202  // moab::BVHTree tree(mbCore);
203  // error = tree.build_tree(melems, &tree_root);MB_CHK_ERR(error);
204 
205  for( size_t ie = 0; ie < selems.size(); ie++ )
206  {
207  moab::EntityHandle selem, leaf;
208  double point[3];
209  selem = selems[ie];
210 
211  // Get the element centroid to be queried
212  error = mbCore->get_coords( &selem, 1, point );MB_CHK_ERR( error );
213 
214  std::vector< moab::EntityHandle > leaf_elems;
215 
216  // Search for the closest source element in the master mesh corresponding
217  // to the target element centroid in the slave mesh
218  error = tree.point_search( point, leaf, treetolerance, btolerance );MB_CHK_ERR( error );
219 
220  // We only care about the dimension that the user specified.
221  // MOAB partitions are ordered by elements anyway.
222  error = mbCore->get_entities_by_dimension( leaf, dimension, leaf_elems, true );MB_CHK_ERR( error );
223 
224  if( leaf != 0 && leaf_elems.size() )
225  {
226 
227  // Now get the master element centroids so that we can compute
228  // the minimum distance to the target point
229  std::vector< double > centroids( leaf_elems.size() * 3 );
230  error = mbCore->get_coords( &leaf_elems[0], leaf_elems.size(), &centroids[0] );MB_CHK_ERR( error );
231 
232  if( !leaf_elems.size() )
233  std::cout << ie << ": "
234  << " No leaf elements found." << std::endl;
235  // else std::cout << ie << " found " << leaf_elems.size() << " leaves for
236  // current point " << std::endl;
237 
238  double dist = 1e10;
239  int pinelem = -1;
240  for( size_t il = 0; il < leaf_elems.size(); ++il )
241  {
242  const double* centroid = &centroids[il * 3];
243  const double locdist = std::pow( point[0] - centroid[0], 2 ) +
244  std::pow( point[1] - centroid[1], 2 ) +
245  std::pow( point[2] - centroid[2], 2 );
246  if( locdist < dist && locdist < 1.0E-2 )
247  {
248  dist = locdist;
249  pinelem = il;
250 
251 #ifdef VERBOSE
252  int gidMelem;
253  error = mbCore->tag_get_data( gidtag, &leaf_elems[il], 1, &gidMelem );MB_CHK_ERR( error );
254  std::cout << "\t Trial leaf " << il << " set " << gidMelem
255  << " and part = " << get_map_value( mpartvals, gidMelem, -1 )
256  << " with distance = " << locdist << std::endl;
257 #endif
258  }
259  }
260 
261  if( pinelem < 0 )
262  {
263 #ifdef VERBOSE
264  std::cout << ie
265  << ": [Error] - Could not find a minimum distance within the "
266  "leaf nodes."
267  << std::endl;
268 #endif
269  npoints_notfound++;
270  spartvals[defaultpart].insert( selems[ie] );
271  }
272  else
273  {
274  int gidMelem;
275  error = mbCore->tag_get_data( gidtag, &leaf_elems[pinelem], 1, &gidMelem );MB_CHK_ERR( error );
276 
277  int mpartid = get_map_value( mpartvals, gidMelem, -1 );
278  if( mpartid < 0 )
279  std::cout << "[WARNING]: Part number for element " << leaf_elems[pinelem]
280  << " with global ID = " << gidMelem << " not found.\n";
281 
282 #ifdef VERBOSE
283  std::cout << "Adding element " << ie << " set " << mpartid << " with distance = " << dist
284  << std::endl;
285 #endif
286  spartvals[mpartid].insert( selems[ie] );
287  }
288  }
289  else
290  {
291 #ifdef VERBOSE
292  std::cout << "[WARNING]: Adding element " << ie << " to default (part) set " << defaultpart
293  << std::endl;
294 #endif
295  spartvals[defaultpart].insert( selems[ie] );
296  }
297  }
298  error = tree.reset_tree();MB_CHK_ERR( error );
299  }
300  if( npoints_notfound )
301  std::cout << "Could not find " << npoints_notfound
302  << " points in the master mesh. Added to defaultpart set = " << defaultpart << std::endl;
303 
304  if( use_spherical )
305  {
306  error = moab::IntxUtils::ScaleToRadius( mbCore, slavefileset, slave_radius );MB_CHK_ERR( error );
307  }
308 
309  error = mbCore->delete_entities( &masterfileset, 1 );MB_CHK_ERR( error );
310  // Find parallel partition sets in the slave mesh - and delete it since we are going to
311  // overwrite the sets
312  if( !keepsparts )
313  {
314  std::cout << "Deleting " << ssets.size() << " sets in the slave mesh" << std::endl;
315  error = mbCore->remove_entities( slavefileset, ssets );MB_CHK_ERR( error );
316  ssets.clear();
317  }
318 
319  size_t ntotslave_elems = 0, ntotslave_parts = 0;
320  for( std::map< int, moab::Range >::iterator it = spartvals.begin(); it != spartvals.end(); ++it )
321  {
322  int partID = it->first;
323  moab::EntityHandle pset;
324  error = mbCore->create_meshset( moab::MESHSET_SET, pset );MB_CHK_ERR( error );
325  error = mbCore->add_entities( pset, it->second );MB_CHK_ERR( error );
326  error = mbCore->add_parent_child( slavefileset, pset );MB_CHK_ERR( error );
327 
328 #ifdef VERBOSE
329  std::cout << "Slave Part " << partID << " has " << it->second.size() << " elements." << std::endl;
330 #endif
331  ntotslave_elems += it->second.size();
332  ntotslave_parts++;
333 
334  if( keepsparts )
335  {
336  error = mbCore->tag_set_data( sparttag, &pset, 1, &partID );MB_CHK_ERR( error );
337  }
338  else
339  {
340  error = mbCore->tag_set_data( parttag, &pset, 1, &partID );MB_CHK_ERR( error );
341  }
342  }
343  std::cout << "Slave mesh: given " << selems.size() << " elements, and assigned " << ntotslave_elems
344  << " elements to " << ntotslave_parts << " parts.\n";
345  assert( ntotslave_elems == selems.size() );
346 
347  // mbCore->print_database();
348 
349  // Write the re-partitioned slave mesh to disk
350  if( nprocs == 1 )
351  {
352  error = mbCore->write_file( "slavemesh.vtk", "VTK", NULL, &slavefileset, 1 );MB_CHK_ERR( error );
353  }
354  error = mbCore->write_file( outfile.c_str(), NULL, write_options.c_str(), &slavefileset, 1 );MB_CHK_ERR( error );
355  // error = mbCore->write_file(outfile.c_str(), NULL,
356  // write_options.c_str());MB_CHK_ERR(error);
357  }
358  catch( std::exception& e )
359  {
360  std::cout << " exception caught during tree initialization " << e.what() << std::endl;
361  }
362  delete mbCore;
363 
364 #ifdef MOAB_HAVE_MPI
365  MPI_Finalize();
366 #endif
367  exit( 0 );
368 }

References moab::Core::add_entities(), moab::Core::add_parent_child(), ProgOptions::addOpt(), moab::AdaptiveKDTree::build_tree(), moab::Range::clear(), moab::Core::create_meshset(), moab::Core::delete_entities(), moab::error(), ErrorCode, moab::Core::get_coords(), moab::Core::get_entities_by_dimension(), moab::Core::get_entities_by_type_and_tag(), get_map_value(), moab::Core::globalId_tag(), moab::Core::load_file(), MB_CHK_ERR, MB_TAG_CREAT, MB_TAG_SPARSE, MB_TYPE_INTEGER, MBENTITYSET, moab::Range::merge(), MESHSET_SET, MESHSET_TRACK_OWNER, ProgOptions::parseCommandLine(), moab::AdaptiveKDTree::point_search(), ProgOptions::printHelp(), moab::Core::remove_entities(), moab::AdaptiveKDTree::reset_tree(), moab::IntxUtils::ScaleToRadius(), moab::Range::size(), moab::Core::tag_get_data(), moab::Core::tag_get_handle(), moab::Core::tag_set_data(), moab::tolerance, moab::Interface::UNION, and moab::Core::write_file().