MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
submesh_utils.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "submesh_utils.hpp"
13#include "ncsubmesh.hpp"
14#include "submesh.hpp"
15#include "pncsubmesh.hpp"
16#include "psubmesh.hpp"
17
18#include <numeric>
19
20namespace mfem
21{
22namespace SubMeshUtils
23{
24int UniqueIndexGenerator::Get(int i, bool &new_index)
25{
26 auto f = idx.find(i);
27 if (f == idx.end())
28 {
29 idx[i] = counter;
30 new_index = true;
31 return counter++;
32 }
33 else
34 {
35 new_index = false;
36 return (*f).second;
37 }
38}
39
40template <typename ElementT>
41bool ElementHasAttribute(const ElementT &el, const Array<int> &attributes)
42{
43 for (int a = 0; a < attributes.Size(); a++)
44 {
45 if (el.GetAttribute() == attributes[a])
46 {
47 return true;
48 }
49 }
50 return false;
51}
52
53std::tuple< Array<int>, Array<int> >
54AddElementsToMesh(const Mesh& parent,
55 Mesh& mesh,
56 const Array<int> &attributes,
57 bool from_boundary)
58{
59 UniqueIndexGenerator vertex_ids;
60 Array<int> parent_vertex_ids, parent_element_ids;
61 Array<int> vert, submesh_vert;
62
63 const int ne = from_boundary ? parent.GetNBE() : parent.GetNE();
64 for (int i = 0; i < ne; i++)
65 {
66 const Element *pel = from_boundary ?
67 parent.GetBdrElement(i) : parent.GetElement(i);
68 if (!HasAttribute(*pel, attributes)) { continue; }
69 pel->GetVertices(vert);
70 submesh_vert.SetSize(vert.Size());
71 for (int iv = 0; iv < vert.Size(); iv++)
72 {
73 bool new_vertex;
74 int mesh_vertex_id = vert[iv];
75 int submesh_vertex_id = vertex_ids.Get(mesh_vertex_id, new_vertex);
76 if (new_vertex)
77 {
78 mesh.AddVertex(parent.GetVertex(mesh_vertex_id));
79 parent_vertex_ids.Append(mesh_vertex_id);
80 }
81 submesh_vert[iv] = submesh_vertex_id;
82 }
83 Element *el = mesh.NewElement(from_boundary ?
84 parent.GetBdrElementType(i) : parent.GetElementType(i));
85 el->SetVertices(submesh_vert);
86 el->SetAttribute(pel->GetAttribute());
87 mesh.AddElement(el);
88 parent_element_ids.Append(i);
89 }
90 return {parent_vertex_ids, parent_element_ids};
91}
92
94 const FiniteElementSpace& parentfes,
95 const SubMesh::From& from,
96 const Array<int>& parent_element_ids,
97 Array<int>& vdof_to_vdof_map)
98{
99 auto *m = subfes.GetMesh();
100 vdof_to_vdof_map.SetSize(subfes.GetVSize());
101 const int vdim = parentfes.GetVDim();
102
104 DenseMatrix T;
105 Array<int> z1;
106 for (int i = 0; i < m->GetNE(); i++)
107 {
108 Array<int> parent_vdofs;
109 if (from == SubMesh::From::Domain)
110 {
111 parentfes.GetElementVDofs(parent_element_ids[i], parent_vdofs);
112 }
113 else if (from == SubMesh::From::Boundary)
114 {
115 if (parentfes.IsDGSpace())
116 {
117 MFEM_ASSERT(static_cast<const L2_FECollection*>
118 (parentfes.FEColl())->GetBasisType() == BasisType::GaussLobatto,
119 "Only BasisType::GaussLobatto is supported for L2 spaces");
120
121 auto pm = parentfes.GetMesh();
122
123 const Geometry::Type face_geom =
124 pm->GetBdrElementGeometry(parent_element_ids[i]);
125 int face_info, parent_volel_id;
126 pm->GetBdrElementAdjacentElement(
127 parent_element_ids[i], parent_volel_id, face_info);
128 face_info = Mesh::EncodeFaceInfo(
131 face_geom, Mesh::DecodeFaceInfoOrientation(face_info)));
132 pm->GetLocalFaceTransformation(
133 pm->GetBdrElementType(parent_element_ids[i]),
134 pm->GetElementType(parent_volel_id),
135 Tr.Transf,
136 face_info);
137
138 const FiniteElement *face_el =
139 parentfes.GetTraceElement(parent_volel_id, face_geom);
140 MFEM_VERIFY(dynamic_cast<const NodalFiniteElement*>(face_el),
141 "Nodal Finite Element is required");
142
143 face_el->GetTransferMatrix(*parentfes.GetFE(parent_volel_id),
144 Tr.Transf,
145 T);
146
147 parentfes.GetElementVDofs(parent_volel_id, z1);
148
149 parent_vdofs.SetSize(vdim * T.Height());
150 for (int j = 0; j < T.Height(); j++)
151 {
152 for (int k = 0; k < T.Width(); k++)
153 {
154 if (T(j, k) != 0.0)
155 {
156 for (int vd=0; vd<vdim; vd++)
157 {
158 int sub_vdof = j + T.Height() * vd;
159 int parent_vdof = k + T.Width() * vd;
160 parent_vdofs[sub_vdof] =
161 z1[static_cast<int>(parent_vdof)];
162 }
163 }
164 }
165 }
166 }
167 else
168 {
169 parentfes.GetBdrElementVDofs(parent_element_ids[i], parent_vdofs);
170 }
171 }
172 else
173 {
174 MFEM_ABORT("SubMesh::From type unknown");
175 }
176
177 Array<int> sub_vdofs;
178 subfes.GetElementVDofs(i, sub_vdofs);
179
180 MFEM_ASSERT(parent_vdofs.Size() == sub_vdofs.Size(),
181 "elem " << i << ' ' << parent_vdofs.Size() << ' ' << sub_vdofs.Size());
182 for (int j = 0; j < parent_vdofs.Size(); j++)
183 {
184 real_t sub_sign = 1.0;
185 int sub_vdof = subfes.DecodeDof(sub_vdofs[j], sub_sign);
186
187 real_t parent_sign = 1.0;
188 int parent_vdof = parentfes.DecodeDof(parent_vdofs[j], parent_sign);
189
190 vdof_to_vdof_map[sub_vdof] =
191 (sub_sign * parent_sign > 0.0) ? parent_vdof : (-1-parent_vdof);
192 }
193 }
194
195#ifdef MFEM_DEBUG
196 auto tmp = vdof_to_vdof_map;
197 tmp.Sort();
198 tmp.Unique();
199
200 if (tmp.Size() != vdof_to_vdof_map.Size())
201 {
202 std::stringstream msg;
203 for (int i = 0; i < vdof_to_vdof_map.Size(); i++)
204 for (int j = i + 1; j < vdof_to_vdof_map.Size(); j++)
205 {
206 auto x = vdof_to_vdof_map[i];
207 auto y = vdof_to_vdof_map[j];
208 if (x == y)
209 {
210 msg << "i " << i << " (" << x << ") j " << j << " (" << y << ")\n";
211 }
212 }
213 MFEM_ABORT("vdof_to_vdof_map should be 1 to 1:\n" << msg.str());
214 }
215#endif
216
217}
218
219Array<int> BuildFaceMap(const Mesh& pm, const Mesh& sm,
220 const Array<int> &parent_element_ids)
221{
222 // TODO: Check if parent is really a parent of mesh
223
224 Array<int> pfids(sm.GetNumFaces());
225 pfids = -1;
226 for (int i = 0; i < sm.GetNE(); i++)
227 {
228 int peid = parent_element_ids[i];
229
230 Array<int> sel_faces, pel_faces, o;
231 if (pm.Dimension() == 2)
232 {
233 sm.GetElementEdges(i, sel_faces, o);
234 pm.GetElementEdges(peid, pel_faces, o);
235 }
236 else
237 {
238 sm.GetElementFaces(i, sel_faces, o);
239 pm.GetElementFaces(peid, pel_faces, o);
240 }
241
242 MFEM_ASSERT(sel_faces.Size() == pel_faces.Size(), "internal error");
243
244 for (int j = 0; j < sel_faces.Size(); j++)
245 {
246 if (pfids[sel_faces[j]] != -1)
247 {
248 MFEM_ASSERT(pfids[sel_faces[j]] == pel_faces[j], "internal error");
249 }
250 pfids[sel_faces[j]] = pel_faces[j];
251 }
252 }
253 return pfids;
254}
255
256template <typename SubMeshT>
257void AddBoundaryElements(SubMeshT &mesh,
258 const std::unordered_map<int,int> &lface_to_boundary_attribute)
259{
260 const int num_codim_1 = mesh.GetNumFaces();
261
262 if (mesh.Dimension() == 3)
263 {
264 // In 3D we check for `bel_to_edge`. It shouldn't have been set
265 // previously.
266 mesh.DeleteBoundaryElementToEdge();
267 }
268 int NumOfBdrElements = 0;
269 for (int i = 0; i < num_codim_1; i++)
270 {
271 if (mesh.GetFaceInformation(i).IsBoundary())
272 {
273 NumOfBdrElements++;
274 }
275 }
276
277 Array<Element *> boundary;
278 Array<int> be_to_face;
279 boundary.Reserve(NumOfBdrElements);
280 be_to_face.Reserve(NumOfBdrElements);
281
282 const auto &parent = *mesh.GetParent();
283 const auto &parent_face_ids = mesh.GetParentFaceIDMap();
284 const auto &parent_edge_ids = mesh.GetParentEdgeIDMap();
285 const auto &parent_vertex_ids = mesh.GetParentVertexIDMap();
286 const auto &parent_face_to_be = parent.GetFaceToBdrElMap();
287 const auto &face_to_be = mesh.GetFaceToBdrElMap();
288 int max_bdr_attr = parent.bdr_attributes.Max();
289 for (int i = 0; i < num_codim_1; i++)
290 {
291 auto pfid = [&](int i)
292 {
293 switch (mesh.Dimension())
294 {
295 case 3: return parent_face_ids[i];
296 case 2: return parent_edge_ids[i];
297 case 1: return parent_vertex_ids[i];
298 }
299 MFEM_ABORT("!");
300 return -1;
301 };
302 if (mesh.GetFaceInformation(i).IsBoundary()
303 && (face_to_be.IsEmpty() || face_to_be[i] == -1))
304 {
305 auto * be = mesh.GetFace(i)->Duplicate(&mesh);
306
307 if (mesh.GetFrom() == SubMesh::From::Domain && mesh.Dimension() >= 2)
308 {
309 int pbeid = parent_face_to_be[pfid(i)];
310 if (pbeid != -1)
311 {
312 be->SetAttribute(parent.GetBdrAttribute(pbeid));
313 }
314 else
315 {
316 auto ghost_attr = lface_to_boundary_attribute.find(pfid(i));
317 int battr = ghost_attr != lface_to_boundary_attribute.end() ?
318 ghost_attr->second : max_bdr_attr + 1;
319 be->SetAttribute(battr);
320 }
321 }
322 else
323 {
324 auto ghost_attr = lface_to_boundary_attribute.find(pfid(i));
325 int battr = ghost_attr != lface_to_boundary_attribute.end() ?
326 ghost_attr->second : max_bdr_attr + 1;
327 be->SetAttribute(battr);
328 }
329 be_to_face.Append(i);
330 boundary.Append(be);
331 }
332 }
333
334 if (mesh.GetFrom() == SubMesh::From::Domain && mesh.Dimension() >= 2)
335 {
336 // Search for and count interior boundary elements
337 int InteriorBdrElems = 0;
338 for (int i=0; i<parent.GetNBE(); i++)
339 {
340 const int parentFaceIdx = parent.GetBdrElementFaceIndex(i);
341 const int submeshFaceIdx =
342 mesh.Dimension() == 3 ?
343 mesh.GetSubMeshFaceFromParent(parentFaceIdx) :
344 mesh.GetSubMeshEdgeFromParent(parentFaceIdx);
345
346 if (submeshFaceIdx == -1) { continue; }
347 if (mesh.GetFaceInformation(submeshFaceIdx).IsBoundary()) { continue; }
348 InteriorBdrElems++;
349 }
350
351 if (InteriorBdrElems > 0)
352 {
353 NumOfBdrElements += InteriorBdrElems;
354 boundary.Reserve(NumOfBdrElements);
355 be_to_face.Reserve(NumOfBdrElements);
356
357 // Search for and transfer interior boundary elements
358 for (int i = 0; i < parent.GetNBE(); i++)
359 {
360 const int parentFaceIdx = parent.GetBdrElementFaceIndex(i);
361 const int submeshFaceIdx =
362 mesh.Dimension() == 3 ?
363 mesh.GetSubMeshFaceFromParent(parentFaceIdx) :
364 mesh.GetSubMeshEdgeFromParent(parentFaceIdx);
365
366 if (submeshFaceIdx == -1) { continue; }
367 if (mesh.GetFaceInformation(submeshFaceIdx).IsBoundary())
368 { continue; }
369
370 auto * be = mesh.GetFace(submeshFaceIdx)->Duplicate(&mesh);
371 be->SetAttribute(parent.GetBdrAttribute(i));
372 boundary.Append(be);
373 be_to_face.Append(submeshFaceIdx);
374 }
375 }
376 }
377 mesh.AddBdrElements(boundary, be_to_face);
378}
379
380// Explicit instantiations
381template void AddBoundaryElements(SubMesh &mesh,
382 const std::unordered_map<int,int> &);
383
384#ifdef MFEM_USE_MPI
385template void AddBoundaryElements(ParSubMesh &mesh,
386 const std::unordered_map<int,int> &);
387#endif
388
389namespace
390{
391/**
392 * @brief Helper class for storing and comparing arrays of face nodes.
393 * @details The comparison operator uses the sorted nodes and a lexicographic
394 * compare so that two different orientations of the same set of nodes will be
395 * identical. The actual nodes are stored unsorted as the ordering is important
396 * for constructing the leaf-root relations.
397 */
398struct FaceNodes
399{
400 std::array<int, NCMesh::MaxFaceNodes> nodes;
401 bool operator<(FaceNodes t2) const
402 {
403 std::array<int, NCMesh::MaxFaceNodes> t1 = nodes;
404 std::sort(t1.begin(), t1.end());
405 std::sort(t2.nodes.begin(), t2.nodes.end());
406 return std::lexicographical_compare(t1.begin(), t1.end(),
407 t2.nodes.begin(), t2.nodes.end());
408 };
409};
410
411/**
412 * @brief Establish the Geometry::Type from an array of nodes
413 *
414 * @param nodes
415 * @return Geometry::Type
416 */
417Geometry::Type FaceGeomFromNodes(const std::array<int, NCMesh::MaxFaceNodes>
418 &nodes)
419{
420 if (nodes[3] == -1) { return Geometry::Type::TRIANGLE; }
421 if (nodes[0] == nodes[1] && nodes[2] == nodes[3]) { return Geometry::Type::SEGMENT; }
423};
424
425} // namespace
426
427template<typename NCSubMeshT>
428void ConstructFaceTree(NCSubMeshT &submesh, const Array<int> &attributes)
429{
430 // Convenience references to avoid `submesh.` repeatedly.
431 auto &parent_node_ids = submesh.parent_node_ids_;
432 auto &parent_element_ids = submesh.parent_element_ids_;
433 auto &parent_to_submesh_node_ids = submesh.parent_to_submesh_node_ids_;
434 auto &parent_to_submesh_element_ids = submesh.parent_to_submesh_element_ids_;
435 const auto &parent = *submesh.GetParent();
436
437 // Collect parent vertex nodes to add in sequence. Map from parent nodes to
438 // the new element in the ncsubmesh.
439 UniqueIndexGenerator node_ids;
440 std::map<FaceNodes, int> pnodes_new_elem;
441 std::set<int> new_nodes;
442 parent_to_submesh_element_ids.reserve(parent.GetNumFaces());
443 parent_element_ids.Reserve(parent.GetNumFaces());
444 // Base class cast then const cast because GetFaceList uses just in time
445 // construction.
446 const auto &face_list = const_cast<NCMesh&>(static_cast<const NCMesh&>
447 (parent)).GetFaceList();
448 // Double indexing loop because begin() and end() do not align with index 0
449 // and size-1.
450 for (int i = 0, ipe = 0; ipe < parent.GetNumFaces(); i++)
451 {
452 const auto &face = parent.GetFace(i);
453 if (face.Unused()) { continue; }
454 ipe++; // actual possible parent element.
455 if (!HasAttribute(face, attributes)
456 || face_list.GetMeshIdType(face.index) == NCMesh::NCList::MeshIdType::MASTER
457 ) { continue; }
458
459 FaceNodes fn{submesh.parent_->FindFaceNodes(face)};
460 if (pnodes_new_elem.find(fn) != pnodes_new_elem.end()) { continue; }
461
462 // TODO: Internal nc submesh can be constructed and solved on, but the
463 // transfer to the parent mesh can be erroneous, this is likely due to not
464 // treating the changing orientation of internal faces for ncmesh within
465 // the ptransfermap.
466 MFEM_ASSERT(face.elem[0] < 0 || face.elem[1] < 0,
467 "Internal nonconforming boundaries are not reliably supported yet.");
468 auto face_geom = FaceGeomFromNodes(fn.nodes);
469 int new_elem_id = submesh.AddElement(face_geom, face.attribute);
470
471 // Rank needs to be established by presence (or lack of) in the submesh.
472 submesh.elements[new_elem_id].rank = [&parent, &face]()
473 {
474 auto rank0 = face.elem[0] >= 0 ? parent.GetElement(face.elem[0]).rank : -1;
475 auto rank1 = face.elem[1] >= 0 ? parent.GetElement(face.elem[1]).rank : -1;
476 if (rank0 < 0) { return rank1; }
477 if (rank1 < 0) { return rank0; }
478 return rank0 < rank1 ? rank0 : rank1;
479 }();
480 pnodes_new_elem[fn] = new_elem_id;
481 parent_element_ids.Append(i);
482 parent_to_submesh_element_ids[i] = new_elem_id;
483
484 // Copy in the parent nodes. These will be relabeled once the tree is
485 // built.
486 std::copy(fn.nodes.begin(), fn.nodes.end(), submesh.elements[new_elem_id].node);
487 for (auto x : fn.nodes)
488 if (x != -1)
489 {
490 new_nodes.insert(x);
491 }
492 auto &gi = submesh.GI[face_geom];
493 gi.InitGeom(face_geom);
494 for (int e = 0; e < gi.ne; e++)
495 {
496 new_nodes.insert(submesh.ParentNodes().FindId(fn.nodes[gi.edges[e][0]],
497 fn.nodes[gi.edges[e][1]]));
498 }
499
500 /*
501 - Check not top level face
502 - Check for parent of the newly entered element
503 - if not present, add in
504 - if present but different order and this path is non-ambiguous,
505 reorder so consistent with child elements.
506 - Set .parent in the newly entered element
507 Break if top level face or joined existing branch (without reordering).
508
509 child element indices will be set afterwards because the orientation can change
510 during traversal.
511 */
512 bool root_path_is_ambiguous=false;
513 bool fix_parent = false, tri_face = (face_geom == Geometry::TRIANGLE);
514 while (true)
515 {
516 int child = submesh.parent_->ParentFaceNodes(fn.nodes);
517 if (tri_face && child == 3)
518 {
519 // Traversing a central triangle face involves flipping the face orientation.
520 // Do not use this pathway for reordering any parent face's nodes.
521 root_path_is_ambiguous = true;
522 }
523
524 if (child == -1) // A root face
525 {
526 submesh.elements[new_elem_id].parent = -1;
527 break;
528 }
529 auto pelem = pnodes_new_elem.find(fn);
530 bool new_parent = pelem == pnodes_new_elem.end();
531 if (new_parent)
532 {
533 // Add in this parent
534 int pelem_id = submesh.AddElement(FaceGeomFromNodes(fn.nodes), face.attribute);
535 pelem = pnodes_new_elem.emplace(fn, pelem_id).first;
536 auto parent_face_id = submesh.ParentFaces().FindId(fn.nodes[0], fn.nodes[1],
537 fn.nodes[2],
538 fn.nodes[3]);
539 parent_element_ids.Append(parent_face_id);
540 }
541 else
542 {
543 // There are two scenarios where the parent nodes should be
544 // rearranged:
545 // 1. The found face is a slave, then the master might have been
546 // added in reverse orientation
547 // 2. The parent face was added from the central face of a triangle,
548 // the orientation of the parent face is only fixed relative to
549 // the outer child faces not the interior. If either of these
550 // scenarios, and there's a mismatch, then reorder the parent and
551 // all ancestors if necessary.
552 if (!root_path_is_ambiguous &&
553 !std::equal(fn.nodes.begin(), fn.nodes.end(), pelem->first.nodes.begin()))
554 {
555 fix_parent = true;
556 auto pelem_id = pelem->second;
557 MFEM_ASSERT(!submesh.elements[pelem_id].IsLeaf(), pelem_id);
558
559 // Re-key the map, the existing entry is inconsistent with the tree.
560 pnodes_new_elem.erase(pelem->first);
561 pelem = pnodes_new_elem.emplace(fn, pelem_id).first;
562 }
563 }
564 // Ensure parent element is marked as non-leaf, and attach to the child.
565 submesh.elements[pelem->second].ref_type = submesh.Dim == 2 ? Refinement::XY :
567 submesh.elements[new_elem_id].parent = pelem->second;
568
569 // If this was neither new nor a fixed parent, the higher levels of the
570 // tree have been built, otherwise we recurse up the tree to add more parents, or
571 // to potentially fix any ambiguously added FaceNodes.
572 if (!new_parent && !fix_parent) { break; }
573
574 new_elem_id = pelem->second;
575 }
576 }
577 parent_element_ids.ShrinkToFit();
578 MFEM_ASSERT(parent_element_ids.Size() == submesh.elements.Size(),
579 parent_element_ids.Size() << ' ' << submesh.elements.Size());
580
581 // All elements have been added, with their parents, and the nodal orientation of parents is
582 // consistent with children, but the children indices have not been marked. Traverse the
583 // tree from root to leaf to fill the child arrays.
584 for (const auto & fn_elem : pnodes_new_elem)
585 {
586 auto fn = fn_elem.first;
587 const auto &child_elem = submesh.elements[fn_elem.second];
588 if (child_elem.parent == -1) { continue; }
589 int child = submesh.parent_->ParentFaceNodes(fn.nodes);
590 MFEM_ASSERT(pnodes_new_elem[fn] == child_elem.parent,
591 pnodes_new_elem[fn] << ' ' << child_elem.parent);
592 MFEM_ASSERT(submesh.elements[child_elem.parent].ref_type != char(0),
593 int(submesh.elements[child_elem.parent].ref_type));
594 submesh.elements[child_elem.parent].child[child] = fn_elem.second;
595 }
596
597 /*
598 All elements have been added into the tree but a) The nodes are all from
599 the parent ncmesh b) The nodes do not know their parents c) The element
600 ordering is wrong, root elements are not first d) The parent and child
601 element numbers reflect the incorrect ordering
602
603 1. Add in nodes in the same order from the parent ncmesh
604 2. Compute reordering of elements with parent elements first, that is
605 stable across processors.
606 */
607 // Build an inverse (and consecutive) map.
608 Array<FaceNodes> new_elem_to_parent_face_nodes(static_cast<int>
609 (pnodes_new_elem.size()));
610 for (const auto &kv : pnodes_new_elem)
611 {
612 new_elem_to_parent_face_nodes[kv.second] = kv.first;
613 }
614 pnodes_new_elem.clear(); // no longer needed
615
616 // Add new nodes preserving parent mesh ordering
617 parent_node_ids.Reserve(static_cast<int>(new_nodes.size()));
618 parent_to_submesh_node_ids.reserve(new_nodes.size());
619 for (auto n : new_nodes)
620 {
621 bool new_node;
622 auto new_node_id = node_ids.Get(n, new_node);
623 MFEM_ASSERT(new_node, "!");
624 submesh.nodes.Alloc(new_node_id, new_node_id, new_node_id);
625 parent_node_ids.Append(n);
626 parent_to_submesh_node_ids[n] = new_node_id;
627 }
628 parent_node_ids.ShrinkToFit();
629 new_nodes.clear(); // not needed any more.
630
631 // Comparator for deciding order of elements. Building the ordering from the
632 // parent ncmesh ensures the root ordering is common across ranks.
633 auto comp_elements = [&](int l, int r)
634 {
635 const auto &elem_l = submesh.elements[l];
636 const auto &elem_r = submesh.elements[r];
637 if (elem_l.parent == elem_r.parent)
638 {
639 const auto &fnl = new_elem_to_parent_face_nodes[l].nodes;
640 const auto &fnr = new_elem_to_parent_face_nodes[r].nodes;
641 return std::lexicographical_compare(fnl.begin(), fnl.end(), fnr.begin(),
642 fnr.end());
643 }
644 else
645 {
646 return elem_l.parent < elem_r.parent;
647 }
648 };
649 Array<int> indices(submesh.elements.Size());
650 auto parental_sorted = [&]()
651 {
652 std::iota(indices.begin(), indices.end(), 0);
653 return std::is_sorted(indices.begin(), indices.end(), comp_elements);
654 };
655
656 Array<int> new_to_old(submesh.elements.Size()),
657 old_to_new(submesh.elements.Size());
658 while (!parental_sorted())
659 {
660 // Stably reorder elements in order of refinement, and by parental nodes
661 // within a nuclear family.
662 new_to_old.SetSize(submesh.elements.Size()),
663 old_to_new.SetSize(submesh.elements.Size());
664 std::iota(new_to_old.begin(), new_to_old.end(), 0);
665 std::stable_sort(new_to_old.begin(), new_to_old.end(), comp_elements);
666 // Build the inverse relation for converting the old elements to new
667 for (int i = 0; i < submesh.elements.Size(); i++)
668 {
669 old_to_new[new_to_old[i]] = i;
670 }
671
672 // Permute whilst reordering new_to_old. Avoids unnecessary copies.
673 Permute(std::move(new_to_old), submesh.elements, parent_element_ids,
674 new_elem_to_parent_face_nodes);
675 parent_to_submesh_element_ids.clear();
676 for (int i = 0; i < parent_element_ids.Size(); i++)
677 {
678 if (parent_element_ids[i] == -1) {continue;}
679 parent_to_submesh_element_ids[parent_element_ids[i]] = i;
680 }
681
682 // Apply the new ordering to child and parent elements
683 for (auto &elem : submesh.elements)
684 {
685 if (!elem.IsLeaf())
686 {
687 // Parent rank is minimum of child ranks.
688 elem.rank = std::numeric_limits<int>::max();
689 for (int c = 0; c < NCMesh::MaxElemChildren && elem.child[c] >= 0; c++)
690 {
691 elem.child[c] = old_to_new[elem.child[c]];
692 elem.rank = std::min(elem.rank, submesh.elements[elem.child[c]].rank);
693 }
694 }
695 elem.parent = elem.parent == -1 ? -1 : old_to_new[elem.parent];
696 }
697 }
698
699 // Apply new node ordering to relations, and sign in on edges/vertices
700 for (auto &elem : submesh.elements)
701 {
702 if (elem.IsLeaf())
703 {
704 bool new_id;
705 auto &gi = submesh.GI[elem.Geom()];
706 gi.InitGeom(elem.Geom());
707 for (int e = 0; e < gi.ne; e++)
708 {
709 const int pid = submesh.ParentNodes().FindId(
710 elem.node[gi.edges[e][0]], elem.node[gi.edges[e][1]]);
711 MFEM_ASSERT(pid >= 0,
712 elem.node[gi.edges[e][0]] << ' ' << elem.node[gi.edges[e][1]]);
713 auto submesh_node_id = node_ids.Get(pid, new_id);
714 MFEM_ASSERT(!new_id, "!");
715 submesh.nodes[submesh_node_id].edge_refc++;
716 }
717 for (int n = 0; n < gi.nv; n++)
718 {
719 MFEM_ASSERT(parent_to_submesh_node_ids.find(elem.node[n]) !=
720 parent_to_submesh_node_ids.end(), "!");
721 elem.node[n] = parent_to_submesh_node_ids[elem.node[n]];
722 submesh.nodes[elem.node[n]].vert_refc++;
723 }
724 // Register faces
725 for (int f = 0; f < gi.nf; f++)
726 {
727 auto *face = submesh.faces.Get(
728 elem.node[gi.faces[f][0]],
729 elem.node[gi.faces[f][1]],
730 elem.node[gi.faces[f][2]],
731 elem.node[gi.faces[f][3]]);
732 face->attribute = -1;
733 face->index = -1;
734 }
735 }
736 }
737}
738
739// Explicit instantiations
740template void ConstructFaceTree(NCSubMesh &submesh,
741 const Array<int> &attributes);
742#ifdef MFEM_USE_MPI
743template void ConstructFaceTree(ParNCSubMesh &submesh,
744 const Array<int> &attributes);
745#endif
746
747template <typename NCSubMeshT>
748void ConstructVolumeTree(NCSubMeshT &submesh, const Array<int> &attributes)
749{
750 // Convenience references to avoid `submesh.` repeatedly.
751 auto &parent_node_ids = submesh.parent_node_ids_;
752 auto &parent_element_ids = submesh.parent_element_ids_;
753 auto &parent_to_submesh_node_ids = submesh.parent_to_submesh_node_ids_;
754 auto &parent_to_submesh_element_ids = submesh.parent_to_submesh_element_ids_;
755 const auto &parent = *submesh.GetParent();
756
757 UniqueIndexGenerator node_ids;
758 parent_to_submesh_element_ids.reserve(parent.GetNumElements());
759 std::set<int> new_nodes;
760 for (int ipe = 0; ipe < parent.GetNumElements(); ipe++)
761 {
762 const auto& pe = parent.GetElement(ipe);
763 if (!HasAttribute(pe, attributes)) { continue; }
764 const int elem_id = submesh.AddElement(pe);
765 auto &el = submesh.elements[elem_id];
766 parent_element_ids.Append(ipe); // submesh -> parent
767 parent_to_submesh_element_ids[ipe] = elem_id; // parent -> submesh
768 if (!pe.IsLeaf()) { continue; }
769 const auto gi = submesh.GI[pe.Geom()];
770 for (int n = 0; n < gi.nv; n++)
771 {
772 new_nodes.insert(el.node[n]);
773 }
774 for (int e = 0; e < gi.ne; e++)
775 {
776 new_nodes.insert(submesh.ParentNodes().FindId(el.node[gi.edges[e][0]],
777 el.node[gi.edges[e][1]]));
778 }
779 }
780
781 parent_node_ids.Reserve(static_cast<int>(new_nodes.size()));
782 parent_to_submesh_node_ids.reserve(new_nodes.size());
783 for (const auto &n : new_nodes)
784 {
785 bool new_node;
786 auto new_node_id = node_ids.Get(n, new_node);
787 MFEM_ASSERT(new_node, "!");
788 submesh.nodes.Alloc(new_node_id, new_node_id, new_node_id);
789 parent_node_ids.Append(n);
790 parent_to_submesh_node_ids[n] = new_node_id;
791 }
792
793 // Loop over elements and reference edges and faces (creating any nodes on
794 // first encounter).
795 for (auto &el : submesh.elements)
796 {
797 if (el.IsLeaf())
798 {
799 const auto gi = submesh.GI[el.Geom()];
800 bool new_id = false;
801
802 for (int n = 0; n < gi.nv; n++)
803 {
804 // Relabel nodes from parent to submesh.
805 el.node[n] = node_ids.Get(el.node[n], new_id);
806 MFEM_ASSERT(new_id == false, "Should not be new.");
807 submesh.nodes[el.node[n]].vert_refc++;
808 }
809 for (int e = 0; e < gi.ne; e++)
810 {
811 const int pid = submesh.ParentNodes().FindId(
812 parent_node_ids[el.node[gi.edges[e][0]]],
813 parent_node_ids[el.node[gi.edges[e][1]]]);
814 MFEM_ASSERT(pid >= 0, "Edge not found");
815 auto submesh_node_id = node_ids.Get(pid, new_id);
816 MFEM_ASSERT(new_id == false, "Should not be new.");
817 submesh.nodes[submesh_node_id].edge_refc++; // Register the edge
818 }
819 for (int f = 0; f < gi.nf; f++)
820 {
821 const int *fv = gi.faces[f];
822 const int pid = submesh.ParentFaces().FindId(
823 parent_node_ids[el.node[fv[0]]],
824 parent_node_ids[el.node[fv[1]]],
825 parent_node_ids[el.node[fv[2]]],
826 el.node[fv[3]] >= 0 ? parent_node_ids[el.node[fv[3]]]: - 1);
827 MFEM_ASSERT(pid >= 0, "Face not found");
828 const int id = submesh.faces.GetId(
829 el.node[fv[0]], el.node[fv[1]], el.node[fv[2]], el.node[fv[3]]);
830 submesh.faces[id].attribute = submesh.ParentFaces()[pid].attribute;
831 }
832 }
833 else
834 {
835 // All elements have been collected, remap the child ids.
836 for (int i = 0; i < NCMesh::MaxElemChildren && el.child[i] >= 0; i++)
837 {
838 el.child[i] = parent_to_submesh_element_ids[el.child[i]];
839 }
840 }
841 el.parent = el.parent < 0 ? el.parent
842 : parent_to_submesh_element_ids.at(el.parent);
843 }
844}
845
846// Explicit instantiations
847template void ConstructVolumeTree(NCSubMesh &submesh,
848 const Array<int> &attributes);
849#ifdef MFEM_USE_MPI
850template void ConstructVolumeTree(ParNCSubMesh &submesh,
851 const Array<int> &attributes);
852#endif
853} // namespace SubMeshUtils
854} // namespace mfem
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
Definition array.hpp:312
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
Definition array.hpp:184
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:840
int Size() const
Return the logical size of the array.
Definition array.hpp:166
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:912
T * end()
STL-like end. Returns pointer after the last element of the array.
Definition array.hpp:369
T * begin()
STL-like begin. Returns pointer to the first element of the array.
Definition array.hpp:366
@ GaussLobatto
Closed type.
Definition fe_base.hpp:36
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Abstract data type element.
Definition element.hpp:29
virtual void GetVertices(Array< int > &v) const =0
Get the indices defining the vertices.
void SetAttribute(const int attr)
Set element's attribute.
Definition element.hpp:61
int GetAttribute() const
Return element's attribute.
Definition element.hpp:58
virtual void SetVertices(const Array< int > &v)=0
Set the indices defining the vertices.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:208
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:304
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3824
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element 'i' to the given 'geom_type'.
Definition fespace.cpp:3942
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:856
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:644
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:826
int GetVDim() const
Returns the vector dimension of the finite element space.
Definition fespace.hpp:819
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition fespace.hpp:1503
static int DecodeDof(int dof)
Helper to return the DOF associated with a sign encoded DOF.
Definition fespace.hpp:1155
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i'th boundary element. The returned indices are offsets int...
Definition fespace.cpp:319
Abstract class for all finite elements.
Definition fe_base.hpp:247
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
Definition fe_base.cpp:129
static int GetInverseOrientation(Type geom_type, int orientation)
Return the inverse of the given orientation for the specified geometry type.
Definition geom.cpp:264
IsoparametricTransformation Transf
Definition eltrans.hpp:733
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:350
int GetBasisType() const
Definition fe_coll.hpp:394
Mesh data type.
Definition mesh.hpp:65
Element * NewElement(int geom)
Definition mesh.cpp:4818
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition mesh.cpp:7962
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition mesh.cpp:7967
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition mesh.cpp:6846
Geometry::Type GetBdrElementGeometry(int i) const
Definition mesh.hpp:1547
const Element * GetElement(int i) const
Return pointer to the i'th element object.
Definition mesh.hpp:1434
static int EncodeFaceInfo(int local_face_index, int orientation)
Given local_face_index and orientation, return the corresponding encoded "face info int".
Definition mesh.hpp:2177
int AddVertex(real_t x, real_t y=0.0, real_t z=0.0)
Definition mesh.cpp:2000
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1377
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1306
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
Definition mesh.hpp:1449
int AddElement(Element *elem)
Definition mesh.cpp:2363
static int DecodeFaceInfoLocalIndex(int info)
Given a "face info int", return the local face index.
Definition mesh.hpp:2173
void GetElementFaces(int i, Array< int > &faces, Array< int > &ori) const
Return the indices and the orientations of all faces of element i.
Definition mesh.cpp:7835
static int DecodeFaceInfoOrientation(int info)
Given a "face info int", return the face orientation.
Definition mesh.hpp:2170
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1380
void GetElementEdges(int i, Array< int > &edges, Array< int > &cor) const
Return the indices and the orientations of all edges of element i.
Definition mesh.cpp:7588
const real_t * GetVertex(int i) const
Return pointer to vertex i's coordinates.
Definition mesh.hpp:1416
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
Definition ncmesh.hpp:190
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
Definition ncmesh.hpp:369
static constexpr int MaxElemChildren
Number of children an element can have.
Definition ncmesh.hpp:558
Class representing a Nonconformal SubMesh. This is only used by SubMesh.
Definition ncsubmesh.hpp:28
Class for standard nodal finite elements.
Definition fe_base.hpp:724
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
Class representing a Parallel Nonconformal SubMesh. This is only used by ParSubMesh.
Subdomain representation of a topological parent in another ParMesh.
Definition psubmesh.hpp:54
Subdomain representation of a topological parent in another Mesh.
Definition submesh.hpp:44
From
Indicator from which part of the parent Mesh the SubMesh is created.
Definition submesh.hpp:49
real_t a
Definition lissajous.cpp:41
std::tuple< Array< int >, Array< int > > AddElementsToMesh(const Mesh &parent, Mesh &mesh, const Array< int > &attributes, bool from_boundary)
Given a Mesh parent and another Mesh mesh using the list of attributes in attributes,...
void ConstructFaceTree(NCSubMeshT &submesh, const Array< int > &attributes)
Construct a nonconformal mesh (serial or parallel) for a surface submesh, from an existing nonconform...
void BuildVdofToVdofMap(const FiniteElementSpace &subfes, const FiniteElementSpace &parentfes, const SubMesh::From &from, const Array< int > &parent_element_ids, Array< int > &vdof_to_vdof_map)
Build the vdof to vdof mapping between two FiniteElementSpace objects.
bool HasAttribute(const T &el, const Array< int > &attributes)
Helper for checking if an object's attributes match a list.
bool ElementHasAttribute(const ElementT &el, const Array< int > &attributes)
void Permute(const Array< int > &indices, T1 &t1, T2 &t2, T3 &t3)
Apply permutation to a container type.
void ConstructVolumeTree(NCSubMeshT &submesh, const Array< int > &attributes)
Construct a nonconformal mesh (serial or parallel) for a volume submesh, from an existing nonconforma...
Array< int > BuildFaceMap(const Mesh &pm, const Mesh &sm, const Array< int > &parent_element_ids)
Given two meshes that have a parent to SubMesh relationship create a face map, using a SubMesh to par...
void AddBoundaryElements(SubMeshT &mesh, const std::unordered_map< int, int > &lface_to_boundary_attribute)
Add boundary elements to the SubMesh.
float real_t
Definition config.hpp:46
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
bool operator<(const Pair< A, B > &p, const Pair< A, B > &q)
Comparison operator for class Pair, based on the first element only.
Convenience object to create unique indices.
std::unordered_map< int, int > idx
int Get(int i, bool &new_index)
Returns the unique index from an index set.
std::array< int, NCMesh::MaxFaceNodes > nodes