21 #include "libmesh/parmetis_partitioner.h" 
   24 #include "libmesh/libmesh_config.h" 
   25 #include "libmesh/mesh_base.h" 
   26 #include "libmesh/mesh_serializer.h" 
   27 #include "libmesh/mesh_tools.h" 
   28 #include "libmesh/mesh_communication.h" 
   29 #include "libmesh/metis_partitioner.h" 
   30 #include "libmesh/parallel_only.h" 
   31 #include "libmesh/libmesh_logging.h" 
   32 #include "libmesh/elem.h" 
   33 #include "libmesh/parmetis_helper.h" 
   40 #ifdef LIBMESH_HAVE_PARMETIS 
   48 #     include "libmesh/ignore_warnings.h" 
   49 #     include "parmetis.h" 
   50 #     include "libmesh/restore_warnings.h" 
   58 #include <unordered_map> 
   66 #ifdef LIBMESH_HAVE_PARMETIS 
   73 #ifdef LIBMESH_HAVE_PARMETIS 
   74   :  _pmetis(libmesh_make_unique<ParmetisHelper>())
 
   82 #ifdef LIBMESH_HAVE_PARMETIS
 
   83   , _pmetis(libmesh_make_unique<ParmetisHelper>(*(other._pmetis)))
 
   90 ParmetisPartitioner::~ParmetisPartitioner() = 
default;
 
   94 void ParmetisPartitioner::_do_partition (MeshBase & 
mesh,
 
   95                                          const unsigned int n_sbdmns)
 
   97   this->_do_repartition (
mesh, n_sbdmns);
 
  103                                            const unsigned int n_sbdmns)
 
  106   libmesh_parallel_only(
mesh.comm());
 
  114       this->single_partition(
mesh);
 
  118   libmesh_assert_greater (n_sbdmns, 0);
 
  121 #ifndef LIBMESH_HAVE_PARMETIS 
  124   libMesh::out << 
"ERROR: The library has been built without" << std::endl
 
  125                << 
"Parmetis support.  Using a Metis"          << std::endl
 
  126                << 
"partitioner instead!"                      << std::endl;);
 
  134                       mesh.active_elements_end(), n_sbdmns);
 
  140   if (
mesh.n_processors() == 1)
 
  150                           mesh.active_elements_end(), n_sbdmns);
 
  154   LOG_SCOPE(
"repartition()", 
"ParmetisPartitioner");
 
  160   this->build_graph (
mesh);
 
  167     bool ready_for_parmetis = 
true;
 
  168     for (
const auto & nelem : _n_active_elem_on_proc)
 
  170         ready_for_parmetis = 
false;
 
  172     std::size_t my_adjacency = _pmetis->adjncy.size();
 
  173     mesh.comm().min(my_adjacency);
 
  175       ready_for_parmetis = 
false;
 
  180     if (!ready_for_parmetis)
 
  192   std::vector<Parmetis::idx_t> vsize(_pmetis->vwgt.size(), 1);
 
  193   Parmetis::real_t itr = 1000000.0;
 
  194   MPI_Comm mpi_comm = 
mesh.comm().get();
 
  199   Parmetis::ParMETIS_V3_AdaptiveRepart(_pmetis->vtxdist.empty() ? nullptr : _pmetis->vtxdist.data(),
 
  200                                        _pmetis->xadj.empty()    ? nullptr : _pmetis->xadj.data(),
 
  201                                        _pmetis->adjncy.empty()  ? nullptr : _pmetis->adjncy.data(),
 
  202                                        _pmetis->vwgt.empty()    ? nullptr : _pmetis->vwgt.data(),
 
  203                                        vsize.empty()            ? nullptr : vsize.data(),
 
  209                                        _pmetis->tpwgts.empty()  ? nullptr : _pmetis->tpwgts.data(),
 
  210                                        _pmetis->ubvec.empty()   ? nullptr : _pmetis->ubvec.data(),
 
  212                                        _pmetis->options.data(),
 
  214                                        _pmetis->part.empty()    ? nullptr : reinterpret_cast<Parmetis::idx_t *>(_pmetis->part.data()),
 
  218   this->assign_partitioning (
mesh, _pmetis->part);
 
  220 #endif // #ifndef LIBMESH_HAVE_PARMETIS ... else ... 
  227 #ifdef LIBMESH_HAVE_PARMETIS 
  230                                       const unsigned int n_sbdmns)
 
  232   LOG_SCOPE(
"initialize()", 
"ParmetisPartitioner");
 
  236   _pmetis->wgtflag = 2;                                      
 
  238   _pmetis->numflag = 0;                                      
 
  239   _pmetis->nparts  = static_cast<Parmetis::idx_t>(n_sbdmns); 
 
  240   _pmetis->edgecut = 0;                                      
 
  244   _pmetis->vtxdist.assign (
mesh.n_processors()+1, 0);
 
  245   _pmetis->tpwgts.assign  (_pmetis->nparts, 1./_pmetis->nparts);
 
  246   _pmetis->ubvec.assign   (_pmetis->ncon, 1.05);
 
  247   _pmetis->part.assign    (n_active_local_elem, 0);
 
  248   _pmetis->options.resize (5);
 
  249   _pmetis->vwgt.resize    (n_active_local_elem);
 
  252   _pmetis->options[0] = 1;  
 
  253   _pmetis->options[1] = 0;  
 
  254   _pmetis->options[2] = 15; 
 
  255   _pmetis->options[3] = 2;  
 
  269   _find_global_index_by_pid_map(
mesh);
 
  281   libmesh_assert_equal_to (_pmetis->vtxdist.size(),
 
  282                            cast_int<std::size_t>(
mesh.n_processors()+1));
 
  283   libmesh_assert_equal_to (_pmetis->vtxdist[0], 0);
 
  287       _pmetis->vtxdist[pid+1] = _pmetis->vtxdist[pid] + _n_active_elem_on_proc[pid];
 
  288       n_active_elem += _n_active_elem_on_proc[pid];
 
  290   libmesh_assert_equal_to (_pmetis->vtxdist.back(), static_cast<Parmetis::idx_t>(n_active_elem));
 
  295   std::unordered_map<dof_id_type, dof_id_type> global_index_map;
 
  298     std::vector<dof_id_type> global_index;
 
  313           const Elem * elem = *it;
 
  317           libmesh_assert_less (cnt, global_index.size());
 
  318           libmesh_assert_less (global_index[cnt], n_active_elem);
 
  320           global_index_map.insert(std::make_pair(elem->
id(), global_index[cnt++]));
 
  324     libmesh_assert_less_equal (global_index_map.size(), n_active_elem);
 
  325     libmesh_assert_less_equal (_global_index_by_pid_map.size(), n_active_elem);
 
  331     if (global_index_map.size() != _global_index_by_pid_map.size())
 
  332       libmesh_error_msg(
"ERROR:  ParmetisPartitioner cannot handle unpartitioned objects!");
 
  339     std::vector<dof_id_type> subdomain_bounds(
mesh.n_processors());
 
  341     const dof_id_type first_local_elem = _pmetis->vtxdist[
mesh.processor_id()];
 
  348         if (pid < static_cast<unsigned int>(_pmetis->nparts))
 
  350             tgt_subdomain_size = n_active_elem/std::min
 
  351               (cast_int<Parmetis::idx_t>(
mesh.n_processors()), _pmetis->nparts);
 
  353             if (pid < n_active_elem%_pmetis->nparts)
 
  354               tgt_subdomain_size++;
 
  357           subdomain_bounds[0] = tgt_subdomain_size;
 
  359           subdomain_bounds[pid] = subdomain_bounds[pid-1] + tgt_subdomain_size;
 
  362     libmesh_assert_equal_to (subdomain_bounds.back(), n_active_elem);
 
  364     for (
const auto & elem : 
mesh.active_local_element_ptr_range())
 
  368           _global_index_by_pid_map[elem->id()];
 
  369         libmesh_assert_less (global_index_by_pid, n_active_elem);
 
  372           global_index_by_pid - first_local_elem;
 
  374         libmesh_assert_less (local_index, n_active_local_elem);
 
  375         libmesh_assert_less (local_index, _pmetis->vwgt.size());
 
  378         _pmetis->vwgt[local_index] = elem->n_nodes();
 
  383           global_index_map[elem->id()];
 
  385         libmesh_assert_less (global_index, subdomain_bounds.back());
 
  387         const unsigned int subdomain_id =
 
  388           cast_int<unsigned int>
 
  390                          std::lower_bound(subdomain_bounds.begin(),
 
  391                                           subdomain_bounds.end(),
 
  393         libmesh_assert_less (subdomain_id, _pmetis->nparts);
 
  394         libmesh_assert_less (local_index, _pmetis->part.size());
 
  396         _pmetis->part[local_index] = subdomain_id;
 
  405   LOG_SCOPE(
"build_graph()", 
"ParmetisPartitioner");
 
  412   Partitioner::build_graph(
mesh);
 
  416   for (
auto & row: _dual_graph)
 
  417    graph_size += cast_int<dof_id_type>(row.size());
 
  420   _pmetis->xadj.clear();
 
  421   _pmetis->xadj.reserve (n_active_local_elem + 1);
 
  422   _pmetis->adjncy.clear();
 
  423   _pmetis->adjncy.reserve (graph_size);
 
  425   for (
auto & graph_row : _dual_graph)
 
  427       _pmetis->xadj.push_back(cast_int<int>(_pmetis->adjncy.size()));
 
  428       _pmetis->adjncy.insert(_pmetis->adjncy.end(),
 
  434   _pmetis->xadj.push_back(cast_int<int>(_pmetis->adjncy.size()));
 
  436   libmesh_assert_equal_to (_pmetis->xadj.size(), n_active_local_elem+1);
 
  437   libmesh_assert_equal_to (_pmetis->adjncy.size(), graph_size);
 
  440 #endif // #ifdef LIBMESH_HAVE_PARMETIS