Line data Source code
1 : // The libMesh Finite Element Library.
2 : // Copyright (C) 2002-2026 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/elem_range.h"
28 : #include "libmesh/mesh_tools.h" // For n_levels
29 : #include "libmesh/parallel.h"
30 : #include "libmesh/remote_elem.h"
31 : #include "libmesh/namebased_io.h"
32 : #include "libmesh/partitioner.h"
33 : #include "libmesh/enum_order.h"
34 : #include "libmesh/mesh_communication.h"
35 : #include "libmesh/enum_to_string.h"
36 : #include "libmesh/mesh_serializer.h"
37 : #include "libmesh/utility.h"
38 :
39 : #ifdef LIBMESH_HAVE_NANOFLANN
40 : #include "libmesh/nanoflann.hpp"
41 : #endif
42 :
43 : // C++ includes
44 : #include <algorithm> // std::all_of
45 : #include <atomic>
46 : #include <fstream>
47 : #include <iomanip>
48 : #include <map>
49 : #include <sstream>
50 : #include <unordered_map>
51 :
52 : // for disjoint neighbors
53 : #include "libmesh/periodic_boundaries.h"
54 : #include "libmesh/periodic_boundary.h"
55 :
56 : namespace {
57 :
58 : using namespace libMesh;
59 :
60 : // Helper functions for all_second_order, all_complete_order
61 :
62 : std::map<std::vector<dof_id_type>, Node *>::iterator
63 33790129 : map_hi_order_node(unsigned int hon,
64 : const Elem & hi_elem,
65 : std::map<std::vector<dof_id_type>, Node *> & adj_vertices_to_ho_nodes)
66 : {
67 : /*
68 : * form a vector that will hold the node id's of
69 : * the vertices that are adjacent to the nth
70 : * higher-order node.
71 : */
72 : const unsigned int n_adjacent_vertices =
73 33790129 : hi_elem.n_second_order_adjacent_vertices(hon);
74 :
75 33790129 : std::vector<dof_id_type> adjacent_vertices_ids(n_adjacent_vertices);
76 :
77 113280998 : for (unsigned int v=0; v<n_adjacent_vertices; v++)
78 81535323 : adjacent_vertices_ids[v] =
79 79490869 : hi_elem.node_id( hi_elem.second_order_adjacent_vertex(hon,v) );
80 :
81 : /*
82 : * \p adjacent_vertices_ids is now in order of the current
83 : * side. sort it, so that comparisons with the
84 : * \p adjacent_vertices_ids created through other elements'
85 : * sides can match
86 : */
87 33790129 : std::sort(adjacent_vertices_ids.begin(),
88 : adjacent_vertices_ids.end());
89 :
90 : // Does this set of vertices already have a mid-node added? If not
91 : // we'll want to add it.
92 35522517 : return adj_vertices_to_ho_nodes.try_emplace(adjacent_vertices_ids, nullptr).first;
93 : }
94 :
95 3625402 : void transfer_elem(Elem & lo_elem,
96 : std::unique_ptr<Elem> hi_elem,
97 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
98 : unique_id_type max_unique_id,
99 : unique_id_type max_new_nodes_per_elem,
100 : #endif
101 : UnstructuredMesh & mesh,
102 : std::map<std::vector<dof_id_type>, Node *> & adj_vertices_to_ho_nodes,
103 : std::unordered_map<Elem *, std::vector<Elem *>> & exterior_children_of)
104 : {
105 92944 : libmesh_assert_equal_to (lo_elem.n_vertices(), hi_elem->n_vertices());
106 :
107 185888 : const processor_id_type my_pid = mesh.processor_id();
108 3625402 : const processor_id_type lo_pid = lo_elem.processor_id();
109 :
110 : /*
111 : * Now handle the additional higher-order nodes. This
112 : * is simply handled through a map that remembers
113 : * the already-added nodes. This map maps the global
114 : * ids of the vertices (that uniquely define this
115 : * higher-order node) to the new node.
116 : * Notation: hon = high-order node
117 : */
118 3625402 : const unsigned int hon_begin = lo_elem.n_nodes();
119 3625402 : const unsigned int hon_end = hi_elem->n_nodes();
120 :
121 37129968 : for (unsigned int hon=hon_begin; hon<hon_end; hon++)
122 : {
123 33504566 : auto pos = map_hi_order_node(hon, *hi_elem, adj_vertices_to_ho_nodes);
124 :
125 : // no, not added yet
126 33504566 : if (!pos->second)
127 : {
128 305940 : const auto & adjacent_vertices_ids = pos->first;
129 :
130 : /*
131 : * for this set of vertices, there is no
132 : * second_order node yet. Add it.
133 : *
134 : * compute the location of the new node as
135 : * the average over the adjacent vertices.
136 : */
137 305940 : Point new_location = 0;
138 42355866 : for (dof_id_type vertex_id : adjacent_vertices_ids)
139 30419800 : new_location += mesh.point(vertex_id);
140 :
141 12242006 : new_location /= static_cast<Real>(adjacent_vertices_ids.size());
142 :
143 : /* Add the new point to the mesh.
144 : *
145 : * If we are on a serialized mesh, then we're doing this
146 : * all in sync, and the node processor_id will be
147 : * consistent between processors.
148 : *
149 : * If we are on a distributed mesh, we can fix
150 : * inconsistent processor ids later, but only if every
151 : * processor gives new nodes a *locally* consistent
152 : * processor id, so we'll give the new node the
153 : * processor id of an adjacent element for now and then
154 : * we'll update that later if appropriate.
155 : */
156 : Node * hi_node = mesh.add_point
157 11936066 : (new_location, DofObject::invalid_id, lo_pid);
158 :
159 : /* Come up with a unique unique_id for a potentially new
160 : * node. On a distributed mesh we don't yet know what
161 : * processor_id will definitely own it, so we can't let
162 : * the pid determine the unique_id. But we're not
163 : * adding unpartitioned nodes in sync, so we can't let
164 : * the mesh autodetermine a unique_id for a new
165 : * unpartitioned node either. So we have to pick unique
166 : * unique_id values manually.
167 : *
168 : * We don't have to pick the *same* unique_id value as
169 : * will be picked on other processors, though; we'll
170 : * sync up each node later. We just need to make sure
171 : * we don't duplicate any unique_id that might be chosen
172 : * by the same process elsewhere.
173 : */
174 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
175 11936066 : unique_id_type new_unique_id = max_unique_id +
176 12242006 : max_new_nodes_per_elem * lo_elem.id() +
177 11936066 : hon - hon_begin;
178 :
179 305940 : hi_node->set_unique_id(new_unique_id);
180 : #endif
181 :
182 : /*
183 : * insert the new node with its defining vertex
184 : * set into the map, and relocate pos to this
185 : * new entry, so that the hi_elem can use
186 : * \p pos for inserting the node
187 : */
188 11936066 : pos->second = hi_node;
189 :
190 11936066 : hi_elem->set_node(hon, hi_node);
191 : }
192 : // yes, already added.
193 : else
194 : {
195 546492 : Node * hi_node = pos->second;
196 546492 : libmesh_assert(hi_node);
197 546492 : libmesh_assert_equal_to(mesh.node_ptr(hi_node->id()), hi_node);
198 :
199 21568500 : hi_elem->set_node(hon, hi_node);
200 :
201 : // We need to ensure that the processor who should own a
202 : // node *knows* they own the node. And because
203 : // Node::choose_processor_id() may depend on Node id,
204 : // which may not yet be authoritative, we still have to
205 : // use a dumb-but-id-independent partitioning heuristic.
206 : processor_id_type chosen_pid =
207 21568500 : std::min (hi_node->processor_id(), lo_pid);
208 :
209 : // Plus, if we just discovered that we own this node,
210 : // then on a distributed mesh we need to make sure to
211 : // give it a valid id, not just a placeholder id!
212 21568500 : if (!mesh.is_replicated() &&
213 21568500 : hi_node->processor_id() != my_pid &&
214 : chosen_pid == my_pid)
215 4130 : mesh.own_node(*hi_node);
216 :
217 21568500 : hi_node->processor_id() = chosen_pid;
218 : }
219 : }
220 :
221 : /*
222 : * find_neighbors relies on remote_elem neighbor links being
223 : * properly maintained. Our own code here relies on ordinary
224 : * neighbor links being properly maintained, so let's just keep
225 : * everything up to date.
226 : */
227 18500576 : for (auto s : lo_elem.side_index_range())
228 : {
229 14875174 : Elem * neigh = lo_elem.neighbor_ptr(s);
230 14875174 : if (!neigh)
231 14167291 : continue;
232 :
233 342207 : if (neigh != remote_elem)
234 : {
235 : // We don't support AMR even outside our own range yet.
236 18422 : libmesh_assert_equal_to (neigh->level(), 0);
237 :
238 327168 : const unsigned int ns = neigh->which_neighbor_am_i(&lo_elem);
239 18422 : libmesh_assert_not_equal_to(ns, libMesh::invalid_uint);
240 :
241 36844 : neigh->set_neighbor(ns, hi_elem.get());
242 : }
243 :
244 37164 : hi_elem->set_neighbor(s, neigh);
245 : }
246 :
247 : /**
248 : * If the old element has an interior_parent(), transfer it to the
249 : * new element ... and if the interior_parent itself might be
250 : * getting upgraded, make sure we later consider the new element to
251 : * be its exterior child, not the old element.
252 : */
253 3625402 : Elem * interior_p = lo_elem.interior_parent();
254 3625402 : if (interior_p)
255 0 : hi_elem->set_interior_parent(interior_p);
256 :
257 3625402 : if (auto parent_exterior_it = exterior_children_of.find(interior_p);
258 92944 : parent_exterior_it != exterior_children_of.end())
259 : {
260 0 : auto & exteriors = parent_exterior_it->second;
261 0 : for (std::size_t i : index_range(exteriors))
262 0 : if (exteriors[i] == &lo_elem)
263 : {
264 0 : exteriors[i] = hi_elem.get();
265 0 : break;
266 : }
267 : }
268 :
269 : /**
270 : * If we had interior_parent() links to the old element, transfer
271 : * them to the new element.
272 : */
273 3718346 : if (auto exterior_it = exterior_children_of.find(&lo_elem);
274 92944 : exterior_it != exterior_children_of.end())
275 : {
276 3625402 : for (Elem * exterior_elem : exterior_it->second)
277 : {
278 0 : libmesh_assert(exterior_elem->interior_parent() == &lo_elem);
279 0 : exterior_elem->set_interior_parent(hi_elem.get());
280 : }
281 : }
282 :
283 : /**
284 : * If the old element had any boundary conditions they
285 : * should be transferred to the second-order element. The old
286 : * boundary conditions will be removed from the BoundaryInfo
287 : * data structure by insert_elem.
288 : *
289 : * Also, prepare_for_use() will reconstruct most of our neighbor
290 : * links, but if we have any remote_elem links in a distributed
291 : * mesh, they need to be preserved. We do that in the same loop
292 : * here.
293 : */
294 92944 : mesh.get_boundary_info().copy_boundary_ids
295 3625402 : (mesh.get_boundary_info(), &lo_elem, hi_elem.get());
296 :
297 : /*
298 : * The new second-order element is ready.
299 : * Inserting it into the mesh will replace and delete
300 : * the first-order element.
301 : */
302 185888 : hi_elem->set_id(lo_elem.id());
303 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
304 185888 : hi_elem->set_unique_id(lo_elem.unique_id());
305 : #endif
306 :
307 3625402 : const unsigned int nei = lo_elem.n_extra_integers();
308 3625402 : hi_elem->add_extra_integers(nei);
309 3625721 : for (unsigned int i=0; i != nei; ++i)
310 319 : hi_elem->set_extra_integer(i, lo_elem.get_extra_integer(i));
311 :
312 3532458 : hi_elem->inherit_data_from(lo_elem);
313 :
314 3718346 : mesh.insert_elem(std::move(hi_elem));
315 3625402 : }
316 :
317 :
318 : template <typename ElemTypeConverter>
319 : void
320 41430 : all_increased_order_range (UnstructuredMesh & mesh,
321 : const SimpleRange<MeshBase::element_iterator> & range,
322 : const unsigned int max_new_nodes_per_elem,
323 : const ElemTypeConverter & elem_type_converter)
324 : {
325 : // This function must be run on all processors at once
326 1182 : timpi_parallel_only(mesh.comm());
327 :
328 : /*
329 : * The maximum number of new higher-order nodes we might be adding,
330 : * for use when picking unique unique_id values later. This variable
331 : * is not used unless unique ids are enabled, so libmesh_ignore() it
332 : * to avoid warnings in that case.
333 : */
334 1182 : libmesh_ignore(max_new_nodes_per_elem);
335 :
336 : /*
337 : * The mesh should at least be consistent enough for us to add new
338 : * nodes consistently.
339 : */
340 1182 : libmesh_assert(mesh.comm().verify(mesh.n_elem()));
341 1182 : libmesh_assert(mesh.comm().verify(mesh.max_elem_id()));
342 :
343 : /*
344 : * If the mesh is empty then we have nothing to do
345 : */
346 41430 : if (!mesh.n_elem())
347 3984 : return;
348 :
349 : // If every element in the range _on every proc_ is already of the
350 : // requested higher order then we have nothing to do. However, if
351 : // any proc has some lower-order elements in the range, then _all_
352 : // processors need to continue this function because it is
353 : // parallel_only().
354 : //
355 : // Note: std::all_of() returns true for an empty range, which can
356 : // happen for example in the DistributedMesh case when there are
357 : // more processors than elements. In the case of an empty range we
358 : // therefore set already_second_order to true on that proc.
359 101568 : auto is_higher_order = [&elem_type_converter](const Elem * elem) {
360 39669 : ElemType old_type = elem->type();
361 2048 : ElemType new_type = elem_type_converter(old_type);
362 39669 : return old_type == new_type;
363 : };
364 :
365 41430 : bool already_higher_order =
366 81678 : std::all_of(range.begin(), range.end(), is_higher_order);
367 :
368 : // Check with other processors and possibly return early
369 41430 : mesh.comm().min(already_higher_order);
370 41430 : if (already_higher_order)
371 116 : return;
372 :
373 : /*
374 : * this map helps in identifying higher order
375 : * nodes. Namely, a higher-order node:
376 : * - edge node
377 : * - face node
378 : * - bubble node
379 : * is uniquely defined through a set of adjacent
380 : * vertices. This set of adjacent vertices is
381 : * used to identify already added higher-order
382 : * nodes. We are safe to use node id's since we
383 : * make sure that these are correctly numbered.
384 : *
385 : * We lazily use an ordered map here to avoid having to implement a
386 : * good hash for vector<dof_id_type>
387 : */
388 2132 : std::map<std::vector<dof_id_type>, Node *> adj_vertices_to_ho_nodes;
389 :
390 : /*
391 : * This map helps us reset any interior_parent() values from the
392 : * lower order element to its higher order replacement. Unlike with
393 : * neighbor pointers, we don't have backlinks here, so we have to
394 : * iterate over the mesh to track forward links.
395 : */
396 2132 : std::unordered_map<Elem *, std::vector<Elem *>> exterior_children_of;
397 :
398 : /*
399 : * max_new_nodes_per_elem is the maximum number of new higher order
400 : * nodes we might be adding, for use when picking unique unique_id
401 : * values later. This variable is not used unless unique ids are
402 : * enabled.
403 : */
404 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
405 37446 : unique_id_type max_unique_id = mesh.parallel_max_unique_id();
406 : #endif
407 :
408 : /**
409 : * On distributed meshes we currently only support unpartitioned
410 : * meshes (where we'll add every node in sync) or
411 : * completely-partitioned meshes (where we'll sync nodes later);
412 : * let's keep track to make sure we're not in any in-between state.
413 : */
414 37446 : dof_id_type n_unpartitioned_elem = 0,
415 37446 : n_partitioned_elem = 0;
416 :
417 : /**
418 : * Loop over the elements in the given range. If any are
419 : * already at higher than first-order, track their higher-order
420 : * nodes in case we need them for neighboring elements later.
421 : *
422 : * In this way we can use this method to "fix up" a mesh which has
423 : * otherwise inconsistent neighbor pairs of lower and higher order
424 : * geometric elements.
425 : *
426 : * If any elements are not at the desired order yet, we need to
427 : * check their neighbors and even their edge neighbors for higher
428 : * order; we may need to share elements with a neighbor not in the
429 : * range.
430 : */
431 7249428 : auto track_if_necessary = [&adj_vertices_to_ho_nodes,
432 : &exterior_children_of,
433 189752 : &elem_type_converter](Elem * elem) {
434 3662755 : if (elem && elem != remote_elem)
435 : {
436 3662755 : if (elem->default_order() != FIRST)
437 305187 : for (unsigned int hon : make_range(elem->n_vertices(), elem->n_nodes()))
438 : {
439 285563 : auto pos = map_hi_order_node(hon, *elem, adj_vertices_to_ho_nodes);
440 299325 : pos->second = elem->node_ptr(hon);
441 : }
442 :
443 3662755 : const ElemType old_type = elem->type();
444 127378 : const ElemType new_type = elem_type_converter(old_type);
445 3662755 : if (old_type != new_type)
446 3643131 : exterior_children_of.emplace(elem, std::vector<Elem *>());
447 : }
448 : };
449 :
450 : // If we're in the common case then just track everything; otherwise
451 : // find point neighbors to track
452 113180 : if (range.begin() == mesh.elements_begin() &&
453 137106 : range.end() == mesh.elements_end())
454 : {
455 7088814 : for (auto & elem : range)
456 3620410 : track_if_necessary(elem);
457 : }
458 : else
459 : {
460 240 : GhostingFunctor::map_type point_neighbor_elements;
461 :
462 240 : GhostPointNeighbors point_neighbor_finder(mesh);
463 11868 : point_neighbor_finder(range.begin(), range.end(),
464 : mesh.n_processors(),
465 : point_neighbor_elements);
466 :
467 46381 : for (auto & [elem, coupling_map] : point_neighbor_elements)
468 : {
469 2168 : libmesh_ignore(coupling_map);
470 42345 : track_if_necessary(const_cast<Elem *>(elem));
471 : }
472 : }
473 :
474 : /**
475 : * Loop over all mesh elements to look for interior_parent links we
476 : * need to upgrade later.
477 : */
478 7341408 : for (auto & elem : mesh.element_ptr_range())
479 3810122 : if (auto exterior_map_it = exterior_children_of.find(elem->interior_parent());
480 96994 : exterior_map_it != exterior_children_of.end())
481 0 : exterior_map_it->second.push_back(elem);
482 :
483 : /**
484 : * Loop over the low-ordered elements in the _elements vector.
485 : * First make sure they _are_ indeed low-order, and then replace
486 : * them with an equivalent second-order element. Don't
487 : * forget to delete the low-order element, or else it will leak!
488 : */
489 7103300 : for (auto & lo_elem : range)
490 : {
491 : // Now we can skip the elements in the range that are already
492 : // higher-order.
493 3625893 : const ElemType old_type = lo_elem->type();
494 125006 : const ElemType new_type = elem_type_converter(old_type);
495 :
496 3625893 : if (old_type == new_type)
497 491 : continue;
498 :
499 : // this does _not_ work for refined elements
500 92944 : libmesh_assert_equal_to (lo_elem->level(), 0);
501 :
502 3625402 : if (lo_elem->processor_id() == DofObject::invalid_processor_id)
503 3554790 : ++n_unpartitioned_elem;
504 : else
505 70612 : ++n_partitioned_elem;
506 :
507 : /*
508 : * Build the higher-order equivalent; add to
509 : * the new_elements list.
510 : */
511 3625402 : auto ho_elem = Elem::build (new_type);
512 :
513 92944 : libmesh_assert_equal_to (lo_elem->n_vertices(), ho_elem->n_vertices());
514 :
515 : /*
516 : * By definition the initial nodes of the lower and higher order
517 : * element are identically numbered. Transfer these.
518 : */
519 18597808 : for (unsigned int v=0, lnn=lo_elem->n_nodes(); v < lnn; v++)
520 15361716 : ho_elem->set_node(v, lo_elem->node_ptr(v));
521 :
522 3811290 : transfer_elem(*lo_elem, std::move(ho_elem),
523 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
524 : max_unique_id, max_new_nodes_per_elem,
525 : #endif
526 : mesh, adj_vertices_to_ho_nodes,
527 : exterior_children_of);
528 : } // end for (auto & lo_elem : range)
529 :
530 : // we can clear the map at this point.
531 1066 : adj_vertices_to_ho_nodes.clear();
532 :
533 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
534 37446 : const unique_id_type new_max_unique_id = max_unique_id +
535 37446 : max_new_nodes_per_elem * mesh.n_elem();
536 37446 : mesh.set_next_unique_id(new_max_unique_id);
537 : #endif
538 :
539 : // On a DistributedMesh our ghost node processor ids may be bad,
540 : // the ids of nodes touching remote elements may be inconsistent,
541 : // unique_ids of newly added non-local nodes remain unset, and our
542 : // partitioning of new nodes may not be well balanced.
543 : //
544 : // make_nodes_parallel_consistent() will fix all this.
545 37446 : if (!mesh.is_replicated())
546 : {
547 32354 : dof_id_type max_unpartitioned_elem = n_unpartitioned_elem;
548 32354 : mesh.comm().max(max_unpartitioned_elem);
549 32354 : if (max_unpartitioned_elem)
550 : {
551 : // We'd better be effectively serialized here. In theory we
552 : // could support more complicated cases but in practice we
553 : // only support "completely partitioned" and/or "serialized"
554 50948 : if (!mesh.comm().verify(n_unpartitioned_elem) ||
555 50970 : !mesh.comm().verify(n_partitioned_elem) ||
556 25485 : !mesh.is_serial())
557 0 : libmesh_not_implemented();
558 : }
559 : else
560 : {
561 6869 : MeshCommunication().make_nodes_parallel_consistent (mesh);
562 : }
563 : }
564 :
565 : // renumber nodes, repartition nodes, etc. We may no longer need a
566 : // find_neighbors() here since we're keeping neighbor links intact
567 : // ourselves, *except* that if we're not already prepared we may
568 : // have user code that was expecting this call to prepare neighbors.
569 2132 : const bool old_find_neighbors = mesh.allow_find_neighbors();
570 37446 : if (mesh.is_prepared())
571 244 : mesh.allow_find_neighbors(false);
572 37446 : mesh.prepare_for_use();
573 1066 : mesh.allow_find_neighbors(old_find_neighbors);
574 : }
575 :
576 :
577 : } // anonymous namespace
578 :
579 :
580 : namespace libMesh
581 : {
582 :
583 : // This class adapts a vector of Nodes (represented by a pair of a Point and a dof_id_type)
584 : // for use in a nanoflann KD-Tree
585 :
586 0 : class VectorOfNodesAdaptor
587 : {
588 : private:
589 : const std::vector<std::pair<Point, dof_id_type>> _nodes;
590 :
591 : public:
592 0 : VectorOfNodesAdaptor(const std::vector<std::pair<Point, dof_id_type>> & nodes) :
593 0 : _nodes(nodes)
594 0 : {}
595 :
596 : /**
597 : * Must return the number of data points
598 : */
599 0 : inline size_t kdtree_get_point_count() const { return _nodes.size(); }
600 :
601 : /**
602 : * \returns The dim'th component of the idx'th point in the class:
603 : * Since this is inlined and the "dim" argument is typically an immediate value, the
604 : * "if's" are actually solved at compile time.
605 : */
606 0 : inline Real kdtree_get_pt(const size_t idx, int dim) const
607 : {
608 0 : libmesh_assert_less (idx, _nodes.size());
609 0 : libmesh_assert_less (dim, 3);
610 :
611 0 : const Point & p(_nodes[idx].first);
612 :
613 0 : if (dim==0) return p(0);
614 0 : if (dim==1) return p(1);
615 0 : return p(2);
616 : }
617 :
618 : /*
619 : * Optional bounding-box computation
620 : */
621 : template <class BBOX>
622 0 : bool kdtree_get_bbox(BBOX & /* bb */) const { return false; }
623 : };
624 :
625 :
626 : // ------------------------------------------------------------
627 : // UnstructuredMesh class member functions
628 331691 : UnstructuredMesh::UnstructuredMesh (const Parallel::Communicator & comm_in,
629 331691 : unsigned char d) :
630 331691 : MeshBase (comm_in,d)
631 : {
632 9630 : libmesh_assert (libMesh::initialized());
633 331691 : }
634 :
635 :
636 :
637 24392 : UnstructuredMesh::UnstructuredMesh (const MeshBase & other_mesh) :
638 24392 : MeshBase (other_mesh)
639 : {
640 804 : libmesh_assert (libMesh::initialized());
641 24392 : }
642 :
643 :
644 :
645 26167 : void UnstructuredMesh::copy_nodes_and_elements(const MeshBase & other_mesh,
646 : const bool skip_find_neighbors,
647 : dof_id_type element_id_offset,
648 : dof_id_type node_id_offset,
649 : unique_id_type
650 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
651 : unique_id_offset
652 : #endif
653 : ,
654 : std::unordered_map<subdomain_id_type, subdomain_id_type> *
655 : id_remapping,
656 : const bool skip_preparation)
657 : {
658 1708 : LOG_SCOPE("copy_nodes_and_elements()", "UnstructuredMesh");
659 :
660 : // If we're asked to skip all preparation, we should be skipping
661 : // find_neighbors specifically.
662 854 : libmesh_assert(!skip_preparation || skip_find_neighbors);
663 :
664 : std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
665 27859 : extra_int_maps = this->merge_extra_integer_names(other_mesh);
666 :
667 26167 : const unsigned int n_old_node_ints = extra_int_maps.second.size(),
668 26167 : n_new_node_ints = _node_integer_names.size(),
669 26167 : n_old_elem_ints = extra_int_maps.first.size(),
670 26167 : n_new_elem_ints = _elem_integer_names.size();
671 :
672 : // If we are partitioned into fewer parts than the incoming mesh has
673 : // processors to handle, then we need to "wrap" the other Mesh's
674 : // processor ids to fit within our range. This can happen, for
675 : // example, while stitching meshes with small numbers of elements in
676 : // parallel...
677 1692 : bool wrap_proc_ids = (this->n_processors() <
678 1692 : other_mesh.n_partitions());
679 :
680 : // We're assuming the other mesh has proper element number ordering,
681 : // so that we add parents before their children, and that the other
682 : // mesh is consistently partitioned. We're not assuming that node
683 : // proc ids are topologically consistent, so we don't just
684 : // libmesh_assert_valid_procids.
685 : #ifdef DEBUG
686 854 : MeshTools::libmesh_assert_valid_amr_elem_ids(other_mesh);
687 854 : MeshTools::libmesh_assert_parallel_consistent_procids<Node>(other_mesh);
688 : #endif
689 :
690 : //Copy in Nodes
691 : {
692 : //Preallocate Memory if necessary
693 26167 : this->reserve_nodes(other_mesh.n_nodes());
694 :
695 7672928 : for (const auto & oldn : other_mesh.node_ptr_range())
696 : {
697 : processor_id_type added_pid = cast_int<processor_id_type>
698 3995224 : (wrap_proc_ids ? oldn->processor_id() % this->n_processors() : oldn->processor_id());
699 :
700 : // Add new nodes in old node Point locations
701 : Node * newn =
702 4551940 : this->add_point(*oldn,
703 3995224 : oldn->id() + node_id_offset,
704 369000 : added_pid);
705 :
706 3995224 : newn->add_extra_integers(n_new_node_ints);
707 4220072 : for (unsigned int i = 0; i != n_old_node_ints; ++i)
708 236796 : newn->set_extra_integer(extra_int_maps.second[i],
709 224848 : oldn->get_extra_integer(i));
710 :
711 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
712 3995224 : newn->set_unique_id(oldn->unique_id() + unique_id_offset);
713 : #endif
714 24475 : }
715 : }
716 :
717 : //Copy in Elements
718 : {
719 : //Preallocate Memory if necessary
720 26167 : this->reserve_elem(other_mesh.n_elem());
721 :
722 : // Declare a map linking old and new elements, needed to copy the neighbor lists
723 : typedef std::unordered_map<const Elem *, Elem *> map_type;
724 1708 : map_type old_elems_to_new_elems, ip_map;
725 :
726 : // Loop over the elements
727 19570998 : for (const auto & old : other_mesh.element_ptr_range())
728 : {
729 : // Build a new element
730 10389951 : Elem * newparent = old->parent() ?
731 240715 : this->elem_ptr(old->parent()->id() + element_id_offset) :
732 323214 : nullptr;
733 10400775 : auto el = old->disconnected_clone();
734 635604 : el->set_parent(newparent);
735 :
736 10077561 : subdomain_id_type sbd_id = old->subdomain_id();
737 10077561 : if (id_remapping)
738 : {
739 538 : auto remapping_it = id_remapping->find(sbd_id);
740 19099 : if (remapping_it != id_remapping->end())
741 568 : sbd_id = remapping_it->second;
742 : }
743 10077561 : el->subdomain_id() = sbd_id;
744 :
745 : // Hold off on trying to set the interior parent because we may actually
746 : // add lower dimensional elements before their interior parents
747 10077561 : if (old->interior_parent())
748 688 : ip_map[old] = el.get();
749 :
750 : #ifdef LIBMESH_ENABLE_AMR
751 10077561 : if (old->has_children())
752 392533 : for (unsigned int c = 0, nc = old->n_children(); c != nc; ++c)
753 332628 : if (old->child_ptr(c) == remote_elem)
754 67489 : el->add_child(const_cast<RemoteElem *>(remote_elem), c);
755 :
756 : //Create the parent's child pointers if necessary
757 10077561 : if (newparent)
758 : {
759 265119 : unsigned int oldc = old->parent()->which_child_am_i(old);
760 240715 : newparent->add_child(el.get(), oldc);
761 : }
762 :
763 : // Copy the refinement flags
764 10077561 : el->set_refinement_flag(old->refinement_flag());
765 :
766 : // Use hack_p_level since we may not have sibling elements
767 : // added yet
768 635604 : el->hack_p_level(old->p_level());
769 :
770 635604 : el->set_p_refinement_flag(old->p_refinement_flag());
771 : #endif // #ifdef LIBMESH_ENABLE_AMR
772 :
773 : //Assign all the nodes
774 53027341 : for (auto i : el->node_index_range())
775 44413088 : el->set_node(i,
776 42949780 : this->node_ptr(old->node_id(i) + node_id_offset));
777 :
778 : // And start it off with the same processor id (mod _n_parts).
779 10077561 : el->processor_id() = cast_int<processor_id_type>
780 10077561 : (wrap_proc_ids ? old->processor_id() % this->n_processors() : old->processor_id());
781 :
782 : // Give it the same element and unique ids
783 10077561 : el->set_id(old->id() + element_id_offset);
784 :
785 10077561 : el->add_extra_integers(n_new_elem_ints);
786 10094814 : for (unsigned int i = 0; i != n_old_elem_ints; ++i)
787 18225 : el->set_extra_integer(extra_int_maps.first[i],
788 17253 : old->get_extra_integer(i));
789 :
790 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
791 10077561 : el->set_unique_id(old->unique_id() + unique_id_offset);
792 : #endif
793 :
794 : //Hold onto it
795 10077561 : if (!skip_find_neighbors)
796 : {
797 220630 : for (auto s : old->side_index_range())
798 182136 : if (old->neighbor_ptr(s) == remote_elem)
799 256 : el->set_neighbor(s, const_cast<RemoteElem *>(remote_elem));
800 46942 : this->add_elem(std::move(el));
801 : }
802 : else
803 : {
804 10344417 : Elem * new_el = this->add_elem(std::move(el));
805 10033435 : old_elems_to_new_elems[old] = new_el;
806 : }
807 9466432 : }
808 :
809 : // If the other_mesh had some interior parents, we may need to
810 : // copy those pointers (if they're to elements in a third mesh),
811 : // or create new equivalent pointers (if they're to elements we
812 : // just copied), or scream and die (if the other mesh had interior
813 : // parents from a third mesh but we already have interior parents
814 : // that aren't to that same third mesh.
815 26167 : if (!ip_map.empty())
816 : {
817 156 : std::atomic<bool> existing_interior_parents{false};
818 :
819 : Threads::parallel_for
820 156 : (this->element_stored_range(),
821 468 : [&existing_interior_parents](const ElemRange & range)
822 : {
823 2750 : for (Elem * elem : range)
824 2594 : if (elem->interior_parent())
825 : {
826 0 : existing_interior_parents = true;
827 0 : break;
828 : }
829 156 : });
830 :
831 : MeshBase * other_interior_mesh =
832 0 : const_cast<MeshBase *>(&other_mesh.interior_mesh());
833 :
834 : // If we don't already have interior parents, then we can just
835 : // use whatever interior_mesh we need for the incoming
836 : // elements.
837 156 : if (!existing_interior_parents)
838 : {
839 156 : if (other_interior_mesh == &other_mesh)
840 0 : this->set_interior_mesh(*this);
841 : else
842 0 : this->set_interior_mesh(*other_interior_mesh);
843 : }
844 :
845 156 : if (other_interior_mesh == &other_mesh &&
846 156 : _interior_mesh == this)
847 844 : for (auto & elem_pair : ip_map)
848 688 : elem_pair.second->set_interior_parent(
849 688 : this->elem_ptr(elem_pair.first->interior_parent()->id() + element_id_offset));
850 0 : else if (other_interior_mesh == _interior_mesh)
851 0 : for (auto & elem_pair : ip_map)
852 : {
853 0 : Elem * ip = const_cast<Elem *>(elem_pair.first->interior_parent());
854 0 : libmesh_assert(ip == remote_elem ||
855 : ip == other_interior_mesh->elem_ptr(ip->id()));
856 0 : elem_pair.second->set_interior_parent(ip);
857 : }
858 : else
859 0 : libmesh_error_msg("Cannot copy boundary elements between meshes with different interior meshes");
860 : }
861 :
862 : // Loop (again) over the elements to fill in the neighbors
863 26167 : if (skip_find_neighbors)
864 : {
865 25883 : old_elems_to_new_elems[remote_elem] = const_cast<RemoteElem*>(remote_elem);
866 :
867 19485002 : for (const auto & old_elem : other_mesh.element_ptr_range())
868 : {
869 10033435 : Elem * new_elem = old_elems_to_new_elems[old_elem];
870 50350488 : for (auto s : old_elem->side_index_range())
871 : {
872 41241311 : const Elem * old_neighbor = old_elem->neighbor_ptr(s);
873 40006071 : Elem * new_neighbor = old_elems_to_new_elems[old_neighbor];
874 2513718 : new_elem->set_neighbor(s, new_neighbor);
875 : }
876 24207 : }
877 : }
878 : }
879 :
880 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
881 : // We set the unique ids of nodes after adding them to the mesh such that our value of
882 : // _next_unique_id may be wrong. So we amend that here
883 26167 : this->set_next_unique_id(other_mesh.parallel_max_unique_id() + unique_id_offset + 1);
884 : #endif
885 :
886 : // Finally, partially prepare the new Mesh for use, if that isn't
887 : // being skipped.
888 : // Even the default behavior here is for backwards compatibility,
889 : // and we don't want to prepare everything.
890 :
891 26167 : if (!skip_preparation)
892 : {
893 : // Keep the same numbering and partitioning and distribution
894 : // status for now, but save our original policies to restore
895 : // later.
896 100 : const bool allowed_renumbering = this->allow_renumbering();
897 100 : const bool allowed_find_neighbors = this->allow_find_neighbors();
898 100 : const bool allowed_elem_removal = this->allow_remote_element_removal();
899 100 : const bool allowed_detect_detect_interior_parents = this->allow_detect_interior_parents();
900 50 : this->allow_renumbering(false);
901 50 : this->allow_remote_element_removal(false);
902 50 : this->allow_find_neighbors(!skip_find_neighbors);
903 100 : this->allow_detect_interior_parents(other_mesh.allow_detect_interior_parents());
904 :
905 : // We should generally be able to skip *all* partitioning here
906 : // because we're only adding one already-consistent mesh to
907 : // another.
908 100 : const bool skipped_partitioning = this->skip_partitioning();
909 50 : this->skip_partitioning(true);
910 :
911 1775 : const bool was_prepared = this->is_prepared();
912 1775 : this->prepare_for_use();
913 :
914 : //But in the long term, don't change our policies.
915 50 : this->allow_find_neighbors(allowed_find_neighbors);
916 50 : this->allow_renumbering(allowed_renumbering);
917 50 : this->allow_remote_element_removal(allowed_elem_removal);
918 50 : this->skip_partitioning(skipped_partitioning);
919 50 : this->allow_detect_interior_parents(allowed_detect_detect_interior_parents);
920 :
921 : // That prepare_for_use() call marked us as prepared, but we
922 : // specifically avoided some important preparation, so we might not
923 : // actually be prepared now.
924 1783 : if (skip_find_neighbors ||
925 1775 : !was_prepared || !other_mesh.is_prepared())
926 1491 : this->unset_is_prepared();
927 : }
928 :
929 : // In general we've just invalidated just about everything, and we'd
930 : // like to unset_is_prepared(), but specific use cases might know a
931 : // priori that they're still partitioned well, or that they've
932 : // copied in a disjoint mesh component and don't need new neighbor
933 : // pointers, or that they're not adding anything that would change
934 : // cached subdomain/element/boundary sets, etc., so we'll rely on
935 : // users of the "advanced" skip_preparation option to also set what
936 : // preparation they still need.
937 :
938 : // else
939 : // this->unset_is_prepared();
940 26167 : }
941 :
942 :
943 :
944 356083 : UnstructuredMesh::~UnstructuredMesh ()
945 : {
946 : // this->clear (); // Nothing to clear at this level
947 :
948 10434 : libmesh_exceptionless_assert (!libMesh::closed());
949 356083 : }
950 :
951 :
952 :
953 :
954 :
955 533211 : void UnstructuredMesh::find_neighbors (const bool reset_remote_elements,
956 : const bool reset_current_list,
957 : const bool assert_valid)
958 : {
959 : // We might actually want to run this on an empty mesh
960 : // (e.g. the boundary mesh for a nonexistent bcid!)
961 : // libmesh_assert_not_equal_to (this->n_nodes(), 0);
962 : // libmesh_assert_not_equal_to (this->n_elem(), 0);
963 :
964 : // This function must be run on all processors at once
965 13198 : parallel_object_only();
966 :
967 13198 : LOG_SCOPE("find_neighbors()", "Mesh");
968 :
969 : //TODO:[BSK] This should be removed later?!
970 533211 : if (reset_current_list)
971 : Threads::parallel_for
972 532797 : (this->element_stored_range(),
973 1549372 : [reset_remote_elements](const ElemRange & range)
974 : {
975 69835198 : for (Elem * e : range)
976 346322674 : for (auto s : e->side_index_range())
977 282655299 : if (e->neighbor_ptr(s) != remote_elem || reset_remote_elements)
978 5674571 : e->set_neighbor(s, nullptr);
979 533181 : });
980 :
981 : // Find neighboring elements by first finding elements
982 : // with identical side keys and then check to see if they
983 : // are neighbors
984 : {
985 : // data structures -- Use the hash_multimap if available
986 : typedef dof_id_type key_type;
987 : typedef std::pair<Elem *, unsigned char> val_type;
988 : typedef std::unordered_multimap<key_type, val_type> map_type;
989 :
990 : // A map from side keys to corresponding elements & side numbers
991 26396 : map_type side_to_elem_map;
992 :
993 : // Pull objects out of the loop to reduce heap operations
994 533211 : std::unique_ptr<Elem> my_side, their_side;
995 :
996 136836384 : for (const auto & element : this->element_ptr_range())
997 : {
998 347734679 : for (auto ms : element->side_index_range())
999 : {
1000 405125905 : next_side:
1001 : // If we haven't yet found a neighbor on this side, try.
1002 : // Even if we think our neighbor is remote, that
1003 : // information may be out of date.
1004 415917722 : if (element->neighbor_ptr(ms) == nullptr ||
1005 130101147 : element->neighbor_ptr(ms) == remote_elem)
1006 : {
1007 : // Get the key for the side of this element. Use the
1008 : // low_order_key so we can find neighbors in
1009 : // mixed-order meshes if necessary.
1010 277118874 : const dof_id_type key = element->low_order_key(ms);
1011 :
1012 : // Look for elements that have an identical side key
1013 5682756 : auto bounds = side_to_elem_map.equal_range(key);
1014 :
1015 : // May be multiple keys, check all the possible
1016 : // elements which _might_ be neighbors.
1017 277118874 : if (bounds.first != bounds.second)
1018 : {
1019 : // Get the side for this element
1020 128183258 : element->side_ptr(my_side, ms);
1021 :
1022 : // Look at all the entries with an equivalent key
1023 128380268 : while (bounds.first != bounds.second)
1024 : {
1025 : // Get the potential element
1026 128297618 : Elem * neighbor = bounds.first->second.first;
1027 :
1028 : // Get the side for the neighboring element
1029 128297618 : const unsigned int ns = bounds.first->second.second;
1030 128297618 : neighbor->side_ptr(their_side, ns);
1031 : //libmesh_assert(my_side.get());
1032 : //libmesh_assert(their_side.get());
1033 :
1034 : // If found a match with my side
1035 : //
1036 : // In 1D, since parents and children have an
1037 : // equal side (i.e. a node) we need to check
1038 : // for matching level() to avoid setting our
1039 : // neighbor pointer to any of our neighbor's
1040 : // descendants.
1041 256595236 : if ((*my_side == *their_side) &&
1042 128297618 : (element->level() == neighbor->level()))
1043 : {
1044 : // So share a side. Is this a mixed pair
1045 : // of subactive and active/ancestor
1046 : // elements?
1047 : // If not, then we're neighbors.
1048 : // If so, then the subactive's neighbor is
1049 :
1050 130669066 : if (element->subactive() ==
1051 128100608 : neighbor->subactive())
1052 : {
1053 : // an element is only subactive if it has
1054 : // been coarsened but not deleted
1055 5143670 : element->set_neighbor (ms,neighbor);
1056 127926903 : neighbor->set_neighbor(ns,element);
1057 : }
1058 173705 : else if (element->subactive())
1059 : {
1060 5096 : element->set_neighbor(ms,neighbor);
1061 : }
1062 94537 : else if (neighbor->subactive())
1063 : {
1064 9368 : neighbor->set_neighbor(ns,element);
1065 : }
1066 2589676 : side_to_elem_map.erase (bounds.first);
1067 :
1068 : // get out of this nested crap
1069 128100608 : goto next_side;
1070 : }
1071 :
1072 2936 : ++bounds.first;
1073 : }
1074 : }
1075 :
1076 : // didn't find a match...
1077 : // Build the map entry for this element
1078 : side_to_elem_map.emplace
1079 149018266 : (key, std::make_pair(element, cast_int<unsigned char>(ms)));
1080 : }
1081 : }
1082 506831 : }
1083 506831 : }
1084 :
1085 : #ifdef LIBMESH_ENABLE_PERIODIC
1086 : // Get the disjoint neighbor boundary pairs object (from periodic BCs)
1087 533211 : auto * db = this->get_disjoint_neighbor_boundary_pairs();
1088 :
1089 533211 : if (db)
1090 : {
1091 : // Obtain a point locator
1092 1045 : std::unique_ptr<PointLocatorBase> point_locator = this->sub_point_locator();
1093 :
1094 7130 : for (const auto & element : this->element_ptr_range())
1095 : {
1096 13062 : for (auto ms : element->side_index_range())
1097 : {
1098 : // Skip if this side already has a valid neighbor (including remote neighbors)
1099 10660 : if (element->neighbor_ptr(ms) != nullptr &&
1100 2003 : element->neighbor_ptr(ms) != remote_elem)
1101 1959 : continue;
1102 :
1103 32183 : for (const auto & [id, boundary_ptr] : *db)
1104 : {
1105 23778 : if (!this->get_boundary_info().has_boundary_id(element, ms, id))
1106 22343 : continue;
1107 :
1108 : unsigned int neigh_side;
1109 : const Elem * neigh =
1110 1463 : db->neighbor(id, *point_locator, element, ms, &neigh_side);
1111 :
1112 1435 : if (neigh && neigh != remote_elem && neigh != element)
1113 : {
1114 1435 : auto neigh_changeable = this->elem_ptr(neigh->id());
1115 1435 : element->set_neighbor(ms, neigh_changeable);
1116 1435 : neigh_changeable->set_neighbor(neigh_side, element);
1117 : }
1118 : }
1119 : }
1120 985 : }
1121 985 : }
1122 : #endif // LIBMESH_ENABLE_PERIODIC
1123 :
1124 : #ifdef LIBMESH_ENABLE_AMR
1125 :
1126 : /**
1127 : * Here we look at all of the child elements which
1128 : * don't already have valid neighbors.
1129 : *
1130 : * If a child element has a nullptr neighbor it is
1131 : * either because it is on the boundary or because
1132 : * its neighbor is at a different level. In the
1133 : * latter case we must get the neighbor from the
1134 : * parent.
1135 : *
1136 : * If a child element has a remote_elem neighbor
1137 : * on a boundary it shares with its parent, that
1138 : * info may have become out-dated through coarsening
1139 : * of the neighbor's parent. In this case, if the
1140 : * parent's neighbor is active then the child should
1141 : * share it.
1142 : *
1143 : * Furthermore, that neighbor better be active,
1144 : * otherwise we missed a child somewhere.
1145 : *
1146 : *
1147 : * We also need to look through children ordered by increasing
1148 : * refinement level in order to add new interior_parent() links in
1149 : * boundary elements which have just been generated by refinement,
1150 : * and fix links in boundary elements whose previous
1151 : * interior_parent() has just been coarsened away.
1152 : */
1153 533211 : const unsigned int n_levels = MeshTools::n_levels(*this);
1154 738167 : for (unsigned int level = 1; level < n_levels; ++level)
1155 : {
1156 1228830 : for (auto & current_elem : as_range(level_elements_begin(level),
1157 92204748 : level_elements_end(level)))
1158 : {
1159 823374 : libmesh_assert(current_elem);
1160 46108605 : Elem * parent = current_elem->parent();
1161 823374 : libmesh_assert(parent);
1162 46108605 : const unsigned int my_child_num = parent->which_child_am_i(current_elem);
1163 :
1164 228078507 : for (auto s : current_elem->side_index_range())
1165 : {
1166 187357638 : if (current_elem->neighbor_ptr(s) == nullptr ||
1167 173707577 : (current_elem->neighbor_ptr(s) == remote_elem &&
1168 939507 : parent->is_child_on_side(my_child_num, s)))
1169 : {
1170 463256 : Elem * neigh = parent->neighbor_ptr(s);
1171 :
1172 : // If neigh was refined and had non-subactive children
1173 : // made remote earlier, then our current elem should
1174 : // actually have one of those remote children as a
1175 : // neighbor
1176 12617970 : if (neigh &&
1177 3753967 : (neigh->ancestor() ||
1178 : // If neigh has subactive children which should have
1179 : // matched as neighbors of the current element but
1180 : // did not, then those likewise must be remote
1181 : // children.
1182 3617419 : (current_elem->subactive() && neigh->has_children() &&
1183 151 : (neigh->level()+1) == current_elem->level())))
1184 : {
1185 : #ifdef DEBUG
1186 : // Let's make sure that "had children made remote"
1187 : // situation is actually the case
1188 0 : libmesh_assert(neigh->has_children());
1189 0 : bool neigh_has_remote_children = false;
1190 0 : for (auto & child : neigh->child_ref_range())
1191 0 : if (&child == remote_elem)
1192 0 : neigh_has_remote_children = true;
1193 0 : libmesh_assert(neigh_has_remote_children);
1194 :
1195 : // And let's double-check that we don't have
1196 : // a remote_elem neighboring an active local element
1197 0 : if (current_elem->active())
1198 0 : libmesh_assert_not_equal_to (current_elem->processor_id(),
1199 : this->processor_id());
1200 : #endif // DEBUG
1201 28928 : neigh = const_cast<RemoteElem *>(remote_elem);
1202 : }
1203 : // If neigh and current_elem are more than one level
1204 : // apart, figuring out whether we have a remote
1205 : // neighbor here becomes much harder.
1206 8950897 : else if (neigh && (current_elem->subactive() &&
1207 7900 : neigh->has_children()))
1208 : {
1209 : // Find the deepest descendant of neigh which
1210 : // we could consider for a neighbor. If we run
1211 : // out of neigh children, then that's our
1212 : // neighbor. If we find a potential neighbor
1213 : // with remote_children and we don't find any
1214 : // potential neighbors among its non-remote
1215 : // children, then our neighbor must be remote.
1216 0 : while (neigh != remote_elem &&
1217 0 : neigh->has_children())
1218 : {
1219 0 : bool found_neigh = false;
1220 0 : for (unsigned int c = 0, nc = neigh->n_children();
1221 0 : !found_neigh && c != nc; ++c)
1222 : {
1223 0 : Elem * child = neigh->child_ptr(c);
1224 0 : if (child == remote_elem)
1225 0 : continue;
1226 0 : for (auto ncn : child->neighbor_ptr_range())
1227 : {
1228 0 : if (ncn != remote_elem &&
1229 0 : ncn->is_ancestor_of(current_elem))
1230 : {
1231 0 : neigh = ncn;
1232 0 : found_neigh = true;
1233 0 : break;
1234 : }
1235 : }
1236 : }
1237 0 : if (!found_neigh)
1238 0 : neigh = const_cast<RemoteElem *>(remote_elem);
1239 : }
1240 : }
1241 8971925 : current_elem->set_neighbor(s, neigh);
1242 : #ifdef DEBUG
1243 231628 : if (neigh != nullptr && neigh != remote_elem)
1244 : // We ignore subactive elements here because
1245 : // we don't care about neighbors of subactive element.
1246 107152 : if ((!neigh->active()) && (!current_elem->subactive()))
1247 : {
1248 0 : libMesh::err << "On processor " << this->processor_id()
1249 0 : << std::endl;
1250 0 : libMesh::err << "Bad element ID = " << current_elem->id()
1251 0 : << ", Side " << s << ", Bad neighbor ID = " << neigh->id() << std::endl;
1252 0 : libMesh::err << "Bad element proc_ID = " << current_elem->processor_id()
1253 0 : << ", Bad neighbor proc_ID = " << neigh->processor_id() << std::endl;
1254 0 : libMesh::err << "Bad element size = " << current_elem->hmin()
1255 0 : << ", Bad neighbor size = " << neigh->hmin() << std::endl;
1256 0 : libMesh::err << "Bad element center = " << current_elem->vertex_average()
1257 0 : << ", Bad neighbor center = " << neigh->vertex_average() << std::endl;
1258 0 : libMesh::err << "ERROR: "
1259 0 : << (current_elem->active()?"Active":"Ancestor")
1260 0 : << " Element at level "
1261 0 : << current_elem->level() << std::endl;
1262 0 : libMesh::err << "with "
1263 0 : << (parent->active()?"active":
1264 0 : (parent->subactive()?"subactive":"ancestor"))
1265 0 : << " parent share "
1266 0 : << (neigh->subactive()?"subactive":"ancestor")
1267 0 : << " neighbor at level " << neigh->level()
1268 0 : << std::endl;
1269 0 : NameBasedIO(*this).write ("bad_mesh.gmv");
1270 0 : libmesh_error_msg("Problematic mesh written to bad_mesh.gmv.");
1271 : }
1272 : #endif // DEBUG
1273 : }
1274 : }
1275 :
1276 : // We can skip to the next element if we're full-dimension
1277 : // and therefore don't have any interior parents
1278 46108605 : if (current_elem->dim() >= LIBMESH_DIM)
1279 4930718 : continue;
1280 :
1281 : // We have no interior parents unless we can find one later
1282 41602693 : current_elem->set_interior_parent(nullptr);
1283 :
1284 41602693 : Elem * pip = parent->interior_parent();
1285 :
1286 41602693 : if (!pip)
1287 40955992 : continue;
1288 :
1289 : // If there's no interior_parent children, whether due to a
1290 : // remote element or a non-conformity, then there's no
1291 : // children to search.
1292 23111 : if (pip == remote_elem || pip->active())
1293 : {
1294 2052 : current_elem->set_interior_parent(pip);
1295 2052 : continue;
1296 : }
1297 :
1298 : // For node comparisons we'll need a sensible tolerance
1299 21059 : Real node_tolerance = current_elem->hmin() * TOLERANCE;
1300 :
1301 : // Otherwise our interior_parent should be a child of our
1302 : // parent's interior_parent.
1303 70139 : for (auto & child : pip->child_ref_range())
1304 : {
1305 : // If we have a remote_elem, that might be our
1306 : // interior_parent. We'll set it provisionally now and
1307 : // keep trying to find something better.
1308 69001 : if (&child == remote_elem)
1309 : {
1310 : current_elem->set_interior_parent
1311 4703 : (const_cast<RemoteElem *>(remote_elem));
1312 4703 : continue;
1313 : }
1314 :
1315 3008 : bool child_contains_our_nodes = true;
1316 130673 : for (auto & n : current_elem->node_ref_range())
1317 : {
1318 5220 : bool child_contains_this_node = false;
1319 711982 : for (auto & cn : child.node_ref_range())
1320 667605 : if (cn.absolute_fuzzy_equals
1321 636405 : (n, node_tolerance))
1322 : {
1323 3212 : child_contains_this_node = true;
1324 3212 : break;
1325 : }
1326 107744 : if (!child_contains_this_node)
1327 : {
1328 2008 : child_contains_our_nodes = false;
1329 2008 : break;
1330 : }
1331 : }
1332 64298 : if (child_contains_our_nodes)
1333 : {
1334 19921 : current_elem->set_interior_parent(&child);
1335 1000 : break;
1336 : }
1337 : }
1338 :
1339 : // We should have found *some* interior_parent at this
1340 : // point, whether semilocal or remote.
1341 1000 : libmesh_assert(current_elem->interior_parent());
1342 196044 : }
1343 : }
1344 : #endif // AMR
1345 :
1346 : #ifdef DEBUG
1347 13198 : if (assert_valid)
1348 : {
1349 13198 : MeshTools::libmesh_assert_valid_neighbors(*this,
1350 13198 : !reset_remote_elements);
1351 13198 : MeshTools::libmesh_assert_valid_amr_interior_parents(*this);
1352 : }
1353 : #else
1354 : libmesh_ignore(assert_valid);
1355 : #endif
1356 :
1357 533211 : this->_preparation.has_neighbor_ptrs = true;
1358 533211 : }
1359 :
1360 :
1361 :
1362 5634 : void UnstructuredMesh::read (const std::string & name,
1363 : void *,
1364 : bool skip_renumber_nodes_and_elements,
1365 : bool skip_find_neighbors,
1366 : bool skip_detect_interior_parents)
1367 : {
1368 : // Set the skip_renumber_nodes_and_elements flag on all processors
1369 : // if necessary.
1370 : // This ensures that renumber_nodes_and_elements is *not* called
1371 : // during prepare_for_use() for certain types of mesh files.
1372 : // This is required in cases where there is an associated solution
1373 : // file which expects a certain ordering of the nodes.
1374 5634 : if (Utility::ends_with(name, ".gmv"))
1375 0 : this->allow_renumbering(false);
1376 :
1377 5634 : NameBasedIO(*this).read(name);
1378 :
1379 5634 : if (skip_renumber_nodes_and_elements)
1380 : {
1381 : // Use MeshBase::allow_renumbering() yourself instead.
1382 : libmesh_deprecated();
1383 0 : this->allow_renumbering(false);
1384 : }
1385 :
1386 : // Done reading the mesh. Now prepare it for use.
1387 348 : const bool old_allow_find_neighbors = this->allow_find_neighbors();
1388 348 : const bool old_allow_detect_interior_parents = this->allow_detect_interior_parents();
1389 :
1390 174 : this->allow_find_neighbors(!skip_find_neighbors);
1391 174 : this->allow_detect_interior_parents(!skip_detect_interior_parents);
1392 :
1393 5634 : this->prepare_for_use();
1394 :
1395 174 : this->allow_find_neighbors(old_allow_find_neighbors);
1396 174 : this->allow_detect_interior_parents(old_allow_detect_interior_parents);
1397 5634 : }
1398 :
1399 :
1400 :
1401 2938 : void UnstructuredMesh::write (const std::string & name) const
1402 : {
1403 598 : LOG_SCOPE("write()", "Mesh");
1404 :
1405 3536 : NameBasedIO(*this).write(name);
1406 2938 : }
1407 :
1408 :
1409 :
1410 0 : void UnstructuredMesh::write (const std::string & name,
1411 : const std::vector<Number> & v,
1412 : const std::vector<std::string> & vn) const
1413 : {
1414 0 : LOG_SCOPE("write()", "Mesh");
1415 :
1416 0 : NameBasedIO(*this).write_nodal_data(name, v, vn);
1417 0 : }
1418 :
1419 :
1420 :
1421 :
1422 :
1423 0 : void UnstructuredMesh::create_pid_mesh(UnstructuredMesh & pid_mesh,
1424 : const processor_id_type pid) const
1425 : {
1426 :
1427 : // Issue a warning if the number the number of processors
1428 : // currently available is less that that requested for
1429 : // partitioning. This is not necessarily an error since
1430 : // you may run on one processor and still partition the
1431 : // mesh into several partitions.
1432 : #ifdef DEBUG
1433 0 : if (this->n_processors() < pid)
1434 : {
1435 0 : libMesh::out << "WARNING: You are creating a "
1436 0 : << "mesh for a processor id (="
1437 0 : << pid
1438 0 : << ") greater than "
1439 0 : << "the number of processors available for "
1440 0 : << "the calculation. (="
1441 0 : << this->n_processors()
1442 0 : << ")."
1443 0 : << std::endl;
1444 : }
1445 : #endif
1446 :
1447 0 : this->create_submesh (pid_mesh,
1448 0 : this->active_pid_elements_begin(pid),
1449 0 : this->active_pid_elements_end(pid));
1450 0 : }
1451 :
1452 :
1453 :
1454 :
1455 :
1456 :
1457 :
1458 0 : void UnstructuredMesh::create_submesh (UnstructuredMesh & new_mesh,
1459 : const const_element_iterator & it,
1460 : const const_element_iterator & it_end) const
1461 : {
1462 : // Just in case the subdomain_mesh already has some information
1463 : // in it, get rid of it.
1464 0 : new_mesh.clear();
1465 :
1466 : // If we're not serial, our submesh isn't either.
1467 : // There are no remote elements to delete on an empty mesh, but
1468 : // calling the method to do so marks the mesh as parallel.
1469 0 : if (!this->is_serial())
1470 0 : new_mesh.delete_remote_elements();
1471 :
1472 : // Fail if (*this == new_mesh), we cannot create a submesh inside ourself!
1473 : // This may happen if the user accidentally passes the original mesh into
1474 : // this function! We will check this by making sure we did not just
1475 : // clear ourself.
1476 0 : libmesh_assert_not_equal_to (this->n_nodes(), 0);
1477 0 : libmesh_assert_not_equal_to (this->n_elem(), 0);
1478 :
1479 : // Container to catch boundary IDs handed back by BoundaryInfo
1480 0 : std::vector<boundary_id_type> bc_ids;
1481 :
1482 : // Put any extra integers on the new mesh too
1483 0 : new_mesh.merge_extra_integer_names(*this);
1484 0 : const unsigned int n_node_ints = _node_integer_names.size();
1485 :
1486 0 : for (const auto & old_elem : as_range(it, it_end))
1487 : {
1488 : // Add an equivalent element type to the new_mesh.
1489 : // disconnected_clone() copies ids, extra element integers, etc.
1490 0 : auto uelem = old_elem->disconnected_clone();
1491 0 : Elem * new_elem = new_mesh.add_elem(std::move(uelem));
1492 0 : libmesh_assert(new_elem);
1493 :
1494 : // Loop over the nodes on this element.
1495 0 : for (auto n : old_elem->node_index_range())
1496 : {
1497 0 : const dof_id_type this_node_id = old_elem->node_id(n);
1498 :
1499 : // Add this node to the new mesh if it's not there already
1500 0 : if (!new_mesh.query_node_ptr(this_node_id))
1501 : {
1502 : Node * newn =
1503 0 : new_mesh.add_point (old_elem->point(n),
1504 : this_node_id,
1505 0 : old_elem->node_ptr(n)->processor_id());
1506 :
1507 0 : newn->add_extra_integers(n_node_ints);
1508 0 : for (unsigned int i = 0; i != n_node_ints; ++i)
1509 0 : newn->set_extra_integer(i, old_elem->node_ptr(n)->get_extra_integer(i));
1510 :
1511 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
1512 0 : newn->set_unique_id(old_elem->node_ptr(n)->unique_id());
1513 : #endif
1514 : }
1515 :
1516 : // Define this element's connectivity on the new mesh
1517 0 : new_elem->set_node(n, new_mesh.node_ptr(this_node_id));
1518 : }
1519 :
1520 : // Maybe add boundary conditions for this element
1521 0 : for (auto s : old_elem->side_index_range())
1522 : {
1523 0 : this->get_boundary_info().boundary_ids(old_elem, s, bc_ids);
1524 0 : new_mesh.get_boundary_info().add_side (new_elem, s, bc_ids);
1525 : }
1526 0 : } // end loop over elements
1527 :
1528 : // Prepare the new_mesh for use
1529 0 : new_mesh.prepare_for_use();
1530 0 : }
1531 :
1532 :
1533 :
1534 : #ifdef LIBMESH_ENABLE_AMR
1535 25826 : bool UnstructuredMesh::contract ()
1536 : {
1537 858 : LOG_SCOPE ("contract()", "Mesh");
1538 :
1539 : // Flag indicating if this call actually changes the mesh
1540 858 : bool mesh_changed = false;
1541 :
1542 : #ifdef DEBUG
1543 545312 : for (const auto & elem : this->element_ptr_range())
1544 544454 : libmesh_assert(elem->active() || elem->subactive() || elem->ancestor());
1545 : #endif
1546 :
1547 : // Loop over the elements.
1548 20050034 : for (auto & elem : this->element_ptr_range())
1549 : {
1550 : // Delete all the subactive ones
1551 10544074 : if (elem->subactive())
1552 : {
1553 : // No level-0 element should be subactive.
1554 : // Note that we CAN'T test elem->level(), as that
1555 : // touches elem->parent()->dim(), and elem->parent()
1556 : // might have already been deleted!
1557 71776 : libmesh_assert(elem->parent());
1558 :
1559 : // Delete the element
1560 : // This just sets a pointer to nullptr, and doesn't
1561 : // invalidate any iterators
1562 1498518 : this->delete_elem(elem);
1563 :
1564 : // the mesh has certainly changed
1565 71776 : mesh_changed = true;
1566 : }
1567 : else
1568 : {
1569 : // Compress all the active ones
1570 472678 : if (elem->active())
1571 6776797 : elem->contract();
1572 : else
1573 113362 : libmesh_assert (elem->ancestor());
1574 : }
1575 24110 : }
1576 :
1577 : // Strip any newly-created nullptr voids out of the element array
1578 25826 : this->renumber_nodes_and_elements();
1579 :
1580 : // FIXME: Need to understand why deleting subactive children
1581 : // invalidates the point locator. For now we will clear it explicitly
1582 25826 : this->clear_point_locator();
1583 :
1584 : // Allow our GhostingFunctor objects to reinit if necessary.
1585 27788 : for (auto & gf : as_range(this->ghosting_functors_begin(),
1586 116257 : this->ghosting_functors_end()))
1587 : {
1588 2820 : libmesh_assert(gf);
1589 86753 : gf->mesh_reinit();
1590 : }
1591 :
1592 26684 : return mesh_changed;
1593 : }
1594 : #endif // #ifdef LIBMESH_ENABLE_AMR
1595 :
1596 :
1597 :
1598 11259 : void UnstructuredMesh::all_first_order ()
1599 : {
1600 784 : LOG_SCOPE("all_first_order()", "Mesh");
1601 :
1602 : /**
1603 : * Prepare to identify (and then delete) a bunch of no-longer-used nodes.
1604 : */
1605 11651 : std::vector<bool> node_touched_by_me(this->max_node_id(), false);
1606 :
1607 : // Loop over the high-ordered elements.
1608 : // First make sure they _are_ indeed high-order, and then replace
1609 : // them with an equivalent first-order element.
1610 481878 : for (auto & so_elem : element_ptr_range())
1611 : {
1612 25848 : libmesh_assert(so_elem);
1613 :
1614 : /*
1615 : * build the first-order equivalent, add to
1616 : * the new_elements list.
1617 : */
1618 : auto lo_elem = Elem::build
1619 : (Elem::first_order_equivalent_type
1620 281572 : (so_elem->type()), so_elem->parent());
1621 :
1622 255724 : const unsigned short n_sides = so_elem->n_sides();
1623 :
1624 1220978 : for (unsigned short s=0; s != n_sides; ++s)
1625 1063058 : if (so_elem->neighbor_ptr(s) == remote_elem)
1626 0 : lo_elem->set_neighbor(s, const_cast<RemoteElem *>(remote_elem));
1627 :
1628 : #ifdef LIBMESH_ENABLE_AMR
1629 : /*
1630 : * Reset the parent links of any child elements
1631 : */
1632 255724 : if (so_elem->has_children())
1633 376697 : for (unsigned int c = 0, nc = so_elem->n_children(); c != nc; ++c)
1634 : {
1635 295596 : Elem * child = so_elem->child_ptr(c);
1636 295596 : if (child != remote_elem)
1637 48208 : child->set_parent(lo_elem.get());
1638 295596 : lo_elem->add_child(child, c);
1639 : }
1640 :
1641 : /*
1642 : * Reset the child link of any parent element
1643 : */
1644 281572 : if (so_elem->parent())
1645 : {
1646 : unsigned int c =
1647 231063 : so_elem->parent()->which_child_am_i(so_elem);
1648 255167 : lo_elem->parent()->replace_child(lo_elem.get(), c);
1649 : }
1650 :
1651 : /*
1652 : * Copy as much data to the new element as makes sense
1653 : */
1654 281572 : lo_elem->set_p_level(so_elem->p_level());
1655 255724 : lo_elem->set_refinement_flag(so_elem->refinement_flag());
1656 51696 : lo_elem->set_p_refinement_flag(so_elem->p_refinement_flag());
1657 : #endif
1658 :
1659 25848 : libmesh_assert_equal_to (lo_elem->n_vertices(), so_elem->n_vertices());
1660 :
1661 : /*
1662 : * By definition the vertices of the linear and
1663 : * second order element are identically numbered.
1664 : * transfer these.
1665 : */
1666 1229918 : for (unsigned int v=0, snv=so_elem->n_vertices(); v < snv; v++)
1667 : {
1668 1072478 : lo_elem->set_node(v, so_elem->node_ptr(v));
1669 294852 : node_touched_by_me[lo_elem->node_id(v)] = true;
1670 : }
1671 :
1672 : /*
1673 : * find_neighbors relies on remote_elem neighbor links being
1674 : * properly maintained.
1675 : */
1676 1220978 : for (unsigned short s=0; s != n_sides; s++)
1677 : {
1678 1063058 : if (so_elem->neighbor_ptr(s) == remote_elem)
1679 0 : lo_elem->set_neighbor(s, const_cast<RemoteElem*>(remote_elem));
1680 : }
1681 :
1682 : /**
1683 : * If the second order element had any boundary conditions they
1684 : * should be transferred to the first-order element. The old
1685 : * boundary conditions will be removed from the BoundaryInfo
1686 : * data structure by insert_elem.
1687 : */
1688 25848 : this->get_boundary_info().copy_boundary_ids
1689 255724 : (this->get_boundary_info(), so_elem, lo_elem.get());
1690 :
1691 : /*
1692 : * The new first-order element is ready.
1693 : * Inserting it into the mesh will replace and delete
1694 : * the second-order element.
1695 : */
1696 255724 : lo_elem->set_id(so_elem->id());
1697 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
1698 51696 : lo_elem->set_unique_id(so_elem->unique_id());
1699 : #endif
1700 :
1701 255724 : const unsigned int nei = so_elem->n_extra_integers();
1702 255724 : lo_elem->add_extra_integers(nei);
1703 258551 : for (unsigned int i=0; i != nei; ++i)
1704 2827 : lo_elem->set_extra_integer(i, so_elem->get_extra_integer(i));
1705 :
1706 255724 : lo_elem->inherit_data_from(*so_elem);
1707 :
1708 307420 : this->insert_elem(std::move(lo_elem));
1709 214503 : }
1710 :
1711 : // Deleting nodes does not invalidate iterators, so this is safe.
1712 1694776 : for (const auto & node : this->node_ptr_range())
1713 1011185 : if (!node_touched_by_me[node->id()])
1714 628888 : this->delete_node(node);
1715 :
1716 : // If crazy people applied boundary info to non-vertices and then
1717 : // deleted those non-vertices, we should make sure their boundary id
1718 : // caches are correct.
1719 11259 : this->get_boundary_info().regenerate_id_sets();
1720 :
1721 : // On hanging nodes that used to also be second order nodes, we
1722 : // might now have an invalid nodal processor_id()
1723 11259 : Partitioner::set_node_processor_ids(*this);
1724 :
1725 : // delete or renumber nodes if desired
1726 11259 : this->prepare_for_use();
1727 11259 : }
1728 :
1729 :
1730 :
1731 : void
1732 21134 : UnstructuredMesh::all_second_order_range (const SimpleRange<element_iterator> & range,
1733 : const bool full_ordered)
1734 : {
1735 614 : LOG_SCOPE("all_second_order_range()", "Mesh");
1736 :
1737 : /*
1738 : * The maximum number of new second order nodes we might be adding,
1739 : * for use when picking unique unique_id values later. This variable
1740 : * is not used unless unique ids are enabled.
1741 : */
1742 : unsigned int max_new_nodes_per_elem;
1743 :
1744 : /*
1745 : * For speed-up of the \p add_point() method, we
1746 : * can reserve memory. Guess the number of additional
1747 : * nodes based on the element spatial dimensions and the
1748 : * total number of nodes in the mesh as an upper bound.
1749 : */
1750 21134 : switch (this->mesh_dimension())
1751 : {
1752 923 : case 1:
1753 : /*
1754 : * in 1D, there can only be order-increase from Edge2
1755 : * to Edge3. Something like 1/2 of n_nodes() have
1756 : * to be added
1757 : */
1758 26 : max_new_nodes_per_elem = 3 - 2;
1759 1846 : this->reserve_nodes(static_cast<unsigned int>
1760 923 : (1.5*static_cast<double>(this->n_nodes())));
1761 897 : break;
1762 :
1763 3586 : case 2:
1764 : /*
1765 : * in 2D, either refine from Tri3 to Tri6 (double the nodes)
1766 : * or from Quad4 to Quad8 (again, double) or Quad9 (2.25 that much)
1767 : */
1768 116 : max_new_nodes_per_elem = 9 - 4;
1769 7172 : this->reserve_nodes(static_cast<unsigned int>
1770 3586 : (2*static_cast<double>(this->n_nodes())));
1771 3470 : break;
1772 :
1773 :
1774 16625 : case 3:
1775 : /*
1776 : * in 3D, either refine from Tet4 to Tet10 (factor = 2.5) up to
1777 : * Hex8 to Hex27 (something > 3). Since in 3D there _are_ already
1778 : * quite some nodes, and since we do not want to overburden the memory by
1779 : * a too conservative guess, use the lower bound
1780 : */
1781 472 : max_new_nodes_per_elem = 27 - 8;
1782 33250 : this->reserve_nodes(static_cast<unsigned int>
1783 16625 : (2.5*static_cast<double>(this->n_nodes())));
1784 16153 : break;
1785 :
1786 0 : default:
1787 : // Hm?
1788 0 : libmesh_error_msg("Unknown mesh dimension " << this->mesh_dimension());
1789 : }
1790 :
1791 : // All the real work is done in the helper function
1792 21134 : all_increased_order_range(*this, range, max_new_nodes_per_elem,
1793 85186 : [full_ordered](ElemType t) {
1794 2074815 : return Elem::second_order_equivalent_type(t, full_ordered);
1795 : });
1796 21134 : }
1797 :
1798 :
1799 :
1800 20296 : void UnstructuredMesh::all_complete_order_range(const SimpleRange<element_iterator> & range)
1801 : {
1802 568 : LOG_SCOPE("all_complete_order()", "Mesh");
1803 :
1804 : /*
1805 : * The maximum number of new higher-order nodes we might be adding,
1806 : * for use when picking unique unique_id values later. This variable
1807 : * is not used unless unique ids are enabled.
1808 : */
1809 : unsigned int max_new_nodes_per_elem;
1810 :
1811 : /*
1812 : * for speed-up of the \p add_point() method, we
1813 : * can reserve memory. Guess the number of additional
1814 : * nodes based on the element spatial dimensions and the
1815 : * total number of nodes in the mesh as an upper bound.
1816 : */
1817 20296 : switch (this->mesh_dimension())
1818 : {
1819 0 : case 1:
1820 : /*
1821 : * in 1D, there can only be order-increase from Edge2
1822 : * to Edge3. Something like 1/2 of n_nodes() have
1823 : * to be added
1824 : */
1825 0 : max_new_nodes_per_elem = 3 - 2;
1826 0 : this->reserve_nodes(static_cast<unsigned int>
1827 0 : (1.5*static_cast<double>(this->n_nodes())));
1828 0 : break;
1829 :
1830 1633 : case 2:
1831 : /*
1832 : * in 2D, we typically refine from Tri6 to Tri7 (1.1667 times
1833 : * the nodes) but might refine from Quad4 to Quad9
1834 : * (2.25 times the nodes)
1835 : */
1836 46 : max_new_nodes_per_elem = 9 - 4;
1837 3266 : this->reserve_nodes(static_cast<unsigned int>
1838 1633 : (2*static_cast<double>(this->n_nodes())));
1839 1587 : break;
1840 :
1841 :
1842 18663 : case 3:
1843 : /*
1844 : * in 3D, we typically refine from Tet10 to Tet14 (factor = 1.4)
1845 : * but may go Hex8 to Hex27 (something > 3). Since in 3D there
1846 : * _are_ already quite some nodes, and since we do not want to
1847 : * overburden the memory by a too conservative guess, use a
1848 : * moderate bound
1849 : */
1850 522 : max_new_nodes_per_elem = 27 - 8;
1851 37326 : this->reserve_nodes(static_cast<unsigned int>
1852 18663 : (2.5*static_cast<double>(this->n_nodes())));
1853 18141 : break;
1854 :
1855 0 : default:
1856 : // Hm?
1857 0 : libmesh_error_msg("Unknown mesh dimension " << this->mesh_dimension());
1858 : }
1859 :
1860 : // All the real work is done in the helper function
1861 20296 : all_increased_order_range(*this, range, max_new_nodes_per_elem,
1862 143032 : [](ElemType t) {
1863 5253502 : return Elem::complete_order_equivalent_type(t);
1864 : });
1865 20296 : }
1866 :
1867 :
1868 : std::size_t
1869 1633 : UnstructuredMesh::stitch_meshes (const MeshBase & other_mesh,
1870 : boundary_id_type this_mesh_boundary_id,
1871 : boundary_id_type other_mesh_boundary_id,
1872 : Real tol,
1873 : bool clear_stitched_boundary_ids,
1874 : bool verbose,
1875 : bool use_binary_search,
1876 : bool enforce_all_nodes_match_on_boundaries,
1877 : bool merge_boundary_nodes_all_or_nothing,
1878 : bool remap_subdomain_ids,
1879 : bool prepare_after_stitching)
1880 : {
1881 92 : LOG_SCOPE("stitch_meshes()", "UnstructuredMesh");
1882 1633 : return stitching_helper(&other_mesh,
1883 : this_mesh_boundary_id,
1884 : other_mesh_boundary_id,
1885 : tol,
1886 : clear_stitched_boundary_ids,
1887 : verbose,
1888 : use_binary_search,
1889 : enforce_all_nodes_match_on_boundaries,
1890 : true,
1891 : merge_boundary_nodes_all_or_nothing,
1892 : remap_subdomain_ids,
1893 1533 : prepare_after_stitching);
1894 : }
1895 :
1896 :
1897 : std::size_t
1898 213 : UnstructuredMesh::stitch_surfaces (boundary_id_type boundary_id_1,
1899 : boundary_id_type boundary_id_2,
1900 : Real tol,
1901 : bool clear_stitched_boundary_ids,
1902 : bool verbose,
1903 : bool use_binary_search,
1904 : bool enforce_all_nodes_match_on_boundaries,
1905 : bool merge_boundary_nodes_all_or_nothing,
1906 : bool prepare_after_stitching)
1907 :
1908 : {
1909 213 : return stitching_helper(nullptr,
1910 : boundary_id_1,
1911 : boundary_id_2,
1912 : tol,
1913 : clear_stitched_boundary_ids,
1914 : verbose,
1915 : use_binary_search,
1916 : enforce_all_nodes_match_on_boundaries,
1917 : /* skip_find_neighbors = */ true,
1918 : merge_boundary_nodes_all_or_nothing,
1919 : /* remap_subdomain_ids = */ false,
1920 213 : prepare_after_stitching);
1921 : }
1922 :
1923 :
1924 : std::size_t
1925 1846 : UnstructuredMesh::stitching_helper (const MeshBase * other_mesh,
1926 : boundary_id_type this_mesh_boundary_id,
1927 : boundary_id_type other_mesh_boundary_id,
1928 : Real tol,
1929 : bool clear_stitched_boundary_ids,
1930 : bool verbose,
1931 : bool use_binary_search,
1932 : bool enforce_all_nodes_match_on_boundaries,
1933 : bool skip_find_neighbors,
1934 : bool merge_boundary_nodes_all_or_nothing,
1935 : bool remap_subdomain_ids,
1936 : bool prepare_after_stitching)
1937 : {
1938 : #ifdef DEBUG
1939 : // We rely on neighbor links here
1940 52 : MeshTools::libmesh_assert_valid_neighbors(*this);
1941 : #endif
1942 :
1943 52 : bool is_valid_disjoint_pair_to_stitch = false;
1944 :
1945 : #ifdef LIBMESH_ENABLE_PERIODIC
1946 1846 : auto * this_db = this->get_disjoint_neighbor_boundary_pairs();
1947 1846 : auto * other_db = (other_mesh ? other_mesh->get_disjoint_neighbor_boundary_pairs() : nullptr);
1948 : const bool have_disc_bdys =
1949 1846 : (this_db && !this_db->empty()) || (other_db && !other_db->empty());
1950 :
1951 52 : if (have_disc_bdys)
1952 : {
1953 10 : const boundary_id_type a = this_mesh_boundary_id;
1954 10 : const boundary_id_type b = other_mesh_boundary_id;
1955 :
1956 40 : auto get_pb = [](const PeriodicBoundaries * db, boundary_id_type id)
1957 : {
1958 539 : return db ? db->boundary(id) : nullptr;
1959 : };
1960 :
1961 : // this mesh
1962 355 : const auto * pb_this_a = get_pb(this_db, a);
1963 355 : const auto * pb_this_b = get_pb(this_db, b);
1964 10 : const bool in_this =
1965 355 : (pb_this_a && pb_this_a->pairedboundary == b) ||
1966 0 : (pb_this_b && pb_this_b->pairedboundary == a);
1967 :
1968 : // other mesh
1969 345 : const auto * pb_other_b = get_pb(other_db, b);
1970 10 : const auto * pb_other_a = get_pb(other_db, a);
1971 10 : const bool in_other =
1972 355 : (pb_other_b && pb_other_b->pairedboundary == a) ||
1973 0 : (pb_other_a && pb_other_a->pairedboundary == b);
1974 :
1975 : // Conflict conditions:
1976 : // Case 1: On "this" mesh, a or b exist but are not paired,
1977 : // while the other mesh pairs them.
1978 355 : if (!in_this && (pb_this_a || pb_this_b) && in_other)
1979 146 : libmesh_error_msg("Disjoint neighbor boundary pairing mismatch: on 'this' mesh, "
1980 : "boundary (" << a << " or " << b
1981 : << ") exists but is not paired; on 'other' mesh the pair is present.");
1982 :
1983 : // Case 2: On "other" mesh, a or b exist but are not paired,
1984 : // while this mesh pairs them.
1985 284 : if (!in_other && (pb_other_a || pb_other_b) && in_this)
1986 0 : libmesh_error_msg("Disjoint neighbor boundary pairing mismatch: on 'other' mesh, "
1987 : "boundary (" << a << " or " << b
1988 : << ") exists but is not paired; on 'this' mesh the pair is present.");
1989 :
1990 : // Legal conditions: either side has a correct pairing
1991 284 : if (in_this || in_other)
1992 6 : is_valid_disjoint_pair_to_stitch = true;
1993 : }
1994 : #endif // LIBMESH_ENABLE_PERIODIC
1995 :
1996 : // We can't even afford any unset neighbor links here.
1997 1775 : if (!this->is_prepared())
1998 0 : this->find_neighbors();
1999 :
2000 : // FIXME: make distributed mesh support efficient.
2001 : // Yes, we currently suck.
2002 1875 : MeshSerializer serialize(*this);
2003 :
2004 : // *Badly*.
2005 1725 : std::unique_ptr<MeshSerializer> serialize_other;
2006 1775 : if (other_mesh)
2007 : serialize_other = std::make_unique<MeshSerializer>
2008 3036 : (*const_cast<MeshBase *>(other_mesh));
2009 :
2010 102 : 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
2011 100 : std::map<dof_id_type, std::vector<dof_id_type>> node_to_elems_map;
2012 :
2013 : typedef dof_id_type key_type;
2014 : typedef std::pair<const Elem *, unsigned char> val_type;
2015 : typedef std::pair<key_type, val_type> key_val_pair;
2016 : typedef std::unordered_multimap<key_type, val_type> map_type;
2017 : // Mapping between all side keys in this mesh and elements+side numbers relevant to the boundary in this mesh as well.
2018 100 : map_type side_to_elem_map;
2019 :
2020 : // If there is only one mesh (i.e. other_mesh == nullptr), then loop over this mesh twice
2021 1775 : if (!other_mesh)
2022 : {
2023 6 : other_mesh = this;
2024 : }
2025 :
2026 1775 : if ((this_mesh_boundary_id != BoundaryInfo::invalid_id) &&
2027 50 : (other_mesh_boundary_id != BoundaryInfo::invalid_id))
2028 : {
2029 100 : LOG_SCOPE("stitch_meshes node merging", "UnstructuredMesh");
2030 :
2031 : // While finding nodes on the boundary, also find the minimum edge length
2032 : // of all faces on both boundaries. This will later be used in relative
2033 : // distance checks when stitching nodes.
2034 1775 : Real h_min = std::numeric_limits<Real>::max();
2035 50 : bool h_min_updated = false;
2036 :
2037 : // Loop below fills in these sets for the two meshes.
2038 100 : std::set<dof_id_type> this_boundary_node_ids, other_boundary_node_ids;
2039 :
2040 : // Pull objects out of the loop to reduce heap operations
2041 1775 : std::unique_ptr<const Elem> side;
2042 :
2043 : {
2044 : // Make temporary fixed-size arrays for loop
2045 1775 : boundary_id_type id_array[2] = {this_mesh_boundary_id, other_mesh_boundary_id};
2046 1775 : std::set<dof_id_type> * set_array[2] = {&this_boundary_node_ids, &other_boundary_node_ids};
2047 1775 : const MeshBase * mesh_array[2] = {this, other_mesh};
2048 :
2049 5325 : for (unsigned i=0; i<2; ++i)
2050 : {
2051 : // First we deal with node boundary IDs. We only enter
2052 : // this loop if we have at least one nodeset. Note that we
2053 : // do not attempt to make an h_min determination here.
2054 : // The h_min determination is done while looping over the
2055 : // Elems and checking their sides and edges for boundary
2056 : // information, below.
2057 3650 : if (mesh_array[i]->get_boundary_info().n_nodeset_conds() > 0)
2058 : {
2059 : // build_node_list() returns a vector of (node-id, bc-id) tuples
2060 324270 : for (const auto & t : mesh_array[i]->get_boundary_info().build_node_list())
2061 : {
2062 321204 : boundary_id_type node_bc_id = std::get<1>(t);
2063 321204 : if (node_bc_id == id_array[i])
2064 : {
2065 61415 : dof_id_type this_node_id = std::get<0>(t);
2066 61415 : set_array[i]->insert( this_node_id );
2067 : }
2068 : }
2069 : }
2070 :
2071 : // Container to catch boundary IDs passed back from BoundaryInfo.
2072 200 : std::vector<boundary_id_type> bc_ids;
2073 :
2074 : // Pointers to boundary NodeElems encountered while looping over the entire Mesh
2075 : // and checking side and edge boundary ids. The Nodes associated with NodeElems
2076 : // may be in a boundary nodeset, but not connected to any other Elems. In this
2077 : // case, we also consider the "minimum node separation distance" amongst all
2078 : // NodeElems when determining the relevant h_min value for this mesh.
2079 200 : std::vector<const Elem *> boundary_node_elems;
2080 :
2081 86902 : for (auto & el : mesh_array[i]->element_ptr_range())
2082 : {
2083 : // Now check whether elem has a face on the specified boundary
2084 267337 : for (auto side_id : el->side_index_range())
2085 : {
2086 : bool should_stitch_this_side =
2087 231450 : (el->neighbor_ptr(side_id) == nullptr) ||
2088 1420 : (is_valid_disjoint_pair_to_stitch &&
2089 1460 : mesh_array[i]->get_boundary_info().has_boundary_id(el, side_id, id_array[i]));
2090 :
2091 6340 : if (should_stitch_this_side)
2092 : {
2093 : // Get *all* boundary IDs on this side, not just the first one!
2094 99718 : mesh_array[i]->get_boundary_info().boundary_ids (el, side_id, bc_ids);
2095 :
2096 99718 : if (std::find(bc_ids.begin(), bc_ids.end(), id_array[i]) != bc_ids.end())
2097 : {
2098 16472 : el->build_side_ptr(side, side_id);
2099 112535 : for (auto & n : side->node_ref_range())
2100 96063 : set_array[i]->insert(n.id());
2101 :
2102 16472 : h_min = std::min(h_min, side->hmin());
2103 464 : h_min_updated = true;
2104 :
2105 : // This side is on the boundary, add its information to side_to_elem
2106 16472 : if (skip_find_neighbors && (i==0))
2107 : {
2108 8236 : key_type key = el->low_order_key(side_id);
2109 232 : val_type val;
2110 8236 : val.first = el;
2111 232 : val.second = cast_int<unsigned char>(side_id);
2112 :
2113 8236 : key_val_pair kvp;
2114 8236 : kvp.first = key;
2115 232 : kvp.second = val;
2116 232 : side_to_elem_map.insert (kvp);
2117 : }
2118 : }
2119 :
2120 : // Also, check the edges on this side. We don't have to worry about
2121 : // updating neighbor info in this case since elements don't store
2122 : // neighbor info on edges.
2123 1211294 : for (auto edge_id : el->edge_index_range())
2124 : {
2125 1111576 : if (el->is_edge_on_side(edge_id, side_id))
2126 : {
2127 : // Get *all* boundary IDs on this edge, not just the first one!
2128 368490 : mesh_array[i]->get_boundary_info().edge_boundary_ids (el, edge_id, bc_ids);
2129 :
2130 368490 : if (std::find(bc_ids.begin(), bc_ids.end(), id_array[i]) != bc_ids.end())
2131 : {
2132 0 : std::unique_ptr<const Elem> edge (el->build_edge_ptr(edge_id));
2133 0 : for (auto & n : edge->node_ref_range())
2134 0 : set_array[i]->insert( n.id() );
2135 :
2136 0 : h_min = std::min(h_min, edge->hmin());
2137 0 : h_min_updated = true;
2138 0 : }
2139 : }
2140 : } // end for (edge_id)
2141 : } // end if (should_stitch_this_side)
2142 : } // end for (side_id)
2143 :
2144 : // Alternatively, is this a boundary NodeElem? If so,
2145 : // add it to a list of NodeElems that will later be
2146 : // used to set h_min based on the minimum node
2147 : // separation distance between all pairs of boundary
2148 : // NodeElems.
2149 41109 : if (el->type() == NODEELEM)
2150 : {
2151 2700 : mesh_array[i]->get_boundary_info().boundary_ids(el->node_ptr(0), bc_ids);
2152 2628 : if (std::find(bc_ids.begin(), bc_ids.end(), id_array[i]) != bc_ids.end())
2153 : {
2154 2556 : boundary_node_elems.push_back(el);
2155 :
2156 : // Debugging:
2157 : // libMesh::out << "Elem " << el->id() << " is a NodeElem on boundary " << id_array[i] << std::endl;
2158 : }
2159 : } // end if (el->type() == NODEELEM)
2160 3350 : } // end for (el)
2161 :
2162 : // Compute the minimum node separation distance amongst
2163 : // all boundary NodeElem pairs.
2164 : {
2165 200 : const auto N = boundary_node_elems.size();
2166 6106 : for (auto node_elem_i : make_range(N))
2167 47286 : for (auto node_elem_j : make_range(node_elem_i+1, N))
2168 : {
2169 : Real node_sep =
2170 47250 : (boundary_node_elems[node_elem_i]->point(0) - boundary_node_elems[node_elem_j]->point(0)).norm();
2171 :
2172 : // We only want to consider non-coincident
2173 : // boundary NodeElem pairs when determining the
2174 : // minimum node separation distance.
2175 44730 : if (node_sep > 0.)
2176 : {
2177 44730 : h_min = std::min(h_min, node_sep);
2178 1260 : h_min_updated = true;
2179 : }
2180 : } // end for (node_elem_j)
2181 : } // end minimum NodeElem separation scope
2182 : } // end for (i)
2183 : } // end scope
2184 :
2185 1775 : if (verbose)
2186 : {
2187 14 : libMesh::out << "In UnstructuredMesh::stitch_meshes:\n"
2188 14 : << "This mesh has " << this_boundary_node_ids.size()
2189 14 : << " nodes on boundary `"
2190 497 : << this->get_boundary_info().get_sideset_name(this_mesh_boundary_id)
2191 14 : << "' (" << this_mesh_boundary_id << ").\n"
2192 14 : << "Other mesh has " << other_boundary_node_ids.size()
2193 14 : << " nodes on boundary `"
2194 497 : << other_mesh->get_boundary_info().get_sideset_name(other_mesh_boundary_id)
2195 14 : << "' (" << other_mesh_boundary_id << ").\n";
2196 :
2197 497 : if (h_min_updated)
2198 : {
2199 28 : libMesh::out << "Minimum edge length on both surfaces is " << h_min << ".\n";
2200 : }
2201 : else
2202 : {
2203 0 : libMesh::out << "No minimum edge length determined on specified surfaces." << std::endl;
2204 : }
2205 : }
2206 :
2207 : // At this point, if h_min==0 it means that there were at least two coincident
2208 : // nodes on the surfaces being stitched, and we don't currently support that case.
2209 : // (It might be possible to support, but getting it exactly right would be tricky
2210 : // and probably not worth the extra complications to the "normal" case.)
2211 1775 : libmesh_error_msg_if(h_min < std::numeric_limits<Real>::epsilon(),
2212 : "Coincident nodes detected on source and/or target "
2213 : "surface, stitching meshes is not possible.");
2214 :
2215 : // We require nanoflann for the "binary search" (really kd-tree)
2216 : // option to work. If it's not available, turn that option off,
2217 : // warn the user, and fall back on the N^2 search algorithm.
2218 : if (use_binary_search)
2219 : {
2220 : #ifndef LIBMESH_HAVE_NANOFLANN
2221 : use_binary_search = false;
2222 : libmesh_warning("The use_binary_search option in the "
2223 : "UnstructuredMesh stitching algorithms requires nanoflann "
2224 : "support. Falling back on N^2 search algorithm.");
2225 : #endif
2226 : }
2227 :
2228 1775 : if (!this_boundary_node_ids.empty())
2229 : {
2230 1775 : if (use_binary_search)
2231 : {
2232 : #ifdef LIBMESH_HAVE_NANOFLANN
2233 : typedef nanoflann::KDTreeSingleIndexAdaptor<nanoflann::L2_Simple_Adaptor<Real, VectorOfNodesAdaptor>,
2234 : VectorOfNodesAdaptor, 3, std::size_t> kd_tree_t;
2235 :
2236 : // Create the dataset needed to build the kd tree with nanoflann
2237 0 : std::vector<std::pair<Point, dof_id_type>> this_mesh_nodes(this_boundary_node_ids.size());
2238 :
2239 0 : for (auto [it, ctr] = std::make_tuple(this_boundary_node_ids.begin(), 0u);
2240 0 : it != this_boundary_node_ids.end(); ++it, ++ctr)
2241 : {
2242 0 : this_mesh_nodes[ctr].first = this->point(*it);
2243 0 : this_mesh_nodes[ctr].second = *it;
2244 : }
2245 :
2246 0 : VectorOfNodesAdaptor vec_nodes_adaptor(this_mesh_nodes);
2247 :
2248 0 : kd_tree_t this_kd_tree(3, vec_nodes_adaptor, 10);
2249 0 : this_kd_tree.buildIndex();
2250 :
2251 : // Storage for nearest neighbor in the loop below
2252 : std::size_t ret_index;
2253 : Real ret_dist_sqr;
2254 :
2255 : // Loop over other mesh. For each node, find its nearest neighbor in this mesh, and fill in the maps.
2256 0 : for (const auto & node_id : other_boundary_node_ids)
2257 : {
2258 0 : const auto & p = other_mesh->point(node_id);
2259 0 : const Real query_pt[] = {p(0), p(1), p(2)};
2260 0 : this_kd_tree.knnSearch(&query_pt[0], 1, &ret_index, &ret_dist_sqr);
2261 :
2262 : // TODO: here we should use the user's specified tolerance
2263 : // and the previously determined value of h_min in the
2264 : // distance comparison, not just TOLERANCE^2.
2265 0 : if (ret_dist_sqr < TOLERANCE*TOLERANCE)
2266 : {
2267 0 : node_to_node_map[this_mesh_nodes[ret_index].second] = node_id;
2268 0 : other_to_this_node_map[node_id] = this_mesh_nodes[ret_index].second;
2269 : }
2270 : }
2271 :
2272 : // If the two maps don't have the same size, it means one
2273 : // node in this mesh is the nearest neighbor of several
2274 : // nodes in other mesh. Since the stitching is ambiguous in
2275 : // this case, we throw an error.
2276 0 : libmesh_error_msg_if(node_to_node_map.size() != other_to_this_node_map.size(),
2277 : "Error: Found multiple matching nodes in stitch_meshes");
2278 : #endif
2279 : }
2280 : else // !use_binary_search
2281 : {
2282 : // In the unlikely event that two meshes composed entirely of
2283 : // NodeElems are being stitched together, we will not have
2284 : // selected a valid h_min value yet, and the distance
2285 : // comparison below will be true for essentially any two
2286 : // nodes. In this case we simply fall back on an absolute
2287 : // distance check.
2288 1775 : if (!h_min_updated)
2289 : {
2290 : libmesh_warning("No valid h_min value was found, falling back on "
2291 : "absolute distance check in the N^2 search algorithm.");
2292 0 : h_min = 1.;
2293 : }
2294 :
2295 : // Otherwise, use a simple N^2 search to find the closest matching points. This can be helpful
2296 : // in the case that we have tolerance issues which cause mismatch between the two surfaces
2297 : // that are being stitched.
2298 32944 : for (const auto & this_node_id : this_boundary_node_ids)
2299 : {
2300 31169 : Node & this_node = this->node_ref(this_node_id);
2301 :
2302 878 : bool found_matching_nodes = false;
2303 :
2304 928538 : for (const auto & other_node_id : other_boundary_node_ids)
2305 : {
2306 897369 : const Node & other_node = other_mesh->node_ref(other_node_id);
2307 :
2308 872091 : Real node_distance = (this_node - other_node).norm();
2309 :
2310 897369 : if (node_distance < tol*h_min)
2311 : {
2312 : // Make sure we didn't already find a matching node!
2313 31169 : libmesh_error_msg_if(found_matching_nodes,
2314 : "Error: Found multiple matching nodes in stitch_meshes");
2315 :
2316 31169 : node_to_node_map[this_node_id] = other_node_id;
2317 31169 : other_to_this_node_map[other_node_id] = this_node_id;
2318 :
2319 878 : found_matching_nodes = true;
2320 : }
2321 : }
2322 : }
2323 : }
2324 : }
2325 :
2326 : // Build up the node_to_elems_map, using only one loop over other_mesh
2327 42830 : for (auto & el : other_mesh->element_ptr_range())
2328 : {
2329 : // For each node on the element, find the corresponding node
2330 : // on "this" Mesh, 'this_node_id', if it exists, and push
2331 : // the current element ID back onto node_to_elems_map[this_node_id].
2332 : // For that we will use the reverse mapping we created at
2333 : // the same time as the forward mapping.
2334 294865 : for (auto & n : el->node_ref_range())
2335 281780 : if (const auto it = other_to_this_node_map.find(/*other_node_id=*/n.id());
2336 7720 : it != other_to_this_node_map.end())
2337 50481 : node_to_elems_map[/*this_node_id=*/it->second].push_back( el->id() );
2338 1675 : }
2339 :
2340 1775 : if (verbose)
2341 : {
2342 14 : libMesh::out << "In UnstructuredMesh::stitch_meshes:\n"
2343 28 : << "Found " << node_to_node_map.size()
2344 14 : << " matching nodes.\n"
2345 14 : << std::endl;
2346 : }
2347 :
2348 1775 : if (enforce_all_nodes_match_on_boundaries)
2349 : {
2350 2 : std::size_t n_matching_nodes = node_to_node_map.size();
2351 2 : std::size_t this_mesh_n_nodes = this_boundary_node_ids.size();
2352 2 : std::size_t other_mesh_n_nodes = other_boundary_node_ids.size();
2353 71 : libmesh_error_msg_if((n_matching_nodes != this_mesh_n_nodes) || (n_matching_nodes != other_mesh_n_nodes),
2354 : "Error: We expected the number of nodes to match.");
2355 : }
2356 :
2357 1775 : if (merge_boundary_nodes_all_or_nothing)
2358 : {
2359 2 : std::size_t n_matching_nodes = node_to_node_map.size();
2360 2 : std::size_t this_mesh_n_nodes = this_boundary_node_ids.size();
2361 2 : std::size_t other_mesh_n_nodes = other_boundary_node_ids.size();
2362 71 : if ((n_matching_nodes != this_mesh_n_nodes) || (n_matching_nodes != other_mesh_n_nodes))
2363 : {
2364 0 : if (verbose)
2365 : {
2366 : libMesh::out << "Skipping node merging in "
2367 : "UnstructuredMesh::stitch_meshes because not "
2368 0 : "all boundary nodes were matched."
2369 0 : << std::endl;
2370 : }
2371 0 : node_to_node_map.clear();
2372 0 : other_to_this_node_map.clear();
2373 0 : node_to_elems_map.clear();
2374 : }
2375 100 : }
2376 3350 : }
2377 : else
2378 : {
2379 0 : if (verbose)
2380 : {
2381 0 : libMesh::out << "Skip node merging in UnstructuredMesh::stitch_meshes:" << std::endl;
2382 : }
2383 : }
2384 :
2385 1775 : dof_id_type node_delta = this->max_node_id();
2386 1775 : dof_id_type elem_delta = this->max_elem_id();
2387 :
2388 : unique_id_type unique_delta =
2389 : #ifdef LIBMESH_ENABLE_UNIQUE_ID
2390 1775 : this->parallel_max_unique_id();
2391 : #else
2392 : 0;
2393 : #endif
2394 :
2395 : // If other_mesh != nullptr, then we have to do a bunch of work
2396 : // in order to copy it to this mesh
2397 1775 : if (this!=other_mesh)
2398 : {
2399 88 : LOG_SCOPE("stitch_meshes copying", "UnstructuredMesh");
2400 :
2401 : #ifdef LIBMESH_ENABLE_PERIODIC
2402 : // Copy disjoint neighbor boundary pairs (PeriodicBoundary objects)
2403 : // from `other_mesh` to `this` mesh
2404 1562 : if (other_db && !other_db->empty())
2405 : {
2406 213 : for (const auto & [bdy_id, pb_ptr] : *other_db)
2407 : {
2408 4 : const auto & pb = *pb_ptr;
2409 142 : const boundary_id_type a = pb.myboundary;
2410 142 : const boundary_id_type b = pb.pairedboundary;
2411 :
2412 142 : if (this_db)
2413 : {
2414 : // Skip if identical pair already exists
2415 142 : if (const auto * existing_pb = this_db->boundary(a))
2416 71 : if ((existing_pb->myboundary == a && existing_pb->pairedboundary == b) ||
2417 0 : (existing_pb->myboundary == b && existing_pb->pairedboundary == a))
2418 69 : continue;
2419 :
2420 : // If both boundary ids exist on this mesh but aren't paired here, refuse to create a new pair
2421 2 : const auto & bdy_ids = this->get_boundary_info().get_boundary_ids();
2422 4 : const bool a_exists = bdy_ids.count(a);
2423 4 : const bool b_exists = bdy_ids.count(b);
2424 : // If a and b already exist on `this`, we should be screaming and dying
2425 : // unless they already have PeriodicBoundary objects connecting them too
2426 71 : if (a_exists && b_exists && !this_db->boundary(a))
2427 0 : libmesh_error_msg("Conflict: boundaries " << a << " and " << b
2428 : << " already exist on this mesh but are not paired.");
2429 : }
2430 :
2431 140 : this->add_disjoint_neighbor_boundary_pairs(a, b, pb.get_corresponding_pos(Point(0.0,0.0,0.0)));
2432 : }
2433 : }
2434 : #endif // LIBMESH_ENABLE_PERIODIC
2435 :
2436 :
2437 : // Increment the node_to_node_map and node_to_elems_map
2438 : // to account for id offsets
2439 32305 : for (auto & pr : node_to_node_map)
2440 30743 : pr.second += node_delta;
2441 :
2442 32305 : for (auto & pr : node_to_elems_map)
2443 80798 : for (auto & entry : pr.second)
2444 50055 : entry += elem_delta;
2445 :
2446 : // We run into problems when the libMesh subdomain standard (the
2447 : // id defines the subdomain; the name was an afterthought) and
2448 : // the MOOSE standard (the name defines the subdomain; the id
2449 : // might be autogenerated) clash.
2450 : //
2451 : // Subdomain ids with the same name in both meshes are surely
2452 : // meant to represent the same subdomain. We can just merge
2453 : // them.
2454 : //
2455 : // Subdomain ids which don't have a name in either mesh are
2456 : // almost surely meant to represent the same subdomain. We'll
2457 : // just merge them.
2458 : //
2459 : // Subdomain ids with different names in different meshes, or
2460 : // names with different ids in different meshes, are trickier.
2461 : // For backwards compatibility we default to the old "just copy
2462 : // all the subdomain ids over" behavior, but if requested we'll
2463 : // remap any ids that appear to be clear conflicts, and we'll
2464 : // scream and die if we see any ids that are ambiguous due to
2465 : // being named in one mesh but not the other.
2466 88 : std::unordered_map<subdomain_id_type, subdomain_id_type> id_remapping;
2467 1562 : if (remap_subdomain_ids)
2468 : {
2469 4 : const auto & this_map = this->get_subdomain_name_map();
2470 4 : const auto & other_map = other_mesh->get_subdomain_name_map();
2471 8 : std::unordered_map<std::string, subdomain_id_type> other_map_reversed;
2472 284 : for (auto & [sid, sname] : other_map)
2473 4 : other_map_reversed.emplace(sname, sid);
2474 :
2475 8 : std::unordered_map<std::string, subdomain_id_type> this_map_reversed;
2476 213 : for (auto & [sid, sname] : this_map)
2477 2 : this_map_reversed.emplace(sname, sid);
2478 :
2479 : // We don't require either mesh to be prepared, but that
2480 : // means we need to check for subdomains manually.
2481 284 : auto get_subdomains = [](const MeshBase & mesh) {
2482 8 : std::set<subdomain_id_type> all_subdomains;
2483 4976 : for (auto & el : mesh.element_ptr_range())
2484 2540 : all_subdomains.insert(el->subdomain_id());
2485 284 : return all_subdomains;
2486 : };
2487 :
2488 146 : const auto this_subdomains = get_subdomains(*this);
2489 146 : const auto other_subdomains = get_subdomains(*other_mesh);
2490 :
2491 213 : for (auto & [sid, sname] : this_map)
2492 : {
2493 : // The same name with the same id means we're fine. The
2494 : // same name with another id means we remap their id to
2495 : // ours
2496 2 : if (const auto other_reverse_it = other_map_reversed.find(sname);
2497 71 : other_reverse_it != other_map_reversed.end() && other_reverse_it->second != sid)
2498 71 : id_remapping[other_reverse_it->second] = sid;
2499 :
2500 : // The same id with a different name, we'll get to
2501 : // later. The same id without any name means we don't
2502 : // know what the user wants.
2503 2 : if (other_subdomains.count(sid) && !other_map.count(sid))
2504 0 : libmesh_error_msg("Can't safely stitch with a mesh sharing subdomain id "
2505 : << sid << " but not subdomain name " << sname);
2506 : }
2507 :
2508 142 : subdomain_id_type next_free_id = 0;
2509 : // We might try to stitch empty meshes ...
2510 142 : if (!this_subdomains.empty())
2511 142 : next_free_id = *this_subdomains.rbegin() + 1;
2512 142 : if (!other_subdomains.empty())
2513 142 : next_free_id =
2514 142 : std::max(next_free_id,
2515 : cast_int<subdomain_id_type>
2516 215 : (*other_subdomains.rbegin() + 1));
2517 :
2518 213 : for (auto & [sid, sname] : other_map)
2519 : {
2520 : // At this point we've figured out any remapping
2521 : // necessary for an sname that we share. And we don't
2522 : // need to remap any sid we don't share.
2523 8 : if (!this_map_reversed.count(sname))
2524 : {
2525 : // But if we don't have this sname and we do have this
2526 : // sid then we can't just merge into that.
2527 2 : if (this_subdomains.count(sid))
2528 : {
2529 : // If we have this sid with no name, we don't
2530 : // know what the user wants.
2531 2 : if (!this_map.count(sid))
2532 211 : libmesh_error_msg("Can't safely stitch with a mesh sharing subdomain id "
2533 : << sid << " but under subdomain name " << sname);
2534 :
2535 : // We have this sid under a different name, so
2536 : // we just need to give the other elements a new
2537 : // id.
2538 :
2539 : // Users might have done crazy things with id
2540 : // choice so let's make sure they didn't get too
2541 : // crazy.
2542 0 : libmesh_error_msg_if ((!this_subdomains.empty() &&
2543 : next_free_id < *this_subdomains.rbegin()) ||
2544 : (!other_subdomains.empty() &&
2545 : next_free_id < *other_subdomains.rbegin()),
2546 : "Subdomain id overflow");
2547 :
2548 0 : id_remapping[sid] = next_free_id++;
2549 0 : this->subdomain_name(next_free_id) = sname;
2550 : }
2551 : // If we don't have this subdomain id, well, we're
2552 : // about to, so we should have its name too.
2553 : else
2554 0 : this->subdomain_name(sid) = sname;
2555 : }
2556 : }
2557 : }
2558 :
2559 : // Copy mesh data. If we skip the call to find_neighbors(), the lists
2560 : // of neighbors will be copied verbatim from the other mesh
2561 1491 : this->copy_nodes_and_elements(*other_mesh, skip_find_neighbors,
2562 : elem_delta, node_delta,
2563 84 : unique_delta, &id_remapping);
2564 :
2565 : // Copy BoundaryInfo from other_mesh too. We do this via the
2566 : // list APIs rather than element-by-element for speed.
2567 42 : BoundaryInfo & boundary = this->get_boundary_info();
2568 42 : const BoundaryInfo & other_boundary = other_mesh->get_boundary_info();
2569 :
2570 158869 : for (const auto & t : other_boundary.build_node_list())
2571 157336 : boundary.add_node(std::get<0>(t) + node_delta,
2572 157336 : std::get<1>(t));
2573 :
2574 42855 : for (const auto & t : other_boundary.build_side_list())
2575 42486 : boundary.add_side(std::get<0>(t) + elem_delta,
2576 41322 : std::get<1>(t),
2577 41322 : std::get<2>(t));
2578 :
2579 1533 : for (const auto & t : other_boundary.build_edge_list())
2580 0 : boundary.add_edge(std::get<0>(t) + elem_delta,
2581 0 : std::get<1>(t),
2582 0 : std::get<2>(t));
2583 :
2584 1533 : for (const auto & t : other_boundary.build_shellface_list())
2585 0 : boundary.add_shellface(std::get<0>(t) + elem_delta,
2586 0 : std::get<1>(t),
2587 0 : std::get<2>(t));
2588 :
2589 42 : const auto & other_ns_id_to_name = other_boundary.get_nodeset_name_map();
2590 42 : auto & ns_id_to_name = boundary.set_nodeset_name_map();
2591 1491 : ns_id_to_name.insert(other_ns_id_to_name.begin(), other_ns_id_to_name.end());
2592 :
2593 42 : const auto & other_ss_id_to_name = other_boundary.get_sideset_name_map();
2594 42 : auto & ss_id_to_name = boundary.set_sideset_name_map();
2595 1491 : ss_id_to_name.insert(other_ss_id_to_name.begin(), other_ss_id_to_name.end());
2596 :
2597 42 : const auto & other_es_id_to_name = other_boundary.get_edgeset_name_map();
2598 42 : auto & es_id_to_name = boundary.set_edgeset_name_map();
2599 1491 : es_id_to_name.insert(other_es_id_to_name.begin(), other_es_id_to_name.end());
2600 :
2601 : // Merge other_mesh's elemset information with ours. Throw an
2602 : // error if this and other_mesh have overlapping elemset codes
2603 : // that refer to different elemset ids.
2604 1533 : std::vector<dof_id_type> this_elemset_codes = this->get_elemset_codes();
2605 84 : MeshBase::elemset_type this_id_set_to_fill, other_id_set_to_fill;
2606 1817 : for (const auto & elemset_code : other_mesh->get_elemset_codes())
2607 : {
2608 : // Get the elemset ids for this elemset_code on other_mesh
2609 284 : other_mesh->get_elemsets(elemset_code, other_id_set_to_fill);
2610 :
2611 : // Check that this elemset code does not already exist
2612 : // in this mesh, or if it does, that it has the same elemset
2613 : // ids associated with it.
2614 : //
2615 : // Note: get_elemset_codes() is guaranteed to return a
2616 : // sorted vector, so we can binary search in it.
2617 268 : auto it = Utility::binary_find(this_elemset_codes.begin(),
2618 : this_elemset_codes.end(),
2619 16 : elemset_code);
2620 :
2621 284 : if (it != this_elemset_codes.end())
2622 : {
2623 : // This mesh has the same elemset code. Does it refer to
2624 : // the same elemset ids?
2625 0 : this->get_elemsets(elemset_code, this_id_set_to_fill);
2626 :
2627 : // Throw an error if they don't match, otherwise we
2628 : // don't need to do anything
2629 0 : libmesh_error_msg_if(other_id_set_to_fill != this_id_set_to_fill,
2630 : "Attempted to stitch together meshes with conflicting elemset codes.");
2631 : }
2632 : else
2633 : {
2634 : // Add other_mesh's elemset code to this mesh
2635 560 : this->add_elemset_code(elemset_code, other_id_set_to_fill);
2636 : }
2637 : }
2638 :
2639 : } // end if (other_mesh)
2640 :
2641 : // Finally, we need to "merge" the overlapping nodes
2642 : // We do this by iterating over node_to_elems_map and updating
2643 : // the elements so that they "point" to the nodes that came
2644 : // from this mesh, rather than from other_mesh.
2645 : // Then we iterate over node_to_node_map and delete the
2646 : // duplicate nodes that came from other_mesh.
2647 :
2648 : {
2649 96 : LOG_SCOPE("stitch_meshes node updates", "UnstructuredMesh");
2650 :
2651 : // Container to catch boundary IDs passed back from BoundaryInfo.
2652 96 : std::vector<boundary_id_type> bc_ids;
2653 :
2654 32234 : for (const auto & [target_node_id, elem_vec] : node_to_elems_map)
2655 : {
2656 30530 : dof_id_type other_node_id = node_to_node_map[target_node_id];
2657 30530 : Node & target_node = this->node_ref(target_node_id);
2658 :
2659 1720 : std::size_t n_elems = elem_vec.size();
2660 79875 : for (std::size_t i=0; i<n_elems; i++)
2661 : {
2662 49345 : dof_id_type elem_id = elem_vec[i];
2663 49345 : Elem * el = this->elem_ptr(elem_id);
2664 :
2665 : // find the local node index that we want to update
2666 47955 : unsigned int local_node_index = el->local_node(other_node_id);
2667 1390 : libmesh_assert_not_equal_to(local_node_index, libMesh::invalid_uint);
2668 :
2669 : // We also need to copy over the nodeset info here,
2670 : // because the node will get deleted below
2671 50735 : this->get_boundary_info().boundary_ids(el->node_ptr(local_node_index), bc_ids);
2672 49345 : el->set_node(local_node_index, &target_node);
2673 49345 : this->get_boundary_info().add_node(&target_node, bc_ids);
2674 : }
2675 : }
2676 : }
2677 :
2678 : {
2679 96 : LOG_SCOPE("stitch_meshes node deletion", "UnstructuredMesh");
2680 32234 : for (const auto & [other_node_id, this_node_id] : node_to_node_map)
2681 : {
2682 : // In the case that this==other_mesh, the two nodes might be the same (e.g. if
2683 : // we're stitching a "sliver"), hence we need to skip node deletion in that case.
2684 30530 : if ((this == other_mesh) && (this_node_id == other_node_id))
2685 0 : continue;
2686 :
2687 30530 : this->delete_node( this->node_ptr(this_node_id) );
2688 : }
2689 : }
2690 :
2691 : // If find_neighbors() wasn't called in prepare_for_use(), we need to
2692 : // manually loop once more over all elements adjacent to the stitched boundary
2693 : // and fix their lists of neighbors.
2694 : // This is done according to the following steps:
2695 : // 1. Loop over all copied elements adjacent to the boundary using node_to_elems_map (trying to avoid duplicates)
2696 : // 2. Look at all their sides with a nullptr neighbor and update them using side_to_elem_map if necessary
2697 : // 3. Update the corresponding side in side_to_elem_map as well
2698 1704 : if (skip_find_neighbors)
2699 : {
2700 96 : LOG_SCOPE("stitch_meshes neighbor fixes", "UnstructuredMesh");
2701 :
2702 : // Pull objects out of the loop to reduce heap operations
2703 1704 : std::unique_ptr<const Elem> my_side, their_side;
2704 :
2705 96 : std::set<dof_id_type> fixed_elems;
2706 32234 : for (const auto & pr : node_to_elems_map)
2707 : {
2708 1720 : std::size_t n_elems = pr.second.size();
2709 79875 : for (std::size_t i=0; i<n_elems; i++)
2710 : {
2711 50735 : dof_id_type elem_id = pr.second[i];
2712 1390 : if (!fixed_elems.count(elem_id))
2713 : {
2714 10508 : Elem * el = this->elem_ptr(elem_id);
2715 10212 : fixed_elems.insert(elem_id);
2716 57226 : for (auto s : el->side_index_range())
2717 : {
2718 46718 : bool has_real_neighbor = (el->neighbor_ptr(s) != nullptr);
2719 47570 : bool has_disdjoint_neighbor = is_valid_disjoint_pair_to_stitch &&
2720 852 : (this->get_boundary_info().has_boundary_id(el, s, this_mesh_boundary_id)
2721 852 : || this->get_boundary_info().has_boundary_id(el, s, other_mesh_boundary_id));
2722 :
2723 46718 : if (!has_real_neighbor || has_disdjoint_neighbor)
2724 : {
2725 21016 : key_type key = el->low_order_key(s);
2726 592 : auto bounds = side_to_elem_map.equal_range(key);
2727 :
2728 21016 : if (bounds.first != bounds.second)
2729 : {
2730 : // Get the side for this element
2731 7952 : el->side_ptr(my_side, s);
2732 :
2733 : // Look at all the entries with an equivalent key
2734 7952 : while (bounds.first != bounds.second)
2735 : {
2736 : // Get the potential element
2737 7952 : Elem * neighbor = const_cast<Elem *>(bounds.first->second.first);
2738 :
2739 : // Get the side for the neighboring element
2740 7952 : const unsigned int ns = bounds.first->second.second;
2741 7952 : neighbor->side_ptr(their_side, ns);
2742 : //libmesh_assert(my_side.get());
2743 : //libmesh_assert(their_side.get());
2744 :
2745 : // If found a match with my side
2746 : //
2747 : // We need special tests here for 1D:
2748 : // since parents and children have an equal
2749 : // side (i.e. a node), we need to check
2750 : // ns != ms, and we also check level() to
2751 : // avoid setting our neighbor pointer to
2752 : // any of our neighbor's descendants
2753 15680 : if ((*my_side == *their_side) &&
2754 15904 : (el->level() == neighbor->level()) &&
2755 7952 : ((el->dim() != 1) || (ns != s)))
2756 : {
2757 : // So share a side. Is this a mixed pair
2758 : // of subactive and active/ancestor
2759 : // elements?
2760 : // If not, then we're neighbors.
2761 : // If so, then the subactive's neighbor is
2762 :
2763 8176 : if (el->subactive() ==
2764 7952 : neighbor->subactive())
2765 : {
2766 : // an element is only subactive if it has
2767 : // been coarsened but not deleted
2768 448 : el->set_neighbor (s,neighbor);
2769 448 : neighbor->set_neighbor(ns,el);
2770 : }
2771 0 : else if (el->subactive())
2772 : {
2773 0 : el->set_neighbor(s,neighbor);
2774 : }
2775 0 : else if (neighbor->subactive())
2776 : {
2777 0 : neighbor->set_neighbor(ns,el);
2778 : }
2779 : // It's OK to invalidate the
2780 : // bounds.first iterator here,
2781 : // as we are immediately going
2782 : // to break out of this while
2783 : // loop. bounds.first will
2784 : // therefore not be used for
2785 : // anything else.
2786 224 : side_to_elem_map.erase (bounds.first);
2787 224 : break;
2788 : }
2789 :
2790 0 : ++bounds.first;
2791 : }
2792 : }
2793 : }
2794 : }
2795 : }
2796 : }
2797 : }
2798 1608 : }
2799 :
2800 : #ifdef LIBMESH_ENABLE_PERIODIC
2801 : // Remove only the disjoint pair that was actually stitched.
2802 : // Safe because `is_valid_disjoint_pair_to_stitch` is true
2803 : // only if this exact (a,b) pair exists in the registry.
2804 : // Other disjoint pairs remain untouched.
2805 1704 : if (is_valid_disjoint_pair_to_stitch)
2806 213 : this->remove_disjoint_boundary_pair(this_mesh_boundary_id, other_mesh_boundary_id);
2807 : #endif
2808 :
2809 1704 : if (prepare_after_stitching)
2810 : {
2811 : // We set our new neighbor pointers already
2812 92 : const bool old_allow_find_neighbors = this->allow_find_neighbors();
2813 46 : this->allow_find_neighbors(!skip_find_neighbors);
2814 :
2815 : // We haven't newly remoted any elements
2816 92 : const bool old_allow_remote_element_removal = this->allow_remote_element_removal();
2817 46 : this->allow_remote_element_removal(false);
2818 :
2819 1633 : this->prepare_for_use();
2820 :
2821 46 : this->allow_find_neighbors(old_allow_find_neighbors);
2822 46 : this->allow_remote_element_removal(old_allow_remote_element_removal);
2823 : }
2824 :
2825 : // After the stitching, we may want to clear boundary IDs from element
2826 : // faces that are now internal to the mesh
2827 1704 : if (clear_stitched_boundary_ids)
2828 : {
2829 92 : LOG_SCOPE("stitch_meshes clear bcids", "UnstructuredMesh");
2830 :
2831 1633 : this->get_boundary_info().clear_stitched_boundary_side_ids(
2832 : this_mesh_boundary_id, other_mesh_boundary_id, /*clear_nodeset_data=*/true);
2833 : }
2834 :
2835 : // Return the number of nodes which were merged.
2836 1752 : return node_to_node_map.size();
2837 1742 : }
2838 :
2839 :
2840 : } // namespace libMesh
|