void MeshTools::Subdivision::add_boundary_ghosts(MeshBase & mesh) { static const Real tol = 1e-5; // add the mirrored ghost elements (without using iterators, because the mesh is modified in the course) std::vector<Tri3Subdivision *> ghost_elems; std::vector<Node *> ghost_nodes; const unsigned int n_elem = mesh.n_elem(); for (unsigned int eid = 0; eid < n_elem; ++eid) { Elem * elem = mesh.elem(eid); libmesh_assert_equal_to(elem->type(), TRI3SUBDIVISION); // If the triangle happens to be in a corner (two boundary // edges), we perform a counter-clockwise loop by mirroring the // previous triangle until we come back to the original // triangle. This prevents degenerated triangles in the mesh // corners and guarantees that the node in the middle of the // loop is of valence=6. for (unsigned int i = 0; i < elem->n_sides(); ++i) { libmesh_assert_not_equal_to(elem->neighbor(i), elem); if (elem->neighbor(i) == libmesh_nullptr && elem->neighbor(next[i]) == libmesh_nullptr) { Elem * nelem = elem; unsigned int k = i; for (unsigned int l=0;l<4;l++) { // this is the vertex to be mirrored Point point = nelem->point(k) + nelem->point(next[k]) - nelem->point(prev[k]); // Check if the proposed vertex doesn't coincide // with one of the existing vertices. This is // necessary because for some triangulations, it can // happen that two mirrored ghost vertices coincide, // which would then lead to a zero size ghost // element below. Node * node = libmesh_nullptr; for (unsigned int j = 0; j < ghost_nodes.size(); ++j) { if ((*ghost_nodes[j] - point).size() < tol * (elem->point(k) - point).size()) { node = ghost_nodes[j]; break; } } // add the new vertex only if no other is nearby if (node == libmesh_nullptr) { node = mesh.add_point(point); ghost_nodes.push_back(node); } Tri3Subdivision * newelem = new Tri3Subdivision(); // add the first new ghost element to the list just as in the non-corner case if (l == 0) ghost_elems.push_back(newelem); newelem->set_node(0) = nelem->get_node(next[k]); newelem->set_node(1) = nelem->get_node(k); newelem->set_node(2) = node; newelem->set_neighbor(0, nelem); newelem->set_ghost(true); if (l>0) newelem->set_neighbor(2, libmesh_nullptr); nelem->set_neighbor(k, newelem); mesh.add_elem(newelem); mesh.get_boundary_info().add_node(nelem->get_node(k), 1); mesh.get_boundary_info().add_node(nelem->get_node(next[k]), 1); mesh.get_boundary_info().add_node(nelem->get_node(prev[k]), 1); mesh.get_boundary_info().add_node(node, 1); nelem = newelem; k = 2 ; } Tri3Subdivision * newelem = new Tri3Subdivision(); newelem->set_node(0) = elem->get_node(next[i]); newelem->set_node(1) = nelem->get_node(2); newelem->set_node(2) = elem->get_node(prev[i]); newelem->set_neighbor(0, nelem); nelem->set_neighbor(2, newelem); newelem->set_ghost(true); newelem->set_neighbor(2, elem); elem->set_neighbor(next[i],newelem); mesh.add_elem(newelem); break; } } for (unsigned int i = 0; i < elem->n_sides(); ++i) { libmesh_assert_not_equal_to(elem->neighbor(i), elem); if (elem->neighbor(i) == libmesh_nullptr) { // this is the vertex to be mirrored Point point = elem->point(i) + elem->point(next[i]) - elem->point(prev[i]); // Check if the proposed vertex doesn't coincide with // one of the existing vertices. This is necessary // because for some triangulations, it can happen that // two mirrored ghost vertices coincide, which would // then lead to a zero size ghost element below. Node * node = libmesh_nullptr; for (unsigned int j = 0; j < ghost_nodes.size(); ++j) { if ((*ghost_nodes[j] - point).size() < tol * (elem->point(i) - point).size()) { node = ghost_nodes[j]; break; } } // add the new vertex only if no other is nearby if (node == libmesh_nullptr) { node = mesh.add_point(point); ghost_nodes.push_back(node); } Tri3Subdivision * newelem = new Tri3Subdivision(); ghost_elems.push_back(newelem); newelem->set_node(0) = elem->get_node(next[i]); newelem->set_node(1) = elem->get_node(i); newelem->set_node(2) = node; newelem->set_neighbor(0, elem); newelem->set_ghost(true); elem->set_neighbor(i, newelem); mesh.add_elem(newelem); mesh.get_boundary_info().add_node(elem->get_node(i), 1); mesh.get_boundary_info().add_node(elem->get_node(next[i]), 1); mesh.get_boundary_info().add_node(elem->get_node(prev[i]), 1); mesh.get_boundary_info().add_node(node, 1); } } } // add the missing ghost elements (connecting new ghost nodes) std::vector<Tri3Subdivision *> missing_ghost_elems; std::vector<Tri3Subdivision *>::iterator ghost_el = ghost_elems.begin(); const std::vector<Tri3Subdivision *>::iterator end_ghost_el = ghost_elems.end(); for (; ghost_el != end_ghost_el; ++ghost_el) { Tri3Subdivision * elem = *ghost_el; libmesh_assert(elem->is_ghost()); for (unsigned int i = 0; i < elem->n_sides(); ++i) { if (elem->neighbor(i) == libmesh_nullptr && elem->neighbor(prev[i]) != libmesh_nullptr) { // go around counter-clockwise Tri3Subdivision * nb1 = static_cast<Tri3Subdivision *>(elem->neighbor(prev[i])); Tri3Subdivision * nb2 = nb1; unsigned int j = i; unsigned int n_nb = 0; while (nb1 != libmesh_nullptr && nb1->id() != elem->id()) { j = nb1->local_node_number(elem->node(i)); nb2 = nb1; nb1 = static_cast<Tri3Subdivision *>(nb1->neighbor(prev[j])); libmesh_assert(nb1 == libmesh_nullptr || nb1->id() != nb2->id()); n_nb++; } libmesh_assert_not_equal_to(nb2->id(), elem->id()); // Above, we merged coinciding ghost vertices. Therefore, we need // to exclude the case where there is no ghost element to add between // these two (identical) ghost nodes. if (elem->get_node(next[i])->id() == nb2->get_node(prev[j])->id()) break; // If the number of already present neighbors is less than 4, we add another extra element // so that the node in the middle of the loop ends up being of valence=6. // This case usually happens when the middle node corresponds to a corner of the original mesh, // and the extra element below prevents degenerated triangles in the mesh corners. if (n_nb < 4) { // this is the vertex to be mirrored Point point = nb2->point(j) + nb2->point(prev[j]) - nb2->point(next[j]); // Check if the proposed vertex doesn't coincide with one of the existing vertices. // This is necessary because for some triangulations, it can happen that two mirrored // ghost vertices coincide, which would then lead to a zero size ghost element below. Node * node = libmesh_nullptr; for (unsigned int k = 0; k < ghost_nodes.size(); ++k) { if ((*ghost_nodes[k] - point).size() < tol * (nb2->point(j) - point).size()) { node = ghost_nodes[k]; break; } } // add the new vertex only if no other is nearby if (node == libmesh_nullptr) { node = mesh.add_point(point); ghost_nodes.push_back(node); } Tri3Subdivision * newelem = new Tri3Subdivision(); newelem->set_node(0) = nb2->get_node(j); newelem->set_node(1) = nb2->get_node(prev[j]); newelem->set_node(2) = node; newelem->set_neighbor(0, nb2); newelem->set_neighbor(1, libmesh_nullptr); newelem->set_ghost(true); nb2->set_neighbor(prev[j], newelem); mesh.add_elem(newelem); mesh.get_boundary_info().add_node(nb2->get_node(j), 1); mesh.get_boundary_info().add_node(nb2->get_node(prev[j]), 1); mesh.get_boundary_info().add_node(node, 1); nb2 = newelem; j = nb2->local_node_number(elem->node(i)); } Tri3Subdivision * newelem = new Tri3Subdivision(); newelem->set_node(0) = elem->get_node(next[i]); newelem->set_node(1) = elem->get_node(i); newelem->set_node(2) = nb2->get_node(prev[j]); newelem->set_neighbor(0, elem); newelem->set_neighbor(1, nb2); newelem->set_neighbor(2, libmesh_nullptr); newelem->set_ghost(true); elem->set_neighbor(i, newelem); nb2->set_neighbor(prev[j], newelem); missing_ghost_elems.push_back(newelem); break; } } // end side loop } // end ghost element loop // add the missing ghost elements to the mesh std::vector<Tri3Subdivision *>::iterator missing_el = missing_ghost_elems.begin(); const std::vector<Tri3Subdivision *>::iterator end_missing_el = missing_ghost_elems.end(); for (; missing_el != end_missing_el; ++missing_el) mesh.add_elem(*missing_el); }
void UnstructuredMesh::copy_nodes_and_elements(const UnstructuredMesh & other_mesh, const bool skip_find_neighbors) { // We're assuming our subclass data needs no copy libmesh_assert_equal_to (_n_parts, other_mesh._n_parts); libmesh_assert (std::equal(_elem_dims.begin(), _elem_dims.end(), other_mesh._elem_dims.begin())); libmesh_assert_equal_to (_is_prepared, other_mesh._is_prepared); // We're assuming the other mesh has proper element number ordering, // so that we add parents before their children. #ifdef DEBUG MeshTools::libmesh_assert_valid_amr_elem_ids(other_mesh); #endif //Copy in Nodes { //Preallocate Memory if necessary this->reserve_nodes(other_mesh.n_nodes()); const_node_iterator it = other_mesh.nodes_begin(); const_node_iterator end = other_mesh.nodes_end(); for (; it != end; ++it) { const Node * oldn = *it; // Add new nodes in old node Point locations #ifdef LIBMESH_ENABLE_UNIQUE_ID Node *newn = #endif this->add_point(*oldn, oldn->id(), oldn->processor_id()); #ifdef LIBMESH_ENABLE_UNIQUE_ID newn->set_unique_id() = oldn->unique_id(); #endif } } //Copy in Elements { //Preallocate Memory if necessary this->reserve_elem(other_mesh.n_elem()); // Declare a map linking old and new elements, needed to copy the neighbor lists std::map<const Elem *, Elem *> old_elems_to_new_elems; // Loop over the elements MeshBase::const_element_iterator it = other_mesh.elements_begin(); const MeshBase::const_element_iterator end = other_mesh.elements_end(); // FIXME: Where do we set element IDs?? for (; it != end; ++it) { //Look at the old element const Elem * old = *it; //Build a new element Elem * newparent = old->parent() ? this->elem_ptr(old->parent()->id()) : libmesh_nullptr; UniquePtr<Elem> ap = Elem::build(old->type(), newparent); Elem * el = ap.release(); el->subdomain_id() = old->subdomain_id(); for (unsigned int s=0; s != old->n_sides(); ++s) if (old->neighbor_ptr(s) == remote_elem) el->set_neighbor(s, const_cast<RemoteElem *>(remote_elem)); #ifdef LIBMESH_ENABLE_AMR if (old->has_children()) for (unsigned int c=0; c != old->n_children(); ++c) if (old->child_ptr(c) == remote_elem) el->add_child(const_cast<RemoteElem *>(remote_elem), c); //Create the parent's child pointers if necessary if (newparent) { unsigned int oldc = old->parent()->which_child_am_i(old); newparent->add_child(el, oldc); } // Copy the refinement flags el->set_refinement_flag(old->refinement_flag()); // Use hack_p_level since we may not have sibling elements // added yet el->hack_p_level(old->p_level()); el->set_p_refinement_flag(old->p_refinement_flag()); #endif // #ifdef LIBMESH_ENABLE_AMR //Assign all the nodes for(unsigned int i=0;i<el->n_nodes();i++) el->set_node(i) = this->node_ptr(old->node_id(i)); // And start it off in the same subdomain el->processor_id() = old->processor_id(); // Give it the same ids el->set_id(old->id()); #ifdef LIBMESH_ENABLE_UNIQUE_ID el->set_unique_id() = old->unique_id(); #endif //Hold onto it if(!skip_find_neighbors) { this->add_elem(el); } else { Elem * new_el = this->add_elem(el); old_elems_to_new_elems[old] = new_el; } // Add the link between the original element and this copy to the map if(skip_find_neighbors) old_elems_to_new_elems[old] = el; } // Loop (again) over the elements to fill in the neighbors if(skip_find_neighbors) { it = other_mesh.elements_begin(); for (; it != end; ++it) { Elem * old_elem = *it; Elem * new_elem = old_elems_to_new_elems[old_elem]; for (unsigned int s=0; s != old_elem->n_neighbors(); ++s) { const Elem * old_neighbor = old_elem->neighbor_ptr(s); Elem * new_neighbor = old_elems_to_new_elems[old_neighbor]; new_elem->set_neighbor(s, new_neighbor); } } } } //Finally prepare the new Mesh for use. Keep the same numbering and //partitioning but also the same renumbering and partitioning //policies as our source mesh. this->allow_renumbering(false); this->skip_partitioning(true); this->prepare_for_use(false, skip_find_neighbors); this->allow_renumbering(other_mesh.allow_renumbering()); this->skip_partitioning(other_mesh.skip_partitioning()); }
void unpack(std::vector<largest_id_type>::const_iterator in, Elem** out, MeshBase* mesh) { #ifndef NDEBUG const std::vector<largest_id_type>::const_iterator original_in = in; const largest_id_type incoming_header = *in++; libmesh_assert_equal_to (incoming_header, elem_magic_header); #endif // int 0: level const unsigned int level = static_cast<unsigned int>(*in++); #ifdef LIBMESH_ENABLE_AMR // int 1: p level const unsigned int p_level = static_cast<unsigned int>(*in++); // int 2: refinement flag const int rflag = *in++; libmesh_assert_greater_equal (rflag, 0); libmesh_assert_less (rflag, Elem::INVALID_REFINEMENTSTATE); const Elem::RefinementState refinement_flag = static_cast<Elem::RefinementState>(rflag); // int 3: p refinement flag const int pflag = *in++; libmesh_assert_greater_equal (pflag, 0); libmesh_assert_less (pflag, Elem::INVALID_REFINEMENTSTATE); const Elem::RefinementState p_refinement_flag = static_cast<Elem::RefinementState>(pflag); #else in += 3; #endif // LIBMESH_ENABLE_AMR // int 4: element type const int typeint = *in++; libmesh_assert_greater_equal (typeint, 0); libmesh_assert_less (typeint, INVALID_ELEM); const ElemType type = static_cast<ElemType>(typeint); const unsigned int n_nodes = Elem::type_to_n_nodes_map[type]; // int 5: processor id const processor_id_type processor_id = static_cast<processor_id_type>(*in++); libmesh_assert (processor_id < mesh->n_processors() || processor_id == DofObject::invalid_processor_id); // int 6: subdomain id const subdomain_id_type subdomain_id = static_cast<subdomain_id_type>(*in++); // int 7: dof object id const dof_id_type id = static_cast<dof_id_type>(*in++); libmesh_assert_not_equal_to (id, DofObject::invalid_id); #ifdef LIBMESH_ENABLE_UNIQUE_ID // int 8: dof object unique id const unique_id_type unique_id = static_cast<unique_id_type>(*in++); #endif #ifdef LIBMESH_ENABLE_AMR // int 9: parent dof object id const dof_id_type parent_id = static_cast<dof_id_type>(*in++); libmesh_assert (level == 0 || parent_id != DofObject::invalid_id); libmesh_assert (level != 0 || parent_id == DofObject::invalid_id); // int 10: local child id const unsigned int which_child_am_i = static_cast<unsigned int>(*in++); #else in += 2; #endif // LIBMESH_ENABLE_AMR // Make sure we don't miscount above when adding the "magic" header // plus the real data header libmesh_assert_equal_to (in - original_in, header_size + 1); Elem *elem = mesh->query_elem(id); // if we already have this element, make sure its // properties match, and update any missing neighbor // links, but then go on if (elem) { libmesh_assert_equal_to (elem->level(), level); libmesh_assert_equal_to (elem->id(), id); //#ifdef LIBMESH_ENABLE_UNIQUE_ID // No check for unqiue id sanity //#endif libmesh_assert_equal_to (elem->processor_id(), processor_id); libmesh_assert_equal_to (elem->subdomain_id(), subdomain_id); libmesh_assert_equal_to (elem->type(), type); libmesh_assert_equal_to (elem->n_nodes(), n_nodes); #ifndef NDEBUG // All our nodes should be correct for (unsigned int i=0; i != n_nodes; ++i) libmesh_assert(elem->node(i) == static_cast<dof_id_type>(*in++)); #else in += n_nodes; #endif #ifdef LIBMESH_ENABLE_AMR libmesh_assert_equal_to (elem->p_level(), p_level); libmesh_assert_equal_to (elem->refinement_flag(), refinement_flag); libmesh_assert_equal_to (elem->p_refinement_flag(), p_refinement_flag); libmesh_assert (!level || elem->parent() != NULL); libmesh_assert (!level || elem->parent()->id() == parent_id); libmesh_assert (!level || elem->parent()->child(which_child_am_i) == elem); #endif // Our neighbor links should be "close to" correct - we may have // to update them, but we can check for some inconsistencies. for (unsigned int n=0; n != elem->n_neighbors(); ++n) { const dof_id_type neighbor_id = static_cast<dof_id_type>(*in++); // If the sending processor sees a domain boundary here, // we'd better agree. if (neighbor_id == DofObject::invalid_id) { libmesh_assert (!(elem->neighbor(n))); continue; } // If the sending processor has a remote_elem neighbor here, // then all we know is that we'd better *not* have a domain // boundary. if (neighbor_id == remote_elem->id()) { libmesh_assert(elem->neighbor(n)); continue; } Elem *neigh = mesh->query_elem(neighbor_id); // The sending processor sees a neighbor here, so if we // don't have that neighboring element, then we'd better // have a remote_elem signifying that fact. if (!neigh) { libmesh_assert_equal_to (elem->neighbor(n), remote_elem); continue; } // The sending processor has a neighbor here, and we have // that element, but that does *NOT* mean we're already // linking to it. Perhaps we initially received both elem // and neigh from processors on which their mutual link was // remote? libmesh_assert(elem->neighbor(n) == neigh || elem->neighbor(n) == remote_elem); // If the link was originally remote, we should update it, // and make sure the appropriate parts of its family link // back to us. if (elem->neighbor(n) == remote_elem) { elem->set_neighbor(n, neigh); elem->make_links_to_me_local(n); } } // FIXME: We should add some debug mode tests to ensure that the // encoded indexing and boundary conditions are consistent. } else { // We don't already have the element, so we need to create it. // Find the parent if necessary Elem *parent = NULL; #ifdef LIBMESH_ENABLE_AMR // Find a child element's parent if (level > 0) { // Note that we must be very careful to construct the send // connectivity so that parents are encountered before // children. If we get here and can't find the parent that // is a fatal error. parent = mesh->elem(parent_id); } // Or assert that the sending processor sees no parent else libmesh_assert_equal_to (parent_id, static_cast<dof_id_type>(-1)); #else // No non-level-0 elements without AMR libmesh_assert_equal_to (level, 0); #endif elem = Elem::build(type,parent).release(); libmesh_assert (elem); #ifdef LIBMESH_ENABLE_AMR if (level != 0) { // Since this is a newly created element, the parent must // have previously thought of this child as a remote element. libmesh_assert_equal_to (parent->child(which_child_am_i), remote_elem); parent->add_child(elem, which_child_am_i); } // Assign the refinement flags and levels elem->set_p_level(p_level); elem->set_refinement_flag(refinement_flag); elem->set_p_refinement_flag(p_refinement_flag); libmesh_assert_equal_to (elem->level(), level); // If this element definitely should have children, assign // remote_elem to all of them for now, for consistency. Later // unpacked elements may overwrite that. if (!elem->active()) for (unsigned int c=0; c != elem->n_children(); ++c) elem->add_child(const_cast<RemoteElem*>(remote_elem), c); #endif // LIBMESH_ENABLE_AMR // Assign the IDs elem->subdomain_id() = subdomain_id; elem->processor_id() = processor_id; elem->set_id() = id; #ifdef LIBMESH_ENABLE_UNIQUE_ID elem->set_unique_id() = unique_id; #endif // Assign the connectivity libmesh_assert_equal_to (elem->n_nodes(), n_nodes); for (unsigned int n=0; n != n_nodes; n++) elem->set_node(n) = mesh->node_ptr (static_cast<dof_id_type>(*in++)); for (unsigned int n=0; n<elem->n_neighbors(); n++) { const dof_id_type neighbor_id = static_cast<dof_id_type>(*in++); if (neighbor_id == DofObject::invalid_id) continue; // We may be unpacking an element that was a ghost element on the // sender, in which case the element's neighbors may not all be // known by the packed element. We'll have to set such // neighbors to remote_elem ourselves and wait for a later // packed element to give us better information. if (neighbor_id == remote_elem->id()) { elem->set_neighbor(n, const_cast<RemoteElem*>(remote_elem)); continue; } // If we don't have the neighbor element, then it's a // remote_elem until we get it. Elem *neigh = mesh->query_elem(neighbor_id); if (!neigh) { elem->set_neighbor(n, const_cast<RemoteElem*>(remote_elem)); continue; } // If we have the neighbor element, then link to it, and // make sure the appropriate parts of its family link back // to us. elem->set_neighbor(n, neigh); elem->make_links_to_me_local(n); } elem->unpack_indexing(in); } in += elem->packed_indexing_size(); // If this is a coarse element, // add any element side boundary condition ids if (level == 0) for (unsigned int s = 0; s != elem->n_sides(); ++s) { const int num_bcs = *in++; libmesh_assert_greater_equal (num_bcs, 0); for(int bc_it=0; bc_it < num_bcs; bc_it++) mesh->boundary_info->add_side (elem, s, *in++); } // Return the new element *out = elem; }
void UnstructuredMesh::find_neighbors (const bool reset_remote_elements, const bool reset_current_list) { // We might actually want to run this on an empty mesh // (e.g. the boundary mesh for a nonexistant bcid!) // libmesh_assert_not_equal_to (this->n_nodes(), 0); // libmesh_assert_not_equal_to (this->n_elem(), 0); // This function must be run on all processors at once parallel_object_only(); LOG_SCOPE("find_neighbors()", "Mesh"); const element_iterator el_end = this->elements_end(); //TODO:[BSK] This should be removed later?! if (reset_current_list) for (element_iterator el = this->elements_begin(); el != el_end; ++el) { Elem * e = *el; for (unsigned int s=0; s<e->n_neighbors(); s++) if (e->neighbor_ptr(s) != remote_elem || reset_remote_elements) e->set_neighbor(s, libmesh_nullptr); } // Find neighboring elements by first finding elements // with identical side keys and then check to see if they // are neighbors { // data structures -- Use the hash_multimap if available typedef unsigned int key_type; typedef std::pair<Elem *, unsigned char> val_type; typedef std::pair<key_type, val_type> key_val_pair; typedef LIBMESH_BEST_UNORDERED_MULTIMAP<key_type, val_type> map_type; // A map from side keys to corresponding elements & side numbers map_type side_to_elem_map; for (element_iterator el = this->elements_begin(); el != el_end; ++el) { Elem * element = *el; for (unsigned char ms=0; ms<element->n_neighbors(); ms++) { next_side: // If we haven't yet found a neighbor on this side, try. // Even if we think our neighbor is remote, that // information may be out of date. if (element->neighbor_ptr(ms) == libmesh_nullptr || element->neighbor_ptr(ms) == remote_elem) { // Get the key for the side of this element const unsigned int key = element->key(ms); // Look for elements that have an identical side key std::pair <map_type::iterator, map_type::iterator> bounds = side_to_elem_map.equal_range(key); // May be multiple keys, check all the possible // elements which _might_ be neighbors. if (bounds.first != bounds.second) { // Get the side for this element const UniquePtr<Elem> my_side(element->side_ptr(ms)); // Look at all the entries with an equivalent key while (bounds.first != bounds.second) { // Get the potential element Elem * neighbor = bounds.first->second.first; // Get the side for the neighboring element const unsigned int ns = bounds.first->second.second; const UniquePtr<Elem> their_side(neighbor->side_ptr(ns)); //libmesh_assert(my_side.get()); //libmesh_assert(their_side.get()); // If found a match with my side // // We need special tests here for 1D: // since parents and children have an equal // side (i.e. a node), we need to check // ns != ms, and we also check level() to // avoid setting our neighbor pointer to // any of our neighbor's descendants if( (*my_side == *their_side) && (element->level() == neighbor->level()) && ((element->dim() != 1) || (ns != ms)) ) { // So share a side. Is this a mixed pair // of subactive and active/ancestor // elements? // If not, then we're neighbors. // If so, then the subactive's neighbor is if (element->subactive() == neighbor->subactive()) { // an element is only subactive if it has // been coarsened but not deleted element->set_neighbor (ms,neighbor); neighbor->set_neighbor(ns,element); } else if (element->subactive()) { element->set_neighbor(ms,neighbor); } else if (neighbor->subactive()) { neighbor->set_neighbor(ns,element); } side_to_elem_map.erase (bounds.first); // get out of this nested crap goto next_side; } ++bounds.first; } } // didn't find a match... // Build the map entry for this element key_val_pair kvp; kvp.first = key; kvp.second.first = element; kvp.second.second = ms; // use the lower bound as a hint for // where to put it. #if defined(LIBMESH_HAVE_UNORDERED_MAP) || defined(LIBMESH_HAVE_TR1_UNORDERED_MAP) || defined(LIBMESH_HAVE_HASH_MAP) || defined(LIBMESH_HAVE_EXT_HASH_MAP) side_to_elem_map.insert (kvp); #else side_to_elem_map.insert (bounds.first,kvp); #endif } } } } #ifdef LIBMESH_ENABLE_AMR /** * Here we look at all of the child elements which * don't already have valid neighbors. * * If a child element has a NULL neighbor it is * either because it is on the boundary or because * its neighbor is at a different level. In the * latter case we must get the neighbor from the * parent. * * If a child element has a remote_elem neighbor * on a boundary it shares with its parent, that * info may have become out-dated through coarsening * of the neighbor's parent. In this case, if the * parent's neighbor is active then the child should * share it. * * Furthermore, that neighbor better be active, * otherwise we missed a child somewhere. * * * We also need to look through children ordered by increasing * refinement level in order to add new interior_parent() links in * boundary elements which have just been generated by refinement, * and fix links in boundary elements whose previous * interior_parent() has just been coarsened away. */ const unsigned int n_levels = MeshTools::n_levels(*this); for (unsigned int level = 1; level < n_levels; ++level) { element_iterator end = this->level_elements_end(level); for (element_iterator el = this->level_elements_begin(level); el != end; ++el) { Elem * current_elem = *el; libmesh_assert(current_elem); Elem * parent = current_elem->parent(); libmesh_assert(parent); const unsigned int my_child_num = parent->which_child_am_i(current_elem); for (unsigned int s=0; s < current_elem->n_neighbors(); s++) { if (current_elem->neighbor_ptr(s) == libmesh_nullptr || (current_elem->neighbor_ptr(s) == remote_elem && parent->is_child_on_side(my_child_num, s))) { Elem * neigh = parent->neighbor_ptr(s); // If neigh was refined and had non-subactive children // made remote earlier, then a non-subactive elem should // actually have one of those remote children as a // neighbor if (neigh && (neigh->ancestor()) && (!current_elem->subactive())) { #ifdef DEBUG // Let's make sure that "had children made remote" // situation is actually the case libmesh_assert(neigh->has_children()); bool neigh_has_remote_children = false; for (unsigned int c = 0; c != neigh->n_children(); ++c) { if (neigh->child_ptr(c) == remote_elem) neigh_has_remote_children = true; } libmesh_assert(neigh_has_remote_children); // And let's double-check that we don't have // a remote_elem neighboring a local element libmesh_assert_not_equal_to (current_elem->processor_id(), this->processor_id()); #endif // DEBUG neigh = const_cast<RemoteElem *>(remote_elem); } if (!current_elem->subactive()) current_elem->set_neighbor(s, neigh); #ifdef DEBUG if (neigh != libmesh_nullptr && neigh != remote_elem) // We ignore subactive elements here because // we don't care about neighbors of subactive element. if ((!neigh->active()) && (!current_elem->subactive())) { libMesh::err << "On processor " << this->processor_id() << std::endl; libMesh::err << "Bad element ID = " << current_elem->id() << ", Side " << s << ", Bad neighbor ID = " << neigh->id() << std::endl; libMesh::err << "Bad element proc_ID = " << current_elem->processor_id() << ", Bad neighbor proc_ID = " << neigh->processor_id() << std::endl; libMesh::err << "Bad element size = " << current_elem->hmin() << ", Bad neighbor size = " << neigh->hmin() << std::endl; libMesh::err << "Bad element center = " << current_elem->centroid() << ", Bad neighbor center = " << neigh->centroid() << std::endl; libMesh::err << "ERROR: " << (current_elem->active()?"Active":"Ancestor") << " Element at level " << current_elem->level() << std::endl; libMesh::err << "with " << (parent->active()?"active": (parent->subactive()?"subactive":"ancestor")) << " parent share " << (neigh->subactive()?"subactive":"ancestor") << " neighbor at level " << neigh->level() << std::endl; NameBasedIO(*this).write ("bad_mesh.gmv"); libmesh_error_msg("Problematic mesh written to bad_mesh.gmv."); } #endif // DEBUG } } // We can skip to the next element if we're full-dimension // and therefore don't have any interior parents if (current_elem->dim() >= LIBMESH_DIM) continue; // We have no interior parents unless we can find one later current_elem->set_interior_parent(libmesh_nullptr); Elem * pip = parent->interior_parent(); if (!pip) continue; // If there's no interior_parent children, whether due to a // remote element or a non-conformity, then there's no // children to search. if (pip == remote_elem || pip->active()) { current_elem->set_interior_parent(pip); continue; } // For node comparisons we'll need a sensible tolerance Real node_tolerance = current_elem->hmin() * TOLERANCE; // Otherwise our interior_parent should be a child of our // parent's interior_parent. for (unsigned int c=0; c != pip->n_children(); ++c) { Elem * child = pip->child_ptr(c); // If we have a remote_elem, that might be our // interior_parent. We'll set it provisionally now and // keep trying to find something better. if (child == remote_elem) { current_elem->set_interior_parent (const_cast<RemoteElem *>(remote_elem)); continue; } bool child_contains_our_nodes = true; for (unsigned int n=0; n != current_elem->n_nodes(); ++n) { bool child_contains_this_node = false; for (unsigned int cn=0; cn != child->n_nodes(); ++cn) if (child->point(cn).absolute_fuzzy_equals (current_elem->point(n), node_tolerance)) { child_contains_this_node = true; break; } if (!child_contains_this_node) { child_contains_our_nodes = false; break; } } if (child_contains_our_nodes) { current_elem->set_interior_parent(child); break; } } // We should have found *some* interior_parent at this // point, whether semilocal or remote. libmesh_assert(current_elem->interior_parent()); } } #endif // AMR #ifdef DEBUG MeshTools::libmesh_assert_valid_neighbors(*this, !reset_remote_elements); MeshTools::libmesh_assert_valid_amr_interior_parents(*this); #endif }
Elem * Packing<Elem *>::unpack (std::vector<largest_id_type>::const_iterator in, MeshBase * mesh) { #ifndef NDEBUG const std::vector<largest_id_type>::const_iterator original_in = in; const largest_id_type incoming_header = *in++; libmesh_assert_equal_to (incoming_header, elem_magic_header); #endif // int 0: level const unsigned int level = cast_int<unsigned int>(*in++); #ifdef LIBMESH_ENABLE_AMR // int 1: p level const unsigned int p_level = cast_int<unsigned int>(*in++); // int 2: refinement flag and encoded has_children const int rflag = cast_int<int>(*in++); const int invalid_rflag = cast_int<int>(Elem::INVALID_REFINEMENTSTATE); libmesh_assert_greater_equal (rflag, 0); libmesh_assert_less (rflag, invalid_rflag*2+1); const bool has_children = (rflag > invalid_rflag); const Elem::RefinementState refinement_flag = has_children ? cast_int<Elem::RefinementState>(rflag - invalid_rflag - 1) : cast_int<Elem::RefinementState>(rflag); // int 3: p refinement flag const int pflag = cast_int<int>(*in++); libmesh_assert_greater_equal (pflag, 0); libmesh_assert_less (pflag, Elem::INVALID_REFINEMENTSTATE); const Elem::RefinementState p_refinement_flag = cast_int<Elem::RefinementState>(pflag); #else in += 3; #endif // LIBMESH_ENABLE_AMR // int 4: element type const int typeint = cast_int<int>(*in++); libmesh_assert_greater_equal (typeint, 0); libmesh_assert_less (typeint, INVALID_ELEM); const ElemType type = cast_int<ElemType>(typeint); const unsigned int n_nodes = Elem::type_to_n_nodes_map[type]; // int 5: processor id const processor_id_type processor_id = cast_int<processor_id_type>(*in++); libmesh_assert (processor_id < mesh->n_processors() || processor_id == DofObject::invalid_processor_id); // int 6: subdomain id const subdomain_id_type subdomain_id = cast_int<subdomain_id_type>(*in++); // int 7: dof object id const dof_id_type id = cast_int<dof_id_type>(*in++); libmesh_assert_not_equal_to (id, DofObject::invalid_id); #ifdef LIBMESH_ENABLE_UNIQUE_ID // int 8: dof object unique id const unique_id_type unique_id = cast_int<unique_id_type>(*in++); #endif #ifdef LIBMESH_ENABLE_AMR // int 9: parent dof object id. // Note: If level==0, then (*in) == invalid_id. In // this case, the equality check in cast_int<unsigned>(*in) will // never succeed. Therefore, we should only attempt the more // rigorous cast verification in cases where level != 0. const dof_id_type parent_id = (level == 0) ? static_cast<dof_id_type>(*in++) : cast_int<dof_id_type>(*in++); libmesh_assert (level == 0 || parent_id != DofObject::invalid_id); libmesh_assert (level != 0 || parent_id == DofObject::invalid_id); // int 10: local child id // Note: If level==0, then which_child_am_i is not valid, so don't // do the more rigorous cast verification. const unsigned int which_child_am_i = (level == 0) ? static_cast<unsigned int>(*in++) : cast_int<unsigned int>(*in++); #else in += 2; #endif // LIBMESH_ENABLE_AMR const dof_id_type interior_parent_id = static_cast<dof_id_type>(*in++); // Make sure we don't miscount above when adding the "magic" header // plus the real data header libmesh_assert_equal_to (in - original_in, header_size + 1); Elem * elem = mesh->query_elem_ptr(id); // if we already have this element, make sure its // properties match, and update any missing neighbor // links, but then go on if (elem) { libmesh_assert_equal_to (elem->level(), level); libmesh_assert_equal_to (elem->id(), id); //#ifdef LIBMESH_ENABLE_UNIQUE_ID // No check for unique id sanity //#endif libmesh_assert_equal_to (elem->processor_id(), processor_id); libmesh_assert_equal_to (elem->subdomain_id(), subdomain_id); libmesh_assert_equal_to (elem->type(), type); libmesh_assert_equal_to (elem->n_nodes(), n_nodes); #ifndef NDEBUG // All our nodes should be correct for (unsigned int i=0; i != n_nodes; ++i) libmesh_assert(elem->node_id(i) == cast_int<dof_id_type>(*in++)); #else in += n_nodes; #endif #ifdef LIBMESH_ENABLE_AMR libmesh_assert_equal_to (elem->refinement_flag(), refinement_flag); libmesh_assert_equal_to (elem->has_children(), has_children); #ifdef DEBUG if (elem->active()) { libmesh_assert_equal_to (elem->p_level(), p_level); libmesh_assert_equal_to (elem->p_refinement_flag(), p_refinement_flag); } #endif libmesh_assert (!level || elem->parent() != libmesh_nullptr); libmesh_assert (!level || elem->parent()->id() == parent_id); libmesh_assert (!level || elem->parent()->child_ptr(which_child_am_i) == elem); #endif // Our interior_parent link should be "close to" correct - we // may have to update it, but we can check for some // inconsistencies. { // If the sending processor sees no interior_parent here, we'd // better agree. if (interior_parent_id == DofObject::invalid_id) { if (elem->dim() < LIBMESH_DIM) libmesh_assert (!(elem->interior_parent())); } // If the sending processor has a remote_elem interior_parent, // then all we know is that we'd better have *some* // interior_parent else if (interior_parent_id == remote_elem->id()) { libmesh_assert(elem->interior_parent()); } else { Elem * ip = mesh->query_elem_ptr(interior_parent_id); // The sending processor sees an interior parent here, so // if we don't have that interior element, then we'd // better have a remote_elem signifying that fact. if (!ip) libmesh_assert_equal_to (elem->interior_parent(), remote_elem); else { // The sending processor has an interior_parent here, // and we have that element, but that does *NOT* mean // we're already linking to it. Perhaps we initially // received elem from a processor on which the // interior_parent link was remote? libmesh_assert(elem->interior_parent() == ip || elem->interior_parent() == remote_elem); // If the link was originally remote, update it if (elem->interior_parent() == remote_elem) { elem->set_interior_parent(ip); } } } } // Our neighbor links should be "close to" correct - we may have // to update a remote_elem link, and we can check for possible // inconsistencies along the way. // // For subactive elements, we don't bother keeping neighbor // links in good shape, so there's nothing we need to set or can // safely assert here. if (!elem->subactive()) for (auto n : elem->side_index_range()) { const dof_id_type neighbor_id = cast_int<dof_id_type>(*in++); // If the sending processor sees a domain boundary here, // we'd better agree. if (neighbor_id == DofObject::invalid_id) { libmesh_assert (!(elem->neighbor_ptr(n))); continue; } // If the sending processor has a remote_elem neighbor here, // then all we know is that we'd better *not* have a domain // boundary. if (neighbor_id == remote_elem->id()) { libmesh_assert(elem->neighbor_ptr(n)); continue; } Elem * neigh = mesh->query_elem_ptr(neighbor_id); // The sending processor sees a neighbor here, so if we // don't have that neighboring element, then we'd better // have a remote_elem signifying that fact. if (!neigh) { libmesh_assert_equal_to (elem->neighbor_ptr(n), remote_elem); continue; } // The sending processor has a neighbor here, and we have // that element, but that does *NOT* mean we're already // linking to it. Perhaps we initially received both elem // and neigh from processors on which their mutual link was // remote? libmesh_assert(elem->neighbor_ptr(n) == neigh || elem->neighbor_ptr(n) == remote_elem); // If the link was originally remote, we should update it, // and make sure the appropriate parts of its family link // back to us. if (elem->neighbor_ptr(n) == remote_elem) { elem->set_neighbor(n, neigh); elem->make_links_to_me_local(n); } } // Our p level and refinement flags should be "close to" correct // if we're not an active element - we might have a p level // increased or decreased by changes in remote_elem children. // // But if we have remote_elem children, then we shouldn't be // doing a projection on this inactive element on this // processor, so we won't need correct p settings. Couldn't // hurt to update, though. #ifdef LIBMESH_ENABLE_AMR if (elem->processor_id() != mesh->processor_id()) { elem->hack_p_level(p_level); elem->set_p_refinement_flag(p_refinement_flag); } #endif // LIBMESH_ENABLE_AMR // FIXME: We should add some debug mode tests to ensure that the // encoded indexing and boundary conditions are consistent. } else { // We don't already have the element, so we need to create it. // Find the parent if necessary Elem * parent = libmesh_nullptr; #ifdef LIBMESH_ENABLE_AMR // Find a child element's parent if (level > 0) { // Note that we must be very careful to construct the send // connectivity so that parents are encountered before // children. If we get here and can't find the parent that // is a fatal error. parent = mesh->elem_ptr(parent_id); } // Or assert that the sending processor sees no parent else libmesh_assert_equal_to (parent_id, DofObject::invalid_id); #else // No non-level-0 elements without AMR libmesh_assert_equal_to (level, 0); #endif elem = Elem::build(type,parent).release(); libmesh_assert (elem); #ifdef LIBMESH_ENABLE_AMR if (level != 0) { // Since this is a newly created element, the parent must // have previously thought of this child as a remote element. libmesh_assert_equal_to (parent->child_ptr(which_child_am_i), remote_elem); parent->add_child(elem, which_child_am_i); } // Assign the refinement flags and levels elem->set_p_level(p_level); elem->set_refinement_flag(refinement_flag); elem->set_p_refinement_flag(p_refinement_flag); libmesh_assert_equal_to (elem->level(), level); // If this element should have children, assign remote_elem to // all of them for now, for consistency. Later unpacked // elements may overwrite that. if (has_children) { const unsigned int nc = elem->n_children(); for (unsigned int c=0; c != nc; ++c) elem->add_child(const_cast<RemoteElem *>(remote_elem), c); } #endif // LIBMESH_ENABLE_AMR // Assign the IDs elem->subdomain_id() = subdomain_id; elem->processor_id() = processor_id; elem->set_id() = id; #ifdef LIBMESH_ENABLE_UNIQUE_ID elem->set_unique_id() = unique_id; #endif // Assign the connectivity libmesh_assert_equal_to (elem->n_nodes(), n_nodes); for (unsigned int n=0; n != n_nodes; n++) elem->set_node(n) = mesh->node_ptr (cast_int<dof_id_type>(*in++)); // Set interior_parent if found { // We may be unpacking an element that was a ghost element on the // sender, in which case the element's interior_parent may not be // known by the packed element. We'll have to set such // interior_parents to remote_elem ourselves and wait for a // later packed element to give us better information. if (interior_parent_id == remote_elem->id()) { elem->set_interior_parent (const_cast<RemoteElem *>(remote_elem)); } else if (interior_parent_id != DofObject::invalid_id) { // If we don't have the interior parent element, then it's // a remote_elem until we get it. Elem * ip = mesh->query_elem_ptr(interior_parent_id); if (!ip ) elem->set_interior_parent (const_cast<RemoteElem *>(remote_elem)); else elem->set_interior_parent(ip); } } for (auto n : elem->side_index_range()) { const dof_id_type neighbor_id = cast_int<dof_id_type>(*in++); if (neighbor_id == DofObject::invalid_id) continue; // We may be unpacking an element that was a ghost element on the // sender, in which case the element's neighbors may not all be // known by the packed element. We'll have to set such // neighbors to remote_elem ourselves and wait for a later // packed element to give us better information. if (neighbor_id == remote_elem->id()) { elem->set_neighbor(n, const_cast<RemoteElem *>(remote_elem)); continue; } // If we don't have the neighbor element, then it's a // remote_elem until we get it. Elem * neigh = mesh->query_elem_ptr(neighbor_id); if (!neigh) { elem->set_neighbor(n, const_cast<RemoteElem *>(remote_elem)); continue; } // If we have the neighbor element, then link to it, and // make sure the appropriate parts of its family link back // to us. elem->set_neighbor(n, neigh); elem->make_links_to_me_local(n); } elem->unpack_indexing(in); } in += elem->packed_indexing_size(); // If this is a coarse element, // add any element side or edge boundary condition ids if (level == 0) { for (auto s : elem->side_index_range()) { const boundary_id_type num_bcs = cast_int<boundary_id_type>(*in++); for (boundary_id_type bc_it=0; bc_it < num_bcs; bc_it++) mesh->get_boundary_info().add_side (elem, s, cast_int<boundary_id_type>(*in++)); } for (auto e : elem->edge_index_range()) { const boundary_id_type num_bcs = cast_int<boundary_id_type>(*in++); for (boundary_id_type bc_it=0; bc_it < num_bcs; bc_it++) mesh->get_boundary_info().add_edge (elem, e, cast_int<boundary_id_type>(*in++)); } for (unsigned short sf=0; sf != 2; ++sf) { const boundary_id_type num_bcs = cast_int<boundary_id_type>(*in++); for (boundary_id_type bc_it=0; bc_it < num_bcs; bc_it++) mesh->get_boundary_info().add_shellface (elem, sf, cast_int<boundary_id_type>(*in++)); } } // Return the new element return elem; }