21 #include "libmesh/mesh_communication.h" 
   24 #include "libmesh/libmesh_config.h" 
   25 #include "libmesh/libmesh_common.h" 
   26 #include "libmesh/libmesh_logging.h" 
   27 #include "libmesh/mesh_base.h" 
   28 #include "libmesh/mesh_tools.h" 
   29 #include "libmesh/parallel_hilbert.h" 
   30 #include "libmesh/parallel_sort.h" 
   31 #include "libmesh/elem.h" 
   32 #include "libmesh/elem_range.h" 
   33 #include "libmesh/node_range.h" 
   40 #ifdef LIBMESH_HAVE_LIBHILBERT 
   45 #ifdef LIBMESH_HAVE_LIBHILBERT 
   52 void get_hilbert_coords (
const Point & p,
 
   54                          CFixBitVec icoords[3])
 
   56   static const Hilbert::inttype max_inttype = static_cast<Hilbert::inttype>(-1);
 
   59     x = static_cast<long double>
 
   60         ((bbox.first(0) == bbox.second(0)) ? 0. :
 
   61          (p(0)-bbox.first(0))/(bbox.second(0)-bbox.first(0))),
 
   64     y = static_cast<long double>
 
   65         ((bbox.first(1) == bbox.second(1)) ? 0. :
 
   66          (p(1)-bbox.first(1))/(bbox.second(1)-bbox.first(1))),
 
   72     z = static_cast<long double>
 
   73         ((bbox.first(2) == bbox.second(2)) ? 0. :
 
   74          (p(2)-bbox.first(2))/(bbox.second(2)-bbox.first(2)));
 
   80   icoords[0] = static_cast<Hilbert::inttype>(x*max_inttype);
 
   81   icoords[1] = static_cast<Hilbert::inttype>(y*max_inttype);
 
   82   icoords[2] = static_cast<Hilbert::inttype>(z*max_inttype);
 
   88 get_dofobject_key (
const Elem * e,
 
   91   static const unsigned int sizeof_inttype = 
sizeof(Hilbert::inttype);
 
   93   Hilbert::HilbertIndices index;
 
   94   CFixBitVec icoords[3];
 
   95   Hilbert::BitVecType bv;
 
   96   get_hilbert_coords (e->
centroid(), bbox, icoords);
 
   97   Hilbert::coordsToIndex (icoords, 8*sizeof_inttype, 3, bv);
 
  100 #ifdef LIBMESH_ENABLE_UNIQUE_ID 
  101   return std::make_pair(index, e->
unique_id());
 
  111 get_dofobject_key (
const Node * n,
 
  114   static const unsigned int sizeof_inttype = 
sizeof(Hilbert::inttype);
 
  116   Hilbert::HilbertIndices index;
 
  117   CFixBitVec icoords[3];
 
  118   Hilbert::BitVecType bv;
 
  119   get_hilbert_coords (*n, bbox, icoords);
 
  120   Hilbert::coordsToIndex (icoords, 8*sizeof_inttype, 3, bv);
 
  123 #ifdef LIBMESH_ENABLE_UNIQUE_ID 
  124   return std::make_pair(index, n->
unique_id());
 
  131 class ComputeHilbertKeys
 
  135                       std::vector<Parallel::DofObjectKey> & keys) :
 
  144     for (
const auto & node : range)
 
  147         libmesh_assert_less (pos, _keys.size());
 
  148         _keys[pos++] = get_dofobject_key (node, _bbox);
 
  156     for (
const auto & elem : range)
 
  159         libmesh_assert_less (pos, _keys.size());
 
  160         _keys[pos++] = get_dofobject_key (elem, _bbox);
 
  166   std::vector<Parallel::DofObjectKey> & _keys;
 
  177 #if defined(LIBMESH_HAVE_LIBHILBERT) && defined(LIBMESH_HAVE_MPI) 
  180   LOG_SCOPE (
"assign_global_indices()", 
"MeshCommunication");
 
  192   const Parallel::Communicator & communicator (
mesh.
comm());
 
  202   std::vector<Parallel::DofObjectKey>
 
  203     node_keys, elem_keys;
 
  210       node_keys.resize (nr.size());
 
  238       elem_keys.resize (er.size());
 
  280   const std::vector<Parallel::DofObjectKey> & my_node_bin =
 
  287   const std::vector<Parallel::DofObjectKey> & my_elem_bin =
 
  294   std::vector<Parallel::DofObjectKey>
 
  295     node_upper_bounds(communicator.size()),
 
  296     elem_upper_bounds(communicator.size());
 
  299     std::vector<Parallel::DofObjectKey> recvbuf(2*communicator.size());
 
  300     std::vector<unsigned short int> 
 
  301       empty_nodes (communicator.size()),
 
  302       empty_elem  (communicator.size());
 
  303     std::vector<Parallel::DofObjectKey> my_max(2);
 
  305     communicator.allgather (static_cast<unsigned short int>(my_node_bin.empty()), empty_nodes);
 
  306     communicator.allgather (static_cast<unsigned short int>(my_elem_bin.empty()),  empty_elem);
 
  308     if (!my_node_bin.empty()) my_max[0] = my_node_bin.back();
 
  309     if (!my_elem_bin.empty()) my_max[1] = my_elem_bin.back();
 
  311     communicator.allgather (my_max,  
true);
 
  318         node_upper_bounds[p] = my_max[2*p+0];
 
  319         elem_upper_bounds[p] = my_max[2*p+1];
 
  323             if (empty_nodes[p]) node_upper_bounds[p] = node_upper_bounds[p-1];
 
  324             if (empty_elem[p])  elem_upper_bounds[p] = elem_upper_bounds[p-1];
 
  339       std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
 
  343       std::map<dof_id_type, std::vector<dof_id_type>>
 
  351             get_dofobject_key (node, bbox);
 
  353             cast_int<processor_id_type>
 
  355                             std::lower_bound(node_upper_bounds.begin(),
 
  356                                              node_upper_bounds.end(),
 
  359           libmesh_assert_less (pid, communicator.size());
 
  361           requested_ids[pid].push_back(hi);
 
  365       std::vector<dof_id_type> node_bin_sizes(communicator.size());
 
  366       communicator.allgather (static_cast<dof_id_type>(my_node_bin.size()), node_bin_sizes);
 
  371         my_offset += node_bin_sizes[pid];
 
  373       auto gather_functor =
 
  383          const std::vector<Parallel::DofObjectKey> & keys,
 
  384          std::vector<dof_id_type> & global_ids)
 
  387           const std::size_t keys_size = keys.size();
 
  388           global_ids.reserve(keys_size);
 
  389           for (std::size_t 
idx=0; 
idx != keys_size; 
idx++)
 
  392               libmesh_assert_less_equal (hi, node_upper_bounds[communicator.rank()]);
 
  395               std::vector<Parallel::DofObjectKey>::const_iterator pos =
 
  396                 std::lower_bound (my_node_bin.begin(), my_node_bin.end(), hi);
 
  398               libmesh_assert_equal_to (*pos, hi);
 
  402               global_ids.push_back(cast_int<dof_id_type>(
std::distance(my_node_bin.begin(), pos) + my_offset));
 
  406       auto action_functor =
 
  409          const std::vector<Parallel::DofObjectKey> &,
 
  410          const std::vector<dof_id_type> & global_ids)
 
  412           filled_request[pid] = global_ids;
 
  417       Parallel::pull_parallel_vector_data
 
  418         (communicator, requested_ids, gather_functor, action_functor, ex);
 
  423         std::map<dof_id_type, std::vector<dof_id_type>::const_iterator>
 
  425         for (
auto & p : filled_request)
 
  426           next_obj_on_proc[p.first] = p.second.begin();
 
  432               get_dofobject_key (node, bbox);
 
  434               cast_int<processor_id_type>
 
  436                               std::lower_bound(node_upper_bounds.begin(),
 
  437                                                node_upper_bounds.end(),
 
  440             libmesh_assert_less (pid, communicator.size());
 
  443             const dof_id_type global_index = *next_obj_on_proc[pid];
 
  444             libmesh_assert_less (global_index, 
mesh.
n_nodes());
 
  445             node->set_id() = global_index;
 
  447             ++next_obj_on_proc[pid];
 
  456       std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
 
  460       std::map<dof_id_type, std::vector<dof_id_type>>
 
  467             get_dofobject_key (elem, bbox);
 
  469             cast_int<processor_id_type>
 
  471                             std::lower_bound(elem_upper_bounds.begin(),
 
  472                                              elem_upper_bounds.end(),
 
  475           libmesh_assert_less (pid, communicator.size());
 
  477           requested_ids[pid].push_back(hi);
 
  481       std::vector<dof_id_type> elem_bin_sizes(communicator.size());
 
  482       communicator.allgather (static_cast<dof_id_type>(my_elem_bin.size()), elem_bin_sizes);
 
  487         my_offset += elem_bin_sizes[pid];
 
  489       auto gather_functor =
 
  499          const std::vector<Parallel::DofObjectKey> & keys,
 
  500          std::vector<dof_id_type> & global_ids)
 
  503           const std::size_t keys_size = keys.size();
 
  504           global_ids.reserve(keys_size);
 
  505           for (std::size_t 
idx=0; 
idx != keys_size; 
idx++)
 
  508               libmesh_assert_less_equal (hi, elem_upper_bounds[communicator.rank()]);
 
  511               std::vector<Parallel::DofObjectKey>::const_iterator pos =
 
  512                 std::lower_bound (my_elem_bin.begin(), my_elem_bin.end(), hi);
 
  514               libmesh_assert_equal_to (*pos, hi);
 
  518               global_ids.push_back (cast_int<dof_id_type>(
std::distance(my_elem_bin.begin(), pos) + my_offset));
 
  522       auto action_functor =
 
  525          const std::vector<Parallel::DofObjectKey> &,
 
  526          const std::vector<dof_id_type> & global_ids)
 
  528           filled_request[pid] = global_ids;
 
  533       Parallel::pull_parallel_vector_data
 
  534         (communicator, requested_ids, gather_functor, action_functor, ex);
 
  539         std::vector<std::vector<dof_id_type>::const_iterator>
 
  540           next_obj_on_proc; next_obj_on_proc.reserve(communicator.size());
 
  542           next_obj_on_proc.push_back(filled_request[pid].begin());
 
  548               get_dofobject_key (elem, bbox);
 
  550               cast_int<processor_id_type>
 
  552                               std::lower_bound(elem_upper_bounds.begin(),
 
  553                                                elem_upper_bounds.end(),
 
  556             libmesh_assert_less (pid, communicator.size());
 
  559             const dof_id_type global_index = *next_obj_on_proc[pid];
 
  560             libmesh_assert_less (global_index, 
mesh.
n_elem());
 
  561             elem->set_id() = global_index;
 
  563             ++next_obj_on_proc[pid];
 
  569 #else // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 
  573 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 
  576 #if defined(LIBMESH_HAVE_LIBHILBERT) && defined(LIBMESH_HAVE_MPI) 
  579   LOG_SCOPE (
"check_for_duplicate_global_indices()", 
"MeshCommunication");
 
  588   std::vector<Parallel::DofObjectKey>
 
  589     node_keys, elem_keys;
 
  596       node_keys.resize (nr.size());
 
  602       for (std::size_t i = 0; i != node_keys.size(); ++i, ++nodei)
 
  605           for (std::size_t j = 0; j != i; ++j, ++nodej)
 
  607               if (node_keys[i] == node_keys[j])
 
  609                   CFixBitVec icoords[3], jcoords[3];
 
  610                   get_hilbert_coords(**nodej, bbox, jcoords);
 
  612                     "node " << (*nodej)->id() << 
", " <<
 
  613                     *(
Point *)(*nodej) << 
" has HilbertIndices " <<
 
  614                     node_keys[j] << std::endl;
 
  615                   get_hilbert_coords(**nodei, bbox, icoords);
 
  617                     "node " << (*nodei)->id() << 
", " <<
 
  618                     *(
Point *)(*nodei) << 
" has HilbertIndices " <<
 
  619                     node_keys[i] << std::endl;
 
  620                   libmesh_error_msg(
"Error: nodes with duplicate Hilbert keys!");
 
  630       elem_keys.resize (er.size());
 
  637       for (std::size_t i = 0; i != elem_keys.size(); ++i, ++elemi)
 
  640           for (std::size_t j = 0; j != i; ++j, ++elemj)
 
  642               if ((elem_keys[i] == elem_keys[j]) &&
 
  643                   ((*elemi)->level() == (*elemj)->level()))
 
  646                     "level " << (*elemj)->level() << 
" elem\n" <<
 
  647                     (**elemj) << 
" centroid " <<
 
  648                     (*elemj)->centroid() << 
" has HilbertIndices " <<
 
  649                     elem_keys[j] << 
" or " <<
 
  650                     get_dofobject_key((*elemj), bbox) <<
 
  653                     "level " << (*elemi)->level() << 
" elem\n" <<
 
  654                     (**elemi) << 
" centroid " <<
 
  655                     (*elemi)->centroid() << 
" has HilbertIndices " <<
 
  656                     elem_keys[i] << 
" or " <<
 
  657                     get_dofobject_key((*elemi), bbox) <<
 
  659                   libmesh_error_msg(
"Error: level " << (*elemi)->level() << 
" elements with duplicate Hilbert keys!");
 
  666 #else // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 
  670 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 
  672 #if defined(LIBMESH_HAVE_LIBHILBERT) && defined(LIBMESH_HAVE_MPI) 
  673 template <
typename ForwardIterator>
 
  675                                             const ForwardIterator & begin,
 
  676                                             const ForwardIterator & 
end,
 
  677                                             std::unordered_map<dof_id_type, dof_id_type> & index_map)
 const 
  679   LOG_SCOPE (
"find_local_indices()", 
"MeshCommunication");
 
  691   std::map<Parallel::DofObjectKey, dof_id_type> hilbert_keys;
 
  693     LOG_SCOPE(
"local_hilbert_indices", 
"MeshCommunication");
 
  694     for (ForwardIterator it=begin; it!=
end; ++it)
 
  697         hilbert_keys.emplace(hi, (*it)->id());
 
  703     for (
auto key_val : hilbert_keys)
 
  704       index_map[key_val.second] = cnt++;
 
  709 template <
typename ForwardIterator>
 
  712                                              const ForwardIterator & begin,
 
  713                                              const ForwardIterator & 
end,
 
  714                                              std::vector<dof_id_type> & index_map)
 const 
  716   LOG_SCOPE (
"find_global_indices()", 
"MeshCommunication");
 
  729   index_map.reserve(n_objects);
 
  736   std::vector<Parallel::DofObjectKey>
 
  739   sorted_hilbert_keys.reserve(n_objects);
 
  740   hilbert_keys.reserve(n_objects);
 
  742     LOG_SCOPE(
"compute_hilbert_indices()", 
"MeshCommunication");
 
  743     for (ForwardIterator it=begin; it!=
end; ++it)
 
  746         hilbert_keys.push_back(hi);
 
  748         if ((*it)->processor_id() == communicator.rank())
 
  749           sorted_hilbert_keys.push_back(hi);
 
  752         if ((communicator.rank() == 0) &&
 
  754           sorted_hilbert_keys.push_back(hi);
 
  760   START_LOG (
"parallel_sort()", 
"MeshCommunication");
 
  762                                                  sorted_hilbert_keys);
 
  764   STOP_LOG (
"parallel_sort()", 
"MeshCommunication");
 
  765   const std::vector<Parallel::DofObjectKey> & my_bin = sorter.
bin();
 
  768   std::vector<unsigned int> bin_sizes(communicator.size());
 
  769   communicator.allgather (static_cast<unsigned int>(my_bin.size()), bin_sizes);
 
  772   unsigned int my_offset = 0;
 
  774     my_offset += bin_sizes[pid];
 
  778   std::vector<Parallel::DofObjectKey>
 
  782     upper_bounds[0] = my_bin.back();
 
  784   communicator.allgather (upper_bounds,  
true);
 
  790     if (!bin_sizes[p]) upper_bounds[p] = upper_bounds[p-1];
 
  801     std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
 
  804     std::map<processor_id_type, std::vector<dof_id_type>>
 
  808     std::vector<Parallel::DofObjectKey>::const_iterator hi =
 
  809       hilbert_keys.begin();
 
  811     for (ForwardIterator it = begin; it != 
end; ++it)
 
  815         std::vector<Parallel::DofObjectKey>::iterator lb =
 
  816           std::lower_bound(upper_bounds.begin(), upper_bounds.end(),
 
  820           cast_int<processor_id_type>
 
  823         libmesh_assert_less (pid, communicator.size());
 
  825         requested_ids[pid].push_back(*hi);
 
  830         index_map.push_back(pid);
 
  833   auto gather_functor =
 
  844      std::vector<dof_id_type> & global_ids)
 
  850       const std::size_t keys_size = keys.size();
 
  852       global_ids.reserve(keys_size);
 
  853       for (std::size_t 
idx=0; 
idx != keys_size; 
idx++)
 
  856           libmesh_assert_less_equal (hilbert_indices, upper_bounds[communicator.rank()]);
 
  859           std::vector<Parallel::DofObjectKey>::const_iterator pos =
 
  860             std::lower_bound (my_bin.begin(), my_bin.end(), hilbert_indices);
 
  867           if (*pos != hilbert_indices)
 
  872               Hilbert::BitVecType input;
 
  873 #ifdef LIBMESH_ENABLE_UNIQUE_ID 
  874               input = hilbert_indices.first;
 
  876               input = hilbert_indices;
 
  880               std::vector<CBigBitVec> output(3);
 
  883               Hilbert::indexToCoords(output.data(), 8*
sizeof(Hilbert::inttype), 3, input);
 
  889               const Real max_int_as_real =
 
  890                 static_cast<Real>(std::numeric_limits<Hilbert::inttype>::max());
 
  895               Point p_hat(static_cast<Real>(output[0].racks()[0]) / max_int_as_real,
 
  896                           static_cast<Real>(output[1].racks()[0]) / max_int_as_real,
 
  897                           static_cast<Real>(output[2].racks()[0]) / max_int_as_real);
 
  901                 xmin = bbox.first(0),
 
  902                 xmax = bbox.second(0),
 
  903                 ymin = bbox.first(1),
 
  904                 ymax = bbox.second(1),
 
  905                 zmin = bbox.first(2),
 
  906                 zmax = bbox.second(2);
 
  909               Point p(xmin + (xmax-xmin)*p_hat(0),
 
  910                       ymin + (ymax-ymin)*p_hat(1),
 
  911                       zmin + (zmax-zmin)*p_hat(2));
 
  913               libmesh_error_msg(
"Could not find hilbert indices: " 
  915                                 << 
" corresponding to point " << p);
 
  921           global_ids.push_back (cast_int<dof_id_type>(
std::distance(my_bin.begin(), pos) + my_offset));
 
  925   auto action_functor =
 
  928      const std::vector<Parallel::DofObjectKey> &,
 
  929      const std::vector<dof_id_type> & global_ids)
 
  931       filled_request[pid] = global_ids;
 
  935   Parallel::pull_parallel_vector_data
 
  936     (communicator, requested_ids, gather_functor, action_functor, ex);
 
  941       std::vector<std::vector<dof_id_type>::const_iterator>
 
  942         next_obj_on_proc; next_obj_on_proc.reserve(communicator.size());
 
  944         next_obj_on_proc.push_back(filled_request[pid].begin());
 
  947       for (ForwardIterator it = begin; it != 
end; ++it, cnt++)
 
  952           libmesh_assert_less (pid, communicator.size());
 
  955           const dof_id_type global_index = *next_obj_on_proc[pid];
 
  956           index_map[cnt] = global_index;
 
  958           ++next_obj_on_proc[pid];
 
  963   libmesh_assert_equal_to(index_map.size(), n_objects);
 
  965 #else // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 
  966 template <
typename ForwardIterator>
 
  968                                             const ForwardIterator & begin,
 
  969                                             const ForwardIterator & 
end,
 
  970                                             std::unordered_map<dof_id_type, dof_id_type> & index_map)
 const 
  974   for (ForwardIterator it=begin; it!=
end; ++it)
 
  975     index_map[(*it)->id()] = (index++);
 
  978 template <
typename ForwardIterator>
 
  981                                              const ForwardIterator & begin,
 
  982                                              const ForwardIterator & 
end,
 
  983                                              std::vector<dof_id_type> & index_map)
 const 
  987   unsigned int index = 0;
 
  988   for (ForwardIterator it=begin; it!=
end; ++it)
 
  989     index_map.push_back(index++);
 
  991 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 
  996 template void MeshCommunication::find_global_indices<MeshBase::const_node_iterator> (
const Parallel::Communicator &,
 
 1000                                                                                      std::vector<dof_id_type> &) 
const;
 
 1002 template void MeshCommunication::find_global_indices<MeshBase::const_element_iterator> (
const Parallel::Communicator &,
 
 1006                                                                                         std::vector<dof_id_type> &) 
const;
 
 1007 template void MeshCommunication::find_global_indices<MeshBase::node_iterator> (
const Parallel::Communicator &,
 
 1011                                                                                std::vector<dof_id_type> &) 
const;
 
 1013 template void MeshCommunication::find_global_indices<MeshBase::element_iterator> (
const Parallel::Communicator &,
 
 1017                                                                                   std::vector<dof_id_type> &) 
const;
 
 1018 template void MeshCommunication::find_local_indices<MeshBase::const_element_iterator> (
const libMesh::BoundingBox &,
 
 1021                                                                                        std::unordered_map<dof_id_type, dof_id_type> &) 
const;