libMesh
mesh_tools.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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/elem.h"
22 #include "libmesh/elem_range.h"
23 #include "libmesh/libmesh_logging.h"
24 #include "libmesh/mesh_base.h"
25 #include "libmesh/mesh_communication.h"
26 #include "libmesh/mesh_serializer.h"
27 #include "libmesh/mesh_tools.h"
28 #include "libmesh/node_range.h"
29 #include "libmesh/parallel.h"
30 #include "libmesh/parallel_algebra.h"
31 #include "libmesh/parallel_ghost_sync.h"
32 #include "libmesh/sphere.h"
33 #include "libmesh/threads.h"
34 #include "libmesh/enum_to_string.h"
35 #include "libmesh/enum_elem_type.h"
36 #include "libmesh/int_range.h"
37 #include "libmesh/utility.h"
38 #include "libmesh/boundary_info.h"
39 
40 #ifndef NDEBUG
41 # include "libmesh/remote_elem.h"
42 #endif
43 
44 // C++ includes
45 #include <limits>
46 #include <numeric> // for std::accumulate
47 #include <set>
48 #include <unordered_map>
49 #include <unordered_set>
50 
51 
52 
53 // ------------------------------------------------------------
54 // anonymous namespace for helper classes
55 namespace {
56 
57 using namespace libMesh;
58 
66 class SumElemWeight
67 {
68 public:
69  SumElemWeight () :
70  _weight(0)
71  {}
72 
73  SumElemWeight (SumElemWeight &, Threads::split) :
74  _weight(0)
75  {}
76 
77  void operator()(const ConstElemRange & range)
78  {
79  for (const auto & elem : range)
80  _weight += elem->n_nodes();
81  }
82 
83  dof_id_type weight() const
84  { return _weight; }
85 
86  // If we don't have threads we never need a join, and icpc yells a
87  // warning if it sees an anonymous function that's never used
88 #if LIBMESH_USING_THREADS
89  void join (const SumElemWeight & other)
90  { _weight += other.weight(); }
91 #endif
92 
93 private:
94  dof_id_type _weight;
95 };
96 
97 
104 class FindBBox
105 {
106 public:
107  FindBBox () : _bbox()
108  {}
109 
110  FindBBox (FindBBox & other, Threads::split) :
111  _bbox(other._bbox)
112  {}
113 
114  void operator()(const ConstNodeRange & range)
115  {
116  for (const auto & node : range)
117  {
118  libmesh_assert(node);
119  _bbox.union_with(*node);
120  }
121  }
122 
123  void operator()(const ConstElemRange & range)
124  {
125  for (const auto & elem : range)
126  {
127  libmesh_assert(elem);
128  _bbox.union_with(elem->loose_bounding_box());
129  }
130  }
131 
132  Point & min() { return _bbox.min(); }
133 
134  Point & max() { return _bbox.max(); }
135 
136  // If we don't have threads we never need a join, and icpc yells a
137  // warning if it sees an anonymous function that's never used
138 #if LIBMESH_USING_THREADS
139  void join (const FindBBox & other)
140  {
141  _bbox.union_with(other._bbox);
142  }
143 #endif
144 
145  libMesh::BoundingBox & bbox ()
146  {
147  return _bbox;
148  }
149 
150 private:
151  BoundingBox _bbox;
152 };
153 
154 #ifdef DEBUG
155 void assert_semiverify_dofobj(const Parallel::Communicator & communicator,
156  const DofObject * d,
157  unsigned int sysnum = libMesh::invalid_uint)
158 {
159  if (d)
160  {
161  const unsigned int n_sys = d->n_systems();
162 
163  std::vector<unsigned int> n_vars (n_sys, 0);
164  for (unsigned int s = 0; s != n_sys; ++s)
165  if (sysnum == s ||
166  sysnum == libMesh::invalid_uint)
167  n_vars[s] = d->n_vars(s);
168 
169  const unsigned int tot_n_vars =
170  std::accumulate(n_vars.begin(), n_vars.end(), 0);
171 
172  std::vector<unsigned int> n_comp (tot_n_vars, 0);
173  std::vector<dof_id_type> first_dof (tot_n_vars, 0);
174 
175  for (unsigned int s = 0, i=0; s != n_sys; ++s)
176  {
177  if (sysnum != s &&
178  sysnum != libMesh::invalid_uint)
179  continue;
180 
181  for (unsigned int v = 0; v != n_vars[s]; ++v, ++i)
182  {
183  n_comp[i] = d->n_comp(s,v);
184  first_dof[i] = n_comp[i] ? d->dof_number(s,v,0) : DofObject::invalid_id;
185  }
186  }
187 
188  libmesh_assert(communicator.semiverify(&n_sys));
189  libmesh_assert(communicator.semiverify(&n_vars));
190  libmesh_assert(communicator.semiverify(&n_comp));
191  libmesh_assert(communicator.semiverify(&first_dof));
192  }
193  else
194  {
195  const unsigned int * p_ui = nullptr;
196  const std::vector<unsigned int> * p_vui = nullptr;
197  const std::vector<dof_id_type> * p_vdid = nullptr;
198 
199  libmesh_assert(communicator.semiverify(p_ui));
200  libmesh_assert(communicator.semiverify(p_vui));
201  libmesh_assert(communicator.semiverify(p_vui));
202  libmesh_assert(communicator.semiverify(p_vdid));
203  }
204 }
205 
206 
207 
208 #ifdef LIBMESH_ENABLE_UNIQUE_ID
209 void assert_dofobj_unique_id(const Parallel::Communicator & comm,
210  const DofObject * d,
211  const std::unordered_set<unique_id_type> & unique_ids)
212 {
213  // Duplicating some semiverify code here so we can reuse
214  // tempmin,tempmax afterward
215 
216  unique_id_type tempmin, tempmax;
217  if (d)
218  {
219  tempmin = tempmax = d->unique_id();
220  }
221  else
222  {
225  }
226  comm.min(tempmin);
227  comm.max(tempmax);
228  bool invalid = d && ((d->unique_id() != tempmin) ||
229  (d->unique_id() != tempmax));
230  comm.max(invalid);
231 
232  // First verify that everything is in sync
233  libmesh_assert(!invalid);
234 
235  // Then verify that any remote id doesn't duplicate a local one.
236  if (!d && tempmin == tempmax)
237  libmesh_assert(!unique_ids.count(tempmin));
238 }
239 #endif // LIBMESH_ENABLE_UNIQUE_ID
240 #endif // DEBUG
241 
242 void find_nodal_neighbors_helper(const dof_id_type global_id,
243  const std::vector<const Elem *> & node_to_elem_vec,
244  std::vector<const Node *> & neighbors)
245 {
246  // We'll construct a std::set<const Node *> for more efficient
247  // searching while finding the nodal neighbors, and return it to the
248  // user in a std::vector.
249  std::set<const Node *> neighbor_set;
250 
251  // Look through the elements that contain this node
252  // find the local node id... then find the side that
253  // node lives on in the element
254  // next, look for the _other_ node on that side
255  // That other node is a "nodal_neighbor"... save it
256  for (const auto & elem : node_to_elem_vec)
257  {
258  // We only care about active elements...
259  if (elem->active())
260  {
261  // Which local node number is global_id?
262  unsigned local_node_number = elem->local_node(global_id);
263 
264  // Make sure it was found
265  libmesh_assert_not_equal_to(local_node_number, libMesh::invalid_uint);
266 
267  const unsigned short n_edges = elem->n_edges();
268 
269  // If this element has no edges, the edge-based algorithm below doesn't make sense.
270  if (!n_edges)
271  {
272  switch (elem->type())
273  {
274  case EDGE2:
275  {
276  switch (local_node_number)
277  {
278  case 0:
279  // The other node is a nodal neighbor
280  neighbor_set.insert(elem->node_ptr(1));
281  break;
282 
283  case 1:
284  // The other node is a nodal neighbor
285  neighbor_set.insert(elem->node_ptr(0));
286  break;
287 
288  default:
289  libmesh_error_msg("Invalid local node number: " << local_node_number << " found." << std::endl);
290  }
291  break;
292  }
293 
294  case EDGE3:
295  {
296  switch (local_node_number)
297  {
298  // The outside nodes have node 2 as a neighbor
299  case 0:
300  case 1:
301  neighbor_set.insert(elem->node_ptr(2));
302  break;
303 
304  // The middle node has the outer nodes as neighbors
305  case 2:
306  neighbor_set.insert(elem->node_ptr(0));
307  neighbor_set.insert(elem->node_ptr(1));
308  break;
309 
310  default:
311  libmesh_error_msg("Invalid local node number: " << local_node_number << " found." << std::endl);
312  }
313  break;
314  }
315 
316  case EDGE4:
317  {
318  switch (local_node_number)
319  {
320  case 0:
321  // The left-middle node is a nodal neighbor
322  neighbor_set.insert(elem->node_ptr(2));
323  break;
324 
325  case 1:
326  // The right-middle node is a nodal neighbor
327  neighbor_set.insert(elem->node_ptr(3));
328  break;
329 
330  // The left-middle node
331  case 2:
332  neighbor_set.insert(elem->node_ptr(0));
333  neighbor_set.insert(elem->node_ptr(3));
334  break;
335 
336  // The right-middle node
337  case 3:
338  neighbor_set.insert(elem->node_ptr(1));
339  neighbor_set.insert(elem->node_ptr(2));
340  break;
341 
342  default:
343  libmesh_error_msg("Invalid local node number: " << local_node_number << " found." << std::endl);
344  }
345  break;
346  }
347 
348  default:
349  libmesh_error_msg("Unrecognized ElemType: " << Utility::enum_to_string(elem->type()) << std::endl);
350  }
351  }
352 
353  const auto elem_order = Elem::type_to_default_order_map[elem->type()];
354 
355  // Index of the current edge
356  unsigned current_edge = 0;
357 
358  const unsigned short n_nodes = elem->n_nodes();
359 
360  while (current_edge < n_edges)
361  {
362  // Find the edge the node is on
363  bool found_edge = false;
364  for (; current_edge<n_edges; ++current_edge)
365  if (elem->is_node_on_edge(local_node_number, current_edge))
366  {
367  found_edge = true;
368  break;
369  }
370 
371  // Did we find one?
372  if (found_edge)
373  {
374  const Node * node_to_save = nullptr;
375 
376  // Find another node in this element on this edge
377  for (unsigned other_node_this_edge = 0; other_node_this_edge != n_nodes; other_node_this_edge++)
378  {
379  const bool both_vertices = elem->is_vertex(local_node_number) &&
380  elem->is_vertex(other_node_this_edge);
381  if ( elem->is_node_on_edge(other_node_this_edge, current_edge) && // On the current edge
382  elem->node_id(other_node_this_edge) != global_id && // But not the original node
383  // vertex nodes on the same edge of higher order elements are not nodal neighbors
384  (elem_order == 1 || !both_vertices))
385  {
386  // We've found a nodal neighbor! Save a pointer to it..
387  node_to_save = elem->node_ptr(other_node_this_edge);
388 
389  // Make sure we found something
390  libmesh_assert(node_to_save != nullptr);
391 
392  neighbor_set.insert(node_to_save);
393  }
394  }
395  }
396 
397  // Keep looking for edges, node may be on more than one edge
398  current_edge++;
399  }
400  } // if (elem->active())
401  } // for
402 
403  // Assign the entries from the set to the vector. Note: this
404  // replaces any existing contents in neighbors and modifies its size
405  // accordingly.
406  neighbors.assign(neighbor_set.begin(), neighbor_set.end());
407 }
408 
409 }
410 
411 
412 namespace libMesh
413 {
414 
415 // ------------------------------------------------------------
416 // MeshTools functions
417 
418 namespace MeshTools
419 {
420 
422 {
423  if (!mesh.is_serial())
424  {
425  libmesh_parallel_only(mesh.comm());
427  mesh.comm().sum(weight);
428  dof_id_type unpartitioned_weight =
430  return weight + unpartitioned_weight;
431  }
432 
433  SumElemWeight sew;
434 
435  Threads::parallel_reduce (ConstElemRange (mesh.elements_begin(),
436  mesh.elements_end()),
437  sew);
438  return sew.weight();
439 
440 }
441 
442 
443 
445 {
446  SumElemWeight sew;
447 
448  Threads::parallel_reduce (ConstElemRange (mesh.pid_elements_begin(pid),
449  mesh.pid_elements_end(pid)),
450  sew);
451  return sew.weight();
452 }
453 
454 
455 
457  std::vector<std::vector<dof_id_type>> & nodes_to_elem_map)
458 {
459  // A vector indexed over all nodes is too inefficient to use for a
460  // distributed mesh. Use the unordered_map API instead.
461  if (!mesh.is_serial())
462  libmesh_deprecated();
463 
464  nodes_to_elem_map.resize (mesh.max_node_id());
465 
466  for (const auto & elem : mesh.element_ptr_range())
467  for (auto & node : elem->node_ref_range())
468  {
469  libmesh_assert_less (node.id(), nodes_to_elem_map.size());
470  libmesh_assert_less (elem->id(), mesh.n_elem());
471 
472  nodes_to_elem_map[node.id()].push_back(elem->id());
473  }
474 }
475 
476 
477 
479  std::vector<std::vector<const Elem *>> & nodes_to_elem_map)
480 {
481  // A vector indexed over all nodes is too inefficient to use for a
482  // distributed mesh. Use the unordered_map API instead.
483  if (!mesh.is_serial())
484  libmesh_deprecated();
485 
486  nodes_to_elem_map.resize (mesh.max_node_id());
487 
488  for (const auto & elem : mesh.element_ptr_range())
489  for (auto & node : elem->node_ref_range())
490  {
491  libmesh_assert_less (node.id(), nodes_to_elem_map.size());
492 
493  nodes_to_elem_map[node.id()].push_back(elem);
494  }
495 }
496 
497 
498 
500  std::unordered_map<dof_id_type, std::vector<dof_id_type>> & nodes_to_elem_map)
501 {
502  nodes_to_elem_map.clear();
503 
504  for (const auto & elem : mesh.element_ptr_range())
505  for (auto & node : elem->node_ref_range())
506  nodes_to_elem_map[node.id()].push_back(elem->id());
507 }
508 
509 
510 
512  std::unordered_map<dof_id_type, std::vector<const Elem *>> & nodes_to_elem_map)
513 {
514  nodes_to_elem_map.clear();
515 
516  for (const auto & elem : mesh.element_ptr_range())
517  for (auto & node : elem->node_ref_range())
518  nodes_to_elem_map[node.id()].push_back(elem);
519 }
520 
521 
522 
523 std::unordered_set<dof_id_type>
525 {
526  std::unordered_set<dof_id_type> boundary_nodes;
527 
528  // Loop over elements, find those on boundary, and
529  // mark them as true in on_boundary.
530  for (const auto & elem : mesh.active_element_ptr_range())
531  for (auto s : elem->side_index_range())
532  if (elem->neighbor_ptr(s) == nullptr) // on the boundary
533  {
534  auto nodes_on_side = elem->nodes_on_side(s);
535 
536  for (auto & local_id : nodes_on_side)
537  boundary_nodes.insert(elem->node_ptr(local_id)->id());
538  }
539 
540  return boundary_nodes;
541 }
542 
543 std::unordered_set<dof_id_type>
545 {
546  std::unordered_set<dof_id_type> block_boundary_nodes;
547 
548  // Loop over elements, find those on boundary, and
549  // mark them as true in on_boundary.
550  for (const auto & elem : mesh.active_element_ptr_range())
551  for (auto s : elem->side_index_range())
552  if (elem->neighbor_ptr(s) && (elem->neighbor_ptr(s)->subdomain_id() != elem->subdomain_id()))
553  {
554  auto nodes_on_side = elem->nodes_on_side(s);
555 
556  for (auto & local_id : nodes_on_side)
557  block_boundary_nodes.insert(elem->node_ptr(local_id)->id());
558  }
559 
560  return block_boundary_nodes;
561 }
562 
563 
564 
567 {
568  // This function must be run on all processors at once
569  libmesh_parallel_only(mesh.comm());
570 
571  FindBBox find_bbox;
572 
573  // Start with any unpartitioned elements we know about locally
575  mesh.pid_elements_end(DofObject::invalid_processor_id)),
576  find_bbox);
577 
578  // And combine with our local elements
579  find_bbox.bbox().union_with(create_local_bounding_box(mesh));
580 
581  // Compare the bounding boxes across processors
582  mesh.comm().min(find_bbox.min());
583  mesh.comm().max(find_bbox.max());
584 
585  return find_bbox.bbox();
586 }
587 
588 
589 
592 {
593  // This function must be run on all processors at once
594  libmesh_parallel_only(mesh.comm());
595 
596  FindBBox find_bbox;
597 
598  // Start with any unpartitioned nodes we know about locally
600  mesh.pid_nodes_end(DofObject::invalid_processor_id)),
601  find_bbox);
602 
603  // Add our local nodes
604  Threads::parallel_reduce (ConstNodeRange (mesh.local_nodes_begin(),
605  mesh.local_nodes_end()),
606  find_bbox);
607 
608  // Compare the bounding boxes across processors
609  mesh.comm().min(find_bbox.min());
610  mesh.comm().max(find_bbox.max());
611 
612  return find_bbox.bbox();
613 }
614 
615 
616 
617 Sphere
619 {
621 
622  const Real diag = (bbox.second - bbox.first).norm();
623  const Point cent = (bbox.second + bbox.first)/2;
624 
625  return Sphere (cent, .5*diag);
626 }
627 
628 
629 
632 {
633  FindBBox find_bbox;
634 
635  Threads::parallel_reduce (ConstElemRange (mesh.local_elements_begin(),
636  mesh.local_elements_end()),
637  find_bbox);
638 
639  return find_bbox.bbox();
640 }
641 
642 
643 
646  const processor_id_type pid)
647 {
648  // This can only be run in parallel, with consistent arguments.
649  libmesh_parallel_only(mesh.comm());
650  libmesh_assert(mesh.comm().verify(pid));
651 
652  libmesh_assert_less (pid, mesh.n_processors());
653 
654  FindBBox find_bbox;
655 
656  Threads::parallel_reduce (ConstElemRange (mesh.pid_elements_begin(pid),
657  mesh.pid_elements_end(pid)),
658  find_bbox);
659 
660  // Compare the bounding boxes across processors
661  mesh.comm().min(find_bbox.min());
662  mesh.comm().max(find_bbox.max());
663 
664  return find_bbox.bbox();
665 }
666 
667 
668 
669 Sphere
671  const processor_id_type pid)
672 {
673  libMesh::BoundingBox bbox =
675 
676  const Real diag = (bbox.second - bbox.first).norm();
677  const Point cent = (bbox.second + bbox.first)/2;
678 
679  return Sphere (cent, .5*diag);
680 }
681 
682 
683 
686  const subdomain_id_type sid)
687 {
688  // This can only be run in parallel, with consistent arguments.
689  libmesh_parallel_only(mesh.comm());
690  libmesh_assert(mesh.comm().verify(sid));
691 
692  FindBBox find_bbox;
693 
695  (ConstElemRange (mesh.active_local_subdomain_elements_begin(sid),
696  mesh.active_local_subdomain_elements_end(sid)),
697  find_bbox);
698 
699  // Compare the bounding boxes across processors
700  mesh.comm().min(find_bbox.min());
701  mesh.comm().max(find_bbox.max());
702 
703  return find_bbox.bbox();
704 }
705 
706 
707 
708 Sphere
710  const subdomain_id_type sid)
711 {
712  libMesh::BoundingBox bbox =
714 
715  const Real diag = (bbox.second - bbox.first).norm();
716  const Point cent = (bbox.second + bbox.first)/2;
717 
718  return Sphere (cent, .5*diag);
719 }
720 
721 
722 
723 void elem_types (const MeshBase & mesh,
724  std::vector<ElemType> & et)
725 {
726  // Loop over the the elements. If the current element type isn't in
727  // the vector, insert it.
728  for (const auto & elem : mesh.element_ptr_range())
729  if (!std::count(et.begin(), et.end(), elem->type()))
730  et.push_back(elem->type());
731 }
732 
733 
734 
736  const ElemType type)
737 {
738  return static_cast<dof_id_type>(std::distance(mesh.type_elements_begin(type),
739  mesh.type_elements_end (type)));
740 }
741 
742 
743 
745  const ElemType type)
746 {
747  return static_cast<dof_id_type>(std::distance(mesh.active_type_elements_begin(type),
748  mesh.active_type_elements_end (type)));
749 }
750 
752  const ElemType type,
753  const unsigned int level)
754 {
755  dof_id_type cnt = 0;
756 
757  // iterate over the elements of the specified type
758  for (const auto & elem : as_range(mesh.type_elements_begin(type),
759  mesh.type_elements_end(type)))
760  if ((elem->level() == level) && !elem->subactive())
761  cnt++;
762 
763  return cnt;
764 }
765 
766 
767 unsigned int n_active_local_levels(const MeshBase & mesh)
768 {
769  unsigned int nl = 0;
770 
771  for (auto & elem : mesh.active_local_element_ptr_range())
772  nl = std::max(elem->level() + 1, nl);
773 
774  return nl;
775 }
776 
777 
778 
779 unsigned int n_active_levels(const MeshBase & mesh)
780 {
781  libmesh_parallel_only(mesh.comm());
782 
783  unsigned int nl = n_active_local_levels(mesh);
784 
785  for (const auto & elem : as_range(mesh.unpartitioned_elements_begin(),
786  mesh.unpartitioned_elements_end()))
787  if (elem->active())
788  nl = std::max(elem->level() + 1, nl);
789 
790  mesh.comm().max(nl);
791  return nl;
792 }
793 
794 
795 
796 unsigned int n_local_levels(const MeshBase & mesh)
797 {
798  unsigned int nl = 0;
799 
800  for (const auto & elem : as_range(mesh.local_elements_begin(),
801  mesh.local_elements_end()))
802  nl = std::max(elem->level() + 1, nl);
803 
804  return nl;
805 }
806 
807 
808 
809 unsigned int n_levels(const MeshBase & mesh)
810 {
811  libmesh_parallel_only(mesh.comm());
812 
813  unsigned int nl = n_local_levels(mesh);
814 
815  for (const auto & elem : as_range(mesh.unpartitioned_elements_begin(),
816  mesh.unpartitioned_elements_end()))
817  nl = std::max(elem->level() + 1, nl);
818 
819  mesh.comm().max(nl);
820 
821  // n_levels() is only valid and should only be called in cases where
822  // the mesh is validly distributed (or serialized). Let's run an
823  // expensive test in debug mode to make sure this is such a case.
824 #ifdef DEBUG
825  const unsigned int paranoid_nl = paranoid_n_levels(mesh);
826  libmesh_assert_equal_to(nl, paranoid_nl);
827 #endif
828  return nl;
829 }
830 
831 
832 
833 unsigned int paranoid_n_levels(const MeshBase & mesh)
834 {
835  libmesh_parallel_only(mesh.comm());
836 
837  unsigned int nl = 0;
838  for (const auto & elem : mesh.element_ptr_range())
839  nl = std::max(elem->level() + 1, nl);
840 
841  mesh.comm().max(nl);
842  return nl;
843 }
844 
845 
846 
848  Real constraint_tol)
849 {
850  LOG_SCOPE("n_connected_components()", "MeshTools");
851 
852  // Yes, I'm being lazy. This is for mesh analysis before a
853  // simulation, not anything going in any loops.
854  if (!mesh.is_serial_on_zero())
855  libmesh_not_implemented();
856 
857  dof_id_type n_components = 0;
858 
859  if (mesh.processor_id())
860  {
861  mesh.comm().broadcast(n_components);
862  return n_components;
863  }
864 
865  // All nodes in a set here are connected (at least indirectly) to
866  // all other nodes in the same set, but have not yet been discovered
867  // to be connected to nodes in other sets.
868  std::vector<std::unordered_set<const Node *>> components;
869 
870  // With a typical mesh with few components and somewhat-contiguous
871  // ordering, vector performance should be fine. With a mesh with
872  // many components or completely scrambled ordering, performance
873  // can be a disaster.
874  auto find_component = [&components](const Node * n) {
875  std::unordered_set<const Node *> * component = nullptr;
876 
877  for (auto & c: components)
878  if (c.find(n) != c.end())
879  {
880  libmesh_assert(component == nullptr);
881  component = &c;
882  }
883 
884  return component;
885  };
886 
887  auto add_to_component =
888  [&find_component]
889  (std::unordered_set<const Node *> & component, const Node * n) {
890 
891  auto current_component = find_component(n);
892  // We may already know we're in the desired component
893  if (&component == current_component)
894  return;
895 
896  // If we're unknown, we should be in the desired component
897  else if (!current_component)
898  component.insert(n);
899 
900  // If we think we're in another component, it should actually be
901  // part of the desired component
902  else
903  {
904  component.merge(*current_component);
905  libmesh_assert(current_component->empty());
906  }
907  };
908 
909  auto & constraint_rows = mesh.get_constraint_rows();
910 
911  for (const auto & elem : mesh.element_ptr_range())
912  {
913  const Node * first_node = elem->node_ptr(0);
914 
915  auto component = find_component(first_node);
916 
917  // If we didn't find one, make a new one, reusing an existing
918  // slot if possible or growing our vector if necessary
919  if (!component)
920  for (auto & c: components)
921  if (c.empty())
922  component = &c;
923 
924  if (!component)
925  component = &components.emplace_back();
926 
927  for (const Node & n : elem->node_ref_range())
928  {
929  add_to_component(*component, &n);
930 
931  auto it = constraint_rows.find(&n);
932  if (it == constraint_rows.end())
933  continue;
934 
935  for (const auto & [pr, val] : it->second)
936  {
937  // Ignore too-trivial constraint coefficients if
938  // we get a non-default-0 constraint_tol
939  if (std::abs(val) < constraint_tol)
940  continue;
941 
942  const Elem * spline_elem = pr.first;
943  libmesh_assert(spline_elem == mesh.elem_ptr(spline_elem->id()));
944 
945  const Node * spline_node =
946  spline_elem->node_ptr(pr.second);
947 
948  add_to_component(*component, spline_node);
949  }
950  }
951  }
952 
953  for (auto & component : components)
954  if (!component.empty())
955  ++n_components;
956 
957  // We calculated this on proc 0; now let everyone else know too
958  mesh.comm().broadcast(n_components);
959 
960  return n_components;
961 }
962 
963 
964 
966  std::set<dof_id_type> & not_subactive_node_ids)
967 {
968  for (const auto & elem : mesh.element_ptr_range())
969  if (!elem->subactive())
970  for (auto & n : elem->node_ref_range())
971  not_subactive_node_ids.insert(n.id());
972 }
973 
974 
975 
978 {
979  return cast_int<dof_id_type>(std::distance(begin, end));
980 }
981 
982 
983 
985  const MeshBase::const_node_iterator & end)
986 {
987  return cast_int<dof_id_type>(std::distance(begin, end));
988 }
989 
990 
991 
993  unsigned int dim)
994 {
995  libmesh_parallel_only(mesh.comm());
996 
997  if (dim == libMesh::invalid_uint)
998  dim = mesh.mesh_dimension();
999 
1000  Real vol = 0;
1001 
1002  // first my local elements
1003  for (const auto & elem : as_range(mesh.local_elements_begin(),
1004  mesh.local_elements_end()))
1005  if (elem->dim() == dim)
1006  vol += elem->volume();
1007 
1008  // then count any unpartitioned objects, once
1009  if (mesh.processor_id() == 0)
1010  for (const auto & elem : as_range(mesh.unpartitioned_elements_begin(),
1011  mesh.unpartitioned_elements_end()))
1012  if (elem->dim() == dim)
1013  vol += elem->volume();
1014 
1015  mesh.comm().sum(vol);
1016  return vol;
1017 }
1018 
1019 
1020 
1021 unsigned int n_p_levels (const MeshBase & mesh)
1022 {
1023  libmesh_parallel_only(mesh.comm());
1024 
1025  unsigned int max_p_level = 0;
1026 
1027  // first my local elements
1028  for (const auto & elem : as_range(mesh.local_elements_begin(),
1029  mesh.local_elements_end()))
1030  max_p_level = std::max(elem->p_level(), max_p_level);
1031 
1032  // then any unpartitioned objects
1033  for (const auto & elem : as_range(mesh.unpartitioned_elements_begin(),
1034  mesh.unpartitioned_elements_end()))
1035  max_p_level = std::max(elem->p_level(), max_p_level);
1036 
1037  mesh.comm().max(max_p_level);
1038  return max_p_level + 1;
1039 }
1040 
1041 
1042 
1044  const Node & node,
1045  const std::vector<std::vector<const Elem *>> & nodes_to_elem_map,
1046  std::vector<const Node *> & neighbors)
1047 {
1048  find_nodal_neighbors_helper(node.id(), nodes_to_elem_map[node.id()],
1049  neighbors);
1050 }
1051 
1052 
1053 
1055  const Node & node,
1056  const std::unordered_map<dof_id_type, std::vector<const Elem *>> & nodes_to_elem_map,
1057  std::vector<const Node *> & neighbors)
1058 {
1059  const std::vector<const Elem *> node_to_elem_vec =
1060  libmesh_map_find(nodes_to_elem_map, node.id());
1061  find_nodal_neighbors_helper(node.id(), node_to_elem_vec, neighbors);
1062 }
1063 
1064 
1065 
1067  std::map<dof_id_type, std::vector<dof_id_type>> & hanging_nodes)
1068 {
1069  // Loop through all the elements
1070  for (auto & elem : mesh.active_local_element_ptr_range())
1071  if (elem->type() == QUAD4)
1072  for (auto s : elem->side_index_range())
1073  {
1074  // Loop over the sides looking for sides that have hanging nodes
1075  // This code is inspired by compute_proj_constraints()
1076  const Elem * neigh = elem->neighbor_ptr(s);
1077 
1078  // If not a boundary side
1079  if (neigh != nullptr)
1080  {
1081  // Is there a coarser element next to this one?
1082  if (neigh->level() < elem->level())
1083  {
1084  const Elem * ancestor = elem;
1085  while (neigh->level() < ancestor->level())
1086  ancestor = ancestor->parent();
1087  unsigned int s_neigh = neigh->which_neighbor_am_i(ancestor);
1088  libmesh_assert_less (s_neigh, neigh->n_neighbors());
1089 
1090  // Couple of helper uints...
1091  unsigned int local_node1=0;
1092  unsigned int local_node2=0;
1093 
1094  bool found_in_neighbor = false;
1095 
1096  // Find the two vertices that make up this side
1097  while (!elem->is_node_on_side(local_node1++,s)) { }
1098  local_node1--;
1099 
1100  // Start looking for the second one with the next node
1101  local_node2=local_node1+1;
1102 
1103  // Find the other one
1104  while (!elem->is_node_on_side(local_node2++,s)) { }
1105  local_node2--;
1106 
1107  //Pull out their global ids:
1108  dof_id_type node1 = elem->node_id(local_node1);
1109  dof_id_type node2 = elem->node_id(local_node2);
1110 
1111  // Now find which node is present in the neighbor
1112  // FIXME This assumes a level one rule!
1113  // The _other_ one is the hanging node
1114 
1115  // First look for the first one
1116  // FIXME could be streamlined a bit
1117  for (unsigned int n=0;n<neigh->n_sides();n++)
1118  if (neigh->node_id(n) == node1)
1119  found_in_neighbor=true;
1120 
1121  dof_id_type hanging_node=0;
1122 
1123  if (!found_in_neighbor)
1124  hanging_node=node1;
1125  else // If it wasn't node1 then it must be node2!
1126  hanging_node=node2;
1127 
1128  // Reset for reuse
1129  local_node1=0;
1130 
1131  // Find the first node that makes up the side in the neighbor (these should be the parent nodes)
1132  while (!neigh->is_node_on_side(local_node1++,s_neigh)) { }
1133  local_node1--;
1134 
1135  local_node2=local_node1+1;
1136 
1137  // Find the second node...
1138  while (!neigh->is_node_on_side(local_node2++,s_neigh)) { }
1139  local_node2--;
1140 
1141  // Save them if we haven't already found the parents for this one
1142  if (hanging_nodes[hanging_node].size()<2)
1143  {
1144  hanging_nodes[hanging_node].push_back(neigh->node_id(local_node1));
1145  hanging_nodes[hanging_node].push_back(neigh->node_id(local_node2));
1146  }
1147  }
1148  }
1149  }
1150 }
1151 
1152 
1153 
1155 {
1156  std::vector<Elem *> nodeelem_to_delete;
1157 
1158  for (auto & elem : mesh.element_ptr_range())
1159  if (elem->type() == NODEELEM &&
1160  elem->mapping_type() == RATIONAL_BERNSTEIN_MAP)
1161  nodeelem_to_delete.push_back(elem);
1162 
1163  auto & constraint_rows = mesh.get_constraint_rows();
1164 
1165  // All our constraint_rows ought to be for spline constraints we're
1166  // about to get rid of.
1167 #ifndef NDEBUG
1168  for (auto & node_row : constraint_rows)
1169  for (auto pr : node_row.second)
1170  {
1171  const Elem * elem = pr.first.first;
1172  libmesh_assert(elem->type() == NODEELEM);
1174  }
1175 #endif
1176 
1177  constraint_rows.clear();
1178 
1179  for (Elem * elem : nodeelem_to_delete)
1180  {
1181  Node * node = elem->node_ptr(0);
1182  mesh.delete_elem(elem);
1183  mesh.delete_node(node);
1184  }
1185 }
1186 
1187 
1188 
1190 {
1191  LOG_SCOPE("libmesh_assert_valid_is_prepared()", "MeshTools");
1192 
1193  if (!mesh.is_prepared())
1194  return true;
1195 
1196  std::unique_ptr<MeshBase> mesh_clone = mesh.clone();
1197 
1198  // Try preparing (without allowing repartitioning or renumbering, to
1199  // avoid false assertion failures)
1200  bool old_allow_renumbering = mesh_clone->allow_renumbering();
1201  mesh_clone->allow_renumbering(false);
1202  bool old_allow_remote_element_removal =
1203  mesh_clone->allow_remote_element_removal();
1204  bool old_skip_partitioning = mesh_clone->skip_partitioning();
1205  mesh_clone->skip_partitioning(true);
1206  mesh_clone->allow_remote_element_removal(false);
1207  mesh_clone->prepare_for_use();
1208  mesh_clone->allow_renumbering(old_allow_renumbering);
1209  mesh_clone->allow_remote_element_removal(old_allow_remote_element_removal);
1210  mesh_clone->skip_partitioning(old_skip_partitioning);
1211 
1212  return (mesh == *mesh_clone);
1213 }
1214 
1215 
1216 
1217 #ifndef NDEBUG
1218 
1220 {
1221  LOG_SCOPE("libmesh_assert_equal_n_systems()", "MeshTools");
1222 
1223  unsigned int n_sys = libMesh::invalid_uint;
1224 
1225  for (const auto & elem : mesh.element_ptr_range())
1226  {
1227  if (n_sys == libMesh::invalid_uint)
1228  n_sys = elem->n_systems();
1229  else
1230  libmesh_assert_equal_to (elem->n_systems(), n_sys);
1231  }
1232 
1233  for (const auto & node : mesh.node_ptr_range())
1234  {
1235  if (n_sys == libMesh::invalid_uint)
1236  n_sys = node->n_systems();
1237  else
1238  libmesh_assert_equal_to (node->n_systems(), n_sys);
1239  }
1240 }
1241 
1242 
1243 
1244 #ifdef LIBMESH_ENABLE_AMR
1246 {
1247  LOG_SCOPE("libmesh_assert_old_dof_objects()", "MeshTools");
1248 
1249  for (const auto & elem : mesh.element_ptr_range())
1250  {
1251  if (elem->refinement_flag() == Elem::JUST_REFINED ||
1252  elem->refinement_flag() == Elem::INACTIVE)
1253  continue;
1254 
1255  if (elem->has_dofs())
1256  libmesh_assert(elem->get_old_dof_object());
1257 
1258  for (auto & node : elem->node_ref_range())
1259  if (node.has_dofs())
1260  libmesh_assert(node.get_old_dof_object());
1261  }
1262 }
1263 #else
1264 void libmesh_assert_old_dof_objects (const MeshBase &) {}
1265 #endif // LIBMESH_ENABLE_AMR
1266 
1267 
1268 
1270 {
1271  LOG_SCOPE("libmesh_assert_valid_node_pointers()", "MeshTools");
1272 
1273  // Here we specifically do not want "auto &" because we need to
1274  // reseat the (temporary) pointer variable in the loop below,
1275  // without modifying the original.
1276  for (const Elem * elem : mesh.element_ptr_range())
1277  {
1278  libmesh_assert (elem);
1279  while (elem)
1280  {
1281  elem->libmesh_assert_valid_node_pointers();
1282  for (auto n : elem->neighbor_ptr_range())
1283  if (n && n != remote_elem)
1284  n->libmesh_assert_valid_node_pointers();
1285 
1286  libmesh_assert_not_equal_to (elem->parent(), remote_elem);
1287  elem = elem->parent();
1288  }
1289  }
1290 }
1291 
1292 
1293 
1295 {
1296  LOG_SCOPE("libmesh_assert_valid_remote_elems()", "MeshTools");
1297 
1298  for (const auto & elem : as_range(mesh.local_elements_begin(),
1299  mesh.local_elements_end()))
1300  {
1301  libmesh_assert (elem);
1302 
1303  // We currently don't allow active_local_elements to have
1304  // remote_elem neighbors
1305  if (elem->active())
1306  for (auto n : elem->neighbor_ptr_range())
1307  libmesh_assert_not_equal_to (n, remote_elem);
1308 
1309 #ifdef LIBMESH_ENABLE_AMR
1310  const Elem * parent = elem->parent();
1311  if (parent)
1312  libmesh_assert_not_equal_to (parent, remote_elem);
1313 
1314  // We can only be strict about active elements' subactive
1315  // children
1316  if (elem->active() && elem->has_children())
1317  for (auto & child : elem->child_ref_range())
1318  libmesh_assert_not_equal_to (&child, remote_elem);
1319 #endif
1320  }
1321 }
1322 
1323 
1324 
1326 {
1327  LOG_SCOPE("libmesh_assert_valid_elem_ids()", "MeshTools");
1328 
1329  processor_id_type lastprocid = 0;
1330  dof_id_type lastelemid = 0;
1331 
1332  for (const auto & elem : mesh.active_element_ptr_range())
1333  {
1334  libmesh_assert (elem);
1335  processor_id_type elemprocid = elem->processor_id();
1336  dof_id_type elemid = elem->id();
1337 
1338  libmesh_assert_greater_equal (elemid, lastelemid);
1339  libmesh_assert_greater_equal (elemprocid, lastprocid);
1340 
1341  lastelemid = elemid;
1342  lastprocid = elemprocid;
1343  }
1344 }
1345 
1346 
1347 
1349 {
1350  LOG_SCOPE("libmesh_assert_valid_amr_elem_ids()", "MeshTools");
1351 
1352  for (const auto & elem : mesh.element_ptr_range())
1353  {
1354  libmesh_assert (elem);
1355 
1356  const Elem * parent = elem->parent();
1357 
1358  if (parent)
1359  {
1360  libmesh_assert_greater_equal (elem->id(), parent->id());
1361  libmesh_assert_greater_equal (elem->processor_id(), parent->processor_id());
1362  }
1363  }
1364 }
1365 
1366 
1367 
1369 {
1370  LOG_SCOPE("libmesh_assert_valid_amr_interior_parents()", "MeshTools");
1371 
1372  for (const auto & elem : mesh.element_ptr_range())
1373  {
1374  libmesh_assert (elem);
1375 
1376  // We can skip to the next element if we're full-dimension
1377  // and therefore don't have any interior parents
1378  if (elem->dim() >= LIBMESH_DIM)
1379  continue;
1380 
1381  const Elem * ip = elem->interior_parent();
1382 
1383  const Elem * parent = elem->parent();
1384 
1385  if (ip && (ip != remote_elem) && parent)
1386  {
1387  libmesh_assert_equal_to (ip->top_parent(),
1388  elem->top_parent()->interior_parent());
1389 
1390  if (ip->level() == elem->level())
1391  libmesh_assert_equal_to (ip->parent(),
1392  parent->interior_parent());
1393  else
1394  {
1395  libmesh_assert_less (ip->level(), elem->level());
1396  libmesh_assert_equal_to (ip, parent->interior_parent());
1397  }
1398  }
1399  }
1400 }
1401 
1402 
1403 
1404 void libmesh_assert_contiguous_dof_ids(const MeshBase & mesh, unsigned int sysnum)
1405 {
1406  LOG_SCOPE("libmesh_assert_contiguous_dof_ids()", "MeshTools");
1407 
1408  if (mesh.n_processors() == 1)
1409  return;
1410 
1411  libmesh_parallel_only(mesh.comm());
1412 
1413  dof_id_type min_dof_id = std::numeric_limits<dof_id_type>::max(),
1414  max_dof_id = std::numeric_limits<dof_id_type>::min();
1415 
1416  // Figure out what our local dof id range is
1417  for (const auto * node : mesh.local_node_ptr_range())
1418  {
1419  for (auto v : make_range(node->n_vars(sysnum)))
1420  for (auto c : make_range(node->n_comp(sysnum, v)))
1421  {
1422  dof_id_type id = node->dof_number(sysnum, v, c);
1423  min_dof_id = std::min (min_dof_id, id);
1424  max_dof_id = std::max (max_dof_id, id);
1425  }
1426  }
1427 
1428  // Make sure no other processors' ids are inside it
1429  for (const auto * node : mesh.node_ptr_range())
1430  {
1431  if (node->processor_id() == mesh.processor_id())
1432  continue;
1433  for (auto v : make_range(node->n_vars(sysnum)))
1434  for (auto c : make_range(node->n_comp(sysnum, v)))
1435  {
1436  dof_id_type id = node->dof_number(sysnum, v, c);
1437  libmesh_assert (id < min_dof_id ||
1438  id > max_dof_id);
1439  }
1440  }
1441 }
1442 
1443 
1444 
1445 template <>
1447 {
1448  LOG_SCOPE("libmesh_assert_topology_consistent_procids()", "MeshTools");
1449 
1450  // This parameter is not used when !LIBMESH_ENABLE_AMR
1452 
1453  // If we're adaptively refining, check processor ids for consistency
1454  // between parents and children.
1455 #ifdef LIBMESH_ENABLE_AMR
1456 
1457  // Ancestor elements we won't worry about, but subactive and active
1458  // elements ought to have parents with consistent processor ids
1459  for (const auto & elem : mesh.element_ptr_range())
1460  {
1461  libmesh_assert(elem);
1462 
1463  if (!elem->active() && !elem->subactive())
1464  continue;
1465 
1466  const Elem * parent = elem->parent();
1467 
1468  if (parent)
1469  {
1470  libmesh_assert(parent->has_children());
1471  processor_id_type parent_procid = parent->processor_id();
1472  bool matching_child_id = false;
1473  // If we've got a remote_elem then we don't know whether
1474  // it's responsible for the parent's processor id; all
1475  // we can do is assume it is and let its processor fail
1476  // an assert if there's something wrong.
1477  for (auto & child : parent->child_ref_range())
1478  if (&child == remote_elem ||
1479  child.processor_id() == parent_procid)
1480  matching_child_id = true;
1481  libmesh_assert(matching_child_id);
1482  }
1483  }
1484 #endif
1485 }
1486 
1487 
1488 
1489 template <>
1491 {
1492  LOG_SCOPE("libmesh_assert_topology_consistent_procids()", "MeshTools");
1493 
1494  if (mesh.n_processors() == 1)
1495  return;
1496 
1497  libmesh_parallel_only(mesh.comm());
1498 
1499  // We want this test to be valid even when called after nodes have
1500  // been added asynchronously but before they're renumbered.
1501  //
1502  // Plus, some code (looking at you, stitch_meshes) modifies
1503  // DofObject ids without keeping max_elem_id()/max_node_id()
1504  // consistent, but that's done in a safe way for performance
1505  // reasons, so we'll play along and just figure out new max ids
1506  // ourselves.
1507  dof_id_type parallel_max_node_id = 0;
1508  for (const auto & node : mesh.node_ptr_range())
1509  parallel_max_node_id = std::max<dof_id_type>(parallel_max_node_id,
1510  node->id()+1);
1511  mesh.comm().max(parallel_max_node_id);
1512 
1513 
1514  std::vector<bool> node_touched_by_me(parallel_max_node_id, false);
1515 
1516  for (const auto & elem : as_range(mesh.local_elements_begin(),
1517  mesh.local_elements_end()))
1518  {
1519  libmesh_assert (elem);
1520 
1521  for (auto & node : elem->node_ref_range())
1522  {
1523  dof_id_type nodeid = node.id();
1524  node_touched_by_me[nodeid] = true;
1525  }
1526  }
1527  std::vector<bool> node_touched_by_anyone(node_touched_by_me);
1528  mesh.comm().max(node_touched_by_anyone);
1529 
1530  for (const auto & node : mesh.local_node_ptr_range())
1531  {
1532  libmesh_assert(node);
1533  dof_id_type nodeid = node->id();
1534  libmesh_assert(!node_touched_by_anyone[nodeid] ||
1535  node_touched_by_me[nodeid]);
1536  }
1537 }
1538 
1539 
1540 
1542 {
1543  for (const auto & elem : mesh.active_element_ptr_range())
1544  for (auto & node : elem->node_ref_range())
1545  libmesh_assert_equal_to
1546  (node.processor_id(),
1547  node.choose_processor_id(node.processor_id(),
1548  elem->processor_id()));
1549 }
1550 
1551 
1552 
1553 #ifdef LIBMESH_ENABLE_AMR
1555 {
1556  LOG_SCOPE("libmesh_assert_valid_refinement_tree()", "MeshTools");
1557 
1558  for (const auto & elem : mesh.element_ptr_range())
1559  {
1560  libmesh_assert(elem);
1561  if (elem->has_children())
1562  for (auto & child : elem->child_ref_range())
1563  if (&child != remote_elem)
1564  libmesh_assert_equal_to (child.parent(), elem);
1565  if (elem->active())
1566  {
1567  libmesh_assert(!elem->ancestor());
1568  libmesh_assert(!elem->subactive());
1569  }
1570  else if (elem->ancestor())
1571  {
1572  libmesh_assert(!elem->subactive());
1573  }
1574  else
1575  libmesh_assert(elem->subactive());
1576 
1577  if (elem->p_refinement_flag() == Elem::JUST_REFINED)
1578  libmesh_assert_greater(elem->p_level(), 0);
1579  }
1580 }
1581 #else
1583 {
1584 }
1585 #endif // LIBMESH_ENABLE_AMR
1586 
1587 #endif // !NDEBUG
1588 
1589 
1590 
1591 #ifdef DEBUG
1592 
1594  const Elem * bad_elem)
1595 {
1596  for (const auto & elem : mesh.element_ptr_range())
1597  {
1598  libmesh_assert (elem);
1599  libmesh_assert_not_equal_to (elem->parent(), bad_elem);
1600  for (auto n : elem->neighbor_ptr_range())
1601  libmesh_assert_not_equal_to (n, bad_elem);
1602 
1603 #ifdef LIBMESH_ENABLE_AMR
1604  if (elem->has_children())
1605  for (auto & child : elem->child_ref_range())
1606  libmesh_assert_not_equal_to (&child, bad_elem);
1607 #endif
1608  }
1609 }
1610 
1611 
1613 {
1614  LOG_SCOPE("libmesh_assert_equal_points()", "MeshTools");
1615 
1616  dof_id_type pmax_node_id = mesh.max_node_id();
1617  mesh.comm().max(pmax_node_id);
1618 
1619  for (dof_id_type i=0; i != pmax_node_id; ++i)
1620  {
1621  const Point * p = mesh.query_node_ptr(i);
1622 
1624  }
1625 }
1626 
1627 
1629 {
1630  LOG_SCOPE("libmesh_assert_equal_connectivity()", "MeshTools");
1631 
1632  dof_id_type pmax_elem_id = mesh.max_elem_id();
1633  mesh.comm().max(pmax_elem_id);
1634 
1635  for (dof_id_type i=0; i != pmax_elem_id; ++i)
1636  {
1637  const Elem * e = mesh.query_elem_ptr(i);
1638 
1639  std::vector<dof_id_type> nodes;
1640  if (e)
1641  for (auto n : e->node_index_range())
1642  nodes.push_back(e->node_id(n));
1643 
1644  libmesh_assert(mesh.comm().semiverify(e ? &nodes : nullptr));
1645  }
1646 }
1647 
1648 
1650 {
1651  LOG_SCOPE("libmesh_assert_connected_nodes()", "MeshTools");
1652 
1653  std::set<const Node *> used_nodes;
1654 
1655  for (const auto & elem : mesh.element_ptr_range())
1656  {
1657  libmesh_assert (elem);
1658 
1659  for (auto & n : elem->node_ref_range())
1660  used_nodes.insert(&n);
1661  }
1662 
1663  for (const auto & node : mesh.node_ptr_range())
1664  {
1665  libmesh_assert(node);
1666  libmesh_assert(used_nodes.count(node));
1667  }
1668 }
1669 
1670 
1671 
1673 {
1674  libmesh_parallel_only(mesh.comm());
1675 
1676  const auto & constraint_rows = mesh.get_constraint_rows();
1677 
1678  bool have_constraint_rows = !constraint_rows.empty();
1679  mesh.comm().max(have_constraint_rows);
1680  if (!have_constraint_rows)
1681  return;
1682 
1683  for (auto & row : constraint_rows)
1684  {
1685  const Node * node = row.first;
1686  libmesh_assert(node == mesh.node_ptr(node->id()));
1687 
1688  for (auto & pr : row.second)
1689  {
1690  const Elem * spline_elem = pr.first.first;
1691  libmesh_assert(spline_elem == mesh.elem_ptr(spline_elem->id()));
1692  }
1693  }
1694 
1695  dof_id_type pmax_node_id = mesh.max_node_id();
1696  mesh.comm().max(pmax_node_id);
1697 
1698  for (dof_id_type i=0; i != pmax_node_id; ++i)
1699  {
1700  const Node * node = mesh.query_node_ptr(i);
1701 
1702  bool have_constraint = constraint_rows.count(node);
1703 
1704  const std::size_t my_n_constraints = have_constraint ?
1705  libmesh_map_find(constraint_rows, node).size() : std::size_t(-1);
1706  const std::size_t * n_constraints = node ?
1707  &my_n_constraints : nullptr;
1708 
1709  libmesh_assert(mesh.comm().semiverify(n_constraints));
1710  }
1711 }
1712 
1713 
1714 
1716 {
1717  LOG_SCOPE("libmesh_assert_valid_boundary_ids()", "MeshTools");
1718 
1719  if (mesh.n_processors() == 1)
1720  return;
1721 
1722  libmesh_parallel_only(mesh.comm());
1723 
1724  const BoundaryInfo & boundary_info = mesh.get_boundary_info();
1725 
1726  dof_id_type pmax_elem_id = mesh.max_elem_id();
1727  mesh.comm().max(pmax_elem_id);
1728 
1729  for (dof_id_type i=0; i != pmax_elem_id; ++i)
1730  {
1731  const Elem * elem = mesh.query_elem_ptr(i);
1732  const unsigned int my_n_nodes = elem ? elem->n_nodes() : 0;
1733  const unsigned int my_n_edges = elem ? elem->n_edges() : 0;
1734  const unsigned int my_n_sides = elem ? elem->n_sides() : 0;
1735  unsigned int
1736  n_nodes = my_n_nodes,
1737  n_edges = my_n_edges,
1738  n_sides = my_n_sides;
1739 
1740  mesh.comm().max(n_nodes);
1741  mesh.comm().max(n_edges);
1742  mesh.comm().max(n_sides);
1743 
1744  if (elem)
1745  {
1746  libmesh_assert_equal_to(my_n_nodes, n_nodes);
1747  libmesh_assert_equal_to(my_n_edges, n_edges);
1748  libmesh_assert_equal_to(my_n_sides, n_sides);
1749  }
1750 
1751  // Let's test all IDs on the element with one communication
1752  // rather than n_nodes + n_edges + n_sides communications, to
1753  // cut down on latency in dbg modes.
1754  std::vector<boundary_id_type> all_bcids;
1755 
1756  for (unsigned int n=0; n != n_nodes; ++n)
1757  {
1758  std::vector<boundary_id_type> bcids;
1759  if (elem)
1760  {
1761  boundary_info.boundary_ids(elem->node_ptr(n), bcids);
1762 
1763  // Ordering of boundary ids shouldn't matter
1764  std::sort(bcids.begin(), bcids.end());
1765  }
1766  // libmesh_assert(mesh.comm().semiverify (elem ? &bcids : nullptr));
1767 
1768  all_bcids.insert(all_bcids.end(), bcids.begin(),
1769  bcids.end());
1770  // Separator
1771  all_bcids.push_back(BoundaryInfo::invalid_id);
1772  }
1773 
1774  for (unsigned short e=0; e != n_edges; ++e)
1775  {
1776  std::vector<boundary_id_type> bcids;
1777 
1778  if (elem)
1779  {
1780  boundary_info.edge_boundary_ids(elem, e, bcids);
1781 
1782  // Ordering of boundary ids shouldn't matter
1783  std::sort(bcids.begin(), bcids.end());
1784  }
1785 
1786  // libmesh_assert(mesh.comm().semiverify (elem ? &bcids : nullptr));
1787 
1788  all_bcids.insert(all_bcids.end(), bcids.begin(),
1789  bcids.end());
1790  // Separator
1791  all_bcids.push_back(BoundaryInfo::invalid_id);
1792 
1793  if (elem)
1794  {
1795  boundary_info.raw_edge_boundary_ids(elem, e, bcids);
1796 
1797  // Ordering of boundary ids shouldn't matter
1798  std::sort(bcids.begin(), bcids.end());
1799 
1800  all_bcids.insert(all_bcids.end(), bcids.begin(),
1801  bcids.end());
1802  // Separator
1803  all_bcids.push_back(BoundaryInfo::invalid_id);
1804  }
1805 
1806  // libmesh_assert(mesh.comm().semiverify (elem ? &bcids : nullptr));
1807  }
1808 
1809  for (unsigned short s=0; s != n_sides; ++s)
1810  {
1811  std::vector<boundary_id_type> bcids;
1812 
1813  if (elem)
1814  {
1815  boundary_info.boundary_ids(elem, s, bcids);
1816 
1817  // Ordering of boundary ids shouldn't matter
1818  std::sort(bcids.begin(), bcids.end());
1819 
1820  all_bcids.insert(all_bcids.end(), bcids.begin(),
1821  bcids.end());
1822  // Separator
1823  all_bcids.push_back(BoundaryInfo::invalid_id);
1824  }
1825 
1826  // libmesh_assert(mesh.comm().semiverify (elem ? &bcids : nullptr));
1827 
1828  if (elem)
1829  {
1830  boundary_info.raw_boundary_ids(elem, s, bcids);
1831 
1832  // Ordering of boundary ids shouldn't matter
1833  std::sort(bcids.begin(), bcids.end());
1834 
1835  all_bcids.insert(all_bcids.end(), bcids.begin(),
1836  bcids.end());
1837  // Separator
1838  all_bcids.push_back(BoundaryInfo::invalid_id);
1839  }
1840 
1841  // libmesh_assert(mesh.comm().semiverify (elem ? &bcids : nullptr));
1842  }
1843 
1844  for (unsigned short sf=0; sf != 2; ++sf)
1845  {
1846  std::vector<boundary_id_type> bcids;
1847 
1848  if (elem)
1849  {
1850  boundary_info.shellface_boundary_ids(elem, sf, bcids);
1851 
1852  // Ordering of boundary ids shouldn't matter
1853  std::sort(bcids.begin(), bcids.end());
1854 
1855  all_bcids.insert(all_bcids.end(), bcids.begin(),
1856  bcids.end());
1857  // Separator
1858  all_bcids.push_back(BoundaryInfo::invalid_id);
1859  }
1860 
1861  // libmesh_assert(mesh.comm().semiverify (elem ? &bcids : nullptr));
1862 
1863  if (elem)
1864  {
1865  boundary_info.raw_shellface_boundary_ids(elem, sf, bcids);
1866 
1867  // Ordering of boundary ids shouldn't matter
1868  std::sort(bcids.begin(), bcids.end());
1869 
1870  all_bcids.insert(all_bcids.end(), bcids.begin(),
1871  bcids.end());
1872  // Separator
1873  all_bcids.push_back(BoundaryInfo::invalid_id);
1874  }
1875 
1876  // libmesh_assert(mesh.comm().semiverify (elem ? &bcids : nullptr));
1877  }
1878 
1880  (elem ? &all_bcids : nullptr));
1881  }
1882 }
1883 
1884 
1885 void libmesh_assert_valid_dof_ids(const MeshBase & mesh, unsigned int sysnum)
1886 {
1887  LOG_SCOPE("libmesh_assert_valid_dof_ids()", "MeshTools");
1888 
1889  if (mesh.n_processors() == 1)
1890  return;
1891 
1892  libmesh_parallel_only(mesh.comm());
1893 
1894  dof_id_type pmax_elem_id = mesh.max_elem_id();
1895  mesh.comm().max(pmax_elem_id);
1896 
1897  for (dof_id_type i=0; i != pmax_elem_id; ++i)
1898  assert_semiverify_dofobj(mesh.comm(),
1899  mesh.query_elem_ptr(i),
1900  sysnum);
1901 
1902  dof_id_type pmax_node_id = mesh.max_node_id();
1903  mesh.comm().max(pmax_node_id);
1904 
1905  for (dof_id_type i=0; i != pmax_node_id; ++i)
1906  assert_semiverify_dofobj(mesh.comm(),
1907  mesh.query_node_ptr(i),
1908  sysnum);
1909 }
1910 
1911 
1912 #ifdef LIBMESH_ENABLE_UNIQUE_ID
1914 {
1915  LOG_SCOPE("libmesh_assert_valid_unique_ids()", "MeshTools");
1916 
1917  libmesh_parallel_only(mesh.comm());
1918 
1919  // First collect all the unique_ids we can see and make sure there's
1920  // no duplicates
1921  std::unordered_set<unique_id_type> semilocal_unique_ids;
1922 
1923  for (auto const & elem : mesh.active_element_ptr_range())
1924  {
1925  libmesh_assert (!semilocal_unique_ids.count(elem->unique_id()));
1926  semilocal_unique_ids.insert(elem->unique_id());
1927  }
1928 
1929  for (auto const & node : mesh.node_ptr_range())
1930  {
1931  libmesh_assert (!semilocal_unique_ids.count(node->unique_id()));
1932  semilocal_unique_ids.insert(node->unique_id());
1933  }
1934 
1935  // Then make sure elements are all in sync and remote elements don't
1936  // duplicate semilocal
1937 
1938  dof_id_type pmax_elem_id = mesh.max_elem_id();
1939  mesh.comm().max(pmax_elem_id);
1940 
1941  for (dof_id_type i=0; i != pmax_elem_id; ++i)
1942  {
1943  const Elem * elem = mesh.query_elem_ptr(i);
1944  assert_dofobj_unique_id(mesh.comm(), elem, semilocal_unique_ids);
1945  }
1946 
1947  // Then make sure nodes are all in sync and remote elements don't
1948  // duplicate semilocal
1949 
1950  dof_id_type pmax_node_id = mesh.max_node_id();
1951  mesh.comm().max(pmax_node_id);
1952 
1953  for (dof_id_type i=0; i != pmax_node_id; ++i)
1954  {
1955  const Node * node = mesh.query_node_ptr(i);
1956  assert_dofobj_unique_id(mesh.comm(), node, semilocal_unique_ids);
1957  }
1958 }
1959 #endif
1960 
1962 {
1963  libmesh_parallel_only(mesh.comm());
1964 
1965  dof_id_type parallel_max_elem_id = mesh.max_elem_id();
1966  mesh.comm().max(parallel_max_elem_id);
1967 
1968  for (dof_id_type i=0; i != parallel_max_elem_id; ++i)
1969  {
1970  const Elem * elem = mesh.query_elem_ptr(i);
1971  processor_id_type pid =
1973  mesh.comm().min(pid);
1974  libmesh_assert(elem || pid != mesh.processor_id());
1975  }
1976 
1977  dof_id_type parallel_max_node_id = mesh.max_node_id();
1978  mesh.comm().max(parallel_max_node_id);
1979 
1980  for (dof_id_type i=0; i != parallel_max_node_id; ++i)
1981  {
1982  const Node * node = mesh.query_node_ptr(i);
1983  processor_id_type pid =
1985  mesh.comm().min(pid);
1986  libmesh_assert(node || pid != mesh.processor_id());
1987  }
1988 }
1989 
1990 
1992 {
1993  libmesh_parallel_only(mesh.comm());
1994  auto locator = mesh.sub_point_locator();
1995 
1996  dof_id_type parallel_max_elem_id = mesh.max_elem_id();
1997  mesh.comm().max(parallel_max_elem_id);
1998 
1999  for (dof_id_type i=0; i != parallel_max_elem_id; ++i)
2000  {
2001  const Elem * elem = mesh.query_elem_ptr(i);
2002 
2003  const unsigned int my_n_nodes = elem ? elem->n_nodes() : 0;
2004  unsigned int n_nodes = my_n_nodes;
2005  mesh.comm().max(n_nodes);
2006 
2007  if (n_nodes)
2008  libmesh_assert(mesh.comm().semiverify(elem ? &my_n_nodes : nullptr));
2009 
2010  for (unsigned int n=0; n != n_nodes; ++n)
2011  {
2012  const Node * node = elem ? elem->node_ptr(n) : nullptr;
2013  processor_id_type pid =
2015  mesh.comm().min(pid);
2016  libmesh_assert(node || pid != mesh.processor_id());
2017  }
2018  }
2019 }
2020 
2021 
2022 
2023 template <>
2025 {
2026  LOG_SCOPE("libmesh_assert_parallel_consistent_procids()", "MeshTools");
2027 
2028  if (mesh.n_processors() == 1)
2029  return;
2030 
2031  libmesh_parallel_only(mesh.comm());
2032 
2033  // Some code (looking at you, stitch_meshes) modifies DofObject ids
2034  // without keeping max_elem_id()/max_node_id() consistent, but
2035  // that's done in a safe way for performance reasons, so we'll play
2036  // along and just figure out new max ids ourselves.
2037  dof_id_type parallel_max_elem_id = 0;
2038  for (const auto & elem : mesh.element_ptr_range())
2039  parallel_max_elem_id = std::max<dof_id_type>(parallel_max_elem_id,
2040  elem->id()+1);
2041  mesh.comm().max(parallel_max_elem_id);
2042 
2043  // Check processor ids for consistency between processors
2044 
2045  for (dof_id_type i=0; i != parallel_max_elem_id; ++i)
2046  {
2047  const Elem * elem = mesh.query_elem_ptr(i);
2048 
2049  processor_id_type min_id =
2050  elem ? elem->processor_id() :
2051  std::numeric_limits<processor_id_type>::max();
2052  mesh.comm().min(min_id);
2053 
2054  processor_id_type max_id =
2055  elem ? elem->processor_id() :
2056  std::numeric_limits<processor_id_type>::min();
2057  mesh.comm().max(max_id);
2058 
2059  if (elem)
2060  {
2061  libmesh_assert_equal_to (min_id, elem->processor_id());
2062  libmesh_assert_equal_to (max_id, elem->processor_id());
2063  }
2064 
2065  if (min_id == mesh.processor_id())
2066  libmesh_assert(elem);
2067  }
2068 }
2069 
2070 
2071 
2073 {
2074  LOG_SCOPE("libmesh_assert_parallel_consistent_new_node_procids()", "MeshTools");
2075 
2076  if (mesh.n_processors() == 1)
2077  return;
2078 
2079  libmesh_parallel_only(mesh.comm());
2080 
2081  // We want this test to hit every node when called even after nodes
2082  // have been added asynchronously but before everything has been
2083  // renumbered.
2084  dof_id_type parallel_max_elem_id = mesh.max_elem_id();
2085  mesh.comm().max(parallel_max_elem_id);
2086 
2087  std::vector<bool> elem_touched_by_anyone(parallel_max_elem_id, false);
2088 
2089  for (dof_id_type i=0; i != parallel_max_elem_id; ++i)
2090  {
2091  const Elem * elem = mesh.query_elem_ptr(i);
2092 
2093  const unsigned int my_n_nodes = elem ? elem->n_nodes() : 0;
2094  unsigned int n_nodes = my_n_nodes;
2095  mesh.comm().max(n_nodes);
2096 
2097  if (n_nodes)
2098  libmesh_assert(mesh.comm().semiverify(elem ? &my_n_nodes : nullptr));
2099 
2100  for (unsigned int n=0; n != n_nodes; ++n)
2101  {
2102  const Node * node = elem ? elem->node_ptr(n) : nullptr;
2103  const processor_id_type pid = node ? node->processor_id() : 0;
2104  libmesh_assert(mesh.comm().semiverify (node ? &pid : nullptr));
2105  }
2106  }
2107 }
2108 
2109 template <>
2111 {
2112  LOG_SCOPE("libmesh_assert_parallel_consistent_procids()", "MeshTools");
2113 
2114  if (mesh.n_processors() == 1)
2115  return;
2116 
2117  libmesh_parallel_only(mesh.comm());
2118 
2119  // We want this test to be valid even when called even after nodes
2120  // have been added asynchronously but before they're renumbered
2121  //
2122  // Plus, some code (looking at you, stitch_meshes) modifies
2123  // DofObject ids without keeping max_elem_id()/max_node_id()
2124  // consistent, but that's done in a safe way for performance
2125  // reasons, so we'll play along and just figure out new max ids
2126  // ourselves.
2127  dof_id_type parallel_max_node_id = 0;
2128  for (const auto & node : mesh.node_ptr_range())
2129  parallel_max_node_id = std::max<dof_id_type>(parallel_max_node_id,
2130  node->id()+1);
2131  mesh.comm().max(parallel_max_node_id);
2132 
2133  std::vector<bool> node_touched_by_anyone(parallel_max_node_id, false);
2134 
2135  for (const auto & elem : as_range(mesh.local_elements_begin(),
2136  mesh.local_elements_end()))
2137  {
2138  libmesh_assert (elem);
2139 
2140  for (auto & node : elem->node_ref_range())
2141  {
2142  dof_id_type nodeid = node.id();
2143  node_touched_by_anyone[nodeid] = true;
2144  }
2145  }
2146  mesh.comm().max(node_touched_by_anyone);
2147 
2148  // Check processor ids for consistency between processors
2149  // on any node an element touches
2150  for (dof_id_type i=0; i != parallel_max_node_id; ++i)
2151  {
2152  if (!node_touched_by_anyone[i])
2153  continue;
2154 
2155  const Node * node = mesh.query_node_ptr(i);
2156  const processor_id_type pid = node ? node->processor_id() : 0;
2157 
2158  libmesh_assert(mesh.comm().semiverify (node ? &pid : nullptr));
2159  }
2160 }
2161 
2162 
2163 
2164 #ifdef LIBMESH_ENABLE_AMR
2166 {
2167  LOG_SCOPE("libmesh_assert_valid_refinement_flags()", "MeshTools");
2168 
2169  libmesh_parallel_only(mesh.comm());
2170  if (mesh.n_processors() == 1)
2171  return;
2172 
2173  dof_id_type pmax_elem_id = mesh.max_elem_id();
2174  mesh.comm().max(pmax_elem_id);
2175 
2176  std::vector<unsigned char> my_elem_h_state(pmax_elem_id, 255);
2177  std::vector<unsigned char> my_elem_p_state(pmax_elem_id, 255);
2178 
2179  for (const auto & elem : mesh.element_ptr_range())
2180  {
2181  libmesh_assert (elem);
2182  dof_id_type elemid = elem->id();
2183 
2184  my_elem_h_state[elemid] =
2185  static_cast<unsigned char>(elem->refinement_flag());
2186 
2187  my_elem_p_state[elemid] =
2188  static_cast<unsigned char>(elem->p_refinement_flag());
2189  }
2190  std::vector<unsigned char> min_elem_h_state(my_elem_h_state);
2191  mesh.comm().min(min_elem_h_state);
2192 
2193  std::vector<unsigned char> min_elem_p_state(my_elem_p_state);
2194  mesh.comm().min(min_elem_p_state);
2195 
2196  for (dof_id_type i=0; i!= pmax_elem_id; ++i)
2197  {
2198  libmesh_assert(my_elem_h_state[i] == 255 ||
2199  my_elem_h_state[i] == min_elem_h_state[i]);
2200  libmesh_assert(my_elem_p_state[i] == 255 ||
2201  my_elem_p_state[i] == min_elem_p_state[i]);
2202  }
2203 }
2204 #else
2206 {
2207 }
2208 #endif // LIBMESH_ENABLE_AMR
2209 
2210 
2211 
2213  bool assert_valid_remote_elems)
2214 {
2215  LOG_SCOPE("libmesh_assert_valid_neighbors()", "MeshTools");
2216 
2217  for (const auto & elem : mesh.element_ptr_range())
2218  {
2219  libmesh_assert (elem);
2220  elem->libmesh_assert_valid_neighbors();
2221  }
2222 
2223  if (mesh.n_processors() == 1)
2224  return;
2225 
2226  libmesh_parallel_only(mesh.comm());
2227 
2228  dof_id_type pmax_elem_id = mesh.max_elem_id();
2229  mesh.comm().max(pmax_elem_id);
2230 
2231  for (dof_id_type i=0; i != pmax_elem_id; ++i)
2232  {
2233  const Elem * elem = mesh.query_elem_ptr(i);
2234 
2235  const unsigned int my_n_neigh = elem ? elem->n_neighbors() : 0;
2236  unsigned int n_neigh = my_n_neigh;
2237  mesh.comm().max(n_neigh);
2238  if (elem)
2239  libmesh_assert_equal_to (my_n_neigh, n_neigh);
2240 
2241  for (unsigned int n = 0; n != n_neigh; ++n)
2242  {
2243  dof_id_type my_neighbor = DofObject::invalid_id;
2244  dof_id_type * p_my_neighbor = nullptr;
2245 
2246  // If we have a non-remote_elem neighbor link, then we can
2247  // verify it.
2248  if (elem && elem->neighbor_ptr(n) != remote_elem)
2249  {
2250  p_my_neighbor = &my_neighbor;
2251  if (elem->neighbor_ptr(n))
2252  my_neighbor = elem->neighbor_ptr(n)->id();
2253 
2254  // But wait - if we haven't set remote_elem links yet then
2255  // some nullptr links on ghost elements might be
2256  // future-remote_elem links, so we can't verify those.
2257  if (!assert_valid_remote_elems &&
2258  !elem->neighbor_ptr(n) &&
2259  elem->processor_id() != mesh.processor_id())
2260  p_my_neighbor = nullptr;
2261  }
2262  libmesh_assert(mesh.comm().semiverify(p_my_neighbor));
2263  }
2264  }
2265 }
2266 #endif // DEBUG
2267 
2268 
2269 
2270 // Functors for correct_node_proc_ids
2271 namespace {
2272 
2273 typedef std::unordered_map<dof_id_type, processor_id_type> proc_id_map_type;
2274 
2275 struct SyncNodeSet
2276 {
2277  typedef unsigned char datum; // bool but without bit twiddling issues
2278 
2279  SyncNodeSet(std::unordered_set<const Node *> & _set,
2280  MeshBase & _mesh) :
2281  node_set(_set), mesh(_mesh) {}
2282 
2283  std::unordered_set<const Node *> & node_set;
2284 
2286 
2287  // ------------------------------------------------------------
2288  void gather_data (const std::vector<dof_id_type> & ids,
2289  std::vector<datum> & data)
2290  {
2291  // Find whether each requested node belongs in the set
2292  data.resize(ids.size());
2293 
2294  for (auto i : index_range(ids))
2295  {
2296  const dof_id_type id = ids[i];
2297 
2298  // We'd better have every node we're asked for
2299  Node * node = mesh.node_ptr(id);
2300 
2301  // Return if the node is in the set.
2302  data[i] = node_set.count(node);
2303  }
2304  }
2305 
2306  // ------------------------------------------------------------
2307  bool act_on_data (const std::vector<dof_id_type> & ids,
2308  const std::vector<datum> in_set)
2309  {
2310  bool data_changed = false;
2311 
2312  // Add nodes we've been informed of to our own set
2313  for (auto i : index_range(ids))
2314  {
2315  if (in_set[i])
2316  {
2317  Node * node = mesh.node_ptr(ids[i]);
2318  if (!node_set.count(node))
2319  {
2320  node_set.insert(node);
2321  data_changed = true;
2322  }
2323  }
2324  }
2325 
2326  return data_changed;
2327  }
2328 };
2329 
2330 
2331 struct NodesNotInSet
2332 {
2333  NodesNotInSet(const std::unordered_set<const Node *> _set)
2334  : node_set(_set) {}
2335 
2336  bool operator() (const Node * node) const
2337  {
2338  if (node_set.count(node))
2339  return false;
2340  return true;
2341  }
2342 
2343  const std::unordered_set<const Node *> node_set;
2344 };
2345 
2346 
2347 struct SyncProcIdsFromMap
2348 {
2349  typedef processor_id_type datum;
2350 
2351  SyncProcIdsFromMap(const proc_id_map_type & _map,
2352  MeshBase & _mesh) :
2353  new_proc_ids(_map), mesh(_mesh) {}
2354 
2355  const proc_id_map_type & new_proc_ids;
2356 
2357  MeshBase & mesh;
2358 
2359  // ------------------------------------------------------------
2360  void gather_data (const std::vector<dof_id_type> & ids,
2361  std::vector<datum> & data)
2362  {
2363  // Find the new processor id of each requested node
2364  data.resize(ids.size());
2365 
2366  for (auto i : index_range(ids))
2367  {
2368  const dof_id_type id = ids[i];
2369 
2370  // Return the node's new processor id if it has one, or its
2371  // old processor id if not.
2372  if (const auto it = new_proc_ids.find(id);
2373  it != new_proc_ids.end())
2374  data[i] = it->second;
2375  else
2376  {
2377  // We'd better find every node we're asked for
2378  const Node & node = mesh.node_ref(id);
2379  data[i] = node.processor_id();
2380  }
2381  }
2382  }
2383 
2384  // ------------------------------------------------------------
2385  void act_on_data (const std::vector<dof_id_type> & ids,
2386  const std::vector<datum> proc_ids)
2387  {
2388  // Set the node processor ids we've now been informed of
2389  for (auto i : index_range(ids))
2390  {
2391  Node & node = mesh.node_ref(ids[i]);
2392  node.processor_id() = proc_ids[i];
2393  }
2394  }
2395 };
2396 }
2397 
2398 
2399 
2401 {
2402  LOG_SCOPE("correct_node_proc_ids()","MeshTools");
2403 
2404  // This function must be run on all processors at once
2405  libmesh_parallel_only(mesh.comm());
2406 
2407  // We require all processors to agree on nodal processor ids before
2408  // going into this algorithm.
2409 #ifdef DEBUG
2411 #endif
2412 
2413  // If we have any unpartitioned elements at this
2414  // stage there is a problem
2415  libmesh_assert (n_elem(mesh.unpartitioned_elements_begin(),
2416  mesh.unpartitioned_elements_end()) == 0);
2417 
2418  // Fix nodes' processor ids. Coarsening may have left us with nodes
2419  // which are no longer touched by any elements of the same processor
2420  // id, and for DofMap to work we need to fix that.
2421 
2422  // This is harder now that libMesh no longer requires a distributed
2423  // mesh to ghost all nodal neighbors: it is possible for two active
2424  // elements on two different processors to share the same node in
2425  // such a way that neither processor knows the others' element
2426  // exists!
2427 
2428  // While we're at it, if this mesh is configured to allow
2429  // repartitioning, we'll repartition *all* the nodes' processor ids
2430  // using the canonical Node heuristic, to try and improve DoF load
2431  // balancing. But if the mesh is disallowing repartitioning, we
2432  // won't touch processor_id on any node where it's valid, regardless
2433  // of whether or not it's canonical.
2434  bool repartition_all_nodes = !mesh.skip_noncritical_partitioning();
2435  std::unordered_set<const Node *> valid_nodes;
2436 
2437  // If we aren't allowed to repartition, then we're going to leave
2438  // every node we can at its current processor_id, and *only*
2439  // repartition the nodes whose current processor id is incompatible
2440  // with DoFMap (because it doesn't touch an active element, e.g. due
2441  // to coarsening)
2442  if (!repartition_all_nodes)
2443  {
2444  for (const auto & elem : mesh.active_element_ptr_range())
2445  for (const auto & node : elem->node_ref_range())
2446  if (elem->processor_id() == node.processor_id())
2447  valid_nodes.insert(&node);
2448 
2449  SyncNodeSet syncv(valid_nodes, mesh);
2450 
2452  (mesh.comm(), mesh.nodes_begin(), mesh.nodes_end(), syncv);
2453  }
2454 
2455  // We build up a set of compatible processor ids for each node
2456  proc_id_map_type new_proc_ids;
2457 
2458  for (auto & elem : mesh.active_element_ptr_range())
2459  {
2460  processor_id_type pid = elem->processor_id();
2461 
2462  for (auto & node : elem->node_ref_range())
2463  {
2464  const dof_id_type id = node.id();
2465  if (auto it = new_proc_ids.find(id);
2466  it == new_proc_ids.end())
2467  new_proc_ids.emplace(id, pid);
2468  else
2469  it->second = node.choose_processor_id(it->second, pid);
2470  }
2471  }
2472 
2473  // Sort the new pids to push to each processor
2474  std::map<processor_id_type, std::vector<std::pair<dof_id_type, processor_id_type>>>
2475  ids_to_push;
2476 
2477  for (const auto & node : mesh.node_ptr_range())
2478  if (const auto it = std::as_const(new_proc_ids).find(node->id());
2479  it != new_proc_ids.end() && node->processor_id() != DofObject::invalid_processor_id)
2480  ids_to_push[node->processor_id()].emplace_back(node->id(), /*pid=*/it->second);
2481 
2482  auto action_functor =
2483  [& mesh, & new_proc_ids]
2485  const std::vector<std::pair<dof_id_type, processor_id_type>> & data)
2486  {
2487  for (const auto & [id, pid] : data)
2488  {
2489  if (const auto it = new_proc_ids.find(id);
2490  it == new_proc_ids.end())
2491  new_proc_ids.emplace(id, pid);
2492  else
2493  {
2494  const Node & node = mesh.node_ref(id);
2495  it->second = node.choose_processor_id(it->second, pid);
2496  }
2497  }
2498  };
2499 
2500  Parallel::push_parallel_vector_data
2501  (mesh.comm(), ids_to_push, action_functor);
2502 
2503  // Now new_proc_ids is correct for every node we used to own. Let's
2504  // ask every other processor about the nodes they used to own. But
2505  // first we'll need to keep track of which nodes we used to own,
2506  // lest we get them confused with nodes we newly own.
2507  std::unordered_set<Node *> ex_local_nodes;
2508  for (auto & node : mesh.local_node_ptr_range())
2509  if (const auto it = new_proc_ids.find(node->id());
2510  it != new_proc_ids.end() && it->second != mesh.processor_id())
2511  ex_local_nodes.insert(node);
2512 
2513  SyncProcIdsFromMap sync(new_proc_ids, mesh);
2514  if (repartition_all_nodes)
2516  (mesh.comm(), mesh.nodes_begin(), mesh.nodes_end(), sync);
2517  else
2518  {
2519  NodesNotInSet nnis(valid_nodes);
2520 
2522  (mesh.comm(), mesh.nodes_begin(), mesh.nodes_end(), nnis, sync);
2523  }
2524 
2525  // And finally let's update the nodes we used to own.
2526  for (const auto & node : ex_local_nodes)
2527  {
2528  if (valid_nodes.count(node))
2529  continue;
2530 
2531  const dof_id_type id = node->id();
2532  const proc_id_map_type::iterator it = new_proc_ids.find(id);
2533  libmesh_assert(it != new_proc_ids.end());
2534  node->processor_id() = it->second;
2535  }
2536 
2537  // We should still have consistent nodal processor ids coming out of
2538  // this algorithm, but if we're allowed to repartition the mesh then
2539  // they should be canonically correct too.
2540 #ifdef DEBUG
2541  libmesh_assert_valid_procids<Node>(mesh);
2542  //if (repartition_all_nodes)
2543  // libmesh_assert_canonical_node_procids(mesh);
2544 #endif
2545 }
2546 
2547 
2548 
2550 {
2552 }
2553 
2554 } // namespace MeshTools
2555 
2556 } // namespace libMesh
unsigned int n_active_local_levels(const MeshBase &mesh)
Definition: mesh_tools.C:767
static const Order type_to_default_order_map[INVALID_ELEM]
This array maps the integer representation of the ElemType enum to the default approximation order of...
Definition: elem.h:992
void libmesh_assert_equal_points(const MeshBase &mesh)
A function for testing that node locations match across processors.
Definition: mesh_tools.C:1612
ElemType
Defines an enum for geometric element types.
void libmesh_assert_topology_consistent_procids< Node >(const MeshBase &mesh)
Definition: mesh_tools.C:1490
void elem_types(const MeshBase &mesh, std::vector< ElemType > &et)
Fills in a vector of all element types in the mesh.
Definition: mesh_tools.C:723
dof_id_type n_nodes(const MeshBase::const_node_iterator &begin, const MeshBase::const_node_iterator &end)
Count up the number of nodes of a specific type (as defined by an iterator range).
Definition: mesh_tools.C:984
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:1032
const Elem * parent() const
Definition: elem.h:3031
bool is_prepared() const
Definition: mesh_base.h:198
static constexpr processor_id_type invalid_processor_id
An invalid processor_id to distinguish DoFs that have not been assigned to a processor.
Definition: dof_object.h:493
constraint_rows_type & get_constraint_rows()
Constraint rows accessors.
Definition: mesh_base.h:1709
A Node is like a Point, but with more information.
Definition: node.h:52
void libmesh_assert_valid_amr_interior_parents(const MeshBase &mesh)
A function for verifying that any interior_parent pointers on elements are consistent with AMR (paren...
Definition: mesh_tools.C:1368
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
unsigned int n_comp(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:1002
std::unique_ptr< PointLocatorBase > sub_point_locator() const
Definition: mesh_base.C:1675
This class defines a sphere.
Definition: sphere.h:72
void libmesh_assert_valid_remote_elems(const MeshBase &mesh)
A function for verifying that active local elements&#39; neighbors are never remote elements.
Definition: mesh_tools.C:1294
void libmesh_assert_consistent_distributed_nodes(const MeshBase &mesh)
A function for verifying that distribution of nodes is parallel consistent (every processor can see e...
Definition: mesh_tools.C:1991
const Elem * interior_parent() const
Definition: elem.C:1184
void correct_node_proc_ids(MeshBase &)
Changes the processor ids on each node so be the same as the id of the lowest element touching that n...
Definition: mesh_tools.C:2400
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Count up the number of elements of a specific type (as defined by an iterator range).
Definition: mesh_tools.C:976
libMesh::BoundingBox create_bounding_box(const MeshBase &mesh)
Definition: mesh_tools.C:566
The definition of the const_element_iterator struct.
Definition: mesh_base.h:2218
bool valid_is_prepared(const MeshBase &mesh)
A function for testing whether a mesh&#39;s cached is_prepared() setting is not a false positive...
Definition: mesh_tools.C:1189
unsigned int dim
const Elem * top_parent() const
Definition: elem.h:3057
Sphere processor_bounding_sphere(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:670
Dummy "splitting object" used to distinguish splitting constructors from copy constructors.
Definition: threads_none.h:63
void clear_spline_nodes(MeshBase &)
Remove spline node (for IsoGeometric Analysis meshes) elements and nodes and constraints from the mes...
Definition: mesh_tools.C:1154
void libmesh_assert_valid_refinement_tree(const MeshBase &mesh)
A function for verifying that elements on this processor have valid descendants and consistent active...
Definition: mesh_tools.C:1554
void sum(T &r) const
std::unordered_set< dof_id_type > find_block_boundary_nodes(const MeshBase &mesh)
Returns a std::set containing Node IDs for all of the block boundary nodes.
Definition: mesh_tools.C:544
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
void assign_global_indices(MeshBase &) const
This method assigns globally unique, partition-agnostic indices to the nodes and elements in the mesh...
void find_nodal_neighbors(const MeshBase &mesh, const Node &n, const std::vector< std::vector< const Elem *>> &nodes_to_elem_map, std::vector< const Node *> &neighbors)
Given a mesh and a node in the mesh, the vector will be filled with every node directly attached to t...
Definition: mesh_tools.C:1043
void libmesh_assert_valid_constraint_rows(const MeshBase &mesh)
A function for verifying that all mesh constraint rows express relations between nodes and elements t...
Definition: mesh_tools.C:1672
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const =0
unsigned int n_local_levels(const MeshBase &mesh)
Definition: mesh_tools.C:796
void skip_noncritical_partitioning(bool skip)
If true is passed in then the elements on this mesh will no longer be (re)partitioned, and the nodes on this mesh will only be repartitioned if they are found "orphaned" via coarsening or other removal of the last element responsible for their node/element processor id consistency.
Definition: mesh_base.h:1236
unique_id_type unique_id() const
Definition: dof_object.h:844
void libmesh_assert_parallel_consistent_new_node_procids(const MeshBase &mesh)
A function for verifying that processor assignment is parallel consistent (every processor agrees on ...
Definition: mesh_tools.C:2072
const Parallel::Communicator & comm() const
void build_nodes_to_elem_map(const MeshBase &mesh, std::vector< std::vector< dof_id_type >> &nodes_to_elem_map)
After calling this function the input vector nodes_to_elem_map will contain the node to element conne...
Definition: mesh_tools.C:456
The StoredRange class defines a contiguous, divisible set of objects.
Definition: stored_range.h:54
The libMesh namespace provides an interface to certain functionality in the library.
MeshBase & mesh
Definition: mesh_tools.C:2285
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:165
virtual std::unique_ptr< MeshBase > clone() const =0
Virtual "copy constructor".
Real distance(const Point &p)
dof_id_type n_elem_of_type(const MeshBase &mesh, const ElemType type)
Definition: mesh_tools.C:735
processor_id_type choose_processor_id(processor_id_type pid1, processor_id_type pid2) const
Return which of pid1 and pid2 would be preferred by the current load-balancing heuristic applied to t...
Definition: node.C:78
libMesh::BoundingBox create_local_bounding_box(const MeshBase &mesh)
Definition: mesh_tools.C:631
uint8_t processor_id_type
Definition: id_types.h:104
This is the MeshBase class.
Definition: mesh_base.h:75
SimpleRange< ChildRefIter > child_ref_range()
Returns a range with all children of a parent element, usable in range-based for loops.
Definition: elem.h:2347
void libmesh_assert_consistent_distributed(const MeshBase &mesh)
A function for verifying that distribution of dof objects is parallel consistent (every processor can...
Definition: mesh_tools.C:1961
unsigned int paranoid_n_levels(const MeshBase &mesh)
Definition: mesh_tools.C:833
std::unordered_set< dof_id_type > find_boundary_nodes(const MeshBase &mesh)
Calling this function on a 2D mesh will convert all the elements to triangles.
Definition: mesh_tools.C:524
virtual bool is_serial_on_zero() const
Definition: mesh_base.h:218
void libmesh_assert_equal_n_systems(const MeshBase &mesh)
The following functions, only available in builds with NDEBUG undefined, are for asserting internal c...
Definition: mesh_tools.C:1219
void libmesh_assert_valid_dof_ids(const MeshBase &mesh, unsigned int sysnum=libMesh::invalid_uint)
A function for verifying that degree of freedom indexing matches across processors.
Definition: mesh_tools.C:1885
uint8_t processor_id_type
processor_id_type n_processors() const
virtual bool is_serial() const
Definition: mesh_base.h:211
void libmesh_ignore(const Args &...)
const dof_id_type n_nodes
Definition: tecplot_io.C:67
void libmesh_assert_valid_node_pointers(const MeshBase &mesh)
A function for walking across the mesh to try and ferret out invalidated or misassigned pointers...
Definition: mesh_tools.C:1269
static const boundary_id_type invalid_id
Number used for internal use.
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
ElemMappingType mapping_type() const
Definition: elem.h:3121
This is the MeshCommunication class.
void libmesh_assert_canonical_node_procids(const MeshBase &mesh)
A function for verifying that processor assignment of nodes matches the heuristic specified in Node::...
Definition: mesh_tools.C:1541
dof_id_type id() const
Definition: dof_object.h:828
void min(const T &r, T &o, Request &req) const
static constexpr dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:482
unsigned int n_vars
virtual unsigned int n_nodes() const =0
void libmesh_assert_no_links_to_elem(const MeshBase &mesh, const Elem *bad_elem)
The following functions, only available in builds with DEBUG defined, typically have surprisingly slo...
Definition: mesh_tools.C:1593
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:444
unsigned int n_vars(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:967
unsigned int which_neighbor_am_i(const Elem *e) const
This function tells you which neighbor e is.
Definition: elem.h:2920
unsigned int n_systems() const
Definition: dof_object.h:937
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
unsigned int n_levels(const MeshBase &mesh)
Definition: mesh_tools.C:809
virtual const Node * query_node_ptr(const dof_id_type i) const =0
virtual dof_id_type max_elem_id() const =0
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
void libmesh_assert_valid_boundary_ids(const MeshBase &mesh)
A function for verifying that boundary condition ids match across processors.
Definition: mesh_tools.C:1715
virtual void delete_node(Node *n)=0
Removes the Node n from the mesh.
void libmesh_assert_valid_amr_elem_ids(const MeshBase &mesh)
A function for verifying that ids of elements are correctly sorted for AMR (parents have lower ids th...
Definition: mesh_tools.C:1348
void libmesh_assert_parallel_consistent_procids< Elem >(const MeshBase &mesh)
Definition: mesh_tools.C:2024
Sphere subdomain_bounding_sphere(const MeshBase &mesh, const subdomain_id_type sid)
Definition: mesh_tools.C:709
const proc_id_map_type & new_proc_ids
Definition: mesh_tools.C:2355
virtual unsigned int n_edges() const =0
timpi_pure bool semiverify(const T *r) const
auto norm(const T &a) -> decltype(std::abs(a))
Definition: tensor_tools.h:74
void sync_dofobject_data_by_id(const Communicator &comm, const Iterator &range_begin, const Iterator &range_end, SyncFunctor &sync)
Request data about a range of ghost dofobjects uniquely identified by their id.
void globally_renumber_nodes_and_elements(MeshBase &)
There is no reason for a user to ever call this function.
Definition: mesh_tools.C:2549
Real volume(const MeshBase &mesh, unsigned int dim=libMesh::invalid_uint)
Find the total volume of a mesh (interpreting that as area for dim = 2, or total arc length for dim =...
Definition: mesh_tools.C:992
void find_hanging_nodes_and_parents(const MeshBase &mesh, std::map< dof_id_type, std::vector< dof_id_type >> &hanging_nodes)
Given a mesh hanging_nodes will be filled with an associative array keyed off the global id of all th...
Definition: mesh_tools.C:1066
dof_id_type total_weight(const MeshBase &mesh)
Definition: mesh_tools.C:421
void broadcast(T &data, const unsigned int root_id=0, const bool identical_sizes=false) const
void libmesh_assert_parallel_consistent_procids< Node >(const MeshBase &mesh)
Definition: mesh_tools.C:2110
dof_id_type n_connected_components(const MeshBase &mesh, Real constraint_tol=0)
Definition: mesh_tools.C:847
libMesh::BoundingBox create_subdomain_bounding_box(const MeshBase &mesh, const subdomain_id_type sid)
Definition: mesh_tools.C:685
std::string enum_to_string(const T e)
Sphere bounding_sphere(const MeshBase &mesh)
Definition: mesh_tools.C:618
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
virtual const Elem * elem_ptr(const dof_id_type i) const =0
virtual unsigned int n_sides() const =0
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2599
void libmesh_assert_equal_connectivity(const MeshBase &mesh)
A function for testing that element nodal connectivities match across processors. ...
Definition: mesh_tools.C:1628
std::unordered_set< const Node * > & node_set
Definition: mesh_tools.C:2283
unsigned int level() const
Definition: elem.h:3075
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual const Elem * query_elem_ptr(const dof_id_type i) const =0
timpi_pure bool verify(const T &r) const
void libmesh_assert_connected_nodes(const MeshBase &mesh)
A function for verifying that all nodes are connected to at least one element.
Definition: mesh_tools.C:1649
void max(const T &r, T &o, Request &req) const
DIE A HORRIBLE DEATH HERE typedef MPI_Comm communicator
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2508
unsigned int n_neighbors() const
Definition: elem.h:715
void get_not_subactive_node_ids(const MeshBase &mesh, std::set< dof_id_type > &not_subactive_node_ids)
Builds a set of node IDs for nodes which belong to non-subactive elements.
Definition: mesh_tools.C:965
The definition of the const_node_iterator struct.
Definition: mesh_base.h:2269
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
The DofObject defines an abstract base class for objects that have degrees of freedom associated with...
Definition: dof_object.h:54
dof_id_type n_active_elem_of_type(const MeshBase &mesh, const ElemType type)
Definition: mesh_tools.C:744
unsigned int mesh_dimension() const
Definition: mesh_base.C:372
void parallel_reduce(const Range &range, Body &body)
Execute the provided reduction operation in parallel on the specified range.
Definition: threads_none.h:101
IntRange< unsigned short > node_index_range() const
Definition: elem.h:2684
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:596
void libmesh_assert_valid_neighbors(const MeshBase &mesh, bool assert_valid_remote_elems=true)
A function for verifying that neighbor connectivity is correct (each element is a neighbor of or desc...
Definition: mesh_tools.C:2212
libMesh::BoundingBox create_nodal_bounding_box(const MeshBase &mesh)
Definition: mesh_tools.C:591
virtual dof_id_type max_node_id() const =0
virtual dof_id_type n_elem() const =0
virtual const Node * node_ptr(const dof_id_type i) const =0
void libmesh_assert_valid_unique_ids(const MeshBase &mesh)
A function for verifying that unique ids match across processors.
Definition: mesh_tools.C:1913
processor_id_type processor_id() const
libMesh::BoundingBox create_processor_bounding_box(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:645
unsigned int n_p_levels(const MeshBase &mesh)
Definition: mesh_tools.C:1021
processor_id_type processor_id() const
Definition: dof_object.h:905
void libmesh_assert_old_dof_objects(const MeshBase &mesh)
A function for testing that all non-recently-created DofObjects within a mesh have old_dof_object dat...
Definition: mesh_tools.C:1245
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:2476
uint8_t unique_id_type
Definition: id_types.h:86
bool has_children() const
Definition: elem.h:2980
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
void libmesh_assert_contiguous_dof_ids(const MeshBase &mesh, unsigned int sysnum)
A function for verifying that degree of freedom indexes are contiguous on each processor, as is required by libMesh numeric classes.
Definition: mesh_tools.C:1404
dof_id_type n_non_subactive_elem_of_type_at_level(const MeshBase &mesh, const ElemType type, const unsigned int level)
Definition: mesh_tools.C:751
void libmesh_assert_topology_consistent_procids< Elem >(const MeshBase &mesh)
Definition: mesh_tools.C:1446
uint8_t dof_id_type
Definition: id_types.h:67
unsigned int n_active_levels(const MeshBase &mesh)
Definition: mesh_tools.C:779
void libmesh_assert_valid_elem_ids(const MeshBase &mesh)
A function for verifying that ids and processor assignment of elements are correctly sorted (monotone...
Definition: mesh_tools.C:1325
void libmesh_assert_valid_refinement_flags(const MeshBase &mesh)
A function for verifying that refinement flags on elements are consistent between processors...
Definition: mesh_tools.C:2165
const RemoteElem * remote_elem
Definition: remote_elem.C:57