34 #ifdef LIBMESH_HAVE_LIBHILBERT 38 #ifdef LIBMESH_HAVE_LIBHILBERT 45 void get_hilbert_coords (
const Point & p,
47 CFixBitVec icoords[3])
49 static const Hilbert::inttype max_inttype =
static_cast<Hilbert::inttype
>(-1);
52 x = ((bbox.first(0) == bbox.second(0)) ? 0. :
53 (p(0)-bbox.first(0))/(bbox.second(0)-bbox.first(0))),
56 y = ((bbox.first(1) == bbox.second(1)) ? 0. :
57 (p(1)-bbox.first(1))/(bbox.second(1)-bbox.first(1))),
63 z = ((bbox.first(2) == bbox.second(2)) ? 0. :
64 (p(2)-bbox.first(2))/(bbox.second(2)-bbox.first(2)));
70 icoords[0] =
static_cast<Hilbert::inttype
>(x*max_inttype);
71 icoords[1] =
static_cast<Hilbert::inttype
>(y*max_inttype);
72 icoords[2] =
static_cast<Hilbert::inttype
>(z*max_inttype);
78 get_dofobject_key (
const Elem * e,
81 static const unsigned int sizeof_inttype =
sizeof(Hilbert::inttype);
83 Hilbert::HilbertIndices index;
84 CFixBitVec icoords[3];
85 Hilbert::BitVecType bv;
86 get_hilbert_coords (e->
centroid(), bbox, icoords);
87 Hilbert::coordsToIndex (icoords, 8*sizeof_inttype, 3, bv);
90 #ifdef LIBMESH_ENABLE_UNIQUE_ID 91 return std::make_pair(index, e->
unique_id());
101 get_dofobject_key (
const Node * n,
104 static const unsigned int sizeof_inttype =
sizeof(Hilbert::inttype);
106 Hilbert::HilbertIndices index;
107 CFixBitVec icoords[3];
108 Hilbert::BitVecType bv;
109 get_hilbert_coords (*n, bbox, icoords);
110 Hilbert::coordsToIndex (icoords, 8*sizeof_inttype, 3, bv);
113 #ifdef LIBMESH_ENABLE_UNIQUE_ID 114 return std::make_pair(index, n->
unique_id());
121 class ComputeHilbertKeys
125 std::vector<Parallel::DofObjectKey> & keys) :
134 for (
const auto & node : range)
136 libmesh_assert(node);
137 libmesh_assert_less (pos, _keys.size());
138 _keys[pos++] = get_dofobject_key (node, _bbox);
146 for (
const auto & elem : range)
148 libmesh_assert(elem);
149 libmesh_assert_less (pos, _keys.size());
150 _keys[pos++] = get_dofobject_key (elem, _bbox);
156 std::vector<Parallel::DofObjectKey> & _keys;
167 #if defined(LIBMESH_HAVE_LIBHILBERT) && defined(LIBMESH_HAVE_MPI) 170 LOG_SCOPE (
"assign_global_indices()",
"MeshCommunication");
192 std::vector<Parallel::DofObjectKey>
193 node_keys, elem_keys;
200 node_keys.resize (nr.size());
234 elem_keys.resize (er.size());
278 const std::vector<Parallel::DofObjectKey> & my_node_bin =
285 const std::vector<Parallel::DofObjectKey> & my_elem_bin =
292 std::vector<Parallel::DofObjectKey>
297 std::vector<Parallel::DofObjectKey> recvbuf(2*
communicator.size());
298 std::vector<unsigned short int>
301 std::vector<Parallel::DofObjectKey> my_max(2);
303 communicator.allgather (static_cast<unsigned short int>(my_node_bin.empty()), empty_nodes);
304 communicator.allgather (static_cast<unsigned short int>(my_elem_bin.empty()), empty_elem);
306 if (!my_node_bin.empty()) my_max[0] = my_node_bin.back();
307 if (!my_elem_bin.empty()) my_max[1] = my_elem_bin.back();
316 node_upper_bounds[p] = my_max[2*p+0];
317 elem_upper_bounds[p] = my_max[2*p+1];
321 if (empty_nodes[p]) node_upper_bounds[p] = node_upper_bounds[p-1];
322 if (empty_elem[p]) elem_upper_bounds[p] = elem_upper_bounds[p-1];
337 std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
341 std::map<dof_id_type, std::vector<dof_id_type>>
347 libmesh_assert(node);
349 get_dofobject_key (node, bbox);
351 cast_int<processor_id_type>
352 (std::distance (node_upper_bounds.begin(),
353 std::lower_bound(node_upper_bounds.begin(),
354 node_upper_bounds.end(),
359 requested_ids[pid].push_back(hi);
363 std::vector<dof_id_type> node_bin_sizes(
communicator.size());
364 communicator.allgather (static_cast<dof_id_type>(my_node_bin.size()), node_bin_sizes);
369 my_offset += node_bin_sizes[pid];
371 auto gather_functor =
381 const std::vector<Parallel::DofObjectKey> & keys,
382 std::vector<dof_id_type> & global_ids)
385 const std::size_t keys_size = keys.size();
386 global_ids.reserve(keys_size);
387 for (std::size_t
idx=0;
idx != keys_size;
idx++)
390 libmesh_assert_less_equal (hi, node_upper_bounds[
communicator.rank()]);
393 std::vector<Parallel::DofObjectKey>::const_iterator pos =
394 std::lower_bound (my_node_bin.begin(), my_node_bin.end(), hi);
395 libmesh_assert (pos != my_node_bin.end());
396 libmesh_assert_equal_to (*pos, hi);
400 global_ids.push_back(cast_int<dof_id_type>(std::distance(my_node_bin.begin(), pos) + my_offset));
404 auto action_functor =
407 const std::vector<Parallel::DofObjectKey> &,
408 const std::vector<dof_id_type> & global_ids)
410 filled_request[pid] = global_ids;
416 (
communicator, requested_ids, gather_functor, action_functor, ex);
421 std::map<dof_id_type, std::vector<dof_id_type>::const_iterator>
423 for (
auto & p : filled_request)
424 next_obj_on_proc[p.first] = p.second.begin();
428 libmesh_assert(node);
430 get_dofobject_key (node, bbox);
432 cast_int<processor_id_type>
433 (std::distance (node_upper_bounds.begin(),
434 std::lower_bound(node_upper_bounds.begin(),
435 node_upper_bounds.end(),
439 libmesh_assert (next_obj_on_proc[pid] != filled_request[pid].
end());
441 const dof_id_type global_index = *next_obj_on_proc[pid];
442 libmesh_assert_less (global_index,
mesh.
n_nodes());
443 node->set_id() = global_index;
445 ++next_obj_on_proc[pid];
454 std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
458 std::map<dof_id_type, std::vector<dof_id_type>>
463 libmesh_assert(elem);
465 get_dofobject_key (elem, bbox);
467 cast_int<processor_id_type>
468 (std::distance (elem_upper_bounds.begin(),
469 std::lower_bound(elem_upper_bounds.begin(),
470 elem_upper_bounds.end(),
475 requested_ids[pid].push_back(hi);
479 std::vector<dof_id_type> elem_bin_sizes(
communicator.size());
480 communicator.allgather (static_cast<dof_id_type>(my_elem_bin.size()), elem_bin_sizes);
485 my_offset += elem_bin_sizes[pid];
487 auto gather_functor =
497 const std::vector<Parallel::DofObjectKey> & keys,
498 std::vector<dof_id_type> & global_ids)
501 const std::size_t keys_size = keys.size();
502 global_ids.reserve(keys_size);
503 for (std::size_t
idx=0;
idx != keys_size;
idx++)
506 libmesh_assert_less_equal (hi, elem_upper_bounds[
communicator.rank()]);
509 std::vector<Parallel::DofObjectKey>::const_iterator pos =
510 std::lower_bound (my_elem_bin.begin(), my_elem_bin.end(), hi);
511 libmesh_assert (pos != my_elem_bin.end());
512 libmesh_assert_equal_to (*pos, hi);
516 global_ids.push_back (cast_int<dof_id_type>(std::distance(my_elem_bin.begin(), pos) + my_offset));
520 auto action_functor =
523 const std::vector<Parallel::DofObjectKey> &,
524 const std::vector<dof_id_type> & global_ids)
526 filled_request[pid] = global_ids;
532 (
communicator, requested_ids, gather_functor, action_functor, ex);
537 std::vector<std::vector<dof_id_type>::const_iterator>
538 next_obj_on_proc; next_obj_on_proc.reserve(
communicator.size());
540 next_obj_on_proc.push_back(filled_request[pid].begin());
544 libmesh_assert(elem);
546 get_dofobject_key (elem, bbox);
548 cast_int<processor_id_type>
549 (std::distance (elem_upper_bounds.begin(),
550 std::lower_bound(elem_upper_bounds.begin(),
551 elem_upper_bounds.end(),
555 libmesh_assert (next_obj_on_proc[pid] != filled_request[pid].
end());
557 const dof_id_type global_index = *next_obj_on_proc[pid];
558 libmesh_assert_less (global_index,
mesh.
n_elem());
559 elem->set_id() = global_index;
561 ++next_obj_on_proc[pid];
567 #else // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 571 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 574 #if defined(LIBMESH_HAVE_LIBHILBERT) && defined(LIBMESH_HAVE_MPI) 577 LOG_SCOPE (
"check_for_duplicate_global_indices()",
"MeshCommunication");
586 std::vector<Parallel::DofObjectKey>
587 node_keys, elem_keys;
594 node_keys.resize (nr.size());
600 for (std::size_t i = 0; i != node_keys.size(); ++i, ++nodei)
603 for (std::size_t j = 0; j != i; ++j, ++nodej)
605 if (node_keys[i] == node_keys[j])
607 CFixBitVec icoords[3], jcoords[3];
608 get_hilbert_coords(**nodej, bbox, jcoords);
610 "node " << (*nodej)->id() <<
", " <<
611 *(
Point *)(*nodej) <<
" has HilbertIndices " <<
612 node_keys[j] << std::endl;
613 get_hilbert_coords(**nodei, bbox, icoords);
615 "node " << (*nodei)->id() <<
", " <<
616 *(
Point *)(*nodei) <<
" has HilbertIndices " <<
617 node_keys[i] << std::endl;
618 libmesh_error_msg(
"Error: nodes with duplicate Hilbert keys!");
628 elem_keys.resize (er.size());
635 for (std::size_t i = 0; i != elem_keys.size(); ++i, ++elemi)
638 for (std::size_t j = 0; j != i; ++j, ++elemj)
640 if ((elem_keys[i] == elem_keys[j]) &&
641 ((*elemi)->level() == (*elemj)->level()))
644 "level " << (*elemj)->level() <<
" elem\n" <<
645 (**elemj) <<
" centroid " <<
646 (*elemj)->centroid() <<
" has HilbertIndices " <<
647 elem_keys[j] <<
" or " <<
648 get_dofobject_key((*elemj), bbox) <<
651 "level " << (*elemi)->level() <<
" elem\n" <<
652 (**elemi) <<
" centroid " <<
653 (*elemi)->centroid() <<
" has HilbertIndices " <<
654 elem_keys[i] <<
" or " <<
655 get_dofobject_key((*elemi), bbox) <<
657 libmesh_error_msg(
"Error: level " << (*elemi)->level() <<
" elements with duplicate Hilbert keys!");
664 #else // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 668 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 670 #if defined(LIBMESH_HAVE_LIBHILBERT) && defined(LIBMESH_HAVE_MPI) 671 template <
typename ForwardIterator>
673 const ForwardIterator & begin,
674 const ForwardIterator &
end,
675 std::unordered_map<dof_id_type, dof_id_type> & index_map)
const 677 LOG_SCOPE (
"find_local_indices()",
"MeshCommunication");
689 std::map<Parallel::DofObjectKey, dof_id_type> hilbert_keys;
691 LOG_SCOPE(
"local_hilbert_indices",
"MeshCommunication");
692 for (ForwardIterator it=begin; it!=
end; ++it)
695 hilbert_keys.emplace(hi, (*it)->id());
701 for (
auto key_val : hilbert_keys)
702 index_map[key_val.second] = cnt++;
707 template <
typename ForwardIterator>
710 const ForwardIterator & begin,
711 const ForwardIterator &
end,
712 std::vector<dof_id_type> & index_map)
const 714 LOG_SCOPE (
"find_global_indices()",
"MeshCommunication");
726 std::size_t n_objects = std::distance (begin,
end);
727 index_map.reserve(n_objects);
734 std::vector<Parallel::DofObjectKey>
737 sorted_hilbert_keys.reserve(n_objects);
738 hilbert_keys.reserve(n_objects);
740 LOG_SCOPE(
"compute_hilbert_indices()",
"MeshCommunication");
741 for (ForwardIterator it=begin; it!=
end; ++it)
744 hilbert_keys.push_back(hi);
747 sorted_hilbert_keys.push_back(hi);
752 sorted_hilbert_keys.push_back(hi);
758 START_LOG (
"parallel_sort()",
"MeshCommunication");
760 sorted_hilbert_keys);
762 STOP_LOG (
"parallel_sort()",
"MeshCommunication");
763 const std::vector<Parallel::DofObjectKey> & my_bin = sorter.
bin();
766 std::vector<unsigned int> bin_sizes(
communicator.size());
767 communicator.allgather (static_cast<unsigned int>(my_bin.size()), bin_sizes);
770 unsigned int my_offset = 0;
771 for (
unsigned int pid=0; pid<
communicator.rank(); pid++)
772 my_offset += bin_sizes[pid];
776 std::vector<Parallel::DofObjectKey>
780 upper_bounds[0] = my_bin.back();
788 if (!bin_sizes[p]) upper_bounds[p] = upper_bounds[p-1];
799 std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
802 std::map<processor_id_type, std::vector<dof_id_type>>
806 std::vector<Parallel::DofObjectKey>::const_iterator hi =
807 hilbert_keys.begin();
809 for (ForwardIterator it = begin; it !=
end; ++it)
811 libmesh_assert (hi != hilbert_keys.end());
813 std::vector<Parallel::DofObjectKey>::iterator lb =
814 std::lower_bound(upper_bounds.begin(), upper_bounds.end(),
818 cast_int<processor_id_type>
819 (std::distance (upper_bounds.begin(), lb));
823 requested_ids[pid].push_back(*hi);
828 index_map.push_back(pid);
831 auto gather_functor =
842 std::vector<dof_id_type> & global_ids)
848 const std::size_t keys_size = keys.size();
850 global_ids.reserve(keys_size);
851 for (std::size_t
idx=0;
idx != keys_size;
idx++)
854 libmesh_assert_less_equal (hilbert_indices, upper_bounds[
communicator.rank()]);
857 std::vector<Parallel::DofObjectKey>::const_iterator pos =
858 std::lower_bound (my_bin.begin(), my_bin.end(), hilbert_indices);
859 libmesh_assert (pos != my_bin.end());
865 if (*pos != hilbert_indices)
870 Hilbert::BitVecType input;
871 #ifdef LIBMESH_ENABLE_UNIQUE_ID 872 input = hilbert_indices.first;
874 input = hilbert_indices;
878 std::vector<CBigBitVec> output(3);
881 Hilbert::indexToCoords(output.data(), 8*
sizeof(Hilbert::inttype), 3, input);
887 const Real max_int_as_real =
893 Point p_hat(static_cast<Real>(output[0].racks()[0]) / max_int_as_real,
894 static_cast<Real>(output[1].racks()[0]) / max_int_as_real,
895 static_cast<Real>(output[2].racks()[0]) / max_int_as_real);
899 xmin = bbox.first(0),
900 xmax = bbox.second(0),
901 ymin = bbox.first(1),
902 ymax = bbox.second(1),
903 zmin = bbox.first(2),
904 zmax = bbox.second(2);
907 Point p(xmin + (xmax-xmin)*p_hat(0),
908 ymin + (ymax-ymin)*p_hat(1),
909 zmin + (zmax-zmin)*p_hat(2));
911 libmesh_error_msg(
"Could not find hilbert indices: " 913 <<
" corresponding to point " << p);
919 global_ids.push_back (cast_int<dof_id_type>(std::distance(my_bin.begin(), pos) + my_offset));
923 auto action_functor =
926 const std::vector<Parallel::DofObjectKey> &,
927 const std::vector<dof_id_type> & global_ids)
929 filled_request[pid] = global_ids;
934 (
communicator, requested_ids, gather_functor, action_functor, ex);
939 std::vector<std::vector<dof_id_type>::const_iterator>
940 next_obj_on_proc; next_obj_on_proc.reserve(
communicator.size());
942 next_obj_on_proc.push_back(filled_request[pid].begin());
945 for (ForwardIterator it = begin; it !=
end; ++it, cnt++)
951 libmesh_assert (next_obj_on_proc[pid] != filled_request[pid].
end());
953 const dof_id_type global_index = *next_obj_on_proc[pid];
954 index_map[cnt] = global_index;
956 ++next_obj_on_proc[pid];
961 libmesh_assert_equal_to(index_map.size(), n_objects);
963 #else // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 964 template <
typename ForwardIterator>
966 const ForwardIterator & begin,
967 const ForwardIterator &
end,
968 std::unordered_map<dof_id_type, dof_id_type> & index_map)
const 972 for (ForwardIterator it=begin; it!=
end; ++it)
973 index_map[(*it)->id()] = (index++);
976 template <
typename ForwardIterator>
979 const ForwardIterator & begin,
980 const ForwardIterator &
end,
981 std::vector<dof_id_type> & index_map)
const 984 index_map.reserve(std::distance (begin,
end));
985 unsigned int index = 0;
986 for (ForwardIterator it=begin; it!=
end; ++it)
987 index_map.push_back(index++);
989 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI 994 template void MeshCommunication::find_global_indices<MeshBase::const_node_iterator> (
const Parallel::Communicator &,
998 std::vector<dof_id_type> &)
const;
1000 template void MeshCommunication::find_global_indices<MeshBase::const_element_iterator> (
const Parallel::Communicator &,
1004 std::vector<dof_id_type> &)
const;
1005 template void MeshCommunication::find_global_indices<MeshBase::node_iterator> (
const Parallel::Communicator &,
1009 std::vector<dof_id_type> &)
const;
1011 template void MeshCommunication::find_global_indices<MeshBase::element_iterator> (
const Parallel::Communicator &,
1015 std::vector<dof_id_type> &)
const;
1016 template void MeshCommunication::find_local_indices<MeshBase::const_element_iterator> (
const libMesh::BoundingBox &,
1019 std::unordered_map<dof_id_type, dof_id_type> &)
const;
void find_global_indices(const Parallel::Communicator &communicator, const libMesh::BoundingBox &, const ForwardIterator &, const ForwardIterator &, std::vector< dof_id_type > &) const
A geometric point in (x,y,z) space associated with a DOF.
const std::vector< KeyType > & bin()
void parallel_for(const Range &range, const Body &body)
std::pair< Hilbert::HilbertIndices, unique_id_type > DofObjectKey
virtual element_iterator local_elements_begin()=0
std::size_t first_idx() const
The base class for all geometric element types.
uint8_t processor_id_type
void assign_global_indices(MeshBase &) const
unique_id_type unique_id() const
const Parallel::Communicator & comm() const
Utility class for defining generic ranges for threading.
long double max(long double a, double b)
void libmesh_ignore(const Args &...)
void check_for_duplicate_global_indices(MeshBase &) const
virtual SimpleRange< element_iterator > element_ptr_range()=0
virtual node_iterator local_nodes_end()=0
void pull_parallel_vector_data(const Communicator &comm, const MapToVectors &queries, RequestContainer &reqs, GatherFunctor &gather_data, ActionFunctor &act_on_data, const datum *example)
static const processor_id_type invalid_processor_id
virtual SimpleRange< node_iterator > node_ptr_range()=0
void find_local_indices(const libMesh::BoundingBox &, const ForwardIterator &, const ForwardIterator &, std::unordered_map< dof_id_type, dof_id_type > &) const
virtual element_iterator local_elements_end()=0
OStreamProxy err(std::cerr)
virtual node_iterator local_nodes_begin()=0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Object for performing parallel sorts using MPI.
virtual dof_id_type n_elem() const =0
A geometric point in (x,y,z) space.
virtual Point centroid() const
virtual dof_id_type n_nodes() const =0