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