Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
unstructured_mesh.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/boundary_info.h"
22 #include "libmesh/ghosting_functor.h"
23 #include "libmesh/ghost_point_neighbors.h"
24 #include "libmesh/unstructured_mesh.h"
25 #include "libmesh/libmesh_logging.h"
26 #include "libmesh/elem.h"
27 #include "libmesh/mesh_tools.h" // For n_levels
28 #include "libmesh/parallel.h"
29 #include "libmesh/remote_elem.h"
30 #include "libmesh/namebased_io.h"
31 #include "libmesh/partitioner.h"
32 #include "libmesh/enum_order.h"
33 #include "libmesh/mesh_communication.h"
34 #include "libmesh/enum_to_string.h"
35 #include "libmesh/mesh_serializer.h"
36 #include "libmesh/utility.h"
37 
38 #ifdef LIBMESH_HAVE_NANOFLANN
39 #include "libmesh/nanoflann.hpp"
40 #endif
41 
42 // C++ includes
43 #include <algorithm> // std::all_of
44 #include <fstream>
45 #include <iomanip>
46 #include <map>
47 #include <sstream>
48 #include <unordered_map>
49 
50 namespace {
51 
52 using namespace libMesh;
53 
54 // Helper functions for all_second_order, all_complete_order
55 
56 std::map<std::vector<dof_id_type>, Node *>::iterator
57 map_hi_order_node(unsigned int hon,
58  const Elem & hi_elem,
59  std::map<std::vector<dof_id_type>, Node *> & adj_vertices_to_ho_nodes)
60 {
61  /*
62  * form a vector that will hold the node id's of
63  * the vertices that are adjacent to the nth
64  * higher-order node.
65  */
66  const unsigned int n_adjacent_vertices =
68 
69  std::vector<dof_id_type> adjacent_vertices_ids(n_adjacent_vertices);
70 
71  for (unsigned int v=0; v<n_adjacent_vertices; v++)
72  adjacent_vertices_ids[v] =
73  hi_elem.node_id( hi_elem.second_order_adjacent_vertex(hon,v) );
74 
75  /*
76  * \p adjacent_vertices_ids is now in order of the current
77  * side. sort it, so that comparisons with the
78  * \p adjacent_vertices_ids created through other elements'
79  * sides can match
80  */
81  std::sort(adjacent_vertices_ids.begin(),
82  adjacent_vertices_ids.end());
83 
84  // Does this set of vertices already have a mid-node added? If not
85  // we'll want to add it.
86  return adj_vertices_to_ho_nodes.try_emplace(adjacent_vertices_ids, nullptr).first;
87 }
88 
89 void transfer_elem(Elem & lo_elem,
90  std::unique_ptr<Elem> hi_elem,
91 #ifdef LIBMESH_ENABLE_UNIQUE_ID
92  unique_id_type max_unique_id,
93  unique_id_type max_new_nodes_per_elem,
94 #endif
96  std::map<std::vector<dof_id_type>, Node *> & adj_vertices_to_ho_nodes,
97  std::unordered_map<Elem *, std::vector<Elem *>> & exterior_children_of)
98 {
99  libmesh_assert_equal_to (lo_elem.n_vertices(), hi_elem->n_vertices());
100 
101  const processor_id_type my_pid = mesh.processor_id();
102  const processor_id_type lo_pid = lo_elem.processor_id();
103 
104  /*
105  * Now handle the additional higher-order nodes. This
106  * is simply handled through a map that remembers
107  * the already-added nodes. This map maps the global
108  * ids of the vertices (that uniquely define this
109  * higher-order node) to the new node.
110  * Notation: hon = high-order node
111  */
112  const unsigned int hon_begin = lo_elem.n_nodes();
113  const unsigned int hon_end = hi_elem->n_nodes();
114 
115  for (unsigned int hon=hon_begin; hon<hon_end; hon++)
116  {
117  auto pos = map_hi_order_node(hon, *hi_elem, adj_vertices_to_ho_nodes);
118 
119  // no, not added yet
120  if (!pos->second)
121  {
122  const auto & adjacent_vertices_ids = pos->first;
123 
124  /*
125  * for this set of vertices, there is no
126  * second_order node yet. Add it.
127  *
128  * compute the location of the new node as
129  * the average over the adjacent vertices.
130  */
131  Point new_location = 0;
132  for (dof_id_type vertex_id : adjacent_vertices_ids)
133  new_location += mesh.point(vertex_id);
134 
135  new_location /= static_cast<Real>(adjacent_vertices_ids.size());
136 
137  /* Add the new point to the mesh.
138  *
139  * If we are on a serialized mesh, then we're doing this
140  * all in sync, and the node processor_id will be
141  * consistent between processors.
142  *
143  * If we are on a distributed mesh, we can fix
144  * inconsistent processor ids later, but only if every
145  * processor gives new nodes a *locally* consistent
146  * processor id, so we'll give the new node the
147  * processor id of an adjacent element for now and then
148  * we'll update that later if appropriate.
149  */
150  Node * hi_node = mesh.add_point
151  (new_location, DofObject::invalid_id, lo_pid);
152 
153  /* Come up with a unique unique_id for a potentially new
154  * node. On a distributed mesh we don't yet know what
155  * processor_id will definitely own it, so we can't let
156  * the pid determine the unique_id. But we're not
157  * adding unpartitioned nodes in sync, so we can't let
158  * the mesh autodetermine a unique_id for a new
159  * unpartitioned node either. So we have to pick unique
160  * unique_id values manually.
161  *
162  * We don't have to pick the *same* unique_id value as
163  * will be picked on other processors, though; we'll
164  * sync up each node later. We just need to make sure
165  * we don't duplicate any unique_id that might be chosen
166  * by the same process elsewhere.
167  */
168 #ifdef LIBMESH_ENABLE_UNIQUE_ID
169  unique_id_type new_unique_id = max_unique_id +
170  max_new_nodes_per_elem * lo_elem.id() +
171  hon - hon_begin;
172 
173  hi_node->set_unique_id(new_unique_id);
174 #endif
175 
176  /*
177  * insert the new node with its defining vertex
178  * set into the map, and relocate pos to this
179  * new entry, so that the hi_elem can use
180  * \p pos for inserting the node
181  */
182  pos->second = hi_node;
183 
184  hi_elem->set_node(hon) = hi_node;
185  }
186  // yes, already added.
187  else
188  {
189  Node * hi_node = pos->second;
190  libmesh_assert(hi_node);
191  libmesh_assert_equal_to(mesh.node_ptr(hi_node->id()), hi_node);
192 
193  hi_elem->set_node(hon) = hi_node;
194 
195  // We need to ensure that the processor who should own a
196  // node *knows* they own the node. And because
197  // Node::choose_processor_id() may depend on Node id,
198  // which may not yet be authoritative, we still have to
199  // use a dumb-but-id-independent partitioning heuristic.
200  processor_id_type chosen_pid =
201  std::min (hi_node->processor_id(), lo_pid);
202 
203  // Plus, if we just discovered that we own this node,
204  // then on a distributed mesh we need to make sure to
205  // give it a valid id, not just a placeholder id!
206  if (!mesh.is_replicated() &&
207  hi_node->processor_id() != my_pid &&
208  chosen_pid == my_pid)
209  mesh.own_node(*hi_node);
210 
211  hi_node->processor_id() = chosen_pid;
212  }
213  }
214 
215  /*
216  * find_neighbors relies on remote_elem neighbor links being
217  * properly maintained. Our own code here relies on ordinary
218  * neighbor links being properly maintained, so let's just keep
219  * everything up to date.
220  */
221  for (auto s : lo_elem.side_index_range())
222  {
223  Elem * neigh = lo_elem.neighbor_ptr(s);
224  if (!neigh)
225  continue;
226 
227  if (neigh != remote_elem)
228  {
229  // We don't support AMR even outside our own range yet.
230  libmesh_assert_equal_to (neigh->level(), 0);
231 
232  const unsigned int ns = neigh->which_neighbor_am_i(&lo_elem);
233  libmesh_assert_not_equal_to(ns, libMesh::invalid_uint);
234 
235  neigh->set_neighbor(ns, hi_elem.get());
236  }
237 
238  hi_elem->set_neighbor(s, neigh);
239  }
240 
247  Elem * interior_p = lo_elem.interior_parent();
248  if (interior_p)
249  hi_elem->set_interior_parent(interior_p);
250 
251  if (auto parent_exterior_it = exterior_children_of.find(interior_p);
252  parent_exterior_it != exterior_children_of.end())
253  {
254  auto & exteriors = parent_exterior_it->second;
255  for (std::size_t i : index_range(exteriors))
256  if (exteriors[i] == &lo_elem)
257  {
258  exteriors[i] = hi_elem.get();
259  break;
260  }
261  }
262 
267  if (auto exterior_it = exterior_children_of.find(&lo_elem);
268  exterior_it != exterior_children_of.end())
269  {
270  for (Elem * exterior_elem : exterior_it->second)
271  {
272  libmesh_assert(exterior_elem->interior_parent() == &lo_elem);
273  exterior_elem->set_interior_parent(hi_elem.get());
274  }
275  }
276 
289  (mesh.get_boundary_info(), &lo_elem, hi_elem.get());
290 
291  /*
292  * The new second-order element is ready.
293  * Inserting it into the mesh will replace and delete
294  * the first-order element.
295  */
296  hi_elem->set_id(lo_elem.id());
297 #ifdef LIBMESH_ENABLE_UNIQUE_ID
298  hi_elem->set_unique_id(lo_elem.unique_id());
299 #endif
300  hi_elem->processor_id() = lo_pid;
301  hi_elem->subdomain_id() = lo_elem.subdomain_id();
302 
303  const unsigned int nei = lo_elem.n_extra_integers();
304  hi_elem->add_extra_integers(nei);
305  for (unsigned int i=0; i != nei; ++i)
306  hi_elem->set_extra_integer(i, lo_elem.get_extra_integer(i));
307 
308  // This might not help anything but shouldn't hurt.
309  hi_elem->set_mapping_type(lo_elem.mapping_type());
310  hi_elem->set_mapping_data(lo_elem.mapping_data());
311 
312  mesh.insert_elem(std::move(hi_elem));
313 }
314 
315 
316 template <typename ElemTypeConverter>
317 void
318 all_increased_order_range (UnstructuredMesh & mesh,
320  const unsigned int max_new_nodes_per_elem,
321  const ElemTypeConverter & elem_type_converter)
322 {
323  // This function must be run on all processors at once
324  timpi_parallel_only(mesh.comm());
325 
326  /*
327  * The maximum number of new higher-order nodes we might be adding,
328  * for use when picking unique unique_id values later. This variable
329  * is not used unless unique ids are enabled, so libmesh_ignore() it
330  * to avoid warnings in that case.
331  */
332  libmesh_ignore(max_new_nodes_per_elem);
333 
334  /*
335  * The mesh should at least be consistent enough for us to add new
336  * nodes consistently.
337  */
340 
341  /*
342  * If the mesh is empty then we have nothing to do
343  */
344  if (!mesh.n_elem())
345  return;
346 
347  // If every element in the range _on every proc_ is already of the
348  // requested higher order then we have nothing to do. However, if
349  // any proc has some lower-order elements in the range, then _all_
350  // processors need to continue this function because it is
351  // parallel_only().
352  //
353  // Note: std::all_of() returns true for an empty range, which can
354  // happen for example in the DistributedMesh case when there are
355  // more processors than elements. In the case of an empty range we
356  // therefore set already_second_order to true on that proc.
357  auto is_higher_order = [&elem_type_converter](const Elem * elem) {
358  ElemType old_type = elem->type();
359  ElemType new_type = elem_type_converter(old_type);
360  return old_type == new_type;
361  };
362 
363  bool already_higher_order =
364  std::all_of(range.begin(), range.end(), is_higher_order);
365 
366  // Check with other processors and possibly return early
367  mesh.comm().min(already_higher_order);
368  if (already_higher_order)
369  return;
370 
371  /*
372  * this map helps in identifying higher order
373  * nodes. Namely, a higher-order node:
374  * - edge node
375  * - face node
376  * - bubble node
377  * is uniquely defined through a set of adjacent
378  * vertices. This set of adjacent vertices is
379  * used to identify already added higher-order
380  * nodes. We are safe to use node id's since we
381  * make sure that these are correctly numbered.
382  *
383  * We lazily use an ordered map here to avoid having to implement a
384  * good hash for vector<dof_id_type>
385  */
386  std::map<std::vector<dof_id_type>, Node *> adj_vertices_to_ho_nodes;
387 
388  /*
389  * This map helps us reset any interior_parent() values from the
390  * lower order element to its higher order replacement. Unlike with
391  * neighbor pointers, we don't have backlinks here, so we have to
392  * iterate over the mesh to track forward links.
393  */
394  std::unordered_map<Elem *, std::vector<Elem *>> exterior_children_of;
395 
396  /*
397  * max_new_nodes_per_elem is the maximum number of new higher order
398  * nodes we might be adding, for use when picking unique unique_id
399  * values later. This variable is not used unless unique ids are
400  * enabled.
401  */
402 #ifdef LIBMESH_ENABLE_UNIQUE_ID
403  unique_id_type max_unique_id = mesh.parallel_max_unique_id();
404 #endif
405 
412  dof_id_type n_unpartitioned_elem = 0,
413  n_partitioned_elem = 0;
414 
429  auto track_if_necessary = [&adj_vertices_to_ho_nodes,
430  &exterior_children_of,
431  &elem_type_converter](Elem * elem) {
432  if (elem && elem != remote_elem)
433  {
434  if (elem->default_order() != FIRST)
435  for (unsigned int hon : make_range(elem->n_vertices(), elem->n_nodes()))
436  {
437  auto pos = map_hi_order_node(hon, *elem, adj_vertices_to_ho_nodes);
438  pos->second = elem->node_ptr(hon);
439  }
440 
441  const ElemType old_type = elem->type();
442  const ElemType new_type = elem_type_converter(old_type);
443  if (old_type != new_type)
444  exterior_children_of.emplace(elem, std::vector<Elem *>());
445  }
446  };
447 
448  // If we're in the common case then just track everything; otherwise
449  // find point neighbors to track
450  if (range.begin() == mesh.elements_begin() &&
451  range.end() == mesh.elements_end())
452  {
453  for (auto & elem : range)
454  track_if_necessary(elem);
455  }
456  else
457  {
458  GhostingFunctor::map_type point_neighbor_elements;
459 
460  GhostPointNeighbors point_neighbor_finder(mesh);
461  point_neighbor_finder(range.begin(), range.end(),
462  mesh.n_processors(),
463  point_neighbor_elements);
464 
465  for (auto & [elem, coupling_map] : point_neighbor_elements)
466  {
467  libmesh_ignore(coupling_map);
468  track_if_necessary(const_cast<Elem *>(elem));
469  }
470  }
471 
476  for (auto & elem : mesh.element_ptr_range())
477  if (auto exterior_map_it = exterior_children_of.find(elem->interior_parent());
478  exterior_map_it != exterior_children_of.end())
479  exterior_map_it->second.push_back(elem);
480 
487  for (auto & lo_elem : range)
488  {
489  // Now we can skip the elements in the range that are already
490  // higher-order.
491  const ElemType old_type = lo_elem->type();
492  const ElemType new_type = elem_type_converter(old_type);
493 
494  if (old_type == new_type)
495  continue;
496 
497  // this does _not_ work for refined elements
498  libmesh_assert_equal_to (lo_elem->level(), 0);
499 
501  ++n_unpartitioned_elem;
502  else
503  ++n_partitioned_elem;
504 
505  /*
506  * Build the higher-order equivalent; add to
507  * the new_elements list.
508  */
509  auto ho_elem = Elem::build (new_type);
510 
511  libmesh_assert_equal_to (lo_elem->n_vertices(), ho_elem->n_vertices());
512 
513  /*
514  * By definition the initial nodes of the lower and higher order
515  * element are identically numbered. Transfer these.
516  */
517  for (unsigned int v=0, lnn=lo_elem->n_nodes(); v < lnn; v++)
518  ho_elem->set_node(v) = lo_elem->node_ptr(v);
519 
520  transfer_elem(*lo_elem, std::move(ho_elem),
521 #ifdef LIBMESH_ENABLE_UNIQUE_ID
522  max_unique_id, max_new_nodes_per_elem,
523 #endif
524  mesh, adj_vertices_to_ho_nodes,
525  exterior_children_of);
526  } // end for (auto & lo_elem : range)
527 
528  // we can clear the map at this point.
529  adj_vertices_to_ho_nodes.clear();
530 
531 #ifdef LIBMESH_ENABLE_UNIQUE_ID
532  const unique_id_type new_max_unique_id = max_unique_id +
533  max_new_nodes_per_elem * mesh.n_elem();
534  mesh.set_next_unique_id(new_max_unique_id);
535 #endif
536 
537  // On a DistributedMesh our ghost node processor ids may be bad,
538  // the ids of nodes touching remote elements may be inconsistent,
539  // unique_ids of newly added non-local nodes remain unset, and our
540  // partitioning of new nodes may not be well balanced.
541  //
542  // make_nodes_parallel_consistent() will fix all this.
543  if (!mesh.is_replicated())
544  {
545  dof_id_type max_unpartitioned_elem = n_unpartitioned_elem;
546  mesh.comm().max(max_unpartitioned_elem);
547  if (max_unpartitioned_elem)
548  {
549  // We'd better be effectively serialized here. In theory we
550  // could support more complicated cases but in practice we
551  // only support "completely partitioned" and/or "serialized"
552  if (!mesh.comm().verify(n_unpartitioned_elem) ||
553  !mesh.comm().verify(n_partitioned_elem) ||
554  !mesh.is_serial())
555  libmesh_not_implemented();
556  }
557  else
558  {
560  }
561  }
562 
563  // renumber nodes, repartition nodes, etc. We may no longer need a
564  // find_neighbors() here since we're keeping neighbor links intact
565  // ourselves, *except* that if we're not already prepared we may
566  // have user code that was expecting this call to prepare neighbors.
567  const bool old_find_neighbors = mesh.allow_find_neighbors();
568  if (mesh.is_prepared())
569  mesh.allow_find_neighbors(false);
571  mesh.allow_find_neighbors(old_find_neighbors);
572 }
573 
574 
575 } // anonymous namespace
576 
577 
578 namespace libMesh
579 {
580 
581 // This class adapts a vector of Nodes (represented by a pair of a Point and a dof_id_type)
582 // for use in a nanoflann KD-Tree
583 
585 {
586 private:
587  const std::vector<std::pair<Point, dof_id_type>> _nodes;
588 
589 public:
590  VectorOfNodesAdaptor(const std::vector<std::pair<Point, dof_id_type>> & nodes) :
591  _nodes(nodes)
592  {}
593 
597  inline size_t kdtree_get_point_count() const { return _nodes.size(); }
598 
604  inline Real kdtree_get_pt(const size_t idx, int dim) const
605  {
606  libmesh_assert_less (idx, _nodes.size());
607  libmesh_assert_less (dim, 3);
608 
609  const Point & p(_nodes[idx].first);
610 
611  if (dim==0) return p(0);
612  if (dim==1) return p(1);
613  return p(2);
614  }
615 
616  /*
617  * Optional bounding-box computation
618  */
619  template <class BBOX>
620  bool kdtree_get_bbox(BBOX & /* bb */) const { return false; }
621 };
622 
623 
624 // ------------------------------------------------------------
625 // UnstructuredMesh class member functions
627  unsigned char d) :
628  MeshBase (comm_in,d)
629 {
631 }
632 
633 
634 
636  MeshBase (other_mesh)
637 {
639 }
640 
641 
642 
644  const bool skip_find_neighbors,
645  dof_id_type element_id_offset,
646  dof_id_type node_id_offset,
648 #ifdef LIBMESH_ENABLE_UNIQUE_ID
649  unique_id_offset
650 #endif
651  ,
652  std::unordered_map<subdomain_id_type, subdomain_id_type> *
653  id_remapping)
654 {
655  LOG_SCOPE("copy_nodes_and_elements()", "UnstructuredMesh");
656 
657  std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
658  extra_int_maps = this->merge_extra_integer_names(other_mesh);
659 
660  const unsigned int n_old_node_ints = extra_int_maps.second.size(),
661  n_new_node_ints = _node_integer_names.size(),
662  n_old_elem_ints = extra_int_maps.first.size(),
663  n_new_elem_ints = _elem_integer_names.size();
664 
665  // If we are partitioned into fewer parts than the incoming mesh has
666  // processors to handle, then we need to "wrap" the other Mesh's
667  // processor ids to fit within our range. This can happen, for
668  // example, while stitching meshes with small numbers of elements in
669  // parallel...
670  bool wrap_proc_ids = (this->n_processors() <
671  other_mesh.n_partitions());
672 
673  // We're assuming the other mesh has proper element number ordering,
674  // so that we add parents before their children, and that the other
675  // mesh is consistently partitioned.
676 #ifdef DEBUG
678  MeshTools::libmesh_assert_valid_procids<Node>(other_mesh);
679 #endif
680 
681  //Copy in Nodes
682  {
683  //Preallocate Memory if necessary
684  this->reserve_nodes(other_mesh.n_nodes());
685 
686  for (const auto & oldn : other_mesh.node_ptr_range())
687  {
688  processor_id_type added_pid = cast_int<processor_id_type>
689  (wrap_proc_ids ? oldn->processor_id() % this->n_processors() : oldn->processor_id());
690 
691  // Add new nodes in old node Point locations
692  Node * newn =
693  this->add_point(*oldn,
694  oldn->id() + node_id_offset,
695  added_pid);
696 
697  newn->add_extra_integers(n_new_node_ints);
698  for (unsigned int i = 0; i != n_old_node_ints; ++i)
699  newn->set_extra_integer(extra_int_maps.second[i],
700  oldn->get_extra_integer(i));
701 
702 #ifdef LIBMESH_ENABLE_UNIQUE_ID
703  newn->set_unique_id(oldn->unique_id() + unique_id_offset);
704 #endif
705  }
706  }
707 
708  //Copy in Elements
709  {
710  //Preallocate Memory if necessary
711  this->reserve_elem(other_mesh.n_elem());
712 
713  // Declare a map linking old and new elements, needed to copy the neighbor lists
714  typedef std::unordered_map<const Elem *, Elem *> map_type;
715  map_type old_elems_to_new_elems, ip_map;
716 
717  // Loop over the elements
718  for (const auto & old : other_mesh.element_ptr_range())
719  {
720  // Build a new element
721  Elem * newparent = old->parent() ?
722  this->elem_ptr(old->parent()->id() + element_id_offset) :
723  nullptr;
724  auto el = old->disconnected_clone();
725  el->set_parent(newparent);
726 
727  subdomain_id_type sbd_id = old->subdomain_id();
728  if (id_remapping)
729  {
730  auto remapping_it = id_remapping->find(sbd_id);
731  if (remapping_it != id_remapping->end())
732  sbd_id = remapping_it->second;
733  }
734  el->subdomain_id() = sbd_id;
735 
736  // Hold off on trying to set the interior parent because we may actually
737  // add lower dimensional elements before their interior parents
738  if (el->dim() < other_mesh.mesh_dimension() && old->interior_parent())
739  ip_map[old] = el.get();
740 
741 #ifdef LIBMESH_ENABLE_AMR
742  if (old->has_children())
743  for (unsigned int c = 0, nc = old->n_children(); c != nc; ++c)
744  if (old->child_ptr(c) == remote_elem)
745  el->add_child(const_cast<RemoteElem *>(remote_elem), c);
746 
747  //Create the parent's child pointers if necessary
748  if (newparent)
749  {
750  unsigned int oldc = old->parent()->which_child_am_i(old);
751  newparent->add_child(el.get(), oldc);
752  }
753 
754  // Copy the refinement flags
755  el->set_refinement_flag(old->refinement_flag());
756 
757  // Use hack_p_level since we may not have sibling elements
758  // added yet
759  el->hack_p_level(old->p_level());
760 
761  el->set_p_refinement_flag(old->p_refinement_flag());
762 #endif // #ifdef LIBMESH_ENABLE_AMR
763 
764  //Assign all the nodes
765  for (auto i : el->node_index_range())
766  el->set_node(i) =
767  this->node_ptr(old->node_id(i) + node_id_offset);
768 
769  // And start it off with the same processor id (mod _n_parts).
770  el->processor_id() = cast_int<processor_id_type>
771  (wrap_proc_ids ? old->processor_id() % this->n_processors() : old->processor_id());
772 
773  // Give it the same element and unique ids
774  el->set_id(old->id() + element_id_offset);
775 
776  el->add_extra_integers(n_new_elem_ints);
777  for (unsigned int i = 0; i != n_old_elem_ints; ++i)
778  el->set_extra_integer(extra_int_maps.first[i],
779  old->get_extra_integer(i));
780 
781 #ifdef LIBMESH_ENABLE_UNIQUE_ID
782  el->set_unique_id(old->unique_id() + unique_id_offset);
783 #endif
784 
785  //Hold onto it
786  if (!skip_find_neighbors)
787  {
788  for (auto s : old->side_index_range())
789  if (old->neighbor_ptr(s) == remote_elem)
790  el->set_neighbor(s, const_cast<RemoteElem *>(remote_elem));
791  this->add_elem(std::move(el));
792  }
793  else
794  {
795  Elem * new_el = this->add_elem(std::move(el));
796  old_elems_to_new_elems[old] = new_el;
797  }
798  }
799 
800  for (auto & elem_pair : ip_map)
801  elem_pair.second->set_interior_parent(
802  this->elem_ptr(elem_pair.first->interior_parent()->id() + element_id_offset));
803 
804  // Loop (again) over the elements to fill in the neighbors
805  if (skip_find_neighbors)
806  {
807  old_elems_to_new_elems[remote_elem] = const_cast<RemoteElem*>(remote_elem);
808 
809  for (const auto & old_elem : other_mesh.element_ptr_range())
810  {
811  Elem * new_elem = old_elems_to_new_elems[old_elem];
812  for (auto s : old_elem->side_index_range())
813  {
814  const Elem * old_neighbor = old_elem->neighbor_ptr(s);
815  Elem * new_neighbor = old_elems_to_new_elems[old_neighbor];
816  new_elem->set_neighbor(s, new_neighbor);
817  }
818  }
819  }
820  }
821 
822 #ifdef LIBMESH_ENABLE_UNIQUE_ID
823  // We set the unique ids of nodes after adding them to the mesh such that our value of
824  // _next_unique_id may be wrong. So we amend that here
825  this->set_next_unique_id(other_mesh.parallel_max_unique_id() + unique_id_offset + 1);
826 #endif
827 
828  //Finally prepare the new Mesh for use. Keep the same numbering and
829  //partitioning for now.
830  this->allow_renumbering(false);
831  this->allow_remote_element_removal(false);
832  this->allow_find_neighbors(!skip_find_neighbors);
833 
834  // We should generally be able to skip *all* partitioning here
835  // because we're only adding one already-consistent mesh to another.
836  this->skip_partitioning(true);
837 
838  this->prepare_for_use();
839 
840  //But in the long term, use the same renumbering and partitioning
841  //policies as our source mesh.
842  this->allow_find_neighbors(other_mesh.allow_find_neighbors());
843  this->allow_renumbering(other_mesh.allow_renumbering());
845  this->skip_partitioning(other_mesh.skip_partitioning());
846 }
847 
848 
849 
851 {
852  // this->clear (); // Nothing to clear at this level
853 
854  libmesh_exceptionless_assert (!libMesh::closed());
855 }
856 
857 
858 
859 
860 
861 void UnstructuredMesh::find_neighbors (const bool reset_remote_elements,
862  const bool reset_current_list)
863 {
864  // We might actually want to run this on an empty mesh
865  // (e.g. the boundary mesh for a nonexistent bcid!)
866  // libmesh_assert_not_equal_to (this->n_nodes(), 0);
867  // libmesh_assert_not_equal_to (this->n_elem(), 0);
868 
869  // This function must be run on all processors at once
870  parallel_object_only();
871 
872  LOG_SCOPE("find_neighbors()", "Mesh");
873 
874  //TODO:[BSK] This should be removed later?!
875  if (reset_current_list)
876  for (const auto & e : this->element_ptr_range())
877  for (auto s : e->side_index_range())
878  if (e->neighbor_ptr(s) != remote_elem || reset_remote_elements)
879  e->set_neighbor(s, nullptr);
880 
881  // Find neighboring elements by first finding elements
882  // with identical side keys and then check to see if they
883  // are neighbors
884  {
885  // data structures -- Use the hash_multimap if available
886  typedef dof_id_type key_type;
887  typedef std::pair<Elem *, unsigned char> val_type;
888  typedef std::unordered_multimap<key_type, val_type> map_type;
889 
890  // A map from side keys to corresponding elements & side numbers
891  map_type side_to_elem_map;
892 
893  // Pull objects out of the loop to reduce heap operations
894  std::unique_ptr<Elem> my_side, their_side;
895 
896  for (const auto & element : this->element_ptr_range())
897  {
898  for (auto ms : element->side_index_range())
899  {
900  next_side:
901  // If we haven't yet found a neighbor on this side, try.
902  // Even if we think our neighbor is remote, that
903  // information may be out of date.
904  if (element->neighbor_ptr(ms) == nullptr ||
905  element->neighbor_ptr(ms) == remote_elem)
906  {
907  // Get the key for the side of this element. Use the
908  // low_order_key so we can find neighbors in
909  // mixed-order meshes if necessary.
910  const dof_id_type key = element->low_order_key(ms);
911 
912  // Look for elements that have an identical side key
913  auto bounds = side_to_elem_map.equal_range(key);
914 
915  // May be multiple keys, check all the possible
916  // elements which _might_ be neighbors.
917  if (bounds.first != bounds.second)
918  {
919  // Get the side for this element
920  element->side_ptr(my_side, ms);
921 
922  // Look at all the entries with an equivalent key
923  while (bounds.first != bounds.second)
924  {
925  // Get the potential element
926  Elem * neighbor = bounds.first->second.first;
927 
928  // Get the side for the neighboring element
929  const unsigned int ns = bounds.first->second.second;
930  neighbor->side_ptr(their_side, ns);
931  //libmesh_assert(my_side.get());
932  //libmesh_assert(their_side.get());
933 
934  // If found a match with my side
935  //
936  // In 1D, since parents and children have an
937  // equal side (i.e. a node) we need to check
938  // for matching level() to avoid setting our
939  // neighbor pointer to any of our neighbor's
940  // descendants.
941  if ((*my_side == *their_side) &&
942  (element->level() == neighbor->level()))
943  {
944  // So share a side. Is this a mixed pair
945  // of subactive and active/ancestor
946  // elements?
947  // If not, then we're neighbors.
948  // If so, then the subactive's neighbor is
949 
950  if (element->subactive() ==
951  neighbor->subactive())
952  {
953  // an element is only subactive if it has
954  // been coarsened but not deleted
955  element->set_neighbor (ms,neighbor);
956  neighbor->set_neighbor(ns,element);
957  }
958  else if (element->subactive())
959  {
960  element->set_neighbor(ms,neighbor);
961  }
962  else if (neighbor->subactive())
963  {
964  neighbor->set_neighbor(ns,element);
965  }
966  side_to_elem_map.erase (bounds.first);
967 
968  // get out of this nested crap
969  goto next_side;
970  }
971 
972  ++bounds.first;
973  }
974  }
975 
976  // didn't find a match...
977  // Build the map entry for this element
978  side_to_elem_map.emplace
979  (key, std::make_pair(element, cast_int<unsigned char>(ms)));
980  }
981  }
982  }
983  }
984 
985 #ifdef LIBMESH_ENABLE_AMR
986 
1014  const unsigned int n_levels = MeshTools::n_levels(*this);
1015  for (unsigned int level = 1; level < n_levels; ++level)
1016  {
1017  for (auto & current_elem : as_range(level_elements_begin(level),
1018  level_elements_end(level)))
1019  {
1020  libmesh_assert(current_elem);
1021  Elem * parent = current_elem->parent();
1022  libmesh_assert(parent);
1023  const unsigned int my_child_num = parent->which_child_am_i(current_elem);
1024 
1025  for (auto s : current_elem->side_index_range())
1026  {
1027  if (current_elem->neighbor_ptr(s) == nullptr ||
1028  (current_elem->neighbor_ptr(s) == remote_elem &&
1029  parent->is_child_on_side(my_child_num, s)))
1030  {
1031  Elem * neigh = parent->neighbor_ptr(s);
1032 
1033  // If neigh was refined and had non-subactive children
1034  // made remote earlier, then our current elem should
1035  // actually have one of those remote children as a
1036  // neighbor
1037  if (neigh &&
1038  (neigh->ancestor() ||
1039  // If neigh has subactive children which should have
1040  // matched as neighbors of the current element but
1041  // did not, then those likewise must be remote
1042  // children.
1043  (current_elem->subactive() && neigh->has_children() &&
1044  (neigh->level()+1) == current_elem->level())))
1045  {
1046 #ifdef DEBUG
1047  // Let's make sure that "had children made remote"
1048  // situation is actually the case
1049  libmesh_assert(neigh->has_children());
1050  bool neigh_has_remote_children = false;
1051  for (auto & child : neigh->child_ref_range())
1052  if (&child == remote_elem)
1053  neigh_has_remote_children = true;
1054  libmesh_assert(neigh_has_remote_children);
1055 
1056  // And let's double-check that we don't have
1057  // a remote_elem neighboring an active local element
1058  if (current_elem->active())
1059  libmesh_assert_not_equal_to (current_elem->processor_id(),
1060  this->processor_id());
1061 #endif // DEBUG
1062  neigh = const_cast<RemoteElem *>(remote_elem);
1063  }
1064  // If neigh and current_elem are more than one level
1065  // apart, figuring out whether we have a remote
1066  // neighbor here becomes much harder.
1067  else if (neigh && (current_elem->subactive() &&
1068  neigh->has_children()))
1069  {
1070  // Find the deepest descendant of neigh which
1071  // we could consider for a neighbor. If we run
1072  // out of neigh children, then that's our
1073  // neighbor. If we find a potential neighbor
1074  // with remote_children and we don't find any
1075  // potential neighbors among its non-remote
1076  // children, then our neighbor must be remote.
1077  while (neigh != remote_elem &&
1078  neigh->has_children())
1079  {
1080  bool found_neigh = false;
1081  for (unsigned int c = 0, nc = neigh->n_children();
1082  !found_neigh && c != nc; ++c)
1083  {
1084  Elem * child = neigh->child_ptr(c);
1085  if (child == remote_elem)
1086  continue;
1087  for (auto ncn : child->neighbor_ptr_range())
1088  {
1089  if (ncn != remote_elem &&
1090  ncn->is_ancestor_of(current_elem))
1091  {
1092  neigh = ncn;
1093  found_neigh = true;
1094  break;
1095  }
1096  }
1097  }
1098  if (!found_neigh)
1099  neigh = const_cast<RemoteElem *>(remote_elem);
1100  }
1101  }
1102  current_elem->set_neighbor(s, neigh);
1103 #ifdef DEBUG
1104  if (neigh != nullptr && neigh != remote_elem)
1105  // We ignore subactive elements here because
1106  // we don't care about neighbors of subactive element.
1107  if ((!neigh->active()) && (!current_elem->subactive()))
1108  {
1109  libMesh::err << "On processor " << this->processor_id()
1110  << std::endl;
1111  libMesh::err << "Bad element ID = " << current_elem->id()
1112  << ", Side " << s << ", Bad neighbor ID = " << neigh->id() << std::endl;
1113  libMesh::err << "Bad element proc_ID = " << current_elem->processor_id()
1114  << ", Bad neighbor proc_ID = " << neigh->processor_id() << std::endl;
1115  libMesh::err << "Bad element size = " << current_elem->hmin()
1116  << ", Bad neighbor size = " << neigh->hmin() << std::endl;
1117  libMesh::err << "Bad element center = " << current_elem->vertex_average()
1118  << ", Bad neighbor center = " << neigh->vertex_average() << std::endl;
1119  libMesh::err << "ERROR: "
1120  << (current_elem->active()?"Active":"Ancestor")
1121  << " Element at level "
1122  << current_elem->level() << std::endl;
1123  libMesh::err << "with "
1124  << (parent->active()?"active":
1125  (parent->subactive()?"subactive":"ancestor"))
1126  << " parent share "
1127  << (neigh->subactive()?"subactive":"ancestor")
1128  << " neighbor at level " << neigh->level()
1129  << std::endl;
1130  NameBasedIO(*this).write ("bad_mesh.gmv");
1131  libmesh_error_msg("Problematic mesh written to bad_mesh.gmv.");
1132  }
1133 #endif // DEBUG
1134  }
1135  }
1136 
1137  // We can skip to the next element if we're full-dimension
1138  // and therefore don't have any interior parents
1139  if (current_elem->dim() >= LIBMESH_DIM)
1140  continue;
1141 
1142  // We have no interior parents unless we can find one later
1143  current_elem->set_interior_parent(nullptr);
1144 
1145  Elem * pip = parent->interior_parent();
1146 
1147  if (!pip)
1148  continue;
1149 
1150  // If there's no interior_parent children, whether due to a
1151  // remote element or a non-conformity, then there's no
1152  // children to search.
1153  if (pip == remote_elem || pip->active())
1154  {
1155  current_elem->set_interior_parent(pip);
1156  continue;
1157  }
1158 
1159  // For node comparisons we'll need a sensible tolerance
1160  Real node_tolerance = current_elem->hmin() * TOLERANCE;
1161 
1162  // Otherwise our interior_parent should be a child of our
1163  // parent's interior_parent.
1164  for (auto & child : pip->child_ref_range())
1165  {
1166  // If we have a remote_elem, that might be our
1167  // interior_parent. We'll set it provisionally now and
1168  // keep trying to find something better.
1169  if (&child == remote_elem)
1170  {
1171  current_elem->set_interior_parent
1172  (const_cast<RemoteElem *>(remote_elem));
1173  continue;
1174  }
1175 
1176  bool child_contains_our_nodes = true;
1177  for (auto & n : current_elem->node_ref_range())
1178  {
1179  bool child_contains_this_node = false;
1180  for (auto & cn : child.node_ref_range())
1181  if (cn.absolute_fuzzy_equals
1182  (n, node_tolerance))
1183  {
1184  child_contains_this_node = true;
1185  break;
1186  }
1187  if (!child_contains_this_node)
1188  {
1189  child_contains_our_nodes = false;
1190  break;
1191  }
1192  }
1193  if (child_contains_our_nodes)
1194  {
1195  current_elem->set_interior_parent(&child);
1196  break;
1197  }
1198  }
1199 
1200  // We should have found *some* interior_parent at this
1201  // point, whether semilocal or remote.
1202  libmesh_assert(current_elem->interior_parent());
1203  }
1204  }
1205 
1206 #endif // AMR
1207 
1208 
1209 #ifdef DEBUG
1211  !reset_remote_elements);
1213 #endif
1214 }
1215 
1216 
1217 
1218 void UnstructuredMesh::read (const std::string & name,
1219  void *,
1220  bool skip_renumber_nodes_and_elements,
1221  bool skip_find_neighbors)
1222 {
1223  // Set the skip_renumber_nodes_and_elements flag on all processors
1224  // if necessary.
1225  // This ensures that renumber_nodes_and_elements is *not* called
1226  // during prepare_for_use() for certain types of mesh files.
1227  // This is required in cases where there is an associated solution
1228  // file which expects a certain ordering of the nodes.
1229  if (name.rfind(".gmv") == name.size() - 4)
1230  this->allow_renumbering(false);
1231 
1232  NameBasedIO(*this).read(name);
1233 
1234  if (skip_renumber_nodes_and_elements)
1235  {
1236  // Use MeshBase::allow_renumbering() yourself instead.
1237  libmesh_deprecated();
1238  this->allow_renumbering(false);
1239  }
1240 
1241  // Done reading the mesh. Now prepare it for use.
1242  const bool old_allow_find_neighbors = this->allow_find_neighbors();
1243  this->allow_find_neighbors(!skip_find_neighbors);
1244  this->prepare_for_use();
1245  this->allow_find_neighbors(old_allow_find_neighbors);
1246 }
1247 
1248 
1249 
1250 void UnstructuredMesh::write (const std::string & name)
1251 {
1252  LOG_SCOPE("write()", "Mesh");
1253 
1254  NameBasedIO(*this).write(name);
1255 }
1256 
1257 
1258 
1259 void UnstructuredMesh::write (const std::string & name,
1260  const std::vector<Number> & v,
1261  const std::vector<std::string> & vn)
1262 {
1263  LOG_SCOPE("write()", "Mesh");
1264 
1265  NameBasedIO(*this).write_nodal_data(name, v, vn);
1266 }
1267 
1268 
1269 
1270 
1271 
1273  const processor_id_type pid) const
1274 {
1275 
1276  // Issue a warning if the number the number of processors
1277  // currently available is less that that requested for
1278  // partitioning. This is not necessarily an error since
1279  // you may run on one processor and still partition the
1280  // mesh into several partitions.
1281 #ifdef DEBUG
1282  if (this->n_processors() < pid)
1283  {
1284  libMesh::out << "WARNING: You are creating a "
1285  << "mesh for a processor id (="
1286  << pid
1287  << ") greater than "
1288  << "the number of processors available for "
1289  << "the calculation. (="
1290  << this->n_processors()
1291  << ")."
1292  << std::endl;
1293  }
1294 #endif
1295 
1296  this->create_submesh (pid_mesh,
1297  this->active_pid_elements_begin(pid),
1298  this->active_pid_elements_end(pid));
1299 }
1300 
1301 
1302 
1303 
1304 
1305 
1306 
1308  const const_element_iterator & it,
1309  const const_element_iterator & it_end) const
1310 {
1311  // Just in case the subdomain_mesh already has some information
1312  // in it, get rid of it.
1313  new_mesh.clear();
1314 
1315  // If we're not serial, our submesh isn't either.
1316  // There are no remote elements to delete on an empty mesh, but
1317  // calling the method to do so marks the mesh as parallel.
1318  if (!this->is_serial())
1319  new_mesh.delete_remote_elements();
1320 
1321  // Fail if (*this == new_mesh), we cannot create a submesh inside ourself!
1322  // This may happen if the user accidentally passes the original mesh into
1323  // this function! We will check this by making sure we did not just
1324  // clear ourself.
1325  libmesh_assert_not_equal_to (this->n_nodes(), 0);
1326  libmesh_assert_not_equal_to (this->n_elem(), 0);
1327 
1328  // Container to catch boundary IDs handed back by BoundaryInfo
1329  std::vector<boundary_id_type> bc_ids;
1330 
1331  // Put any extra integers on the new mesh too
1332  new_mesh.merge_extra_integer_names(*this);
1333  const unsigned int n_node_ints = _node_integer_names.size();
1334 
1335  for (const auto & old_elem : as_range(it, it_end))
1336  {
1337  // Add an equivalent element type to the new_mesh.
1338  // disconnected_clone() copies ids, extra element integers, etc.
1339  auto uelem = old_elem->disconnected_clone();
1340  Elem * new_elem = new_mesh.add_elem(std::move(uelem));
1341  libmesh_assert(new_elem);
1342 
1343  // Loop over the nodes on this element.
1344  for (auto n : old_elem->node_index_range())
1345  {
1346  const dof_id_type this_node_id = old_elem->node_id(n);
1347 
1348  // Add this node to the new mesh if it's not there already
1349  if (!new_mesh.query_node_ptr(this_node_id))
1350  {
1351  Node * newn =
1352  new_mesh.add_point (old_elem->point(n),
1353  this_node_id,
1354  old_elem->node_ptr(n)->processor_id());
1355 
1356  newn->add_extra_integers(n_node_ints);
1357  for (unsigned int i = 0; i != n_node_ints; ++i)
1358  newn->set_extra_integer(i, old_elem->node_ptr(n)->get_extra_integer(i));
1359 
1360 #ifdef LIBMESH_ENABLE_UNIQUE_ID
1361  newn->set_unique_id(old_elem->node_ptr(n)->unique_id());
1362 #endif
1363  }
1364 
1365  // Define this element's connectivity on the new mesh
1366  new_elem->set_node(n) = new_mesh.node_ptr(this_node_id);
1367  }
1368 
1369  // Maybe add boundary conditions for this element
1370  for (auto s : old_elem->side_index_range())
1371  {
1372  this->get_boundary_info().boundary_ids(old_elem, s, bc_ids);
1373  new_mesh.get_boundary_info().add_side (new_elem, s, bc_ids);
1374  }
1375  } // end loop over elements
1376 
1377  // Prepare the new_mesh for use
1378  new_mesh.prepare_for_use();
1379 }
1380 
1381 
1382 
1383 #ifdef LIBMESH_ENABLE_AMR
1385 {
1386  LOG_SCOPE ("contract()", "Mesh");
1387 
1388  // Flag indicating if this call actually changes the mesh
1389  bool mesh_changed = false;
1390 
1391 #ifdef DEBUG
1392  for (const auto & elem : this->element_ptr_range())
1393  libmesh_assert(elem->active() || elem->subactive() || elem->ancestor());
1394 #endif
1395 
1396  // Loop over the elements.
1397  for (auto & elem : this->element_ptr_range())
1398  {
1399  // Delete all the subactive ones
1400  if (elem->subactive())
1401  {
1402  // No level-0 element should be subactive.
1403  // Note that we CAN'T test elem->level(), as that
1404  // touches elem->parent()->dim(), and elem->parent()
1405  // might have already been deleted!
1406  libmesh_assert(elem->parent());
1407 
1408  // Delete the element
1409  // This just sets a pointer to nullptr, and doesn't
1410  // invalidate any iterators
1411  this->delete_elem(elem);
1412 
1413  // the mesh has certainly changed
1414  mesh_changed = true;
1415  }
1416  else
1417  {
1418  // Compress all the active ones
1419  if (elem->active())
1420  elem->contract();
1421  else
1422  libmesh_assert (elem->ancestor());
1423  }
1424  }
1425 
1426  // Strip any newly-created nullptr voids out of the element array
1428 
1429  // FIXME: Need to understand why deleting subactive children
1430  // invalidates the point locator. For now we will clear it explicitly
1431  this->clear_point_locator();
1432 
1433  // Allow our GhostingFunctor objects to reinit if necessary.
1434  for (auto & gf : as_range(this->ghosting_functors_begin(),
1435  this->ghosting_functors_end()))
1436  {
1437  libmesh_assert(gf);
1438  gf->mesh_reinit();
1439  }
1440 
1441  return mesh_changed;
1442 }
1443 #endif // #ifdef LIBMESH_ENABLE_AMR
1444 
1445 
1446 
1448 {
1449  LOG_SCOPE("all_first_order()", "Mesh");
1450 
1454  std::vector<bool> node_touched_by_me(this->max_node_id(), false);
1455 
1456  // Loop over the high-ordered elements.
1457  // First make sure they _are_ indeed high-order, and then replace
1458  // them with an equivalent first-order element.
1459  for (auto & so_elem : element_ptr_range())
1460  {
1461  libmesh_assert(so_elem);
1462 
1463  /*
1464  * build the first-order equivalent, add to
1465  * the new_elements list.
1466  */
1467  auto lo_elem = Elem::build
1469  (so_elem->type()), so_elem->parent());
1470 
1471  const unsigned short n_sides = so_elem->n_sides();
1472 
1473  for (unsigned short s=0; s != n_sides; ++s)
1474  if (so_elem->neighbor_ptr(s) == remote_elem)
1475  lo_elem->set_neighbor(s, const_cast<RemoteElem *>(remote_elem));
1476 
1477 #ifdef LIBMESH_ENABLE_AMR
1478  /*
1479  * Reset the parent links of any child elements
1480  */
1481  if (so_elem->has_children())
1482  for (unsigned int c = 0, nc = so_elem->n_children(); c != nc; ++c)
1483  {
1484  Elem * child = so_elem->child_ptr(c);
1485  if (child != remote_elem)
1486  child->set_parent(lo_elem.get());
1487  lo_elem->add_child(child, c);
1488  }
1489 
1490  /*
1491  * Reset the child link of any parent element
1492  */
1493  if (so_elem->parent())
1494  {
1495  unsigned int c =
1496  so_elem->parent()->which_child_am_i(so_elem);
1497  lo_elem->parent()->replace_child(lo_elem.get(), c);
1498  }
1499 
1500  /*
1501  * Copy as much data to the new element as makes sense
1502  */
1503  lo_elem->set_p_level(so_elem->p_level());
1504  lo_elem->set_refinement_flag(so_elem->refinement_flag());
1505  lo_elem->set_p_refinement_flag(so_elem->p_refinement_flag());
1506 #endif
1507 
1508  libmesh_assert_equal_to (lo_elem->n_vertices(), so_elem->n_vertices());
1509 
1510  /*
1511  * By definition the vertices of the linear and
1512  * second order element are identically numbered.
1513  * transfer these.
1514  */
1515  for (unsigned int v=0, snv=so_elem->n_vertices(); v < snv; v++)
1516  {
1517  lo_elem->set_node(v) = so_elem->node_ptr(v);
1518  node_touched_by_me[lo_elem->node_id(v)] = true;
1519  }
1520 
1521  /*
1522  * find_neighbors relies on remote_elem neighbor links being
1523  * properly maintained.
1524  */
1525  for (unsigned short s=0; s != n_sides; s++)
1526  {
1527  if (so_elem->neighbor_ptr(s) == remote_elem)
1528  lo_elem->set_neighbor(s, const_cast<RemoteElem*>(remote_elem));
1529  }
1530 
1538  (this->get_boundary_info(), so_elem, lo_elem.get());
1539 
1540  /*
1541  * The new first-order element is ready.
1542  * Inserting it into the mesh will replace and delete
1543  * the second-order element.
1544  */
1545  lo_elem->set_id(so_elem->id());
1546 #ifdef LIBMESH_ENABLE_UNIQUE_ID
1547  lo_elem->set_unique_id(so_elem->unique_id());
1548 #endif
1549  lo_elem->processor_id() = so_elem->processor_id();
1550  lo_elem->subdomain_id() = so_elem->subdomain_id();
1551 
1552  const unsigned int nei = so_elem->n_extra_integers();
1553  lo_elem->add_extra_integers(nei);
1554  for (unsigned int i=0; i != nei; ++i)
1555  lo_elem->set_extra_integer(i, so_elem->get_extra_integer(i));
1556 
1557  // This is probably moot but shouldn't hurt
1558  lo_elem->set_mapping_type(so_elem->mapping_type());
1559  lo_elem->set_mapping_data(so_elem->mapping_data());
1560 
1561  this->insert_elem(std::move(lo_elem));
1562  }
1563 
1564  // Deleting nodes does not invalidate iterators, so this is safe.
1565  for (const auto & node : this->node_ptr_range())
1566  if (!node_touched_by_me[node->id()])
1567  this->delete_node(node);
1568 
1569  // If crazy people applied boundary info to non-vertices and then
1570  // deleted those non-vertices, we should make sure their boundary id
1571  // caches are correct.
1573 
1574  // On hanging nodes that used to also be second order nodes, we
1575  // might now have an invalid nodal processor_id()
1577 
1578  // delete or renumber nodes if desired
1579  this->prepare_for_use();
1580 }
1581 
1582 
1583 
1584 void
1586  const bool full_ordered)
1587 {
1588  LOG_SCOPE("all_second_order_range()", "Mesh");
1589 
1590  /*
1591  * The maximum number of new second order nodes we might be adding,
1592  * for use when picking unique unique_id values later. This variable
1593  * is not used unless unique ids are enabled.
1594  */
1595  unsigned int max_new_nodes_per_elem;
1596 
1597  /*
1598  * For speed-up of the \p add_point() method, we
1599  * can reserve memory. Guess the number of additional
1600  * nodes based on the element spatial dimensions and the
1601  * total number of nodes in the mesh as an upper bound.
1602  */
1603  switch (this->mesh_dimension())
1604  {
1605  case 1:
1606  /*
1607  * in 1D, there can only be order-increase from Edge2
1608  * to Edge3. Something like 1/2 of n_nodes() have
1609  * to be added
1610  */
1611  max_new_nodes_per_elem = 3 - 2;
1612  this->reserve_nodes(static_cast<unsigned int>
1613  (1.5*static_cast<double>(this->n_nodes())));
1614  break;
1615 
1616  case 2:
1617  /*
1618  * in 2D, either refine from Tri3 to Tri6 (double the nodes)
1619  * or from Quad4 to Quad8 (again, double) or Quad9 (2.25 that much)
1620  */
1621  max_new_nodes_per_elem = 9 - 4;
1622  this->reserve_nodes(static_cast<unsigned int>
1623  (2*static_cast<double>(this->n_nodes())));
1624  break;
1625 
1626 
1627  case 3:
1628  /*
1629  * in 3D, either refine from Tet4 to Tet10 (factor = 2.5) up to
1630  * Hex8 to Hex27 (something > 3). Since in 3D there _are_ already
1631  * quite some nodes, and since we do not want to overburden the memory by
1632  * a too conservative guess, use the lower bound
1633  */
1634  max_new_nodes_per_elem = 27 - 8;
1635  this->reserve_nodes(static_cast<unsigned int>
1636  (2.5*static_cast<double>(this->n_nodes())));
1637  break;
1638 
1639  default:
1640  // Hm?
1641  libmesh_error_msg("Unknown mesh dimension " << this->mesh_dimension());
1642  }
1643 
1644  // All the real work is done in the helper function
1645  all_increased_order_range(*this, range, max_new_nodes_per_elem,
1646  [full_ordered](ElemType t) {
1647  return Elem::second_order_equivalent_type(t, full_ordered);
1648  });
1649 }
1650 
1651 
1652 
1654 {
1655  LOG_SCOPE("all_complete_order()", "Mesh");
1656 
1657  /*
1658  * The maximum number of new higher-order nodes we might be adding,
1659  * for use when picking unique unique_id values later. This variable
1660  * is not used unless unique ids are enabled.
1661  */
1662  unsigned int max_new_nodes_per_elem;
1663 
1664  /*
1665  * for speed-up of the \p add_point() method, we
1666  * can reserve memory. Guess the number of additional
1667  * nodes based on the element spatial dimensions and the
1668  * total number of nodes in the mesh as an upper bound.
1669  */
1670  switch (this->mesh_dimension())
1671  {
1672  case 1:
1673  /*
1674  * in 1D, there can only be order-increase from Edge2
1675  * to Edge3. Something like 1/2 of n_nodes() have
1676  * to be added
1677  */
1678  max_new_nodes_per_elem = 3 - 2;
1679  this->reserve_nodes(static_cast<unsigned int>
1680  (1.5*static_cast<double>(this->n_nodes())));
1681  break;
1682 
1683  case 2:
1684  /*
1685  * in 2D, we typically refine from Tri6 to Tri7 (1.1667 times
1686  * the nodes) but might refine from Quad4 to Quad9
1687  * (2.25 times the nodes)
1688  */
1689  max_new_nodes_per_elem = 9 - 4;
1690  this->reserve_nodes(static_cast<unsigned int>
1691  (2*static_cast<double>(this->n_nodes())));
1692  break;
1693 
1694 
1695  case 3:
1696  /*
1697  * in 3D, we typically refine from Tet10 to Tet14 (factor = 1.4)
1698  * but may go Hex8 to Hex27 (something > 3). Since in 3D there
1699  * _are_ already quite some nodes, and since we do not want to
1700  * overburden the memory by a too conservative guess, use a
1701  * moderate bound
1702  */
1703  max_new_nodes_per_elem = 27 - 8;
1704  this->reserve_nodes(static_cast<unsigned int>
1705  (2.5*static_cast<double>(this->n_nodes())));
1706  break;
1707 
1708  default:
1709  // Hm?
1710  libmesh_error_msg("Unknown mesh dimension " << this->mesh_dimension());
1711  }
1712 
1713  // All the real work is done in the helper function
1714  all_increased_order_range(*this, range, max_new_nodes_per_elem,
1715  [](ElemType t) {
1717  });
1718 }
1719 
1720 
1721 std::size_t
1723  boundary_id_type this_mesh_boundary_id,
1724  boundary_id_type other_mesh_boundary_id,
1725  Real tol,
1726  bool clear_stitched_boundary_ids,
1727  bool verbose,
1728  bool use_binary_search,
1729  bool enforce_all_nodes_match_on_boundaries,
1730  bool merge_boundary_nodes_all_or_nothing,
1731  bool remap_subdomain_ids)
1732 {
1733  LOG_SCOPE("stitch_meshes()", "UnstructuredMesh");
1734  return stitching_helper(&other_mesh,
1735  this_mesh_boundary_id,
1736  other_mesh_boundary_id,
1737  tol,
1738  clear_stitched_boundary_ids,
1739  verbose,
1740  use_binary_search,
1741  enforce_all_nodes_match_on_boundaries,
1742  true,
1743  merge_boundary_nodes_all_or_nothing,
1744  remap_subdomain_ids);
1745 }
1746 
1747 
1748 std::size_t
1750  boundary_id_type boundary_id_2,
1751  Real tol,
1752  bool clear_stitched_boundary_ids,
1753  bool verbose,
1754  bool use_binary_search,
1755  bool enforce_all_nodes_match_on_boundaries,
1756  bool merge_boundary_nodes_all_or_nothing)
1757 
1758 {
1759  return stitching_helper(nullptr,
1760  boundary_id_1,
1761  boundary_id_2,
1762  tol,
1763  clear_stitched_boundary_ids,
1764  verbose,
1765  use_binary_search,
1766  enforce_all_nodes_match_on_boundaries,
1767  true,
1768  merge_boundary_nodes_all_or_nothing,
1769  false);
1770 }
1771 
1772 
1773 std::size_t
1775  boundary_id_type this_mesh_boundary_id,
1776  boundary_id_type other_mesh_boundary_id,
1777  Real tol,
1778  bool clear_stitched_boundary_ids,
1779  bool verbose,
1780  bool use_binary_search,
1781  bool enforce_all_nodes_match_on_boundaries,
1782  bool skip_find_neighbors,
1783  bool merge_boundary_nodes_all_or_nothing,
1784  bool remap_subdomain_ids)
1785 {
1786 #ifdef DEBUG
1787  // We rely on neighbor links here
1789 #endif
1790 
1791  // FIXME: make distributed mesh support efficient.
1792  // Yes, we currently suck.
1793  MeshSerializer serialize(*this);
1794 
1795  // *Badly*.
1796  std::unique_ptr<MeshSerializer> serialize_other;
1797  if (other_mesh)
1798  serialize_other = std::make_unique<MeshSerializer>
1799  (*const_cast<MeshBase *>(other_mesh));
1800 
1801  std::map<dof_id_type, dof_id_type> node_to_node_map, other_to_this_node_map; // The second is the inverse map of the first
1802  std::map<dof_id_type, std::vector<dof_id_type>> node_to_elems_map;
1803 
1804  typedef dof_id_type key_type;
1805  typedef std::pair<const Elem *, unsigned char> val_type;
1806  typedef std::pair<key_type, val_type> key_val_pair;
1807  typedef std::unordered_multimap<key_type, val_type> map_type;
1808  // Mapping between all side keys in this mesh and elements+side numbers relevant to the boundary in this mesh as well.
1809  map_type side_to_elem_map;
1810 
1811  // If there is only one mesh (i.e. other_mesh == nullptr), then loop over this mesh twice
1812  if (!other_mesh)
1813  {
1814  other_mesh = this;
1815  }
1816 
1817  if ((this_mesh_boundary_id != BoundaryInfo::invalid_id) &&
1818  (other_mesh_boundary_id != BoundaryInfo::invalid_id))
1819  {
1820  LOG_SCOPE("stitch_meshes node merging", "UnstructuredMesh");
1821 
1822  // While finding nodes on the boundary, also find the minimum edge length
1823  // of all faces on both boundaries. This will later be used in relative
1824  // distance checks when stitching nodes.
1825  Real h_min = std::numeric_limits<Real>::max();
1826  bool h_min_updated = false;
1827 
1828  // Loop below fills in these sets for the two meshes.
1829  std::set<dof_id_type> this_boundary_node_ids, other_boundary_node_ids;
1830 
1831  // Pull objects out of the loop to reduce heap operations
1832  std::unique_ptr<const Elem> side;
1833 
1834  {
1835  // Make temporary fixed-size arrays for loop
1836  boundary_id_type id_array[2] = {this_mesh_boundary_id, other_mesh_boundary_id};
1837  std::set<dof_id_type> * set_array[2] = {&this_boundary_node_ids, &other_boundary_node_ids};
1838  const MeshBase * mesh_array[2] = {this, other_mesh};
1839 
1840  for (unsigned i=0; i<2; ++i)
1841  {
1842  // First we deal with node boundary IDs. We only enter
1843  // this loop if we have at least one nodeset. Note that we
1844  // do not attempt to make an h_min determination here.
1845  // The h_min determination is done while looping over the
1846  // Elems and checking their sides and edges for boundary
1847  // information, below.
1848  if (mesh_array[i]->get_boundary_info().n_nodeset_conds() > 0)
1849  {
1850  // build_node_list() returns a vector of (node-id, bc-id) tuples
1851  for (const auto & t : mesh_array[i]->get_boundary_info().build_node_list())
1852  {
1853  boundary_id_type node_bc_id = std::get<1>(t);
1854  if (node_bc_id == id_array[i])
1855  {
1856  dof_id_type this_node_id = std::get<0>(t);
1857  set_array[i]->insert( this_node_id );
1858  }
1859  }
1860  }
1861 
1862  // Container to catch boundary IDs passed back from BoundaryInfo.
1863  std::vector<boundary_id_type> bc_ids;
1864 
1865  // Pointers to boundary NodeElems encountered while looping over the entire Mesh
1866  // and checking side and edge boundary ids. The Nodes associated with NodeElems
1867  // may be in a boundary nodeset, but not connected to any other Elems. In this
1868  // case, we also consider the "minimum node separation distance" amongst all
1869  // NodeElems when determining the relevant h_min value for this mesh.
1870  std::vector<const Elem *> boundary_node_elems;
1871 
1872  for (auto & el : mesh_array[i]->element_ptr_range())
1873  {
1874  // Now check whether elem has a face on the specified boundary
1875  for (auto side_id : el->side_index_range())
1876  if (el->neighbor_ptr(side_id) == nullptr)
1877  {
1878  // Get *all* boundary IDs on this side, not just the first one!
1879  mesh_array[i]->get_boundary_info().boundary_ids (el, side_id, bc_ids);
1880 
1881  if (std::find(bc_ids.begin(), bc_ids.end(), id_array[i]) != bc_ids.end())
1882  {
1883  el->build_side_ptr(side, side_id);
1884  for (auto & n : side->node_ref_range())
1885  set_array[i]->insert(n.id());
1886 
1887  h_min = std::min(h_min, side->hmin());
1888  h_min_updated = true;
1889 
1890  // This side is on the boundary, add its information to side_to_elem
1891  if (skip_find_neighbors && (i==0))
1892  {
1893  key_type key = el->low_order_key(side_id);
1894  val_type val;
1895  val.first = el;
1896  val.second = cast_int<unsigned char>(side_id);
1897 
1898  key_val_pair kvp;
1899  kvp.first = key;
1900  kvp.second = val;
1901  side_to_elem_map.insert (kvp);
1902  }
1903  }
1904 
1905  // Also, check the edges on this side. We don't have to worry about
1906  // updating neighbor info in this case since elements don't store
1907  // neighbor info on edges.
1908  for (auto edge_id : el->edge_index_range())
1909  {
1910  if (el->is_edge_on_side(edge_id, side_id))
1911  {
1912  // Get *all* boundary IDs on this edge, not just the first one!
1913  mesh_array[i]->get_boundary_info().edge_boundary_ids (el, edge_id, bc_ids);
1914 
1915  if (std::find(bc_ids.begin(), bc_ids.end(), id_array[i]) != bc_ids.end())
1916  {
1917  std::unique_ptr<const Elem> edge (el->build_edge_ptr(edge_id));
1918  for (auto & n : edge->node_ref_range())
1919  set_array[i]->insert( n.id() );
1920 
1921  h_min = std::min(h_min, edge->hmin());
1922  h_min_updated = true;
1923  }
1924  }
1925  } // end for (edge_id)
1926  } // end if (side == nullptr)
1927 
1928  // Alternatively, is this a boundary NodeElem? If so,
1929  // add it to a list of NodeElems that will later be
1930  // used to set h_min based on the minimum node
1931  // separation distance between all pairs of boundary
1932  // NodeElems.
1933  if (el->type() == NODEELEM)
1934  {
1935  mesh_array[i]->get_boundary_info().boundary_ids(el->node_ptr(0), bc_ids);
1936  if (std::find(bc_ids.begin(), bc_ids.end(), id_array[i]) != bc_ids.end())
1937  {
1938  boundary_node_elems.push_back(el);
1939 
1940  // Debugging:
1941  // libMesh::out << "Elem " << el->id() << " is a NodeElem on boundary " << id_array[i] << std::endl;
1942  }
1943  }
1944  } // end for (el)
1945 
1946  // Compute the minimum node separation distance amongst
1947  // all boundary NodeElem pairs.
1948  {
1949  const auto N = boundary_node_elems.size();
1950  for (auto node_elem_i : make_range(N))
1951  for (auto node_elem_j : make_range(node_elem_i+1, N))
1952  {
1953  Real node_sep =
1954  (boundary_node_elems[node_elem_i]->point(0) - boundary_node_elems[node_elem_j]->point(0)).norm();
1955 
1956  // We only want to consider non-coincident
1957  // boundary NodeElem pairs when determining the
1958  // minimum node separation distance.
1959  if (node_sep > 0.)
1960  {
1961  h_min = std::min(h_min, node_sep);
1962  h_min_updated = true;
1963  }
1964  } // end for (node_elem_j)
1965  } // end minimum NodeElem separation scope
1966  } // end for (i)
1967  } // end scope
1968 
1969  if (verbose)
1970  {
1971  libMesh::out << "In UnstructuredMesh::stitch_meshes:\n"
1972  << "This mesh has " << this_boundary_node_ids.size()
1973  << " nodes on boundary `"
1974  << this->get_boundary_info().get_sideset_name(this_mesh_boundary_id)
1975  << "' (" << this_mesh_boundary_id << ").\n"
1976  << "Other mesh has " << other_boundary_node_ids.size()
1977  << " nodes on boundary `"
1978  << this->get_boundary_info().get_sideset_name(other_mesh_boundary_id)
1979  << "' (" << other_mesh_boundary_id << ").\n";
1980 
1981  if (h_min_updated)
1982  {
1983  libMesh::out << "Minimum edge length on both surfaces is " << h_min << ".\n";
1984  }
1985  else
1986  {
1987  libMesh::out << "No minimum edge length determined on specified surfaces." << std::endl;
1988  }
1989  }
1990 
1991  // At this point, if h_min==0 it means that there were at least two coincident
1992  // nodes on the surfaces being stitched, and we don't currently support that case.
1993  // (It might be possible to support, but getting it exactly right would be tricky
1994  // and probably not worth the extra complications to the "normal" case.)
1995  libmesh_error_msg_if(h_min < std::numeric_limits<Real>::epsilon(),
1996  "Coincident nodes detected on source and/or target "
1997  "surface, stitching meshes is not possible.");
1998 
1999  // We require nanoflann for the "binary search" (really kd-tree)
2000  // option to work. If it's not available, turn that option off,
2001  // warn the user, and fall back on the N^2 search algorithm.
2002  if (use_binary_search)
2003  {
2004 #ifndef LIBMESH_HAVE_NANOFLANN
2005  use_binary_search = false;
2006  libmesh_warning("The use_binary_search option in the "
2007  "UnstructuredMesh stitching algorithms requires nanoflann "
2008  "support. Falling back on N^2 search algorithm.");
2009 #endif
2010  }
2011 
2012  if (!this_boundary_node_ids.empty())
2013  {
2014  if (use_binary_search)
2015  {
2016 #ifdef LIBMESH_HAVE_NANOFLANN
2017  typedef nanoflann::KDTreeSingleIndexAdaptor<nanoflann::L2_Simple_Adaptor<Real, VectorOfNodesAdaptor>,
2018  VectorOfNodesAdaptor, 3, std::size_t> kd_tree_t;
2019 
2020  // Create the dataset needed to build the kd tree with nanoflann
2021  std::vector<std::pair<Point, dof_id_type>> this_mesh_nodes(this_boundary_node_ids.size());
2022 
2023  for (auto [it, ctr] = std::make_tuple(this_boundary_node_ids.begin(), 0u);
2024  it != this_boundary_node_ids.end(); ++it, ++ctr)
2025  {
2026  this_mesh_nodes[ctr].first = this->point(*it);
2027  this_mesh_nodes[ctr].second = *it;
2028  }
2029 
2030  VectorOfNodesAdaptor vec_nodes_adaptor(this_mesh_nodes);
2031 
2032  kd_tree_t this_kd_tree(3, vec_nodes_adaptor, 10);
2033  this_kd_tree.buildIndex();
2034 
2035  // Storage for nearest neighbor in the loop below
2036  std::size_t ret_index;
2037  Real ret_dist_sqr;
2038 
2039  // Loop over other mesh. For each node, find its nearest neighbor in this mesh, and fill in the maps.
2040  for (const auto & node_id : other_boundary_node_ids)
2041  {
2042  const auto & p = other_mesh->point(node_id);
2043  const Real query_pt[] = {p(0), p(1), p(2)};
2044  this_kd_tree.knnSearch(&query_pt[0], 1, &ret_index, &ret_dist_sqr);
2045 
2046  // TODO: here we should use the user's specified tolerance
2047  // and the previously determined value of h_min in the
2048  // distance comparison, not just TOLERANCE^2.
2049  if (ret_dist_sqr < TOLERANCE*TOLERANCE)
2050  {
2051  node_to_node_map[this_mesh_nodes[ret_index].second] = node_id;
2052  other_to_this_node_map[node_id] = this_mesh_nodes[ret_index].second;
2053  }
2054  }
2055 
2056  // If the two maps don't have the same size, it means one
2057  // node in this mesh is the nearest neighbor of several
2058  // nodes in other mesh. Since the stitching is ambiguous in
2059  // this case, we throw an error.
2060  libmesh_error_msg_if(node_to_node_map.size() != other_to_this_node_map.size(),
2061  "Error: Found multiple matching nodes in stitch_meshes");
2062 #endif
2063  }
2064  else // !use_binary_search
2065  {
2066  // In the unlikely event that two meshes composed entirely of
2067  // NodeElems are being stitched together, we will not have
2068  // selected a valid h_min value yet, and the distance
2069  // comparison below will be true for essentially any two
2070  // nodes. In this case we simply fall back on an absolute
2071  // distance check.
2072  if (!h_min_updated)
2073  {
2074  libmesh_warning("No valid h_min value was found, falling back on "
2075  "absolute distance check in the N^2 search algorithm.");
2076  h_min = 1.;
2077  }
2078 
2079  // Otherwise, use a simple N^2 search to find the closest matching points. This can be helpful
2080  // in the case that we have tolerance issues which cause mismatch between the two surfaces
2081  // that are being stitched.
2082  for (const auto & this_node_id : this_boundary_node_ids)
2083  {
2084  Node & this_node = this->node_ref(this_node_id);
2085 
2086  bool found_matching_nodes = false;
2087 
2088  for (const auto & other_node_id : other_boundary_node_ids)
2089  {
2090  const Node & other_node = other_mesh->node_ref(other_node_id);
2091 
2092  Real node_distance = (this_node - other_node).norm();
2093 
2094  if (node_distance < tol*h_min)
2095  {
2096  // Make sure we didn't already find a matching node!
2097  libmesh_error_msg_if(found_matching_nodes,
2098  "Error: Found multiple matching nodes in stitch_meshes");
2099 
2100  node_to_node_map[this_node_id] = other_node_id;
2101  other_to_this_node_map[other_node_id] = this_node_id;
2102 
2103  found_matching_nodes = true;
2104  }
2105  }
2106  }
2107  }
2108  }
2109 
2110  // Build up the node_to_elems_map, using only one loop over other_mesh
2111  for (auto & el : other_mesh->element_ptr_range())
2112  {
2113  // For each node on the element, find the corresponding node
2114  // on "this" Mesh, 'this_node_id', if it exists, and push
2115  // the current element ID back onto node_to_elems_map[this_node_id].
2116  // For that we will use the reverse mapping we created at
2117  // the same time as the forward mapping.
2118  for (auto & n : el->node_ref_range())
2119  if (const auto it = other_to_this_node_map.find(/*other_node_id=*/n.id());
2120  it != other_to_this_node_map.end())
2121  node_to_elems_map[/*this_node_id=*/it->second].push_back( el->id() );
2122  }
2123 
2124  if (verbose)
2125  {
2126  libMesh::out << "In UnstructuredMesh::stitch_meshes:\n"
2127  << "Found " << node_to_node_map.size()
2128  << " matching nodes.\n"
2129  << std::endl;
2130  }
2131 
2132  if (enforce_all_nodes_match_on_boundaries)
2133  {
2134  std::size_t n_matching_nodes = node_to_node_map.size();
2135  std::size_t this_mesh_n_nodes = this_boundary_node_ids.size();
2136  std::size_t other_mesh_n_nodes = other_boundary_node_ids.size();
2137  libmesh_error_msg_if((n_matching_nodes != this_mesh_n_nodes) || (n_matching_nodes != other_mesh_n_nodes),
2138  "Error: We expected the number of nodes to match.");
2139  }
2140 
2141  if (merge_boundary_nodes_all_or_nothing)
2142  {
2143  std::size_t n_matching_nodes = node_to_node_map.size();
2144  std::size_t this_mesh_n_nodes = this_boundary_node_ids.size();
2145  std::size_t other_mesh_n_nodes = other_boundary_node_ids.size();
2146  if ((n_matching_nodes != this_mesh_n_nodes) || (n_matching_nodes != other_mesh_n_nodes))
2147  {
2148  if (verbose)
2149  {
2150  libMesh::out << "Skipping node merging in "
2151  "UnstructuredMesh::stitch_meshes because not "
2152  "all boundary nodes were matched."
2153  << std::endl;
2154  }
2155  node_to_node_map.clear();
2156  other_to_this_node_map.clear();
2157  node_to_elems_map.clear();
2158  }
2159  }
2160  }
2161  else
2162  {
2163  if (verbose)
2164  {
2165  libMesh::out << "Skip node merging in UnstructuredMesh::stitch_meshes:" << std::endl;
2166  }
2167  }
2168 
2169  dof_id_type node_delta = this->max_node_id();
2170  dof_id_type elem_delta = this->max_elem_id();
2171 
2172  unique_id_type unique_delta =
2173 #ifdef LIBMESH_ENABLE_UNIQUE_ID
2174  this->parallel_max_unique_id();
2175 #else
2176  0;
2177 #endif
2178 
2179  // If other_mesh != nullptr, then we have to do a bunch of work
2180  // in order to copy it to this mesh
2181  if (this!=other_mesh)
2182  {
2183  LOG_SCOPE("stitch_meshes copying", "UnstructuredMesh");
2184 
2185  // Increment the node_to_node_map and node_to_elems_map
2186  // to account for id offsets
2187  for (auto & pr : node_to_node_map)
2188  pr.second += node_delta;
2189 
2190  for (auto & pr : node_to_elems_map)
2191  for (auto & entry : pr.second)
2192  entry += elem_delta;
2193 
2194  // We run into problems when the libMesh subdomain standard (the
2195  // id defines the subdomain; the name was an afterthought) and
2196  // the MOOSE standard (the name defines the subdomain; the id
2197  // might be autogenerated) clash.
2198  //
2199  // Subdomain ids with the same name in both meshes are surely
2200  // meant to represent the same subdomain. We can just merge
2201  // them.
2202  //
2203  // Subdomain ids which don't have a name in either mesh are
2204  // almost surely meant to represent the same subdomain. We'll
2205  // just merge them.
2206  //
2207  // Subdomain ids with different names in different meshes, or
2208  // names with different ids in different meshes, are trickier.
2209  // For backwards compatibility we default to the old "just copy
2210  // all the subdomain ids over" behavior, but if requested we'll
2211  // remap any ids that appear to be clear conflicts, and we'll
2212  // scream and die if we see any ids that are ambiguous due to
2213  // being named in one mesh but not the other.
2214  std::unordered_map<subdomain_id_type, subdomain_id_type> id_remapping;
2215  if (remap_subdomain_ids)
2216  {
2217  const auto & this_map = this->get_subdomain_name_map();
2218  const auto & other_map = other_mesh->get_subdomain_name_map();
2219  std::unordered_map<std::string, subdomain_id_type> other_map_reversed;
2220  for (auto & [sid, sname] : other_map)
2221  other_map_reversed.emplace(sname, sid);
2222 
2223  std::unordered_map<std::string, subdomain_id_type> this_map_reversed;
2224  for (auto & [sid, sname] : this_map)
2225  this_map_reversed.emplace(sname, sid);
2226 
2227  // We don't require either mesh to be prepared, but that
2228  // means we need to check for subdomains manually.
2229  auto get_subdomains = [](const MeshBase & mesh) {
2230  std::set<subdomain_id_type> all_subdomains;
2231  for (auto & el : mesh.element_ptr_range())
2232  all_subdomains.insert(el->subdomain_id());
2233  return all_subdomains;
2234  };
2235 
2236  const auto this_subdomains = get_subdomains(*this);
2237  const auto other_subdomains = get_subdomains(*other_mesh);
2238 
2239  for (auto & [sid, sname] : this_map)
2240  {
2241  // The same name with the same id means we're fine. The
2242  // same name with another id means we remap their id to
2243  // ours
2244  if (const auto other_reverse_it = other_map_reversed.find(sname);
2245  other_reverse_it != other_map_reversed.end() && other_reverse_it->second != sid)
2246  id_remapping[other_reverse_it->second] = sid;
2247 
2248  // The same id with a different name, we'll get to
2249  // later. The same id without any name means we don't
2250  // know what the user wants.
2251  if (other_subdomains.count(sid) && !other_map.count(sid))
2252  libmesh_error_msg("Can't safely stitch with a mesh sharing subdomain id "
2253  << sid << " but not subdomain name " << sname);
2254  }
2255 
2256  subdomain_id_type next_free_id = 0;
2257  // We might try to stitch empty meshes ...
2258  if (!this_subdomains.empty())
2259  next_free_id = *this_subdomains.rbegin() + 1;
2260  if (!other_subdomains.empty())
2261  next_free_id =
2262  std::max(next_free_id,
2263  cast_int<subdomain_id_type>
2264  (*other_subdomains.rbegin() + 1));
2265 
2266  for (auto & [sid, sname] : other_map)
2267  {
2268  // At this point we've figured out any remapping
2269  // necessary for an sname that we share. And we don't
2270  // need to remap any sid we don't share.
2271  if (!this_map_reversed.count(sname))
2272  {
2273  // But if we don't have this sname and we do have this
2274  // sid then we can't just merge into that.
2275  if (this_subdomains.count(sid))
2276  {
2277  // If we have this sid with no name, we don't
2278  // know what the user wants.
2279  if (!this_map.count(sid))
2280  libmesh_error_msg("Can't safely stitch with a mesh sharing subdomain id "
2281  << sid << " but under subdomain name " << sname);
2282 
2283  // We have this sid under a different name, so
2284  // we just need to give the other elements a new
2285  // id.
2286 
2287  // Users might have done crazy things with id
2288  // choice so let's make sure they didn't get too
2289  // crazy.
2290  libmesh_error_msg_if ((!this_subdomains.empty() &&
2291  next_free_id < *this_subdomains.rbegin()) ||
2292  (!other_subdomains.empty() &&
2293  next_free_id < *other_subdomains.rbegin()),
2294  "Subdomain id overflow");
2295 
2296  id_remapping[sid] = next_free_id++;
2297  this->subdomain_name(next_free_id) = sname;
2298  }
2299  // If we don't have this subdomain id, well, we're
2300  // about to, so we should have its name too.
2301  else
2302  this->subdomain_name(sid) = sname;
2303  }
2304  }
2305  }
2306 
2307  // Copy mesh data. If we skip the call to find_neighbors(), the lists
2308  // of neighbors will be copied verbatim from the other mesh
2309  this->copy_nodes_and_elements(*other_mesh, skip_find_neighbors,
2310  elem_delta, node_delta,
2311  unique_delta, &id_remapping);
2312 
2313  // Copy BoundaryInfo from other_mesh too. We do this via the
2314  // list APIs rather than element-by-element for speed.
2315  BoundaryInfo & boundary = this->get_boundary_info();
2316  const BoundaryInfo & other_boundary = other_mesh->get_boundary_info();
2317 
2318  for (const auto & t : other_boundary.build_node_list())
2319  boundary.add_node(std::get<0>(t) + node_delta,
2320  std::get<1>(t));
2321 
2322  for (const auto & t : other_boundary.build_side_list())
2323  boundary.add_side(std::get<0>(t) + elem_delta,
2324  std::get<1>(t),
2325  std::get<2>(t));
2326 
2327  for (const auto & t : other_boundary.build_edge_list())
2328  boundary.add_edge(std::get<0>(t) + elem_delta,
2329  std::get<1>(t),
2330  std::get<2>(t));
2331 
2332  for (const auto & t : other_boundary.build_shellface_list())
2333  boundary.add_shellface(std::get<0>(t) + elem_delta,
2334  std::get<1>(t),
2335  std::get<2>(t));
2336 
2337  const auto & other_ns_id_to_name = other_boundary.get_nodeset_name_map();
2338  auto & ns_id_to_name = boundary.set_nodeset_name_map();
2339  ns_id_to_name.insert(other_ns_id_to_name.begin(), other_ns_id_to_name.end());
2340 
2341  const auto & other_ss_id_to_name = other_boundary.get_sideset_name_map();
2342  auto & ss_id_to_name = boundary.set_sideset_name_map();
2343  ss_id_to_name.insert(other_ss_id_to_name.begin(), other_ss_id_to_name.end());
2344 
2345  const auto & other_es_id_to_name = other_boundary.get_edgeset_name_map();
2346  auto & es_id_to_name = boundary.set_edgeset_name_map();
2347  es_id_to_name.insert(other_es_id_to_name.begin(), other_es_id_to_name.end());
2348 
2349  // Merge other_mesh's elemset information with ours. Throw an
2350  // error if this and other_mesh have overlapping elemset codes
2351  // that refer to different elemset ids.
2352  std::vector<dof_id_type> this_elemset_codes = this->get_elemset_codes();
2353  MeshBase::elemset_type this_id_set_to_fill, other_id_set_to_fill;
2354  for (const auto & elemset_code : other_mesh->get_elemset_codes())
2355  {
2356  // Get the elemset ids for this elemset_code on other_mesh
2357  other_mesh->get_elemsets(elemset_code, other_id_set_to_fill);
2358 
2359  // Check that this elemset code does not already exist
2360  // in this mesh, or if it does, that it has the same elemset
2361  // ids associated with it.
2362  //
2363  // Note: get_elemset_codes() is guaranteed to return a
2364  // sorted vector, so we can binary search in it.
2365  auto it = Utility::binary_find(this_elemset_codes.begin(),
2366  this_elemset_codes.end(),
2367  elemset_code);
2368 
2369  if (it != this_elemset_codes.end())
2370  {
2371  // This mesh has the same elemset code. Does it refer to
2372  // the same elemset ids?
2373  this->get_elemsets(elemset_code, this_id_set_to_fill);
2374 
2375  // Throw an error if they don't match, otherwise we
2376  // don't need to do anything
2377  libmesh_error_msg_if(other_id_set_to_fill != this_id_set_to_fill,
2378  "Attempted to stitch together meshes with conflicting elemset codes.");
2379  }
2380  else
2381  {
2382  // Add other_mesh's elemset code to this mesh
2383  this->add_elemset_code(elemset_code, other_id_set_to_fill);
2384  }
2385  }
2386  } // end if (other_mesh)
2387 
2388  // Finally, we need to "merge" the overlapping nodes
2389  // We do this by iterating over node_to_elems_map and updating
2390  // the elements so that they "point" to the nodes that came
2391  // from this mesh, rather than from other_mesh.
2392  // Then we iterate over node_to_node_map and delete the
2393  // duplicate nodes that came from other_mesh.
2394 
2395  {
2396  LOG_SCOPE("stitch_meshes node updates", "UnstructuredMesh");
2397 
2398  // Container to catch boundary IDs passed back from BoundaryInfo.
2399  std::vector<boundary_id_type> bc_ids;
2400 
2401  for (const auto & [target_node_id, elem_vec] : node_to_elems_map)
2402  {
2403  dof_id_type other_node_id = node_to_node_map[target_node_id];
2404  Node & target_node = this->node_ref(target_node_id);
2405 
2406  std::size_t n_elems = elem_vec.size();
2407  for (std::size_t i=0; i<n_elems; i++)
2408  {
2409  dof_id_type elem_id = elem_vec[i];
2410  Elem * el = this->elem_ptr(elem_id);
2411 
2412  // find the local node index that we want to update
2413  unsigned int local_node_index = el->local_node(other_node_id);
2414  libmesh_assert_not_equal_to(local_node_index, libMesh::invalid_uint);
2415 
2416  // We also need to copy over the nodeset info here,
2417  // because the node will get deleted below
2418  this->get_boundary_info().boundary_ids(el->node_ptr(local_node_index), bc_ids);
2419  el->set_node(local_node_index) = &target_node;
2420  this->get_boundary_info().add_node(&target_node, bc_ids);
2421  }
2422  }
2423  }
2424 
2425  {
2426  LOG_SCOPE("stitch_meshes node deletion", "UnstructuredMesh");
2427  for (const auto & [other_node_id, this_node_id] : node_to_node_map)
2428  {
2429  // In the case that this==other_mesh, the two nodes might be the same (e.g. if
2430  // we're stitching a "sliver"), hence we need to skip node deletion in that case.
2431  if ((this == other_mesh) && (this_node_id == other_node_id))
2432  continue;
2433 
2434  this->delete_node( this->node_ptr(this_node_id) );
2435  }
2436  }
2437 
2438  // If find_neighbors() wasn't called in prepare_for_use(), we need to
2439  // manually loop once more over all elements adjacent to the stitched boundary
2440  // and fix their lists of neighbors.
2441  // This is done according to the following steps:
2442  // 1. Loop over all copied elements adjacent to the boundary using node_to_elems_map (trying to avoid duplicates)
2443  // 2. Look at all their sides with a nullptr neighbor and update them using side_to_elem_map if necessary
2444  // 3. Update the corresponding side in side_to_elem_map as well
2445  if (skip_find_neighbors)
2446  {
2447  LOG_SCOPE("stitch_meshes neighbor fixes", "UnstructuredMesh");
2448 
2449  // Pull objects out of the loop to reduce heap operations
2450  std::unique_ptr<const Elem> my_side, their_side;
2451 
2452  std::set<dof_id_type> fixed_elems;
2453  for (const auto & pr : node_to_elems_map)
2454  {
2455  std::size_t n_elems = pr.second.size();
2456  for (std::size_t i=0; i<n_elems; i++)
2457  {
2458  dof_id_type elem_id = pr.second[i];
2459  if (!fixed_elems.count(elem_id))
2460  {
2461  Elem * el = this->elem_ptr(elem_id);
2462  fixed_elems.insert(elem_id);
2463  for (auto s : el->side_index_range())
2464  {
2465  if (el->neighbor_ptr(s) == nullptr)
2466  {
2467  key_type key = el->low_order_key(s);
2468  auto bounds = side_to_elem_map.equal_range(key);
2469 
2470  if (bounds.first != bounds.second)
2471  {
2472  // Get the side for this element
2473  el->side_ptr(my_side, s);
2474 
2475  // Look at all the entries with an equivalent key
2476  while (bounds.first != bounds.second)
2477  {
2478  // Get the potential element
2479  Elem * neighbor = const_cast<Elem *>(bounds.first->second.first);
2480 
2481  // Get the side for the neighboring element
2482  const unsigned int ns = bounds.first->second.second;
2483  neighbor->side_ptr(their_side, ns);
2484  //libmesh_assert(my_side.get());
2485  //libmesh_assert(their_side.get());
2486 
2487  // If found a match with my side
2488  //
2489  // We need special tests here for 1D:
2490  // since parents and children have an equal
2491  // side (i.e. a node), we need to check
2492  // ns != ms, and we also check level() to
2493  // avoid setting our neighbor pointer to
2494  // any of our neighbor's descendants
2495  if ((*my_side == *their_side) &&
2496  (el->level() == neighbor->level()) &&
2497  ((el->dim() != 1) || (ns != s)))
2498  {
2499  // So share a side. Is this a mixed pair
2500  // of subactive and active/ancestor
2501  // elements?
2502  // If not, then we're neighbors.
2503  // If so, then the subactive's neighbor is
2504 
2505  if (el->subactive() ==
2506  neighbor->subactive())
2507  {
2508  // an element is only subactive if it has
2509  // been coarsened but not deleted
2510  el->set_neighbor (s,neighbor);
2511  neighbor->set_neighbor(ns,el);
2512  }
2513  else if (el->subactive())
2514  {
2515  el->set_neighbor(s,neighbor);
2516  }
2517  else if (neighbor->subactive())
2518  {
2519  neighbor->set_neighbor(ns,el);
2520  }
2521  // It's OK to invalidate the
2522  // bounds.first iterator here,
2523  // as we are immediately going
2524  // to break out of this while
2525  // loop. bounds.first will
2526  // therefore not be used for
2527  // anything else.
2528  side_to_elem_map.erase (bounds.first);
2529  break;
2530  }
2531 
2532  ++bounds.first;
2533  }
2534  }
2535  }
2536  }
2537  }
2538  }
2539  }
2540  }
2541 
2542  const bool old_allow_find_neighbors = this->allow_find_neighbors();
2543  const bool old_allow_remote_element_removal = this->allow_remote_element_removal();
2544  this->allow_find_neighbors(!skip_find_neighbors);
2545  this->allow_remote_element_removal(false);
2546  this->prepare_for_use();
2547  this->allow_find_neighbors(old_allow_find_neighbors);
2548  this->allow_remote_element_removal(old_allow_remote_element_removal);
2549 
2550  // After the stitching, we may want to clear boundary IDs from element
2551  // faces that are now internal to the mesh
2552  if (clear_stitched_boundary_ids)
2553  {
2554  LOG_SCOPE("stitch_meshes clear bcids", "UnstructuredMesh");
2555 
2557  this_mesh_boundary_id, other_mesh_boundary_id, /*clear_nodeset_data=*/true);
2558  }
2559 
2560  // Return the number of nodes which were merged.
2561  return node_to_node_map.size();
2562 }
2563 
2564 
2565 } // namespace libMesh
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
OStreamProxy err
IndexType begin() const
Definition: simple_range.h:41
void set_p_level(const unsigned int p)
Sets the value of the p-refinement level for the element.
unsigned char mapping_data() const
Definition: elem.h:3086
bool closed()
Checks that the library has been closed.
Definition: libmesh.C:285
ElemType
Defines an enum for geometric element types.
The SimpleRange templated class is intended to make it easy to construct ranges from pairs of iterato...
Definition: simple_range.h:36
This class supports simple reads and writes in any libMesh-supported format, by dispatching to one of...
Definition: namebased_io.h:44
unique_id_type & set_unique_id()
Definition: dof_object.h:858
virtual void reserve_nodes(const dof_id_type nn)=0
Reserves space for a known number of nodes.
std::size_t stitch_surfaces(boundary_id_type boundary_id_1, boundary_id_type boundary_id_2, Real tol=TOLERANCE, bool clear_stitched_boundary_ids=false, bool verbose=true, bool use_binary_search=true, bool enforce_all_nodes_match_on_boundaries=false, bool merge_boundary_nodes_all_or_nothing=false)
Similar to stitch_meshes, except that we stitch two adjacent surfaces within this mesh...
const Elem * parent() const
Definition: elem.h:2980
bool is_prepared() const
Definition: mesh_base.h:193
virtual void all_complete_order_range(const SimpleRange< element_iterator > &range) override
Converts a (conforming, non-refined) mesh with linear elements into a mesh with "complete" order elem...
virtual void copy_nodes_and_elements(const MeshBase &other_mesh, const bool skip_find_neighbors=false, dof_id_type element_id_offset=0, dof_id_type node_id_offset=0, unique_id_type unique_id_offset=0, std::unordered_map< subdomain_id_type, subdomain_id_type > *id_remapping=nullptr)
Deep copy of nodes and elements from another mesh object (used by subclass copy constructors and by m...
In parallel meshes where a ghost element has neighbors which do not exist on the local processor...
Definition: remote_elem.h:59
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2494
A Node is like a Point, but with more information.
Definition: node.h:52
bool is_ancestor_of(const Elem *descendant) const
Definition: elem.h:2959
static ElemType complete_order_equivalent_type(const ElemType et)
Definition: elem.C:3190
virtual unique_id_type parallel_max_unique_id() const =0
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:1355
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:293
void set_parent(Elem *p)
Sets the pointer to the element&#39;s parent.
Definition: elem.h:2996
std::vector< std::string > _elem_integer_names
The array of names for integer data associated with each element in the mesh.
Definition: mesh_base.h:1997
const Elem * interior_parent() const
Definition: elem.C:1057
std::vector< dof_id_type > get_elemset_codes() const
Return a vector of all elemset codes defined on the mesh.
Definition: mesh_base.C:427
virtual void read(const std::string &mesh_file) override
This method implements reading a mesh from a specified file.
Definition: namebased_io.C:69
void set_mapping_type(const ElemMappingType type)
Sets the value of the mapping type for the element.
Definition: elem.h:3078
void allow_renumbering(bool allow)
If false is passed in then this mesh will no longer be renumbered when being prepared for use...
Definition: mesh_base.h:1191
IntRange< unsigned short > side_index_range() const
Definition: elem.h:2623
void skip_partitioning(bool skip)
If true is passed in then nothing on this mesh will be (re)partitioned.
Definition: mesh_base.h:1249
The definition of the const_element_iterator struct.
Definition: mesh_base.h:2180
static constexpr Real TOLERANCE
unsigned int dim
static void set_node_processor_ids(MeshBase &mesh)
This function is called after partitioning to set the processor IDs for the nodes.
Definition: partitioner.C:852
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const =0
size_t kdtree_get_point_count() const
Must return the number of data points.
void prepare_for_use(const bool skip_renumber_nodes_and_elements, const bool skip_find_neighbors)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:741
virtual dof_id_type low_order_key(const unsigned int s) const =0
const boundary_id_type side_id
VectorOfNodesAdaptor(const std::vector< std::pair< Point, dof_id_type >> &nodes)
void add_elemset_code(dof_id_type code, MeshBase::elemset_type id_set)
Tabulate a user-defined "code" for elements which belong to the element sets specified in id_set...
Definition: mesh_base.C:380
const std::map< boundary_id_type, std::string > & get_sideset_name_map() const
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
MeshBase & mesh
void add_child(Elem *elem)
Adds a child pointer to the array of children of this element.
Definition: elem.C:1915
void set_refinement_flag(const RefinementState rflag)
Sets the value of the refinement flag for the element.
Definition: elem.h:3168
unique_id_type unique_id() const
Definition: dof_object.h:844
const Parallel::Communicator & comm() const
virtual unsigned int n_children() const =0
void boundary_ids(const Node *node, std::vector< boundary_id_type > &vec_to_fill) const
Fills a user-provided std::vector with the boundary ids associated with Node node.
std::map< const Elem *, const CouplingMatrix *, CompareDofObjectsByPIDAndThenID > map_type
What elements do we care about and what variables do we care about on each element?
virtual void own_node(Node &)
Takes ownership of node n on this partition of a distributed mesh, by setting n.processor_id() to thi...
Definition: mesh_base.h:735
virtual void all_first_order() override
Converts a mesh with higher-order elements into a mesh with linear elements.
virtual void set_next_unique_id(unique_id_type id)=0
Sets the next available unique id to be used.
virtual void find_neighbors(const bool reset_remote_elements=false, const bool reset_current_list=true) override
Other functions from MeshBase requiring re-definition.
The libMesh namespace provides an interface to certain functionality in the library.
void create_pid_mesh(UnstructuredMesh &pid_mesh, const processor_id_type pid) const
Generates a new mesh containing all the elements which are assigned to processor pid.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:160
std::vector< std::string > _node_integer_names
The array of names for integer data associated with each node in the mesh.
Definition: mesh_base.h:2009
virtual void all_second_order_range(const SimpleRange< element_iterator > &range, const bool full_ordered=true) override
Converts a (conforming, non-refined) mesh with linear elements into a mesh with second-order elements...
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
void set_interior_parent(Elem *p)
Sets the pointer to the element&#39;s interior_parent.
Definition: elem.C:1119
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:2296
dof_id_type & set_id()
Definition: dof_object.h:836
void build_side_list(std::vector< dof_id_type > &element_id_list, std::vector< unsigned short int > &side_list, std::vector< boundary_id_type > &bc_id_list) const
Creates a list of element numbers, sides, and ids for those sides.
bool ancestor() const
Definition: elem.C:1881
void get_elemsets(dof_id_type elemset_code, MeshBase::elemset_type &id_set_to_fill) const
Look up the element sets for a given elemset code and vice-versa.
Definition: mesh_base.C:411
virtual bool contract() override
Delete subactive (i.e.
std::size_t stitching_helper(const MeshBase *other_mesh, boundary_id_type boundary_id_1, boundary_id_type boundary_id_2, Real tol, bool clear_stitched_boundary_ids, bool verbose, bool use_binary_search, bool enforce_all_nodes_match_on_boundaries, bool skip_find_neighbors, bool merge_boundary_nodes_all_or_nothing, bool remap_subdomain_ids)
Helper function for stitch_meshes and stitch_surfaces that does the mesh stitching.
uint8_t processor_id_type
std::map< boundary_id_type, std::string > & set_sideset_name_map()
void replace_child(Elem *elem, unsigned int c)
Replaces the child pointer at the specified index in the child array.
Definition: elem.C:1961
void allow_remote_element_removal(bool allow)
If false is passed in then this mesh will no longer have remote elements deleted when being prepared ...
Definition: mesh_base.h:1207
processor_id_type n_processors() const
virtual bool is_serial() const
Definition: mesh_base.h:206
void libmesh_ignore(const Args &...)
void build_node_list(std::vector< dof_id_type > &node_id_list, std::vector< boundary_id_type > &bc_id_list) const
Creates a list of nodes and ids for those nodes.
void add_node(const Node *node, const boundary_id_type id)
Add Node node with boundary id id to the boundary information data structures.
const std::map< boundary_id_type, std::string > & get_nodeset_name_map() const
int8_t boundary_id_type
Definition: id_types.h:51
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:3070
This is the MeshCommunication class.
dof_id_type id() const
Definition: dof_object.h:828
void min(const T &r, T &o, Request &req) const
const std::map< subdomain_id_type, std::string > & get_subdomain_name_map() const
Definition: mesh_base.h:1689
virtual Real hmin() const
Definition: elem.C:573
std::size_t stitch_meshes(const MeshBase &other_mesh, boundary_id_type this_mesh_boundary, boundary_id_type other_mesh_boundary, Real tol=TOLERANCE, bool clear_stitched_boundary_ids=false, bool verbose=true, bool use_binary_search=true, bool enforce_all_nodes_match_on_boundaries=false, bool merge_boundary_nodes_all_or_nothing=false, bool remap_subdomain_ids=false)
Stitch other_mesh to this mesh so that this mesh is the union of the two meshes.
virtual unsigned int n_nodes() const =0
static const 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
The UnstructuredMesh class is derived from the MeshBase class.
unsigned int which_neighbor_am_i(const Elem *e) const
This function tells you which neighbor e is.
Definition: elem.h:2869
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:321
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:802
virtual const Node * query_node_ptr(const dof_id_type i) const =0
void set_mapping_data(const unsigned char data)
Sets the value of the mapping data for the element.
Definition: elem.h:3094
const std::map< boundary_id_type, std::string > & get_edgeset_name_map() const
virtual dof_id_type max_elem_id() const =0
void clear_point_locator()
Releases the current PointLocator object.
Definition: mesh_base.C:1653
bool allow_find_neighbors() const
Definition: mesh_base.h:1199
virtual void write(const std::string &mesh_file) override
This method implements writing a mesh to a specified file.
Definition: namebased_io.C:303
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:482
virtual void delete_node(Node *n)=0
Removes the Node n from the mesh.
std::pair< std::vector< unsigned int >, std::vector< unsigned int > > merge_extra_integer_names(const MeshBase &other)
Merge extra-integer arrays from an other mesh.
Definition: mesh_base.C:1925
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:1335
std::string & subdomain_name(subdomain_id_type id)
Definition: mesh_base.C:1674
void allow_find_neighbors(bool allow)
If false is passed then this mesh will no longer work to find element neighbors when being prepared f...
Definition: mesh_base.h:1198
virtual std::unique_ptr< Elem > disconnected_clone() const
Definition: elem.C:287
virtual unsigned int n_second_order_adjacent_vertices(const unsigned int n) const
Definition: elem.C:2911
auto norm(const T &a) -> decltype(std::abs(a))
Definition: tensor_tools.h:74
bool allow_remote_element_removal() const
Definition: mesh_base.h:1208
void set_neighbor(const unsigned int i, Elem *n)
Assigns n as the neighbor.
Definition: elem.h:2540
void build_shellface_list(std::vector< dof_id_type > &element_id_list, std::vector< unsigned short int > &shellface_list, std::vector< boundary_id_type > &bc_id_list) const
Creates a list of element numbers, shellfaces, and boundary ids for those shellfaces.
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &) override
This method implements writing a mesh with nodal data to a specified file where the nodal data and va...
Definition: namebased_io.C:464
void regenerate_id_sets()
Clears and regenerates the cached sets of ids.
void prepare_for_use()
Definition: mesh_base.C:776
virtual void clear()
Deletes all the element and node data that is currently stored.
Definition: mesh_base.C:902
static ElemType second_order_equivalent_type(const ElemType et, const bool full_ordered=true)
Definition: elem.C:3006
unsigned int which_child_am_i(const Elem *e) const
Definition: elem.h:3142
bool skip_partitioning() const
Definition: mesh_base.h:1251
unsigned int n_partitions() const
Definition: mesh_base.h:1340
virtual const Elem * elem_ptr(const dof_id_type i) const =0
UnstructuredMesh(const Parallel::Communicator &comm_in, unsigned char dim=1)
Constructor.
Real kdtree_get_pt(const size_t idx, int dim) const
virtual Elem * insert_elem(Elem *e)=0
Insert elem e to the element array, preserving its id and replacing/deleting any existing element wit...
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2520
std::set< GhostingFunctor * >::const_iterator ghosting_functors_begin() const
Beginning of range of ghosting functors.
Definition: mesh_base.h:1286
std::set< elemset_id_type > elemset_type
Typedef for the "set" container used to store elemset ids.
Definition: mesh_base.h:313
unsigned int level() const
Definition: elem.h:3024
virtual unsigned int n_vertices() const =0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< std::pair< Point, dof_id_type > > _nodes
virtual std::unique_ptr< Elem > side_ptr(unsigned int i)=0
virtual ~UnstructuredMesh()
Destructor.
subdomain_id_type subdomain_id() const
Definition: elem.h:2504
timpi_pure bool verify(const T &r) const
void copy_boundary_ids(const BoundaryInfo &old_boundary_info, const Elem *const old_elem, const Elem *const new_elem)
void max(const T &r, T &o, Request &req) const
virtual unsigned short dim() const =0
IndexType end() const
Definition: simple_range.h:43
Temporarily serialize a DistributedMesh for non-distributed-mesh capable code paths.
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2444
std::map< boundary_id_type, std::string > & set_nodeset_name_map()
void build_edge_list(std::vector< dof_id_type > &element_id_list, std::vector< unsigned short int > &edge_list, std::vector< boundary_id_type > &bc_id_list) const
Creates a list of element numbers, edges, and boundary ids for those edges.
void clear_stitched_boundary_side_ids(boundary_id_type sideset_id, boundary_id_type other_sideset_id, bool clear_nodeset_data=false)
Clear sideset information along a stitched mesh interface.
void make_nodes_parallel_consistent(MeshBase &)
Copy processor_ids and ids on ghost nodes from their local processors.
OStreamProxy out
virtual bool is_replicated() const
Definition: mesh_base.h:228
void add_side(const dof_id_type elem, const unsigned short int side, const boundary_id_type id)
Add side side of element number elem with boundary id id to the boundary information data structure...
const std::string & get_sideset_name(boundary_id_type id) const
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const
Definition: elem.C:2919
void add_shellface(const dof_id_type elem, const unsigned short int shellface, const boundary_id_type id)
Add shell face shellface of element number elem with boundary id id to the boundary information data ...
bool subactive() const
Definition: elem.h:2909
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
void add_extra_integers(const unsigned int n_integers)
Assigns a set of extra integers to this DofObject.
Definition: dof_object.C:490
ForwardIterator binary_find(ForwardIterator first, ForwardIterator last, const T &value)
The STL provides std::binary_search() which returns true or false depending on whether the searched-f...
Definition: utility.h:265
unsigned int mesh_dimension() const
Definition: mesh_base.C:354
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:278
unsigned int n_extra_integers() const
Returns how many extra integers are associated to the DofObject.
Definition: dof_object.h:1170
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:591
virtual const Point & point(const dof_id_type i) const =0
unsigned int local_node(const dof_id_type i) const
Definition: elem.h:2424
void edge_boundary_ids(const Elem *const elem, const unsigned short int edge, std::vector< boundary_id_type > &vec_to_fill) const
SimpleRange< NeighborPtrIter > neighbor_ptr_range()
Returns a range with all neighbors of an element, usable in range-based for loops.
Definition: elem.h:3439
void set_p_refinement_flag(const RefinementState pflag)
Sets the value of the p-refinement flag for the element.
Definition: elem.h:3184
bool allow_renumbering() const
Definition: mesh_base.h:1192
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:2199
virtual void delete_remote_elements()
When supported, deletes all nonlocal elements of the mesh except for "ghosts" which touch a local ele...
Definition: mesh_base.h:248
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false) override
Reads the file specified by name.
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
processor_id_type processor_id() const
bool active() const
Definition: elem.h:2891
bool kdtree_get_bbox(BBOX &) const
static ElemType first_order_equivalent_type(const ElemType et)
Definition: elem.C:2937
processor_id_type processor_id() const
Definition: dof_object.h:905
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:2412
virtual void reserve_elem(const dof_id_type ne)=0
Reserves space for a known number of elements.
uint8_t unique_id_type
Definition: id_types.h:86
bool has_children() const
Definition: elem.h:2929
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 set_extra_integer(const unsigned int index, const dof_id_type value)
Sets the value on this object of the extra integer associated with index, which should have been obta...
Definition: dof_object.h:1086
virtual void renumber_nodes_and_elements()=0
After partitioning a mesh it is useful to renumber the nodes and elements so that they lie in contigu...
virtual dof_id_type n_nodes() const =0
This class implements the original default geometry ghosting requirements in libMesh: point neighbors...
void create_submesh(UnstructuredMesh &new_mesh, const const_element_iterator &it, const const_element_iterator &it_end) const
Constructs a mesh called "new_mesh" from the current mesh by iterating over the elements between it a...
Point vertex_average() const
Definition: elem.C:559
virtual void write(const std::string &name) override
Write the file specified by name.
dof_id_type get_extra_integer(const unsigned int index) const
Gets the value on this object of the extra integer associated with index, which should have been obta...
Definition: dof_object.h:1102
const Elem * child_ptr(unsigned int i) const
Definition: elem.h:3113
std::set< GhostingFunctor * >::const_iterator ghosting_functors_end() const
End of range of ghosting functors.
Definition: mesh_base.h:1292
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
A useful inline function which replaces the macros used previously.
uint8_t dof_id_type
Definition: id_types.h:67
std::map< boundary_id_type, std::string > & set_edgeset_name_map()
void add_edge(const dof_id_type elem, const unsigned short int edge, const boundary_id_type id)
Add edge edge of element number elem with boundary id id to the boundary information data structure...
const RemoteElem * remote_elem
Definition: remote_elem.C:57