MOAB: Mesh Oriented datABase  (version 5.3.1)
mbslavepart.cpp
Go to the documentation of this file.
00001 //
00002 // Usage:
00003 // tools/mbslavepart -d 2 -m mpas/x1.2562.grid.h5m -s mpas/x1.10242.grid.h5m -o mpas_slave.h5m -e
00004 // 1e-8 -b 1e-6 -O
00005 //
00006 #include <iostream>
00007 #include <exception>
00008 #include <cmath>
00009 #include <vector>
00010 #include <string>
00011 
00012 #include "moab/ProgOptions.hpp"
00013 #include "moab/Core.hpp"
00014 #include "moab/AdaptiveKDTree.hpp"
00015 #include "moab/BVHTree.hpp"
00016 
00017 #include "moab/IntxMesh/IntxUtils.hpp"
00018 
00019 #ifdef MOAB_HAVE_MPI
00020 #include "moab_mpi.h"
00021 #include "moab/ParallelComm.hpp"
00022 #include "MBParallelConventions.h"
00023 #endif
00024 
00025 using namespace moab;
00026 
00027 // A function to get the non-default value from a std::map
00028 template < typename K, typename V >
00029 static V get_map_value( const std::map< K, V >& m, const K& key, const V& defval )
00030 {
00031     typename std::map< K, V >::const_iterator it = m.find( key );
00032     if( it == m.end() ) { return defval; }
00033     else
00034     {
00035         return it->second;
00036     }
00037 }
00038 
00039 int main( int argc, char* argv[] )
00040 {
00041     int nprocs = 1, dimension = 3;
00042 #ifdef MOAB_HAVE_MPI
00043     int proc_id = 0;
00044     MPI_Init( &argc, &argv );
00045     MPI_Comm_rank( MPI_COMM_WORLD, &proc_id );
00046     MPI_Comm_size( MPI_COMM_WORLD, &nprocs );
00047 #endif
00048 
00049     int defaultpart  = 0;
00050     double tolerance = 1e-2, treetolerance = 1e-13, btolerance = 1e-4;
00051     std::string masterfile, slavefile, outfile( "slavemesh.h5m" );
00052     bool keepsparts    = false;
00053     bool use_spherical = false;
00054     ProgOptions opts;
00055 
00056     opts.addOpt< std::string >( "master,m", "Master mesh filename", &masterfile );
00057     opts.addOpt< std::string >( "slave,s", "Slave mesh filename", &slavefile );
00058     opts.addOpt< std::string >( "output,o", "Output partitioned mesh filename", &outfile );
00059     opts.addOpt< int >( "dim,d", "Dimension of entities to use for partitioning", &dimension );
00060     opts.addOpt< int >( "defaultpart,p", "Default partition number if target element is not found on source grid",
00061                         &defaultpart );
00062     opts.addOpt< double >( "eps,e", "Tolerance for the point search", &tolerance );
00063     opts.addOpt< double >( "beps,b", "Tolerance for the bounding box search", &btolerance );
00064     opts.addOpt< void >( "keep,K",
00065                          "Keep the existing partitions in the slave mesh (use PARALLEL_PARTITION_SLAVE instead)",
00066                          &keepsparts );
00067     opts.addOpt< void >( "spherical", "Hint that the meshes are defined on a spherical surface (Climate problems)",
00068                          &use_spherical );
00069     opts.parseCommandLine( argc, argv );
00070 
00071     if( masterfile.empty() || slavefile.empty() )
00072     {
00073         opts.printHelp();
00074 #ifdef MOAB_HAVE_MPI
00075         MPI_Finalize();
00076 #endif
00077         exit( 1 );
00078     }
00079 
00080     ErrorCode error;
00081     Core* mbCore = new Core();
00082 
00083     // Set the read options for parallel file loading
00084     const std::string partition_set_name = "PARALLEL_PARTITION";
00085     const std::string global_id_name     = "GLOBAL_ID";
00086     std::vector< std::string > read_opts, write_opts;
00087     std::string read_options( "" ), write_options( "" );
00088 
00089     if( nprocs > 1 )
00090     {
00091         read_options  = "PARALLEL=READ_PART;PARTITION=" + partition_set_name + ";PARALLEL_RESOLVE_SHARED_ENTS";
00092         write_options = "PARALLEL=WRITE_PART";
00093     }
00094 
00095     EntityHandle masterfileset, slavefileset;
00096     error = mbCore->create_meshset( moab::MESHSET_TRACK_OWNER | moab::MESHSET_SET, masterfileset );MB_CHK_ERR( error );
00097     error = mbCore->create_meshset( moab::MESHSET_TRACK_OWNER | moab::MESHSET_SET, slavefileset );MB_CHK_ERR( error );
00098 
00099     // Load file
00100     error = mbCore->load_file( masterfile.c_str(), &masterfileset, read_options.c_str() );MB_CHK_ERR( error );
00101     error = mbCore->load_file( slavefile.c_str(), &slavefileset, read_options.c_str() );MB_CHK_ERR( error );
00102     // if (error != MB_SUCCESS && size > 1)
00103     // {
00104     //   std::string newread_options = "PARALLEL=BCAST_DELETE;PARALLEL_RESOLVE_SHARED_ENTS";
00105     //   error = mbCore->load_file(slavefile.c_str(), &slavefileset, newread_options.c_str());
00106     // }
00107     // else MB_CHK_ERR(error);
00108 
00109     Tag gidtag = 0, parttag = 0, sparttag = 0;
00110     int dum_id = -1;
00111     error      = mbCore->tag_get_handle( partition_set_name.c_str(), 1, MB_TYPE_INTEGER, parttag,
00112                                     MB_TAG_SPARSE | MB_TAG_CREAT, &dum_id );MB_CHK_ERR( error );
00113     gidtag = mbCore->globalId_tag();
00114     if( keepsparts )
00115     {
00116         error = mbCore->tag_get_handle( std::string( partition_set_name + "_SLAVE" ).c_str(), 1, MB_TYPE_INTEGER,
00117                                         sparttag, MB_TAG_CREAT | MB_TAG_SPARSE, &dum_id );MB_CHK_ERR( error );
00118     }
00119 
00120     Range melems, msets, selems, ssets;
00121 
00122     // Get the partition sets on the master mesh
00123     std::map< int, int > mpartvals;
00124     error = mbCore->get_entities_by_type_and_tag( masterfileset, MBENTITYSET, &parttag, NULL, 1, msets,
00125                                                   moab::Interface::UNION, true );MB_CHK_ERR( error );
00126     if( msets.size() == 0 )
00127     {
00128         std::cout << "No partition sets found in the master mesh. Quitting..." << std::endl;
00129         exit( 1 );
00130     }
00131 
00132     for( unsigned i = 0; i < msets.size(); ++i )
00133     {
00134         EntityHandle mset = msets[i];
00135 
00136         moab::Range msetelems;
00137         error = mbCore->get_entities_by_dimension( mset, dimension, msetelems );MB_CHK_ERR( error );
00138         melems.merge( msetelems );
00139 
00140         int partID;
00141         error = mbCore->tag_get_data( parttag, &mset, 1, &partID );MB_CHK_ERR( error );
00142 
00143         // Get the global ID and use that as the indicator
00144         std::vector< int > gidMelems( msetelems.size() );
00145         error = mbCore->tag_get_data( gidtag, msetelems, gidMelems.data() );MB_CHK_ERR( error );
00146 
00147         for( unsigned j = 0; j < msetelems.size(); ++j )
00148             mpartvals[gidMelems[j]] = partID;
00149             // mpartvals[msetelems[j]]=partID;
00150 #ifdef VERBOSE
00151         std::cout << "Part " << partID << " has " << msetelems.size() << " elements." << std::endl;
00152 #endif
00153     }
00154 
00155     // Get information about the slave file set
00156     error = mbCore->get_entities_by_type_and_tag( slavefileset, MBENTITYSET, &parttag, NULL, 1, ssets,
00157                                                   moab::Interface::UNION );MB_CHK_ERR( error );
00158     // TODO: expand and add other dimensional elements
00159     error = mbCore->get_entities_by_dimension( slavefileset, dimension, selems );MB_CHK_ERR( error );
00160 
00161     std::cout << "Master (elements, parts) : (" << melems.size() << ", " << msets.size()
00162               << "), Slave (elements, parts) : (" << selems.size() << ", " << ssets.size() << ")" << std::endl;
00163 
00164     double master_radius = 1.0, slave_radius = 1.0;
00165     std::vector< double > mastercoords;
00166     Range masterverts, slaveverts;
00167     {
00168         error = mbCore->get_entities_by_dimension( masterfileset, 0, masterverts );MB_CHK_ERR( error );
00169         error = mbCore->get_entities_by_dimension( slavefileset, 0, slaveverts );MB_CHK_ERR( error );
00170     }
00171     if( use_spherical )
00172     {
00173         double points[6];
00174         EntityHandle mfrontback[2] = { masterverts[0], masterverts[masterverts.size() - 1] };
00175         error                      = mbCore->get_coords( &mfrontback[0], 2, points );MB_CHK_ERR( error );
00176         master_radius = 0.5 * ( std::sqrt( points[0] * points[0] + points[1] * points[1] + points[2] * points[2] ) +
00177                                 std::sqrt( points[3] * points[3] + points[4] * points[4] + points[5] * points[5] ) );
00178         EntityHandle sfrontback[2] = { slaveverts[0], slaveverts[slaveverts.size() - 1] };
00179         error                      = mbCore->get_coords( &sfrontback[0], 2, points );MB_CHK_ERR( error );
00180         slave_radius = 0.5 * ( std::sqrt( points[0] * points[0] + points[1] * points[1] + points[2] * points[2] ) +
00181                                std::sqrt( points[3] * points[3] + points[4] * points[4] + points[5] * points[5] ) );
00182         // Let us rescale both master and slave meshes to a unit sphere
00183         error = moab::IntxUtils::ScaleToRadius( mbCore, masterfileset, 1.0 );MB_CHK_ERR( error );
00184         error = moab::IntxUtils::ScaleToRadius( mbCore, slavefileset, 1.0 );MB_CHK_ERR( error );
00185     }
00186 
00187     try
00188     {
00189         std::map< int, moab::Range > spartvals;
00190         int npoints_notfound = 0;
00191         {
00192             FileOptions fopts( ( use_spherical ? "SPHERICAL" : "" ) );
00193             moab::EntityHandle tree_root = masterfileset;
00194 
00195             moab::AdaptiveKDTree tree( mbCore );
00196             error = tree.build_tree( melems, &tree_root, &fopts );
00197 
00198             // moab::BVHTree tree(mbCore);
00199             // error = tree.build_tree(melems, &tree_root);MB_CHK_ERR(error);
00200 
00201             for( size_t ie = 0; ie < selems.size(); ie++ )
00202             {
00203                 moab::EntityHandle selem, leaf;
00204                 double point[3];
00205                 selem = selems[ie];
00206 
00207                 // Get the element centroid to be queried
00208                 error = mbCore->get_coords( &selem, 1, point );MB_CHK_ERR( error );
00209 
00210                 std::vector< moab::EntityHandle > leaf_elems;
00211 
00212                 // Search for the closest source element in the master mesh corresponding
00213                 // to the target element centroid in the slave mesh
00214                 error = tree.point_search( point, leaf, treetolerance, btolerance );MB_CHK_ERR( error );
00215 
00216                 // We only care about the dimension that the user specified.
00217                 // MOAB partitions are ordered by elements anyway.
00218                 error = mbCore->get_entities_by_dimension( leaf, dimension, leaf_elems, true );MB_CHK_ERR( error );
00219 
00220                 if( leaf != 0 && leaf_elems.size() )
00221                 {
00222 
00223                     // Now get the master element centroids so that we can compute
00224                     // the minimum distance to the target point
00225                     std::vector< double > centroids( leaf_elems.size() * 3 );
00226                     error = mbCore->get_coords( &leaf_elems[0], leaf_elems.size(), &centroids[0] );MB_CHK_ERR( error );
00227 
00228                     if( !leaf_elems.size() )
00229                         std::cout << ie << ": "
00230                                   << " No leaf elements found." << std::endl;
00231                     // else std::cout << ie << " found " << leaf_elems.size() << " leaves for
00232                     // current point " << std::endl;
00233 
00234                     double dist = 1e10;
00235                     int pinelem = -1;
00236                     for( size_t il = 0; il < leaf_elems.size(); ++il )
00237                     {
00238                         const double* centroid = &centroids[il * 3];
00239                         const double locdist   = std::pow( point[0] - centroid[0], 2 ) +
00240                                                std::pow( point[1] - centroid[1], 2 ) +
00241                                                std::pow( point[2] - centroid[2], 2 );
00242                         if( locdist < dist && locdist < 1.0E-2 )
00243                         {
00244                             dist    = locdist;
00245                             pinelem = il;
00246 
00247 #ifdef VERBOSE
00248                             int gidMelem;
00249                             error = mbCore->tag_get_data( gidtag, &leaf_elems[il], 1, &gidMelem );MB_CHK_ERR( error );
00250                             std::cout << "\t Trial leaf " << il << " set " << gidMelem
00251                                       << " and part = " << get_map_value( mpartvals, gidMelem, -1 )
00252                                       << " with distance = " << locdist << std::endl;
00253 #endif
00254                         }
00255                     }
00256 
00257                     if( pinelem < 0 )
00258                     {
00259 #ifdef VERBOSE
00260                         std::cout << ie
00261                                   << ": [Error] - Could not find a minimum distance within the "
00262                                      "leaf nodes."
00263                                   << std::endl;
00264 #endif
00265                         npoints_notfound++;
00266                         spartvals[defaultpart].insert( selems[ie] );
00267                     }
00268                     else
00269                     {
00270                         int gidMelem;
00271                         error = mbCore->tag_get_data( gidtag, &leaf_elems[pinelem], 1, &gidMelem );MB_CHK_ERR( error );
00272 
00273                         int mpartid = get_map_value( mpartvals, gidMelem, -1 );
00274                         if( mpartid < 0 )
00275                             std::cout << "[WARNING]: Part number for element " << leaf_elems[pinelem]
00276                                       << " with global ID = " << gidMelem << " not found.\n";
00277 
00278 #ifdef VERBOSE
00279                         std::cout << "Adding element " << ie << " set " << mpartid << " with distance = " << dist
00280                                   << std::endl;
00281 #endif
00282                         spartvals[mpartid].insert( selems[ie] );
00283                     }
00284                 }
00285                 else
00286                 {
00287 #ifdef VERBOSE
00288                     std::cout << "[WARNING]: Adding element " << ie << " to default (part) set " << defaultpart
00289                               << std::endl;
00290 #endif
00291                     spartvals[defaultpart].insert( selems[ie] );
00292                 }
00293             }
00294             error = tree.reset_tree();MB_CHK_ERR( error );
00295         }
00296         if( npoints_notfound )
00297             std::cout << "Could not find " << npoints_notfound
00298                       << " points in the master mesh. Added to defaultpart set = " << defaultpart << std::endl;
00299 
00300         if( use_spherical )
00301         {
00302             error = moab::IntxUtils::ScaleToRadius( mbCore, slavefileset, slave_radius );MB_CHK_ERR( error );
00303         }
00304 
00305         error = mbCore->delete_entities( &masterfileset, 1 );MB_CHK_ERR( error );
00306         // Find parallel partition sets in the slave mesh - and delete it since we are going to
00307         // overwrite the sets
00308         if( !keepsparts )
00309         {
00310             std::cout << "Deleting " << ssets.size() << " sets in the slave mesh" << std::endl;
00311             error = mbCore->remove_entities( slavefileset, ssets );MB_CHK_ERR( error );
00312             ssets.clear();
00313         }
00314 
00315         size_t ntotslave_elems = 0, ntotslave_parts = 0;
00316         for( std::map< int, moab::Range >::iterator it = spartvals.begin(); it != spartvals.end(); ++it )
00317         {
00318             int partID = it->first;
00319             moab::EntityHandle pset;
00320             error = mbCore->create_meshset( moab::MESHSET_SET, pset );MB_CHK_ERR( error );
00321             error = mbCore->add_entities( pset, it->second );MB_CHK_ERR( error );
00322             error = mbCore->add_parent_child( slavefileset, pset );MB_CHK_ERR( error );
00323 
00324 #ifdef VERBOSE
00325             std::cout << "Slave Part " << partID << " has " << it->second.size() << " elements." << std::endl;
00326 #endif
00327             ntotslave_elems += it->second.size();
00328             ntotslave_parts++;
00329 
00330             if( keepsparts )
00331             {
00332                 error = mbCore->tag_set_data( sparttag, &pset, 1, &partID );MB_CHK_ERR( error );
00333             }
00334             else
00335             {
00336                 error = mbCore->tag_set_data( parttag, &pset, 1, &partID );MB_CHK_ERR( error );
00337             }
00338         }
00339         std::cout << "Slave mesh: given " << selems.size() << " elements, and assigned " << ntotslave_elems
00340                   << " elements to " << ntotslave_parts << " parts.\n";
00341         assert( ntotslave_elems == selems.size() );
00342 
00343         // mbCore->print_database();
00344 
00345         // Write the re-partitioned slave mesh to disk
00346         if( nprocs == 1 )
00347         {
00348             error = mbCore->write_file( "slavemesh.vtk", "VTK", NULL, &slavefileset, 1 );MB_CHK_ERR( error );
00349         }
00350         error = mbCore->write_file( outfile.c_str(), NULL, write_options.c_str(), &slavefileset, 1 );MB_CHK_ERR( error );
00351         // error = mbCore->write_file(outfile.c_str(), NULL,
00352         // write_options.c_str());MB_CHK_ERR(error);
00353     }
00354     catch( std::exception& e )
00355     {
00356         std::cout << " exception caught during tree initialization " << e.what() << std::endl;
00357     }
00358     delete mbCore;
00359 
00360 #ifdef MOAB_HAVE_MPI
00361     MPI_Finalize();
00362 #endif
00363     exit( 0 );
00364 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines