mesh_communication_global_indices.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // Local Includes
21 #include "libmesh/libmesh_config.h"
22 #include "libmesh/libmesh_common.h"
24 #include "libmesh/mesh_base.h"
25 #include "libmesh/mesh_tools.h"
27 #include "libmesh/parallel.h"
29 #include "libmesh/parallel_sort.h"
30 #include "libmesh/parallel_sync.h"
31 #include "libmesh/elem.h"
32 #include "libmesh/elem_range.h"
33 #include "libmesh/node_range.h"
34 #ifdef LIBMESH_HAVE_LIBHILBERT
35 # include "hilbert.h"
36 #endif
37 
38 #ifdef LIBMESH_HAVE_LIBHILBERT
39 namespace { // anonymous namespace for helper functions
40 
41 using namespace libMesh;
42 
43 // Utility function to map (x,y,z) in [bbox.min, bbox.max]^3 into
44 // [0,max_inttype]^3 for computing Hilbert keys
45 void get_hilbert_coords (const Point & p,
46  const libMesh::BoundingBox & bbox,
47  CFixBitVec icoords[3])
48 {
49  static const Hilbert::inttype max_inttype = static_cast<Hilbert::inttype>(-1);
50 
51  const long double // put (x,y,z) in [0,1]^3 (don't divide by 0)
52  x = ((bbox.first(0) == bbox.second(0)) ? 0. :
53  (p(0)-bbox.first(0))/(bbox.second(0)-bbox.first(0))),
54 
55 #if LIBMESH_DIM > 1
56  y = ((bbox.first(1) == bbox.second(1)) ? 0. :
57  (p(1)-bbox.first(1))/(bbox.second(1)-bbox.first(1))),
58 #else
59  y = 0.,
60 #endif
61 
62 #if LIBMESH_DIM > 2
63  z = ((bbox.first(2) == bbox.second(2)) ? 0. :
64  (p(2)-bbox.first(2))/(bbox.second(2)-bbox.first(2)));
65 #else
66  z = 0.;
67 #endif
68 
69  // (icoords) in [0,max_inttype]^3
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);
73 }
74 
75 
76 
78 get_dofobject_key (const Elem * e,
79  const libMesh::BoundingBox & bbox)
80 {
81  static const unsigned int sizeof_inttype = sizeof(Hilbert::inttype);
82 
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);
88  index = bv;
89 
90 #ifdef LIBMESH_ENABLE_UNIQUE_ID
91  return std::make_pair(index, e->unique_id());
92 #else
93  return index;
94 #endif
95 }
96 
97 
98 
99 // Compute the hilbert index
101 get_dofobject_key (const Node * n,
102  const libMesh::BoundingBox & bbox)
103 {
104  static const unsigned int sizeof_inttype = sizeof(Hilbert::inttype);
105 
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);
111  index = bv;
112 
113 #ifdef LIBMESH_ENABLE_UNIQUE_ID
114  return std::make_pair(index, n->unique_id());
115 #else
116  return index;
117 #endif
118 }
119 
120 // Helper class for threaded Hilbert key computation
121 class ComputeHilbertKeys
122 {
123 public:
124  ComputeHilbertKeys (const libMesh::BoundingBox & bbox,
125  std::vector<Parallel::DofObjectKey> & keys) :
126  _bbox(bbox),
127  _keys(keys)
128  {}
129 
130  // computes the hilbert index for a node
131  void operator() (const ConstNodeRange & range) const
132  {
133  std::size_t pos = range.first_idx();
134  for (const auto & node : range)
135  {
136  libmesh_assert(node);
137  libmesh_assert_less (pos, _keys.size());
138  _keys[pos++] = get_dofobject_key (node, _bbox);
139  }
140  }
141 
142  // computes the hilbert index for an element
143  void operator() (const ConstElemRange & range) const
144  {
145  std::size_t pos = range.first_idx();
146  for (const auto & elem : range)
147  {
148  libmesh_assert(elem);
149  libmesh_assert_less (pos, _keys.size());
150  _keys[pos++] = get_dofobject_key (elem, _bbox);
151  }
152  }
153 
154 private:
155  const libMesh::BoundingBox & _bbox;
156  std::vector<Parallel::DofObjectKey> & _keys;
157 };
158 }
159 #endif
160 
161 
162 namespace libMesh
163 {
164 
165 // ------------------------------------------------------------
166 // MeshCommunication class members
167 #if defined(LIBMESH_HAVE_LIBHILBERT) && defined(LIBMESH_HAVE_MPI)
169 {
170  LOG_SCOPE ("assign_global_indices()", "MeshCommunication");
171 
172  // This method determines partition-agnostic global indices
173  // for nodes and elements.
174 
175  // Algorithm:
176  // (1) compute the Hilbert key for each local node/element
177  // (2) perform a parallel sort of the Hilbert key
178  // (3) get the min/max value on each processor
179  // (4) determine the position in the global ranking for
180  // each local object
181 
183 
184  // Global bounding box. We choose the nodal bounding box for
185  // backwards compatibility; the element bounding box may be looser
186  // on curved elements.
187  BoundingBox bbox =
189 
190  //-------------------------------------------------------------
191  // (1) compute Hilbert keys
192  std::vector<Parallel::DofObjectKey>
193  node_keys, elem_keys;
194 
195  {
196  // Nodes first
197  {
200  node_keys.resize (nr.size());
201  Threads::parallel_for (nr, ComputeHilbertKeys (bbox, node_keys));
202 
203  // // It's O(N^2) to check that these keys don't duplicate before the
204  // // sort...
205  // MeshBase::const_node_iterator nodei = mesh.local_nodes_begin();
206  // for (std::size_t i = 0; i != node_keys.size(); ++i, ++nodei)
207  // {
208  // MeshBase::const_node_iterator nodej = mesh.local_nodes_begin();
209  // for (std::size_t j = 0; j != i; ++j, ++nodej)
210  // {
211  // if (node_keys[i] == node_keys[j])
212  // {
213  // CFixBitVec icoords[3], jcoords[3];
214  // get_hilbert_coords(**nodej, bbox, jcoords);
215  // libMesh::err <<
216  // "node " << (*nodej)->id() << ", " <<
217  // *(Point *)(*nodej) << " has HilbertIndices " <<
218  // node_keys[j] << std::endl;
219  // get_hilbert_coords(**nodei, bbox, icoords);
220  // libMesh::err <<
221  // "node " << (*nodei)->id() << ", " <<
222  // *(Point *)(*nodei) << " has HilbertIndices " <<
223  // node_keys[i] << std::endl;
224  // libmesh_error_msg("Error: nodes with duplicate Hilbert keys!");
225  // }
226  // }
227  // }
228  }
229 
230  // Elements next
231  {
234  elem_keys.resize (er.size());
235  Threads::parallel_for (er, ComputeHilbertKeys (bbox, elem_keys));
236 
237  // // For elements, the keys can be (and in the case of TRI, are
238  // // expected to be) duplicates, but only if the elements are at
239  // // different levels
240  // MeshBase::const_element_iterator elemi = mesh.local_elements_begin();
241  // for (std::size_t i = 0; i != elem_keys.size(); ++i, ++elemi)
242  // {
243  // MeshBase::const_element_iterator elemj = mesh.local_elements_begin();
244  // for (std::size_t j = 0; j != i; ++j, ++elemj)
245  // {
246  // if ((elem_keys[i] == elem_keys[j]) &&
247  // ((*elemi)->level() == (*elemj)->level()))
248  // {
249  // libMesh::err <<
250  // "level " << (*elemj)->level() << " elem\n" <<
251  // (**elemj) << " centroid " <<
252  // (*elemj)->centroid() << " has HilbertIndices " <<
253  // elem_keys[j] << " or " <<
254  // get_dofobject_key((*elemj), bbox) <<
255  // std::endl;
256  // libMesh::err <<
257  // "level " << (*elemi)->level() << " elem\n" <<
258  // (**elemi) << " centroid " <<
259  // (*elemi)->centroid() << " has HilbertIndices " <<
260  // elem_keys[i] << " or " <<
261  // get_dofobject_key((*elemi), bbox) <<
262  // std::endl;
263  // libmesh_error_msg("Error: level " << (*elemi)->level() << " elements with duplicate Hilbert keys!");
264  // }
265  // }
266  // }
267  }
268  } // done computing Hilbert keys
269 
270 
271 
272  //-------------------------------------------------------------
273  // (2) parallel sort the Hilbert keys
275  node_keys);
276  node_sorter.sort(); /* done with node_keys */ //node_keys.clear();
277 
278  const std::vector<Parallel::DofObjectKey> & my_node_bin =
279  node_sorter.bin();
280 
282  elem_keys);
283  elem_sorter.sort(); /* done with elem_keys */ //elem_keys.clear();
284 
285  const std::vector<Parallel::DofObjectKey> & my_elem_bin =
286  elem_sorter.bin();
287 
288 
289 
290  //-------------------------------------------------------------
291  // (3) get the max value on each processor
292  std::vector<Parallel::DofObjectKey>
293  node_upper_bounds(communicator.size()),
294  elem_upper_bounds(communicator.size());
295 
296  { // limit scope of temporaries
297  std::vector<Parallel::DofObjectKey> recvbuf(2*communicator.size());
298  std::vector<unsigned short int> /* do not use a vector of bools here since it is not always so! */
299  empty_nodes (communicator.size()),
300  empty_elem (communicator.size());
301  std::vector<Parallel::DofObjectKey> my_max(2);
302 
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);
305 
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();
308 
309  communicator.allgather (my_max, /* identical_buffer_sizes = */ true);
310 
311  // Be careful here. The *_upper_bounds will be used to find the processor
312  // a given object belongs to. So, if a processor contains no objects (possible!)
313  // then copy the bound from the lower processor id.
314  for (processor_id_type p=0; p<communicator.size(); p++)
315  {
316  node_upper_bounds[p] = my_max[2*p+0];
317  elem_upper_bounds[p] = my_max[2*p+1];
318 
319  if (p > 0) // default hilbert index value is the OK upper bound for processor 0.
320  {
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];
323  }
324  }
325  }
326 
327 
328 
329  //-------------------------------------------------------------
330  // (4) determine the position in the global ranking for
331  // each local object
332  {
333  //----------------------------------------------
334  // Nodes first -- all nodes, not just local ones
335  {
336  // Request sets to send to each processor
337  std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
338  requested_ids;
339  // Results to gather from each processor - kept in a map so we
340  // do only one loop over nodes after all receives are done.
341  std::map<dof_id_type, std::vector<dof_id_type>>
342  filled_request;
343 
344  // build up list of requests
345  for (const auto & node : mesh.node_ptr_range())
346  {
347  libmesh_assert(node);
348  const Parallel::DofObjectKey hi =
349  get_dofobject_key (node, bbox);
350  const processor_id_type pid =
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(),
355  hi)));
356 
357  libmesh_assert_less (pid, communicator.size());
358 
359  requested_ids[pid].push_back(hi);
360  }
361 
362  // The number of objects in my_node_bin on each processor
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);
365 
366  // The offset of my first global index
367  dof_id_type my_offset = 0;
368  for (processor_id_type pid=0; pid<communicator.rank(); pid++)
369  my_offset += node_bin_sizes[pid];
370 
371  auto gather_functor =
372  [
373 #ifndef NDEBUG
374  & node_upper_bounds,
375  & communicator,
376 #endif
377  & my_node_bin,
378  my_offset
379  ]
381  const std::vector<Parallel::DofObjectKey> & keys,
382  std::vector<dof_id_type> & global_ids)
383  {
384  // Fill the requests
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++)
388  {
389  const Parallel::DofObjectKey & hi = keys[idx];
390  libmesh_assert_less_equal (hi, node_upper_bounds[communicator.rank()]);
391 
392  // find the requested index in my node bin
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);
397 
398  // Finally, assign the global index based off the position of the index
399  // in my array, properly offset.
400  global_ids.push_back(cast_int<dof_id_type>(std::distance(my_node_bin.begin(), pos) + my_offset));
401  }
402  };
403 
404  auto action_functor =
405  [&filled_request]
406  (processor_id_type pid,
407  const std::vector<Parallel::DofObjectKey> &,
408  const std::vector<dof_id_type> & global_ids)
409  {
410  filled_request[pid] = global_ids;
411  };
412 
413  // Trade requests with other processors
414  const dof_id_type * ex = nullptr;
416  (communicator, requested_ids, gather_functor, action_functor, ex);
417 
418  // We now have all the filled requests, so we can loop through our
419  // nodes once and assign the global index to each one.
420  {
421  std::map<dof_id_type, std::vector<dof_id_type>::const_iterator>
422  next_obj_on_proc;
423  for (auto & p : filled_request)
424  next_obj_on_proc[p.first] = p.second.begin();
425 
426  for (auto & node : mesh.node_ptr_range())
427  {
428  libmesh_assert(node);
429  const Parallel::DofObjectKey hi =
430  get_dofobject_key (node, bbox);
431  const processor_id_type pid =
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(),
436  hi)));
437 
438  libmesh_assert_less (pid, communicator.size());
439  libmesh_assert (next_obj_on_proc[pid] != filled_request[pid].end());
440 
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;
444 
445  ++next_obj_on_proc[pid];
446  }
447  }
448  }
449 
450  //---------------------------------------------------
451  // elements next -- all elements, not just local ones
452  {
453  // Request sets to send to each processor
454  std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
455  requested_ids;
456  // Results to gather from each processor - kept in a map so we
457  // do only one loop over elements after all receives are done.
458  std::map<dof_id_type, std::vector<dof_id_type>>
459  filled_request;
460 
461  for (const auto & elem : mesh.element_ptr_range())
462  {
463  libmesh_assert(elem);
464  const Parallel::DofObjectKey hi =
465  get_dofobject_key (elem, bbox);
466  const processor_id_type pid =
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(),
471  hi)));
472 
473  libmesh_assert_less (pid, communicator.size());
474 
475  requested_ids[pid].push_back(hi);
476  }
477 
478  // The number of objects in my_elem_bin on each processor
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);
481 
482  // The offset of my first global index
483  dof_id_type my_offset = 0;
484  for (processor_id_type pid=0; pid<communicator.rank(); pid++)
485  my_offset += elem_bin_sizes[pid];
486 
487  auto gather_functor =
488  [
489 #ifndef NDEBUG
490  & elem_upper_bounds,
491  & communicator,
492 #endif
493  & my_elem_bin,
494  my_offset
495  ]
497  const std::vector<Parallel::DofObjectKey> & keys,
498  std::vector<dof_id_type> & global_ids)
499  {
500  // Fill the requests
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++)
504  {
505  const Parallel::DofObjectKey & hi = keys[idx];
506  libmesh_assert_less_equal (hi, elem_upper_bounds[communicator.rank()]);
507 
508  // find the requested index in my elem bin
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);
513 
514  // Finally, assign the global index based off the position of the index
515  // in my array, properly offset.
516  global_ids.push_back (cast_int<dof_id_type>(std::distance(my_elem_bin.begin(), pos) + my_offset));
517  }
518  };
519 
520  auto action_functor =
521  [&filled_request]
522  (processor_id_type pid,
523  const std::vector<Parallel::DofObjectKey> &,
524  const std::vector<dof_id_type> & global_ids)
525  {
526  filled_request[pid] = global_ids;
527  };
528 
529  // Trade requests with other processors
530  const dof_id_type * ex = nullptr;
532  (communicator, requested_ids, gather_functor, action_functor, ex);
533 
534  // We now have all the filled requests, so we can loop through our
535  // elements once and assign the global index to each one.
536  {
537  std::vector<std::vector<dof_id_type>::const_iterator>
538  next_obj_on_proc; next_obj_on_proc.reserve(communicator.size());
539  for (processor_id_type pid=0; pid<communicator.size(); pid++)
540  next_obj_on_proc.push_back(filled_request[pid].begin());
541 
542  for (auto & elem : mesh.element_ptr_range())
543  {
544  libmesh_assert(elem);
545  const Parallel::DofObjectKey hi =
546  get_dofobject_key (elem, bbox);
547  const processor_id_type pid =
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(),
552  hi)));
553 
554  libmesh_assert_less (pid, communicator.size());
555  libmesh_assert (next_obj_on_proc[pid] != filled_request[pid].end());
556 
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;
560 
561  ++next_obj_on_proc[pid];
562  }
563  }
564  }
565  }
566 }
567 #else // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI
569 {
570 }
571 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI
572 
573 
574 #if defined(LIBMESH_HAVE_LIBHILBERT) && defined(LIBMESH_HAVE_MPI)
576 {
577  LOG_SCOPE ("check_for_duplicate_global_indices()", "MeshCommunication");
578 
579  // Global bounding box. We choose the nodal bounding box for
580  // backwards compatibility; the element bounding box may be looser
581  // on curved elements.
582  BoundingBox bbox =
584 
585 
586  std::vector<Parallel::DofObjectKey>
587  node_keys, elem_keys;
588 
589  {
590  // Nodes first
591  {
594  node_keys.resize (nr.size());
595  Threads::parallel_for (nr, ComputeHilbertKeys (bbox, node_keys));
596 
597  // It's O(N^2) to check that these keys don't duplicate before the
598  // sort...
600  for (std::size_t i = 0; i != node_keys.size(); ++i, ++nodei)
601  {
603  for (std::size_t j = 0; j != i; ++j, ++nodej)
604  {
605  if (node_keys[i] == node_keys[j])
606  {
607  CFixBitVec icoords[3], jcoords[3];
608  get_hilbert_coords(**nodej, bbox, jcoords);
609  libMesh::err <<
610  "node " << (*nodej)->id() << ", " <<
611  *(Point *)(*nodej) << " has HilbertIndices " <<
612  node_keys[j] << std::endl;
613  get_hilbert_coords(**nodei, bbox, icoords);
614  libMesh::err <<
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!");
619  }
620  }
621  }
622  }
623 
624  // Elements next
625  {
628  elem_keys.resize (er.size());
629  Threads::parallel_for (er, ComputeHilbertKeys (bbox, elem_keys));
630 
631  // For elements, the keys can be (and in the case of TRI, are
632  // expected to be) duplicates, but only if the elements are at
633  // different levels
635  for (std::size_t i = 0; i != elem_keys.size(); ++i, ++elemi)
636  {
638  for (std::size_t j = 0; j != i; ++j, ++elemj)
639  {
640  if ((elem_keys[i] == elem_keys[j]) &&
641  ((*elemi)->level() == (*elemj)->level()))
642  {
643  libMesh::err <<
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) <<
649  std::endl;
650  libMesh::err <<
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) <<
656  std::endl;
657  libmesh_error_msg("Error: level " << (*elemi)->level() << " elements with duplicate Hilbert keys!");
658  }
659  }
660  }
661  }
662  } // done checking Hilbert keys
663 }
664 #else // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI
666 {
667 }
668 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI
669 
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
676 {
677  LOG_SCOPE ("find_local_indices()", "MeshCommunication");
678 
679  // This method determines id-agnostic local indices
680  // for nodes and elements by sorting Hilbert keys.
681 
682  index_map.clear();
683 
684  //-------------------------------------------------------------
685  // (1) compute Hilbert keys
686  // These aren't trivial to compute, and we will need them again.
687  // But the binsort will sort the input vector, trashing the order
688  // that we'd like to rely on. So, two vectors...
689  std::map<Parallel::DofObjectKey, dof_id_type> hilbert_keys;
690  {
691  LOG_SCOPE("local_hilbert_indices", "MeshCommunication");
692  for (ForwardIterator it=begin; it!=end; ++it)
693  {
694  const Parallel::DofObjectKey hi(get_dofobject_key ((*it), bbox));
695  hilbert_keys.emplace(hi, (*it)->id());
696  }
697  }
698 
699  {
700  dof_id_type cnt = 0;
701  for (auto key_val : hilbert_keys)
702  index_map[key_val.second] = cnt++;
703  }
704 }
705 
706 
707 template <typename ForwardIterator>
709  const libMesh::BoundingBox & bbox,
710  const ForwardIterator & begin,
711  const ForwardIterator & end,
712  std::vector<dof_id_type> & index_map) const
713 {
714  LOG_SCOPE ("find_global_indices()", "MeshCommunication");
715 
716  // This method determines partition-agnostic global indices
717  // for nodes and elements.
718 
719  // Algorithm:
720  // (1) compute the Hilbert key for each local node/element
721  // (2) perform a parallel sort of the Hilbert key
722  // (3) get the min/max value on each processor
723  // (4) determine the position in the global ranking for
724  // each local object
725  index_map.clear();
726  std::size_t n_objects = std::distance (begin, end);
727  index_map.reserve(n_objects);
728 
729  //-------------------------------------------------------------
730  // (1) compute Hilbert keys
731  // These aren't trivial to compute, and we will need them again.
732  // But the binsort will sort the input vector, trashing the order
733  // that we'd like to rely on. So, two vectors...
734  std::vector<Parallel::DofObjectKey>
735  sorted_hilbert_keys,
736  hilbert_keys;
737  sorted_hilbert_keys.reserve(n_objects);
738  hilbert_keys.reserve(n_objects);
739  {
740  LOG_SCOPE("compute_hilbert_indices()", "MeshCommunication");
741  for (ForwardIterator it=begin; it!=end; ++it)
742  {
743  const Parallel::DofObjectKey hi(get_dofobject_key (*it, bbox));
744  hilbert_keys.push_back(hi);
745 
746  if ((*it)->processor_id() == communicator.rank())
747  sorted_hilbert_keys.push_back(hi);
748 
749  // someone needs to take care of unpartitioned objects!
750  if ((communicator.rank() == 0) &&
751  ((*it)->processor_id() == DofObject::invalid_processor_id))
752  sorted_hilbert_keys.push_back(hi);
753  }
754  }
755 
756  //-------------------------------------------------------------
757  // (2) parallel sort the Hilbert keys
758  START_LOG ("parallel_sort()", "MeshCommunication");
760  sorted_hilbert_keys);
761  sorter.sort();
762  STOP_LOG ("parallel_sort()", "MeshCommunication");
763  const std::vector<Parallel::DofObjectKey> & my_bin = sorter.bin();
764 
765  // The number of objects in my_bin on each processor
766  std::vector<unsigned int> bin_sizes(communicator.size());
767  communicator.allgather (static_cast<unsigned int>(my_bin.size()), bin_sizes);
768 
769  // The offset of my first global index
770  unsigned int my_offset = 0;
771  for (unsigned int pid=0; pid<communicator.rank(); pid++)
772  my_offset += bin_sizes[pid];
773 
774  //-------------------------------------------------------------
775  // (3) get the max value on each processor
776  std::vector<Parallel::DofObjectKey>
777  upper_bounds(1);
778 
779  if (!my_bin.empty())
780  upper_bounds[0] = my_bin.back();
781 
782  communicator.allgather (upper_bounds, /* identical_buffer_sizes = */ true);
783 
784  // Be careful here. The *_upper_bounds will be used to find the processor
785  // a given object belongs to. So, if a processor contains no objects (possible!)
786  // then copy the bound from the lower processor id.
787  for (unsigned int p=1; p<communicator.size(); p++)
788  if (!bin_sizes[p]) upper_bounds[p] = upper_bounds[p-1];
789 
790 
791  //-------------------------------------------------------------
792  // (4) determine the position in the global ranking for
793  // each local object
794  {
795  //----------------------------------------------
796  // all objects, not just local ones
797 
798  // Request sets to send to each processor
799  std::map<processor_id_type, std::vector<Parallel::DofObjectKey>>
800  requested_ids;
801  // Results to gather from each processor
802  std::map<processor_id_type, std::vector<dof_id_type>>
803  filled_request;
804 
805  // build up list of requests
806  std::vector<Parallel::DofObjectKey>::const_iterator hi =
807  hilbert_keys.begin();
808 
809  for (ForwardIterator it = begin; it != end; ++it)
810  {
811  libmesh_assert (hi != hilbert_keys.end());
812 
813  std::vector<Parallel::DofObjectKey>::iterator lb =
814  std::lower_bound(upper_bounds.begin(), upper_bounds.end(),
815  *hi);
816 
817  const processor_id_type pid =
818  cast_int<processor_id_type>
819  (std::distance (upper_bounds.begin(), lb));
820 
821  libmesh_assert_less (pid, communicator.size());
822 
823  requested_ids[pid].push_back(*hi);
824 
825  ++hi;
826  // go ahead and put pid in index_map, that way we
827  // don't have to repeat the std::lower_bound()
828  index_map.push_back(pid);
829  }
830 
831  auto gather_functor =
832  [
833 #ifndef NDEBUG
834  & upper_bounds,
835  & communicator,
836 #endif
837  & bbox,
838  & my_bin,
839  my_offset
840  ]
841  (processor_id_type, const std::vector<Parallel::DofObjectKey> & keys,
842  std::vector<dof_id_type> & global_ids)
843  {
844  // Ignore unused lambda capture warnings in devel mode
845  libmesh_ignore(bbox);
846 
847  // Fill the requests
848  const std::size_t keys_size = keys.size();
849  global_ids.clear();
850  global_ids.reserve(keys_size);
851  for (std::size_t idx=0; idx != keys_size; idx++)
852  {
853  const Parallel::DofObjectKey & hilbert_indices = keys[idx];
854  libmesh_assert_less_equal (hilbert_indices, upper_bounds[communicator.rank()]);
855 
856  // find the requested index in my node bin
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());
860 #ifdef DEBUG
861  // If we could not find the requested Hilbert index in
862  // my_bin, something went terribly wrong, possibly the
863  // Mesh was displaced differently on different processors,
864  // and therefore the Hilbert indices don't agree.
865  if (*pos != hilbert_indices)
866  {
867  // The input will be hilbert_indices. We convert it
868  // to BitVecType using the operator= provided by the
869  // BitVecType class. BitVecType is a CBigBitVec!
870  Hilbert::BitVecType input;
871 #ifdef LIBMESH_ENABLE_UNIQUE_ID
872  input = hilbert_indices.first;
873 #else
874  input = hilbert_indices;
875 #endif
876 
877  // Get output in a vector of CBigBitVec
878  std::vector<CBigBitVec> output(3);
879 
880  // Call the indexToCoords function
881  Hilbert::indexToCoords(output.data(), 8*sizeof(Hilbert::inttype), 3, input);
882 
883  // The entries in the output racks are integers in the
884  // range [0, Hilbert::inttype::max] which can be
885  // converted to floating point values in [0,1] and
886  // finally to actual values using the bounding box.
887  const Real max_int_as_real =
889 
890  // Get the points in [0,1]^3. The zeroth rack of each entry in
891  // 'output' maps to the normalized x, y, and z locations,
892  // respectively.
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);
896 
897  // Convert the points from [0,1]^3 to their actual (x,y,z) locations
898  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);
905 
906  // Convert the points from [0,1]^3 to their actual (x,y,z) locations
907  Point p(xmin + (xmax-xmin)*p_hat(0),
908  ymin + (ymax-ymin)*p_hat(1),
909  zmin + (zmax-zmin)*p_hat(2));
910 
911  libmesh_error_msg("Could not find hilbert indices: "
912  << hilbert_indices
913  << " corresponding to point " << p);
914  }
915 #endif
916 
917  // Finally, assign the global index based off the position of the index
918  // in my array, properly offset.
919  global_ids.push_back (cast_int<dof_id_type>(std::distance(my_bin.begin(), pos) + my_offset));
920  }
921  };
922 
923  auto action_functor =
924  [&filled_request]
925  (processor_id_type pid,
926  const std::vector<Parallel::DofObjectKey> &,
927  const std::vector<dof_id_type> & global_ids)
928  {
929  filled_request[pid] = global_ids;
930  };
931 
932  const dof_id_type * ex = nullptr;
934  (communicator, requested_ids, gather_functor, action_functor, ex);
935 
936  // We now have all the filled requests, so we can loop through our
937  // nodes once and assign the global index to each one.
938  {
939  std::vector<std::vector<dof_id_type>::const_iterator>
940  next_obj_on_proc; next_obj_on_proc.reserve(communicator.size());
941  for (processor_id_type pid=0; pid<communicator.size(); pid++)
942  next_obj_on_proc.push_back(filled_request[pid].begin());
943 
944  unsigned int cnt=0;
945  for (ForwardIterator it = begin; it != end; ++it, cnt++)
946  {
947  const processor_id_type pid = cast_int<processor_id_type>
948  (index_map[cnt]);
949 
950  libmesh_assert_less (pid, communicator.size());
951  libmesh_assert (next_obj_on_proc[pid] != filled_request[pid].end());
952 
953  const dof_id_type global_index = *next_obj_on_proc[pid];
954  index_map[cnt] = global_index;
955 
956  ++next_obj_on_proc[pid];
957  }
958  }
959  }
960 
961  libmesh_assert_equal_to(index_map.size(), n_objects);
962 }
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
969 {
970  index_map.clear();
971  dof_id_type index = 0;
972  for (ForwardIterator it=begin; it!=end; ++it)
973  index_map[(*it)->id()] = (index++);
974 }
975 
976 template <typename ForwardIterator>
978  const libMesh::BoundingBox &,
979  const ForwardIterator & begin,
980  const ForwardIterator & end,
981  std::vector<dof_id_type> & index_map) const
982 {
983  index_map.clear();
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++);
988 }
989 #endif // LIBMESH_HAVE_LIBHILBERT, LIBMESH_HAVE_MPI
990 
991 
992 
993 //------------------------------------------------------------------
994 template void MeshCommunication::find_global_indices<MeshBase::const_node_iterator> (const Parallel::Communicator &,
995  const libMesh::BoundingBox &,
998  std::vector<dof_id_type> &) const;
999 
1000 template void MeshCommunication::find_global_indices<MeshBase::const_element_iterator> (const Parallel::Communicator &,
1001  const libMesh::BoundingBox &,
1004  std::vector<dof_id_type> &) const;
1005 template void MeshCommunication::find_global_indices<MeshBase::node_iterator> (const Parallel::Communicator &,
1006  const libMesh::BoundingBox &,
1007  const MeshBase::node_iterator &,
1008  const MeshBase::node_iterator &,
1009  std::vector<dof_id_type> &) const;
1010 
1011 template void MeshCommunication::find_global_indices<MeshBase::element_iterator> (const Parallel::Communicator &,
1012  const libMesh::BoundingBox &,
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;
1020 
1021 } // namespace libMesh
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.
Definition: node.h:52
const std::vector< KeyType > & bin()
void parallel_for(const Range &range, const Body &body)
Definition: threads_none.h:73
std::pair< Hilbert::HilbertIndices, unique_id_type > DofObjectKey
MPI_Comm communicator
Definition: communicator.h:57
virtual element_iterator local_elements_begin()=0
std::size_t first_idx() const
Definition: stored_range.h:271
The base class for all geometric element types.
Definition: elem.h:100
MeshBase & mesh
uint8_t processor_id_type
Definition: id_types.h:99
unique_id_type unique_id() const
Definition: dof_object.h:672
const Parallel::Communicator & comm() const
IterBase * end
Utility class for defining generic ranges for threading.
Definition: stored_range.h:52
long double max(long double a, double b)
Base class for Mesh.
Definition: mesh_base.h:77
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
Definition: dof_object.h:358
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.
Definition: parallel_sort.h:53
libMesh::BoundingBox create_nodal_bounding_box(const MeshBase &mesh)
Definition: mesh_tools.C:411
virtual dof_id_type n_elem() const =0
A geometric point in (x,y,z) space.
Definition: point.h:38
virtual Point centroid() const
Definition: elem.C:344
virtual dof_id_type n_nodes() const =0
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
uint8_t dof_id_type
Definition: id_types.h:64