elem.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 00020 // C++ includes 00021 #include <algorithm> // for std::sort 00022 #include <iterator> // for std::ostream_iterator 00023 #include <sstream> 00024 00025 00026 // Local includes 00027 #include "libmesh/elem.h" 00028 #include "libmesh/fe_type.h" 00029 #include "libmesh/fe_interface.h" 00030 #include "libmesh/edge_edge2.h" 00031 #include "libmesh/edge_edge3.h" 00032 #include "libmesh/edge_edge4.h" 00033 #include "libmesh/edge_inf_edge2.h" 00034 #include "libmesh/face_tri3.h" 00035 #include "libmesh/face_tri6.h" 00036 #include "libmesh/face_quad4.h" 00037 #include "libmesh/face_quad8.h" 00038 #include "libmesh/face_quad9.h" 00039 #include "libmesh/face_inf_quad4.h" 00040 #include "libmesh/face_inf_quad6.h" 00041 #include "libmesh/cell_tet4.h" 00042 #include "libmesh/cell_tet10.h" 00043 #include "libmesh/cell_hex8.h" 00044 #include "libmesh/cell_hex20.h" 00045 #include "libmesh/cell_hex27.h" 00046 #include "libmesh/cell_inf_hex8.h" 00047 #include "libmesh/cell_inf_hex16.h" 00048 #include "libmesh/cell_inf_hex18.h" 00049 #include "libmesh/cell_prism6.h" 00050 #include "libmesh/cell_prism15.h" 00051 #include "libmesh/cell_prism18.h" 00052 #include "libmesh/cell_inf_prism6.h" 00053 #include "libmesh/cell_inf_prism12.h" 00054 #include "libmesh/cell_pyramid5.h" 00055 #include "libmesh/fe_base.h" 00056 #include "libmesh/mesh_base.h" 00057 #include "libmesh/quadrature_gauss.h" 00058 #include "libmesh/remote_elem.h" 00059 #include "libmesh/string_to_enum.h" 00060 00061 #ifdef LIBMESH_ENABLE_PERIODIC 00062 #include "libmesh/mesh.h" 00063 #include "libmesh/periodic_boundaries.h" 00064 #include "libmesh/boundary_info.h" 00065 #endif 00066 00067 namespace libMesh 00068 { 00069 00070 // Initialize static member variables 00071 const unsigned int Elem::type_to_n_nodes_map [] = 00072 { 00073 2, // EDGE2 00074 3, // EDGE3 00075 4, // EDGE4 00076 00077 3, // TRI3 00078 6, // TRI6 00079 00080 4, // QUAD4 00081 8, // QUAD8 00082 9, // QUAD9 00083 00084 4, // TET4 00085 10, // TET10 00086 00087 8, // HEX8 00088 20, // HEX20 00089 27, // HEX27 00090 00091 6, // PRISM6 00092 15, // PRISM15 00093 18, // PRISM18 00094 00095 5, // PYRAMID5 00096 00097 2, // INFEDGE2 00098 00099 4, // INFQUAD4 00100 6, // INFQUAD6 00101 00102 8, // INFHEX8 00103 16, // INFHEX16 00104 18, // INFHEX18 00105 00106 6, // INFPRISM6 00107 16, // INFPRISM12 00108 00109 1, // NODEELEM 00110 }; 00111 00112 const unsigned int Elem::type_to_n_sides_map [] = 00113 { 00114 2, // EDGE2 00115 2, // EDGE3 00116 2, // EDGE4 00117 00118 3, // TRI3 00119 3, // TRI6 00120 00121 4, // QUAD4 00122 4, // QUAD8 00123 4, // QUAD9 00124 00125 4, // TET4 00126 4, // TET10 00127 00128 6, // HEX8 00129 6, // HEX20 00130 6, // HEX27 00131 00132 5, // PRISM6 00133 5, // PRISM15 00134 5, // PRISM18 00135 00136 5, // PYRAMID5 00137 00138 2, // INFEDGE2 00139 00140 3, // INFQUAD4 00141 3, // INFQUAD6 00142 00143 5, // INFHEX8 00144 5, // INFHEX16 00145 5, // INFHEX18 00146 00147 4, // INFPRISM6 00148 4, // INFPRISM12 00149 00150 0, // NODEELEM 00151 }; 00152 00153 00154 // ------------------------------------------------------------ 00155 // Elem class member funcions 00156 AutoPtr<Elem> Elem::build(const ElemType type, 00157 Elem* p) 00158 { 00159 Elem* elem = NULL; 00160 00161 switch (type) 00162 { 00163 // 1D elements 00164 case EDGE2: 00165 { 00166 elem = new Edge2(p); 00167 break; 00168 } 00169 case EDGE3: 00170 { 00171 elem = new Edge3(p); 00172 break; 00173 } 00174 case EDGE4: 00175 { 00176 elem = new Edge4(p); 00177 break; 00178 } 00179 00180 00181 00182 // 2D elements 00183 case TRI3: 00184 { 00185 elem = new Tri3(p); 00186 break; 00187 } 00188 case TRI6: 00189 { 00190 elem = new Tri6(p); 00191 break; 00192 } 00193 case QUAD4: 00194 { 00195 elem = new Quad4(p); 00196 break; 00197 } 00198 case QUAD8: 00199 { 00200 elem = new Quad8(p); 00201 break; 00202 } 00203 case QUAD9: 00204 { 00205 elem = new Quad9(p); 00206 break; 00207 } 00208 00209 00210 // 3D elements 00211 case TET4: 00212 { 00213 elem = new Tet4(p); 00214 break; 00215 } 00216 case TET10: 00217 { 00218 elem = new Tet10(p); 00219 break; 00220 } 00221 case HEX8: 00222 { 00223 elem = new Hex8(p); 00224 break; 00225 } 00226 case HEX20: 00227 { 00228 elem = new Hex20(p); 00229 break; 00230 } 00231 case HEX27: 00232 { 00233 elem = new Hex27(p); 00234 break; 00235 } 00236 case PRISM6: 00237 { 00238 elem = new Prism6(p); 00239 break; 00240 } 00241 case PRISM15: 00242 { 00243 elem = new Prism15(p); 00244 break; 00245 } 00246 case PRISM18: 00247 { 00248 elem = new Prism18(p); 00249 break; 00250 } 00251 case PYRAMID5: 00252 { 00253 elem = new Pyramid5(p); 00254 break; 00255 } 00256 00257 00258 00259 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00260 00261 // 1D infinite elements 00262 case INFEDGE2: 00263 { 00264 elem = new InfEdge2(p); 00265 break; 00266 } 00267 00268 00269 // 2D infinite elements 00270 case INFQUAD4: 00271 { 00272 elem = new InfQuad4(p); 00273 break; 00274 } 00275 case INFQUAD6: 00276 { 00277 elem = new InfQuad6(p); 00278 break; 00279 } 00280 00281 00282 // 3D infinite elements 00283 case INFHEX8: 00284 { 00285 elem = new InfHex8(p); 00286 break; 00287 } 00288 case INFHEX16: 00289 { 00290 elem = new InfHex16(p); 00291 break; 00292 } 00293 case INFHEX18: 00294 { 00295 elem = new InfHex18(p); 00296 break; 00297 } 00298 case INFPRISM6: 00299 { 00300 elem = new InfPrism6(p); 00301 break; 00302 } 00303 case INFPRISM12: 00304 { 00305 elem = new InfPrism12(p); 00306 break; 00307 } 00308 00309 #endif 00310 00311 default: 00312 { 00313 libMesh::err << "ERROR: Undefined element type!." << std::endl; 00314 libmesh_error(); 00315 } 00316 } 00317 00318 00319 AutoPtr<Elem> ap(elem); 00320 return ap; 00321 } 00322 00323 00324 00325 Point Elem::centroid() const 00326 { 00327 Point cp; 00328 00329 for (unsigned int n=0; n<this->n_vertices(); n++) 00330 cp.add (this->point(n)); 00331 00332 return (cp /= static_cast<Real>(this->n_vertices())); 00333 } 00334 00335 00336 00337 Real Elem::hmin() const 00338 { 00339 Real h_min=1.e30; 00340 00341 for (unsigned int n_outer=0; n_outer<this->n_vertices(); n_outer++) 00342 for (unsigned int n_inner=n_outer+1; n_inner<this->n_vertices(); n_inner++) 00343 { 00344 const Point diff = (this->point(n_outer) - this->point(n_inner)); 00345 00346 h_min = std::min(h_min,diff.size()); 00347 } 00348 00349 return h_min; 00350 } 00351 00352 00353 00354 Real Elem::hmax() const 00355 { 00356 Real h_max=0; 00357 00358 for (unsigned int n_outer=0; n_outer<this->n_vertices(); n_outer++) 00359 for (unsigned int n_inner=n_outer+1; n_inner<this->n_vertices(); n_inner++) 00360 { 00361 const Point diff = (this->point(n_outer) - this->point(n_inner)); 00362 00363 h_max = std::max(h_max,diff.size()); 00364 } 00365 00366 return h_max; 00367 } 00368 00369 00370 00371 Real Elem::length(const unsigned int n1, 00372 const unsigned int n2) const 00373 { 00374 libmesh_assert_less ( n1, this->n_vertices() ); 00375 libmesh_assert_less ( n2, this->n_vertices() ); 00376 00377 return (this->point(n1) - this->point(n2)).size(); 00378 } 00379 00380 00381 00382 bool Elem::operator == (const Elem& rhs) const 00383 { 00384 00385 // Cast rhs to an Elem* 00386 // const Elem* rhs_elem = dynamic_cast<const Elem*>(&rhs); 00387 const Elem* rhs_elem = &rhs; 00388 00389 // If we cannot cast to an Elem*, rhs must be a Node 00390 // if(rhs_elem == static_cast<const Elem*>(NULL)) 00391 // return false; 00392 00393 // libmesh_assert (n_nodes()); 00394 // libmesh_assert (rhs.n_nodes()); 00395 00396 // // Elements can only be equal if they 00397 // // contain the same number of nodes. 00398 // if (this->n_nodes() == rhs.n_nodes()) 00399 // { 00400 // // Create a set that contains our global 00401 // // node numbers and those of our neighbor. 00402 // // If the set is the same size as the number 00403 // // of nodes in both elements then they must 00404 // // be connected to the same nodes. 00405 // std::set<unsigned int> nodes_set; 00406 00407 // for (unsigned int n=0; n<this->n_nodes(); n++) 00408 // { 00409 // nodes_set.insert(this->node(n)); 00410 // nodes_set.insert(rhs.node(n)); 00411 // } 00412 00413 // // If this passes the elements are connected 00414 // // to the same global nodes 00415 // if (nodes_set.size() == this->n_nodes()) 00416 // return true; 00417 // } 00418 00419 // // If we get here it is because the elements either 00420 // // do not have the same number of nodes or they are 00421 // // connected to different nodes. Either way they 00422 // // are not the same element 00423 // return false; 00424 00425 // Useful typedefs 00426 typedef std::vector<dof_id_type>::iterator iterator; 00427 00428 00429 // Elements can only be equal if they 00430 // contain the same number of nodes. 00431 // However, we will only test the vertices, 00432 // which is sufficient & cheaper 00433 if (this->n_nodes() == rhs_elem->n_nodes()) 00434 { 00435 // The number of nodes in the element 00436 const unsigned int nn = this->n_nodes(); 00437 00438 // Create a vector that contains our global 00439 // node numbers and those of our neighbor. 00440 // If the sorted, unique vector is the same size 00441 // as the number of nodes in both elements then 00442 // they must be connected to the same nodes. 00443 // 00444 // The vector will be no larger than 2*n_nodes(), 00445 // so we might as well reserve the space. 00446 std::vector<dof_id_type> common_nodes; 00447 common_nodes.reserve (2*nn); 00448 00449 // Add the global indices of the nodes 00450 for (unsigned int n=0; n<nn; n++) 00451 { 00452 common_nodes.push_back (this->node(n)); 00453 common_nodes.push_back (rhs_elem->node(n)); 00454 } 00455 00456 // Sort the vector and find out how long 00457 // the sorted vector is. 00458 std::sort (common_nodes.begin(), common_nodes.end()); 00459 00460 iterator new_end = std::unique (common_nodes.begin(), 00461 common_nodes.end()); 00462 00463 const int new_size = libmesh_cast_int<int> 00464 (std::distance (common_nodes.begin(), new_end)); 00465 00466 // If this passes the elements are connected 00467 // to the same global vertex nodes 00468 if (new_size == static_cast<int>(nn)) 00469 return true; 00470 } 00471 00472 // If we get here it is because the elements either 00473 // do not have the same number of nodes or they are 00474 // connected to different nodes. Either way they 00475 // are not the same element 00476 return false; 00477 } 00478 00479 00480 00481 bool Elem::is_semilocal() const 00482 { 00483 std::set<const Elem *> point_neighbors; 00484 00485 this->find_point_neighbors(point_neighbors); 00486 00487 std::set<const Elem*>::const_iterator it = point_neighbors.begin(); 00488 const std::set<const Elem*>::const_iterator end = point_neighbors.end(); 00489 00490 for (; it != end; ++it) 00491 { 00492 const Elem* elem = *it; 00493 if (elem->processor_id() == libMesh::processor_id()) 00494 return true; 00495 } 00496 00497 return false; 00498 } 00499 00500 00501 00502 bool Elem::contains_vertex_of(const Elem *e) const 00503 { 00504 // Our vertices are the first numbered nodes 00505 for (unsigned int n = 0; n != e->n_vertices(); ++n) 00506 if (this->contains_point(e->point(n))) 00507 return true; 00508 return false; 00509 } 00510 00511 00512 00513 bool Elem::contains_edge_of(const Elem *e) const 00514 { 00515 unsigned int num_contained_edges = 0; 00516 00517 // Our vertices are the first numbered nodes 00518 for (unsigned int n = 0; n != e->n_vertices(); ++n) 00519 { 00520 if (this->contains_point(e->point(n))) 00521 { 00522 num_contained_edges++; 00523 if(num_contained_edges>=2) 00524 { 00525 return true; 00526 } 00527 } 00528 } 00529 return false; 00530 } 00531 00532 00533 00534 void Elem::find_point_neighbors(const Point &p, 00535 std::set<const Elem *> &neighbor_set) const 00536 { 00537 libmesh_assert(this->contains_point(p)); 00538 00539 neighbor_set.clear(); 00540 neighbor_set.insert(this); 00541 00542 std::set<const Elem *> untested_set, next_untested_set; 00543 untested_set.insert(this); 00544 00545 while (!untested_set.empty()) 00546 { 00547 // Loop over all the elements in the patch that haven't already 00548 // been tested 00549 std::set<const Elem*>::const_iterator it = untested_set.begin(); 00550 const std::set<const Elem*>::const_iterator end = untested_set.end(); 00551 00552 for (; it != end; ++it) 00553 { 00554 const Elem* elem = *it; 00555 00556 for (unsigned int s=0; s<elem->n_sides(); s++) 00557 { 00558 const Elem* current_neighbor = elem->neighbor(s); 00559 if (current_neighbor && 00560 current_neighbor != remote_elem) // we have a real neighbor on this side 00561 { 00562 if (current_neighbor->active()) // ... if it is active 00563 { 00564 if (current_neighbor->contains_point(p)) // ... and touches p 00565 { 00566 // Make sure we'll test it 00567 if (!neighbor_set.count(current_neighbor)) 00568 next_untested_set.insert (current_neighbor); 00569 00570 // And add it 00571 neighbor_set.insert (current_neighbor); 00572 } 00573 } 00574 #ifdef LIBMESH_ENABLE_AMR 00575 else // ... the neighbor is *not* active, 00576 { // ... so add *all* neighboring 00577 // active children that touch p 00578 std::vector<const Elem*> active_neighbor_children; 00579 00580 current_neighbor->active_family_tree_by_neighbor 00581 (active_neighbor_children, elem); 00582 00583 std::vector<const Elem*>::const_iterator 00584 child_it = active_neighbor_children.begin(); 00585 const std::vector<const Elem*>::const_iterator 00586 child_end = active_neighbor_children.end(); 00587 for (; child_it != child_end; ++child_it) 00588 { 00589 const Elem *current_child = *child_it; 00590 if (current_child->contains_point(p)) 00591 { 00592 // Make sure we'll test it 00593 if (!neighbor_set.count(current_child)) 00594 next_untested_set.insert (current_child); 00595 00596 neighbor_set.insert (current_child); 00597 } 00598 } 00599 } 00600 #endif // #ifdef LIBMESH_ENABLE_AMR 00601 } 00602 } 00603 } 00604 untested_set.swap(next_untested_set); 00605 next_untested_set.clear(); 00606 } 00607 } 00608 00609 00610 00611 void Elem::find_point_neighbors(std::set<const Elem *> &neighbor_set) const 00612 { 00613 neighbor_set.clear(); 00614 neighbor_set.insert(this); 00615 00616 std::set<const Elem *> untested_set, next_untested_set; 00617 untested_set.insert(this); 00618 00619 while (!untested_set.empty()) 00620 { 00621 // Loop over all the elements in the patch that haven't already 00622 // been tested 00623 std::set<const Elem*>::const_iterator it = untested_set.begin(); 00624 const std::set<const Elem*>::const_iterator end = untested_set.end(); 00625 00626 for (; it != end; ++it) 00627 { 00628 const Elem* elem = *it; 00629 00630 for (unsigned int s=0; s<elem->n_sides(); s++) 00631 { 00632 const Elem* current_neighbor = elem->neighbor(s); 00633 if (current_neighbor && 00634 current_neighbor != remote_elem) // we have a real neighbor on this side 00635 { 00636 if (current_neighbor->active()) // ... if it is active 00637 { 00638 if (this->contains_vertex_of(current_neighbor) // ... and touches us 00639 || current_neighbor->contains_vertex_of(this)) 00640 { 00641 // Make sure we'll test it 00642 if (!neighbor_set.count(current_neighbor)) 00643 next_untested_set.insert (current_neighbor); 00644 00645 // And add it 00646 neighbor_set.insert (current_neighbor); 00647 } 00648 } 00649 #ifdef LIBMESH_ENABLE_AMR 00650 else // ... the neighbor is *not* active, 00651 { // ... so add *all* neighboring 00652 // active children 00653 std::vector<const Elem*> active_neighbor_children; 00654 00655 current_neighbor->active_family_tree_by_neighbor 00656 (active_neighbor_children, elem); 00657 00658 std::vector<const Elem*>::const_iterator 00659 child_it = active_neighbor_children.begin(); 00660 const std::vector<const Elem*>::const_iterator 00661 child_end = active_neighbor_children.end(); 00662 for (; child_it != child_end; ++child_it) 00663 { 00664 const Elem *current_child = *child_it; 00665 if (this->contains_vertex_of(current_child) || 00666 (current_child)->contains_vertex_of(this)) 00667 { 00668 // Make sure we'll test it 00669 if (!neighbor_set.count(current_child)) 00670 next_untested_set.insert (current_child); 00671 00672 neighbor_set.insert (current_child); 00673 } 00674 } 00675 } 00676 #endif // #ifdef LIBMESH_ENABLE_AMR 00677 } 00678 } 00679 } 00680 untested_set.swap(next_untested_set); 00681 next_untested_set.clear(); 00682 } 00683 } 00684 00685 00686 00687 void Elem::find_edge_neighbors(const Point& p1, 00688 const Point& p2, 00689 std::set<const Elem *> &neighbor_set) const 00690 { 00691 // Simple but perhaps suboptimal code: find elements containing the 00692 // first point, then winnow this set down by removing elements which 00693 // don't also contain the second point 00694 00695 libmesh_assert(this->contains_point(p2)); 00696 this->find_point_neighbors(p1, neighbor_set); 00697 00698 std::set<const Elem*>::iterator it = neighbor_set.begin(); 00699 const std::set<const Elem*>::iterator end = neighbor_set.end(); 00700 00701 while(it != end) { 00702 std::set<const Elem*>::iterator current = it++; 00703 00704 const Elem* elem = *current; 00705 // This won't invalidate iterator it, because it is already 00706 // pointing to the next element 00707 if (!elem->contains_point(p2)) 00708 neighbor_set.erase(current); 00709 } 00710 } 00711 00712 00713 00714 void Elem::find_edge_neighbors(std::set<const Elem *> &neighbor_set) const 00715 { 00716 neighbor_set.clear(); 00717 neighbor_set.insert(this); 00718 00719 std::set<const Elem *> untested_set, next_untested_set; 00720 untested_set.insert(this); 00721 00722 while (!untested_set.empty()) 00723 { 00724 // Loop over all the elements in the patch that haven't already 00725 // been tested 00726 std::set<const Elem*>::const_iterator it = untested_set.begin(); 00727 const std::set<const Elem*>::const_iterator end = untested_set.end(); 00728 00729 for (; it != end; ++it) 00730 { 00731 const Elem* elem = *it; 00732 00733 for (unsigned int s=0; s<elem->n_sides(); s++) 00734 { 00735 const Elem* current_neighbor = elem->neighbor(s); 00736 if (current_neighbor && 00737 current_neighbor != remote_elem) // we have a real neighbor on this side 00738 { 00739 if (current_neighbor->active()) // ... if it is active 00740 { 00741 if (this->contains_edge_of(current_neighbor) // ... and touches us 00742 || current_neighbor->contains_edge_of(this)) 00743 { 00744 // Make sure we'll test it 00745 if (!neighbor_set.count(current_neighbor)) 00746 next_untested_set.insert (current_neighbor); 00747 00748 // And add it 00749 neighbor_set.insert (current_neighbor); 00750 } 00751 } 00752 #ifdef LIBMESH_ENABLE_AMR 00753 else // ... the neighbor is *not* active, 00754 { // ... so add *all* neighboring 00755 // active children 00756 std::vector<const Elem*> active_neighbor_children; 00757 00758 current_neighbor->active_family_tree_by_neighbor 00759 (active_neighbor_children, elem); 00760 00761 std::vector<const Elem*>::const_iterator 00762 child_it = active_neighbor_children.begin(); 00763 const std::vector<const Elem*>::const_iterator 00764 child_end = active_neighbor_children.end(); 00765 for (; child_it != child_end; ++child_it) 00766 { 00767 const Elem *current_child = *child_it; 00768 if (this->contains_edge_of(*child_it) || 00769 (*child_it)->contains_edge_of(this)) 00770 { 00771 // Make sure we'll test it 00772 if (!neighbor_set.count(current_child)) 00773 next_untested_set.insert (current_child); 00774 00775 neighbor_set.insert (current_child); 00776 } 00777 } 00778 } 00779 #endif // #ifdef LIBMESH_ENABLE_AMR 00780 } 00781 } 00782 } 00783 untested_set.swap(next_untested_set); 00784 next_untested_set.clear(); 00785 } 00786 } 00787 00788 #ifdef LIBMESH_ENABLE_PERIODIC 00789 00790 Elem* Elem::topological_neighbor (const unsigned int i, 00791 MeshBase & mesh, 00792 const PointLocatorBase& point_locator, 00793 const PeriodicBoundaries * pb) 00794 { 00795 libmesh_assert_less (i, this->n_neighbors()); 00796 00797 Elem * neighbor_i = this->neighbor(i); 00798 if (neighbor_i != NULL) 00799 return neighbor_i; 00800 00801 if (pb) 00802 { 00803 // Since the neighbor is NULL it must be on a boundary. We need 00804 // see if this is a periodic boundary in which case it will have a 00805 // topological neighbor 00806 00807 std::vector<boundary_id_type> boundary_ids = mesh.boundary_info->boundary_ids(this, i); 00808 for (std::vector<boundary_id_type>::iterator j = boundary_ids.begin(); j != boundary_ids.end(); ++j) 00809 if (pb->boundary(*j)) 00810 { 00811 // Since the point locator inside of periodic boundaries 00812 // returns a const pointer we will retrieve the proper 00813 // pointer directly from the mesh object. Also since coarse 00814 // elements do not have more refined neighbors we need to make 00815 // sure that we don't return one of these types of neighbors. 00816 neighbor_i = mesh.elem(pb->neighbor(*j, point_locator, this, i)->id()); 00817 if (level() < neighbor_i->level()) 00818 neighbor_i = neighbor_i->parent(); 00819 return neighbor_i; 00820 } 00821 } 00822 00823 return NULL; 00824 } 00825 00826 00827 00828 const Elem* Elem::topological_neighbor (const unsigned int i, 00829 const MeshBase & mesh, 00830 const PointLocatorBase& point_locator, 00831 const PeriodicBoundaries * pb) const 00832 { 00833 libmesh_assert_less (i, this->n_neighbors()); 00834 00835 const Elem * neighbor_i = this->neighbor(i); 00836 if (neighbor_i != NULL) 00837 return neighbor_i; 00838 00839 if (pb) 00840 { 00841 // Since the neighbor is NULL it must be on a boundary. We need 00842 // see if this is a periodic boundary in which case it will have a 00843 // topological neighbor 00844 00845 std::vector<boundary_id_type> boundary_ids = mesh.boundary_info->boundary_ids(this, i); 00846 for (std::vector<boundary_id_type>::iterator j = boundary_ids.begin(); j != boundary_ids.end(); ++j) 00847 if (pb->boundary(*j)) 00848 { 00849 // Since the point locator inside of periodic boundaries 00850 // returns a const pointer we will retrieve the proper 00851 // pointer directly from the mesh object. Also since coarse 00852 // elements do not have more refined neighbors we need to make 00853 // sure that we don't return one of these types of neighbors. 00854 neighbor_i = mesh.elem(pb->neighbor(*j, point_locator, this, i)->id()); 00855 if (level() < neighbor_i->level()) 00856 neighbor_i = neighbor_i->parent(); 00857 return neighbor_i; 00858 } 00859 } 00860 00861 return NULL; 00862 } 00863 00864 00865 bool Elem::has_topological_neighbor (const Elem* elem, 00866 const MeshBase & mesh, 00867 const PointLocatorBase& point_locator, 00868 PeriodicBoundaries * pb) const 00869 { 00870 // First see if this is a normal "interior" neighbor 00871 if (has_neighbor(elem)) 00872 return true; 00873 00874 for (unsigned int n=0; n<this->n_neighbors(); n++) 00875 if (this->topological_neighbor(n, mesh, point_locator, pb)) 00876 return true; 00877 00878 return false; 00879 } 00880 00881 00882 #endif 00883 00884 #ifdef DEBUG 00885 00886 void Elem::libmesh_assert_valid_node_pointers() const 00887 { 00888 libmesh_assert(this->valid_id()); 00889 for (unsigned int n=0; n != this->n_nodes(); ++n) 00890 { 00891 libmesh_assert(this->get_node(n)); 00892 libmesh_assert(this->get_node(n)->valid_id()); 00893 } 00894 } 00895 00896 00897 00898 void Elem::libmesh_assert_valid_neighbors() const 00899 { 00900 for (unsigned int s=0; s<this->n_neighbors(); s++) 00901 { 00902 const Elem *neigh = this->neighbor(s); 00903 00904 // Any element might have a remote neighbor; checking 00905 // to make sure that's not inaccurate is tough. 00906 if (neigh == remote_elem) 00907 continue; 00908 00909 if (neigh) 00910 { 00911 // Only subactive elements have subactive neighbors 00912 libmesh_assert (this->subactive() || !neigh->subactive()); 00913 00914 const Elem *elem = this; 00915 00916 // If we're subactive but our neighbor isn't, its 00917 // return neighbor link will be to our first active 00918 // ancestor OR to our inactive ancestor of the same 00919 // level as neigh, 00920 if (this->subactive() && !neigh->subactive()) 00921 { 00922 for (elem = this; !elem->active(); 00923 elem = elem->parent()) 00924 libmesh_assert(elem); 00925 } 00926 else 00927 { 00928 unsigned int rev = neigh->which_neighbor_am_i(elem); 00929 libmesh_assert_less (rev, neigh->n_neighbors()); 00930 00931 if (this->subactive() && !neigh->subactive()) 00932 { 00933 while (neigh->neighbor(rev) != elem) 00934 { 00935 libmesh_assert(elem->parent()); 00936 elem = elem->parent(); 00937 } 00938 } 00939 else 00940 { 00941 Elem *nn = neigh->neighbor(rev); 00942 libmesh_assert(nn); 00943 00944 for (; elem != nn; elem = elem->parent()) 00945 libmesh_assert(elem); 00946 } 00947 } 00948 } 00949 // If we don't have a neighbor and we're not subactive, our 00950 // ancestors shouldn't have any neighbors in this same 00951 // direction. 00952 else if (!this->subactive()) 00953 { 00954 const Elem *my_parent = this->parent(); 00955 if (my_parent && 00956 // A parent with a different dimension isn't really one of 00957 // our ancestors, it means we're on a boundary mesh and this 00958 // is an interior mesh element for which we're on a side. 00959 // Nothing to test for in that case. 00960 (my_parent->dim() == this->dim())) 00961 libmesh_assert (!my_parent->neighbor(s)); 00962 } 00963 } 00964 } 00965 00966 #endif // DEBUG 00967 00968 00969 00970 void Elem::make_links_to_me_local(unsigned int n) 00971 { 00972 Elem *neigh = this->neighbor(n); 00973 00974 // Don't bother calling this function unless it's necessary 00975 libmesh_assert(neigh); 00976 libmesh_assert(!neigh->is_remote()); 00977 00978 // We never have neighbors more refined than us 00979 libmesh_assert_less_equal (neigh->level(), this->level()); 00980 00981 // We never have subactive neighbors of non subactive elements 00982 libmesh_assert(!neigh->subactive() || this->subactive()); 00983 00984 // If we have a neighbor less refined than us then it must not 00985 // have any more refined active descendants we could have 00986 // pointed to instead. 00987 libmesh_assert(neigh->level() == this->level() || 00988 neigh->active()); 00989 00990 // If neigh is at our level, then its family might have 00991 // remote_elem neighbor links which need to point to us 00992 // instead, but if not, then we're done. 00993 if (neigh->level() != this->level()) 00994 return; 00995 00996 // If neigh is subactive then we're not updating its neighbor links 00997 // FIXME - this needs to change when we start using subactive 00998 // elements for more than just the two-phase 00999 // restriction/prolongation projections. 01000 if (neigh->subactive()) 01001 return; 01002 01003 // What side of neigh are we on? We can't use the usual Elem 01004 // method because we're in the middle of restoring topology 01005 const AutoPtr<Elem> my_side = this->side(n); 01006 unsigned int nn = 0; 01007 for (; nn != neigh->n_sides(); ++nn) 01008 { 01009 const AutoPtr<Elem> neigh_side = neigh->side(nn); 01010 if (*my_side == *neigh_side) 01011 break; 01012 } 01013 01014 // we had better be on *some* side of neigh 01015 libmesh_assert_less (nn, neigh->n_sides()); 01016 01017 // Find any elements that ought to point to elem 01018 std::vector<const Elem*> neigh_family; 01019 #ifdef LIBMESH_ENABLE_AMR 01020 if (this->active()) 01021 neigh->family_tree_by_side(neigh_family, nn); 01022 else 01023 #endif 01024 neigh_family.push_back(neigh); 01025 01026 // And point them to elem 01027 for (unsigned int i = 0; i != neigh_family.size(); ++i) 01028 { 01029 Elem* neigh_family_member = const_cast<Elem*>(neigh_family[i]); 01030 01031 // Ideally, the neighbor link ought to either be correct 01032 // already or ought to be to remote_elem. 01033 // 01034 // However, if we're redistributing a newly created elem, 01035 // after an AMR step but before find_neighbors has fixed up 01036 // neighbor links, we might have an out of date neighbor 01037 // link to elem's parent instead. 01038 #ifdef LIBMESH_ENABLE_AMR 01039 libmesh_assert((neigh_family_member->neighbor(nn) == this) || 01040 (neigh_family_member->neighbor(nn) == remote_elem) 01041 || ((this->refinement_flag() == JUST_REFINED) && 01042 (this->parent() != NULL) && 01043 (neigh_family_member->neighbor(nn) == this->parent()))); 01044 #else 01045 libmesh_assert((neigh_family_member->neighbor(nn) == this) || 01046 (neigh_family_member->neighbor(nn) == remote_elem)); 01047 #endif 01048 01049 neigh_family_member->set_neighbor(nn, this); 01050 } 01051 } 01052 01053 01054 void Elem::make_links_to_me_remote() 01055 { 01056 libmesh_assert_not_equal_to (this, remote_elem); 01057 01058 // We need to have handled any children first 01059 #if defined(LIBMESH_ENABLE_AMR) && defined(DEBUG) 01060 if (this->has_children()) 01061 for (unsigned int c = 0; c != this->n_children(); ++c) 01062 { 01063 Elem *current_child = this->child(c); 01064 libmesh_assert_equal_to (current_child, remote_elem); 01065 } 01066 #endif 01067 01068 // Remotify any neighbor links to non-subactive elements 01069 if (!this->subactive()) 01070 { 01071 for (unsigned int s = 0; s != this->n_sides(); ++s) 01072 { 01073 Elem *neigh = this->neighbor(s); 01074 if (neigh && neigh != remote_elem && !neigh->subactive()) 01075 { 01076 // My neighbor should never be more refined than me; my real 01077 // neighbor would have been its parent in that case. 01078 libmesh_assert_greater_equal (this->level(), neigh->level()); 01079 01080 if (this->level() == neigh->level() && 01081 neigh->has_neighbor(this)) 01082 { 01083 #ifdef LIBMESH_ENABLE_AMR 01084 // My neighbor may have descendants which also consider me a 01085 // neighbor 01086 std::vector<const Elem*> family; 01087 neigh->family_tree_by_neighbor (family, this); 01088 01089 // FIXME - There's a lot of ugly const_casts here; we 01090 // may want to make remote_elem non-const and create 01091 // non-const versions of the family_tree methods 01092 for (unsigned int i=0; i != family.size(); ++i) 01093 { 01094 Elem *n = const_cast<Elem*>(family[i]); 01095 libmesh_assert (n); 01096 if (n == remote_elem) 01097 continue; 01098 unsigned int my_s = n->which_neighbor_am_i(this); 01099 libmesh_assert_less (my_s, n->n_neighbors()); 01100 libmesh_assert_equal_to (n->neighbor(my_s), this); 01101 n->set_neighbor(my_s, const_cast<RemoteElem*>(remote_elem)); 01102 } 01103 #else 01104 unsigned int my_s = neigh->which_neighbor_am_i(this); 01105 libmesh_assert_less (my_s, neigh->n_neighbors()); 01106 libmesh_assert_equal_to (neigh->neighbor(my_s), this); 01107 neigh->set_neighbor(my_s, const_cast<RemoteElem*>(remote_elem)); 01108 #endif 01109 } 01110 #ifdef LIBMESH_ENABLE_AMR 01111 // Even if my neighbor doesn't link back to me, it might 01112 // have subactive descendants which do 01113 else if (neigh->has_children()) 01114 { 01115 // If my neighbor at the same level doesn't have me as a 01116 // neighbor, I must be subactive 01117 libmesh_assert(this->level() > neigh->level() || 01118 this->subactive()); 01119 01120 // My neighbor must have some ancestor of mine as a 01121 // neighbor 01122 Elem *my_ancestor = this->parent(); 01123 libmesh_assert(my_ancestor); 01124 while (!neigh->has_neighbor(my_ancestor)) 01125 { 01126 my_ancestor = my_ancestor->parent(); 01127 libmesh_assert(my_ancestor); 01128 } 01129 01130 // My neighbor may have descendants which consider me a 01131 // neighbor 01132 std::vector<const Elem*> family; 01133 neigh->family_tree_by_subneighbor (family, my_ancestor, this); 01134 01135 // FIXME - There's a lot of ugly const_casts here; we 01136 // may want to make remote_elem non-const and create 01137 // non-const versions of the family_tree methods 01138 for (unsigned int i=0; i != family.size(); ++i) 01139 { 01140 Elem *n = const_cast<Elem*>(family[i]); 01141 libmesh_assert (n); 01142 if (n == remote_elem) 01143 continue; 01144 unsigned int my_s = n->which_neighbor_am_i(this); 01145 libmesh_assert_less (my_s, n->n_neighbors()); 01146 libmesh_assert_equal_to (n->neighbor(my_s), this); 01147 n->set_neighbor(my_s, const_cast<RemoteElem*>(remote_elem)); 01148 } 01149 } 01150 #endif 01151 } 01152 } 01153 } 01154 01155 #ifdef LIBMESH_ENABLE_AMR 01156 // Remotify parent's child link 01157 Elem *my_parent = this->parent(); 01158 if (my_parent && 01159 // As long as it's not already remote 01160 my_parent != remote_elem && 01161 // And it's a real parent, not an interior parent 01162 this->dim() == my_parent->dim()) 01163 { 01164 unsigned int me = my_parent->which_child_am_i(this); 01165 libmesh_assert_equal_to (my_parent->child(me), this); 01166 my_parent->set_child(me, const_cast<RemoteElem*>(remote_elem)); 01167 } 01168 #endif 01169 } 01170 01171 01172 01173 void Elem::write_connectivity (std::ostream& out_stream, 01174 const IOPackage iop) const 01175 { 01176 libmesh_assert (out_stream.good()); 01177 libmesh_assert(_nodes); 01178 libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE); 01179 01180 switch (iop) 01181 { 01182 case TECPLOT: 01183 { 01184 // This connectivity vector will be used repeatedly instead 01185 // of being reconstructed inside the loop. 01186 std::vector<dof_id_type> conn; 01187 for (unsigned int sc=0; sc <this->n_sub_elem(); sc++) 01188 { 01189 this->connectivity(sc, TECPLOT, conn); 01190 01191 std::copy(conn.begin(), 01192 conn.end(), 01193 std::ostream_iterator<dof_id_type>(out_stream, " ")); 01194 01195 out_stream << '\n'; 01196 } 01197 return; 01198 } 01199 01200 case UCD: 01201 { 01202 for (unsigned int i=0; i<this->n_nodes(); i++) 01203 out_stream << this->node(i)+1 << "\t"; 01204 01205 out_stream << '\n'; 01206 return; 01207 } 01208 01209 default: 01210 libmesh_error(); 01211 } 01212 01213 libmesh_error(); 01214 } 01215 01216 01217 // void Elem::write_tecplot_connectivity(std::ostream& out) const 01218 // { 01219 // libmesh_assert (!out.bad()); 01220 // libmesh_assert(_nodes); 01221 01222 // // This connectivity vector will be used repeatedly instead 01223 // // of being reconstructed inside the loop. 01224 // std::vector<dof_id_type> conn; 01225 // for (unsigned int sc=0; sc <this->n_sub_elem(); sc++) 01226 // { 01227 // this->connectivity(sc, TECPLOT, conn); 01228 01229 // std::copy(conn.begin(), 01230 // conn.end(), 01231 // std::ostream_iterator<dof_id_type>(out, " ")); 01232 01233 // out << std::endl; 01234 // } 01235 // } 01236 01237 01238 01239 // void Elem::write_ucd_connectivity(std::ostream &out) const 01240 // { 01241 // libmesh_assert (out); 01242 // libmesh_assert(_nodes); 01243 01244 // for (unsigned int i=0; i<this->n_nodes(); i++) 01245 // out << this->node(i)+1 << "\t"; 01246 01247 // out << std::endl; 01248 // } 01249 01250 01251 01252 Real Elem::quality (const ElemQuality q) const 01253 { 01254 switch (q) 01255 { 01259 default: 01260 { 01261 libmesh_here(); 01262 01263 libMesh::err << "ERROR: unknown quality metric: " 01264 << q 01265 << std::endl 01266 << "Cowardly returning 1." 01267 << std::endl; 01268 01269 return 1.; 01270 } 01271 } 01272 01273 01274 // Will never get here... 01275 libmesh_error(); 01276 return 0.; 01277 } 01278 01279 01280 01281 bool Elem::ancestor() const 01282 { 01283 #ifdef LIBMESH_ENABLE_AMR 01284 01285 if (this->active()) 01286 return false; 01287 01288 if (!this->has_children()) 01289 return false; 01290 if (this->child(0)->active()) 01291 return true; 01292 01293 return this->child(0)->ancestor(); 01294 #else 01295 return false; 01296 #endif 01297 } 01298 01299 01300 01301 #ifdef LIBMESH_ENABLE_AMR 01302 01303 void Elem::add_child (Elem* elem) 01304 { 01305 if(_children == NULL) 01306 { 01307 _children = new Elem*[this->n_children()]; 01308 01309 for (unsigned int c=0; c<this->n_children(); c++) 01310 this->set_child(c, NULL); 01311 } 01312 01313 for (unsigned int c=0; c<this->n_children(); c++) 01314 { 01315 if(this->_children[c] == NULL || this->_children[c] == remote_elem) 01316 { 01317 libmesh_assert_equal_to (this, elem->parent()); 01318 this->set_child(c, elem); 01319 return; 01320 } 01321 } 01322 01323 libMesh::err << "Error: Tried to add a child to an element with full children array" 01324 << std::endl; 01325 libmesh_error(); 01326 } 01327 01328 01329 01330 void Elem::add_child (Elem* elem, unsigned int c) 01331 { 01332 if(!this->has_children()) 01333 { 01334 _children = new Elem*[this->n_children()]; 01335 01336 for (unsigned int i=0; i<this->n_children(); i++) 01337 this->set_child(i, NULL); 01338 } 01339 01340 libmesh_assert (this->_children[c] == NULL || this->child(c) == remote_elem); 01341 libmesh_assert (elem == remote_elem || this == elem->parent()); 01342 01343 this->set_child(c, elem); 01344 } 01345 01346 01347 01348 void Elem::replace_child (Elem* elem, unsigned int c) 01349 { 01350 libmesh_assert(this->has_children()); 01351 01352 libmesh_assert(this->child(c)); 01353 01354 this->set_child(c, elem); 01355 } 01356 01357 01358 01359 bool Elem::is_child_on_edge(const unsigned int libmesh_dbg_var(c), 01360 const unsigned int e) const 01361 { 01362 libmesh_assert_less (c, this->n_children()); 01363 libmesh_assert_less (e, this->n_edges()); 01364 01365 AutoPtr<Elem> my_edge = this->build_edge(e); 01366 AutoPtr<Elem> child_edge = this->build_edge(e); 01367 01368 // We're assuming that an overlapping child edge has the same 01369 // number and orientation as its parent 01370 return (child_edge->node(0) == my_edge->node(0) || 01371 child_edge->node(1) == my_edge->node(1)); 01372 } 01373 01374 01375 void Elem::family_tree (std::vector<const Elem*>& family, 01376 const bool reset) const 01377 { 01378 // The "family tree" doesn't include subactive elements 01379 libmesh_assert(!this->subactive()); 01380 01381 // Clear the vector if the flag reset tells us to. 01382 if (reset) 01383 family.clear(); 01384 01385 // Add this element to the family tree. 01386 family.push_back(this); 01387 01388 // Recurse into the elements children, if it has them. 01389 // Do not clear the vector any more. 01390 if (!this->active()) 01391 for (unsigned int c=0; c<this->n_children(); c++) 01392 if (!this->child(c)->is_remote()) 01393 this->child(c)->family_tree (family, false); 01394 } 01395 01396 01397 01398 void Elem::total_family_tree (std::vector<const Elem*>& family, 01399 const bool reset) const 01400 { 01401 // Clear the vector if the flag reset tells us to. 01402 if (reset) 01403 family.clear(); 01404 01405 // Add this element to the family tree. 01406 family.push_back(this); 01407 01408 // Recurse into the elements children, if it has them. 01409 // Do not clear the vector any more. 01410 if (this->has_children()) 01411 for (unsigned int c=0; c<this->n_children(); c++) 01412 if (!this->child(c)->is_remote()) 01413 this->child(c)->total_family_tree (family, false); 01414 } 01415 01416 01417 01418 void Elem::active_family_tree (std::vector<const Elem*>& active_family, 01419 const bool reset) const 01420 { 01421 // The "family tree" doesn't include subactive elements 01422 libmesh_assert(!this->subactive()); 01423 01424 // Clear the vector if the flag reset tells us to. 01425 if (reset) 01426 active_family.clear(); 01427 01428 // Add this element to the family tree if it is active 01429 if (this->active()) 01430 active_family.push_back(this); 01431 01432 // Otherwise recurse into the element's children. 01433 // Do not clear the vector any more. 01434 else 01435 for (unsigned int c=0; c<this->n_children(); c++) 01436 if (!this->child(c)->is_remote()) 01437 this->child(c)->active_family_tree (active_family, false); 01438 } 01439 01440 01441 01442 void Elem::family_tree_by_side (std::vector<const Elem*>& family, 01443 const unsigned int s, 01444 const bool reset) const 01445 { 01446 // The "family tree" doesn't include subactive elements 01447 libmesh_assert(!this->subactive()); 01448 01449 // Clear the vector if the flag reset tells us to. 01450 if (reset) 01451 family.clear(); 01452 01453 libmesh_assert_less (s, this->n_sides()); 01454 01455 // Add this element to the family tree. 01456 family.push_back(this); 01457 01458 // Recurse into the elements children, if it has them. 01459 // Do not clear the vector any more. 01460 if (!this->active()) 01461 for (unsigned int c=0; c<this->n_children(); c++) 01462 if (!this->child(c)->is_remote() && this->is_child_on_side(c, s)) 01463 this->child(c)->family_tree_by_side (family, s, false); 01464 } 01465 01466 01467 01468 void Elem::active_family_tree_by_side (std::vector<const Elem*>& family, 01469 const unsigned int s, 01470 const bool reset) const 01471 { 01472 // The "family tree" doesn't include subactive elements 01473 libmesh_assert(!this->subactive()); 01474 01475 // Clear the vector if the flag reset tells us to. 01476 if (reset) 01477 family.clear(); 01478 01479 libmesh_assert_less (s, this->n_sides()); 01480 01481 // Add an active element to the family tree. 01482 if (this->active()) 01483 family.push_back(this); 01484 01485 // Or recurse into an ancestor element's children. 01486 // Do not clear the vector any more. 01487 else 01488 for (unsigned int c=0; c<this->n_children(); c++) 01489 if (!this->child(c)->is_remote() && this->is_child_on_side(c, s)) 01490 this->child(c)->active_family_tree_by_side (family, s, false); 01491 } 01492 01493 01494 01495 void Elem::family_tree_by_neighbor (std::vector<const Elem*>& family, 01496 const Elem* neighbor_in, 01497 const bool reset) const 01498 { 01499 // The "family tree" doesn't include subactive elements 01500 libmesh_assert(!this->subactive()); 01501 01502 // Clear the vector if the flag reset tells us to. 01503 if (reset) 01504 family.clear(); 01505 01506 // This only makes sense if we're already a neighbor 01507 libmesh_assert (this->has_neighbor(neighbor_in)); 01508 01509 // Add this element to the family tree. 01510 family.push_back(this); 01511 01512 // Recurse into the elements children, if it's not active. 01513 // Do not clear the vector any more. 01514 if (!this->active()) 01515 for (unsigned int c=0; c<this->n_children(); c++) 01516 { 01517 Elem *current_child = this->child(c); 01518 if (current_child != remote_elem && current_child->has_neighbor(neighbor_in)) 01519 current_child->family_tree_by_neighbor (family, neighbor_in, false); 01520 } 01521 } 01522 01523 01524 01525 void Elem::family_tree_by_subneighbor (std::vector<const Elem*>& family, 01526 const Elem* neighbor_in, 01527 const Elem* subneighbor, 01528 const bool reset) const 01529 { 01530 // The "family tree" doesn't include subactive elements 01531 libmesh_assert(!this->subactive()); 01532 01533 // Clear the vector if the flag reset tells us to. 01534 if (reset) 01535 family.clear(); 01536 01537 // To simplifly this function we need an existing neighbor 01538 libmesh_assert (neighbor_in); 01539 libmesh_assert_not_equal_to (neighbor_in, remote_elem); 01540 libmesh_assert (this->has_neighbor(neighbor_in)); 01541 01542 // This only makes sense if subneighbor descends from neighbor 01543 libmesh_assert (subneighbor); 01544 libmesh_assert_not_equal_to (subneighbor, remote_elem); 01545 libmesh_assert (neighbor_in->is_ancestor_of(subneighbor)); 01546 01547 // Add this element to the family tree if applicable. 01548 if (neighbor_in == subneighbor) 01549 family.push_back(this); 01550 01551 // Recurse into the elements children, if it's not active. 01552 // Do not clear the vector any more. 01553 if (!this->active()) 01554 for (unsigned int c=0; c != this->n_children(); ++c) 01555 { 01556 Elem *current_child = this->child(c); 01557 if (current_child != remote_elem) 01558 for (unsigned int s=0; s != current_child->n_sides(); ++s) 01559 { 01560 Elem *child_neigh = current_child->neighbor(s); 01561 if (child_neigh && 01562 (child_neigh == neighbor_in || 01563 (child_neigh->parent() == neighbor_in && 01564 child_neigh->is_ancestor_of(subneighbor)))) 01565 current_child->family_tree_by_subneighbor (family, child_neigh, 01566 subneighbor, false); 01567 } 01568 } 01569 } 01570 01571 01572 01573 void Elem::active_family_tree_by_neighbor (std::vector<const Elem*>& family, 01574 const Elem* neighbor_in, 01575 const bool reset) const 01576 { 01577 // The "family tree" doesn't include subactive elements 01578 libmesh_assert(!this->subactive()); 01579 01580 // Clear the vector if the flag reset tells us to. 01581 if (reset) 01582 family.clear(); 01583 01584 // This only makes sense if we're already a neighbor 01585 if (this->level() >= neighbor_in->level()) 01586 libmesh_assert (this->has_neighbor(neighbor_in)); 01587 01588 // Add an active element to the family tree. 01589 if (this->active()) 01590 family.push_back(this); 01591 01592 // Or recurse into an ancestor element's children. 01593 // Do not clear the vector any more. 01594 else if (!this->active()) 01595 for (unsigned int c=0; c<this->n_children(); c++) 01596 { 01597 Elem *current_child = this->child(c); 01598 if (current_child != remote_elem && current_child->has_neighbor(neighbor_in)) 01599 current_child->active_family_tree_by_neighbor (family, neighbor_in, false); 01600 } 01601 } 01602 01603 01604 01605 unsigned int Elem::min_p_level_by_neighbor(const Elem* neighbor_in, 01606 unsigned int current_min) const 01607 { 01608 libmesh_assert(!this->subactive()); 01609 libmesh_assert(neighbor_in->active()); 01610 01611 // If we're an active element this is simple 01612 if (this->active()) 01613 return std::min(current_min, this->p_level()); 01614 01615 libmesh_assert(has_neighbor(neighbor_in)); 01616 01617 // The p_level() of an ancestor element is already the minimum 01618 // p_level() of its children - so if that's high enough, we don't 01619 // need to examine any children. 01620 if (current_min <= this->p_level()) 01621 return current_min; 01622 01623 unsigned int min_p_level = current_min; 01624 01625 for (unsigned int c=0; c<this->n_children(); c++) 01626 { 01627 const Elem* const current_child = this->child(c); 01628 if (current_child != remote_elem && current_child->has_neighbor(neighbor_in)) 01629 min_p_level = 01630 current_child->min_p_level_by_neighbor(neighbor_in, 01631 min_p_level); 01632 } 01633 01634 return min_p_level; 01635 } 01636 01637 01638 unsigned int Elem::min_new_p_level_by_neighbor(const Elem* neighbor_in, 01639 unsigned int current_min) const 01640 { 01641 libmesh_assert(!this->subactive()); 01642 libmesh_assert(neighbor_in->active()); 01643 01644 // If we're an active element this is simple 01645 if (this->active()) 01646 { 01647 unsigned int new_p_level = this->p_level(); 01648 if (this->p_refinement_flag() == Elem::REFINE) 01649 new_p_level += 1; 01650 if (this->p_refinement_flag() == Elem::COARSEN) 01651 { 01652 libmesh_assert_greater (new_p_level, 0); 01653 new_p_level -= 1; 01654 } 01655 return std::min(current_min, new_p_level); 01656 } 01657 01658 libmesh_assert(has_neighbor(neighbor_in)); 01659 01660 unsigned int min_p_level = current_min; 01661 01662 for (unsigned int c=0; c<this->n_children(); c++) 01663 { 01664 const Elem* const current_child = this->child(c); 01665 if (current_child && current_child != remote_elem) 01666 if (current_child->has_neighbor(neighbor_in)) 01667 min_p_level = 01668 current_child->min_new_p_level_by_neighbor(neighbor_in, 01669 min_p_level); 01670 } 01671 01672 return min_p_level; 01673 } 01674 01675 #endif // #ifdef LIBMESH_ENABLE_AMR 01676 01677 01678 01679 01680 bool Elem::contains_point (const Point& p, Real tol) const 01681 { 01682 // We currently allow the user to enlarge the bounding box by 01683 // providing a tol > TOLERANCE (so this routine is identical to 01684 // Elem::close_to_point()), but print a warning so that the 01685 // user can eventually switch his code over to calling close_to_point() 01686 // instead, which is intended to be used for this purpose. 01687 if ( tol > TOLERANCE ) 01688 { 01689 libmesh_do_once(libMesh::err 01690 << "WARNING: Resizing bounding box to match user-specified tolerance!\n" 01691 << "In the future, calls to Elem::contains_point() with tol > TOLERANCE\n" 01692 << "will be more optimized, but should not be used\n" 01693 << "to search for points 'close to' elements!\n" 01694 << "Instead, use Elem::close_to_point() for this purpose.\n" 01695 << std::endl;); 01696 return this->point_test(p, tol, tol); 01697 } 01698 else 01699 return this->point_test(p, TOLERANCE, tol); 01700 } 01701 01702 01703 01704 01705 bool Elem::close_to_point (const Point& p, Real tol) const 01706 { 01707 // This test uses the user's passed-in tolerance for the 01708 // bounding box test as well, thereby allowing the routine to 01709 // find points which are not only "in" the element, but also 01710 // "nearby" to within some tolerance. 01711 return this->point_test(p, tol, tol); 01712 } 01713 01714 01715 01716 01717 bool Elem::point_test(const Point& p, Real box_tol, Real map_tol) const 01718 { 01719 libmesh_assert_greater (box_tol, 0.); 01720 libmesh_assert_greater (map_tol, 0.); 01721 01722 // This is a great optimization on first order elements, but it 01723 // could return false negatives on higher orders 01724 if (this->default_order() == FIRST) 01725 { 01726 // Check to make sure the element *could* contain this point, so we 01727 // can avoid an expensive inverse_map call if it doesn't. 01728 bool 01729 #if LIBMESH_DIM > 2 01730 point_above_min_z = false, 01731 point_below_max_z = false, 01732 #endif 01733 #if LIBMESH_DIM > 1 01734 point_above_min_y = false, 01735 point_below_max_y = false, 01736 #endif 01737 point_above_min_x = false, 01738 point_below_max_x = false; 01739 01740 // For relative bounding box checks in physical space 01741 const Real my_hmax = this->hmax(); 01742 01743 for (unsigned int n=0; n != this->n_nodes(); ++n) 01744 { 01745 Point pe = this->point(n); 01746 point_above_min_x = point_above_min_x || (pe(0) - my_hmax*box_tol <= p(0)); 01747 point_below_max_x = point_below_max_x || (pe(0) + my_hmax*box_tol >= p(0)); 01748 #if LIBMESH_DIM > 1 01749 point_above_min_y = point_above_min_y || (pe(1) - my_hmax*box_tol <= p(1)); 01750 point_below_max_y = point_below_max_y || (pe(1) + my_hmax*box_tol >= p(1)); 01751 #endif 01752 #if LIBMESH_DIM > 2 01753 point_above_min_z = point_above_min_z || (pe(2) - my_hmax*box_tol <= p(2)); 01754 point_below_max_z = point_below_max_z || (pe(2) + my_hmax*box_tol >= p(2)); 01755 #endif 01756 } 01757 01758 if ( 01759 #if LIBMESH_DIM > 2 01760 !point_above_min_z || 01761 !point_below_max_z || 01762 #endif 01763 #if LIBMESH_DIM > 1 01764 !point_above_min_y || 01765 !point_below_max_y || 01766 #endif 01767 !point_above_min_x || 01768 !point_below_max_x) 01769 return false; 01770 } 01771 01772 // Declare a basic FEType. Will be a Lagrange 01773 // element by default. 01774 FEType fe_type(this->default_order()); 01775 01776 // To be on the safe side, we converge the inverse_map() iteration 01777 // to a slightly tighter tolerance than that requested by the 01778 // user... 01779 const Point mapped_point = FEInterface::inverse_map(this->dim(), 01780 fe_type, 01781 this, 01782 p, 01783 0.1*map_tol, // <- this is |dx| tolerance, the Newton residual should be ~ |dx|^2 01784 /*secure=*/ false); 01785 01786 // Check that the refspace point maps back to p! This is only necessary 01787 // for 1D and 2D elements, 3D elements always live in 3D. 01788 // 01789 // TODO: The contains_point() function could most likely be implemented 01790 // more efficiently in the element sub-classes themselves, at least for 01791 // the linear element types. 01792 if (this->dim() < 3) 01793 { 01794 Point xyz = FEInterface::map(this->dim(), 01795 fe_type, 01796 this, 01797 mapped_point); 01798 01799 // Compute the distance between the original point and the re-mapped point. 01800 // They should be in the same place. 01801 Real dist = (xyz - p).size(); 01802 01803 01804 // If dist is larger than some fraction of the tolerance, then return false. 01805 // This can happen when e.g. a 2D element is living in 3D, and 01806 // FEInterface::inverse_map() maps p onto the projection of the element, 01807 // effectively "tricking" FEInterface::on_reference_element(). 01808 if (dist > this->hmax() * map_tol) 01809 return false; 01810 } 01811 01812 01813 01814 return FEInterface::on_reference_element(mapped_point, this->type(), map_tol); 01815 } 01816 01817 01818 01819 01820 void Elem::print_info (std::ostream& os) const 01821 { 01822 os << this->get_info() 01823 << std::endl; 01824 } 01825 01826 01827 01828 std::string Elem::get_info () const 01829 { 01830 std::ostringstream oss; 01831 01832 oss << " Elem Information" << '\n' 01833 << " id()="; 01834 01835 if (this->valid_id()) 01836 oss << this->id(); 01837 else 01838 oss << "invalid"; 01839 01840 oss << ", processor_id()=" << this->processor_id() << '\n'; 01841 01842 oss << " type()=" << Utility::enum_to_string(this->type()) << '\n' 01843 << " dim()=" << this->dim() << '\n' 01844 << " n_nodes()=" << this->n_nodes() << '\n'; 01845 01846 for (unsigned int n=0; n != this->n_nodes(); ++n) 01847 oss << " " << n << *this->get_node(n); 01848 01849 oss << " n_sides()=" << this->n_sides() << '\n'; 01850 01851 for (unsigned int s=0; s != this->n_sides(); ++s) 01852 { 01853 oss << " neighbor(" << s << ")="; 01854 if (this->neighbor(s)) 01855 oss << this->neighbor(s)->id() << '\n'; 01856 else 01857 oss << "NULL\n"; 01858 } 01859 01860 oss << " hmin()=" << this->hmin() 01861 << ", hmax()=" << this->hmax() << '\n' 01862 << " volume()=" << this->volume() << '\n' 01863 << " active()=" << this->active() 01864 << ", ancestor()=" << this->ancestor() 01865 << ", subactive()=" << this->subactive() 01866 << ", has_children()=" << this->has_children() << '\n' 01867 << " parent()="; 01868 if (this->parent()) 01869 oss << this->parent()->id() << '\n'; 01870 else 01871 oss << "NULL\n"; 01872 oss << " level()=" << this->level() 01873 << ", p_level()=" << this->p_level() << '\n' 01874 #ifdef LIBMESH_ENABLE_AMR 01875 << " refinement_flag()=" << Utility::enum_to_string(this->refinement_flag()) << '\n' 01876 << " p_refinement_flag()=" << Utility::enum_to_string(this->p_refinement_flag()) << '\n' 01877 #endif 01878 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 01879 << " infinite()=" << this->infinite() << '\n'; 01880 if (this->infinite()) 01881 oss << " origin()=" << this->origin() << '\n' 01882 #endif 01883 ; 01884 01885 oss << " DoFs="; 01886 for (unsigned int s=0; s != this->n_systems(); ++s) 01887 for (unsigned int v=0; v != this->n_vars(s); ++v) 01888 for (unsigned int c=0; c != this->n_comp(s,v); ++c) 01889 oss << '(' << s << '/' << v << '/' << this->dof_number(s,v,c) << ") "; 01890 01891 01892 return oss.str(); 01893 } 01894 01895 01896 01897 void Elem::nullify_neighbors () 01898 { 01899 // Tell any of my neighbors about my death... 01900 // Looks strange, huh? 01901 for (unsigned int n=0; n<this->n_neighbors(); n++) 01902 { 01903 Elem* current_neighbor = this->neighbor(n); 01904 if (current_neighbor && current_neighbor != remote_elem) 01905 { 01906 // Note: it is possible that I see the neighbor 01907 // (which is coarser than me) 01908 // but they don't see me, so avoid that case. 01909 if (current_neighbor->level() == this->level()) 01910 { 01911 const unsigned int w_n_a_i = current_neighbor->which_neighbor_am_i(this); 01912 libmesh_assert_less (w_n_a_i, current_neighbor->n_neighbors()); 01913 current_neighbor->set_neighbor(w_n_a_i, NULL); 01914 this->set_neighbor(n, NULL); 01915 } 01916 } 01917 } 01918 } 01919 01920 01921 01922 01923 unsigned int Elem::n_second_order_adjacent_vertices (const unsigned int) const 01924 { 01925 // for linear elements, always return 0 01926 return 0; 01927 } 01928 01929 01930 01931 unsigned short int Elem::second_order_adjacent_vertex (const unsigned int, 01932 const unsigned int) const 01933 { 01934 // for linear elements, always return 0 01935 return 0; 01936 } 01937 01938 01939 01940 std::pair<unsigned short int, unsigned short int> 01941 Elem::second_order_child_vertex (const unsigned int) const 01942 { 01943 // for linear elements, always return 0 01944 return std::pair<unsigned short int, unsigned short int>(0,0); 01945 } 01946 01947 01948 01949 ElemType Elem::first_order_equivalent_type (const ElemType et) 01950 { 01951 switch (et) 01952 { 01953 case EDGE2: 01954 case EDGE3: 01955 case EDGE4: 01956 return EDGE2; 01957 case TRI3: 01958 case TRI6: 01959 return TRI3; 01960 case QUAD4: 01961 case QUAD8: 01962 case QUAD9: 01963 return QUAD4; 01964 case TET4: 01965 case TET10: 01966 return TET4; 01967 case HEX8: 01968 case HEX27: 01969 case HEX20: 01970 return HEX8; 01971 case PRISM6: 01972 case PRISM15: 01973 case PRISM18: 01974 return PRISM6; 01975 01976 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 01977 01978 case INFQUAD4: 01979 case INFQUAD6: 01980 return INFQUAD4; 01981 case INFHEX8: 01982 case INFHEX16: 01983 case INFHEX18: 01984 return INFHEX8; 01985 case INFPRISM6: 01986 case INFPRISM12: 01987 return INFPRISM6; 01988 01989 #endif 01990 01991 default: 01992 // unknown element 01993 return INVALID_ELEM; 01994 } 01995 } 01996 01997 01998 01999 ElemType Elem::second_order_equivalent_type (const ElemType et, 02000 const bool full_ordered) 02001 { 02002 /* for second-order elements, always return \p INVALID_ELEM 02003 * since second-order elements should not be converted 02004 * into something else. Only linear elements should 02005 * return something sensible here 02006 */ 02007 switch (et) 02008 { 02009 case EDGE2: 02010 { 02011 // full_ordered not relevant 02012 return EDGE3; 02013 } 02014 02015 case TRI3: 02016 { 02017 // full_ordered not relevant 02018 return TRI6; 02019 } 02020 02021 case QUAD4: 02022 { 02023 if (full_ordered) 02024 return QUAD9; 02025 else 02026 return QUAD8; 02027 } 02028 02029 case TET4: 02030 { 02031 // full_ordered not relevant 02032 return TET10; 02033 } 02034 02035 case HEX8: 02036 { 02037 // see below how this correlates with INFHEX8 02038 if (full_ordered) 02039 return HEX27; 02040 else 02041 return HEX20; 02042 } 02043 02044 case PRISM6: 02045 { 02046 if (full_ordered) 02047 return PRISM18; 02048 else 02049 return PRISM15; 02050 } 02051 02052 case PYRAMID5: 02053 { 02054 // libmesh_error(); 02055 return INVALID_ELEM; 02056 } 02057 02058 02059 02060 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 02061 02062 // infinite elements 02063 case INFEDGE2: 02064 { 02065 // libmesh_error(); 02066 return INVALID_ELEM; 02067 } 02068 02069 case INFQUAD4: 02070 { 02071 // full_ordered not relevant 02072 return INFQUAD6; 02073 } 02074 02075 case INFHEX8: 02076 { 02077 /* 02078 * Note that this matches with \p Hex8: 02079 * For full-ordered, \p InfHex18 and \p Hex27 02080 * belong together, and for not full-ordered, 02081 * \p InfHex16 and \p Hex20 belong together. 02082 */ 02083 if (full_ordered) 02084 return INFHEX18; 02085 else 02086 return INFHEX16; 02087 } 02088 02089 case INFPRISM6: 02090 { 02091 // full_ordered not relevant 02092 return INFPRISM12; 02093 } 02094 02095 #endif 02096 02097 02098 default: 02099 { 02100 // second-order element 02101 return INVALID_ELEM; 02102 } 02103 } 02104 } 02105 02106 02107 02108 Elem::side_iterator Elem::boundary_sides_begin() 02109 { 02110 Predicates::BoundarySide<SideIter> bsp; 02111 return side_iterator(this->_first_side(), this->_last_side(), bsp); 02112 } 02113 02114 02115 02116 02117 Elem::side_iterator Elem::boundary_sides_end() 02118 { 02119 Predicates::BoundarySide<SideIter> bsp; 02120 return side_iterator(this->_last_side(), this->_last_side(), bsp); 02121 } 02122 02123 02124 02125 02126 Real Elem::volume () const 02127 { 02128 // The default implementation builds a finite element of the correct 02129 // order and sums up the JxW contributions. This can be expensive, 02130 // so the various element types can overload this method and compute 02131 // the volume more efficiently. 02132 FEType fe_type (this->default_order() , LAGRANGE); 02133 02134 AutoPtr<FEBase> fe (FEBase::build(this->dim(), 02135 fe_type)); 02136 02137 const std::vector<Real>& JxW = fe->get_JxW(); 02138 02139 // The default quadrature rule should integrate the mass matrix, 02140 // thus it should be plenty to compute the area 02141 QGauss qrule (this->dim(), fe_type.default_quadrature_order()); 02142 02143 fe->attach_quadrature_rule(&qrule); 02144 02145 fe->reinit(this); 02146 02147 Real vol=0.; 02148 for (unsigned int qp=0; qp<qrule.n_points(); ++qp) 02149 vol += JxW[qp]; 02150 02151 return vol; 02152 02153 } 02154 02155 02156 // ------------------------------------------------------------ 02157 // Elem::PackedElem static data 02158 const unsigned int Elem::PackedElem::header_size = 10; 02159 02160 02161 // Elem::PackedElem member funcions 02162 void Elem::PackedElem::pack (std::vector<int> &conn, const Elem* elem) 02163 { 02164 libmesh_assert(elem); 02165 02166 // we can do at least this good. note that hopefully in general 02167 // the user will already have reserved the full space, which will render 02168 // this redundant 02169 conn.reserve (conn.size() + elem->packed_size()); 02170 02171 #ifdef LIBMESH_ENABLE_AMR 02172 conn.push_back (static_cast<int>(elem->level())); 02173 conn.push_back (static_cast<int>(elem->p_level())); 02174 conn.push_back (static_cast<int>(elem->refinement_flag())); 02175 conn.push_back (static_cast<int>(elem->p_refinement_flag())); 02176 #else 02177 conn.push_back (0); 02178 conn.push_back (0); 02179 conn.push_back (0); 02180 conn.push_back (0); 02181 #endif 02182 conn.push_back (static_cast<int>(elem->type())); 02183 conn.push_back (elem->processor_id()); 02184 conn.push_back (elem->subdomain_id()); 02185 conn.push_back (elem->id()); 02186 02187 #ifdef LIBMESH_ENABLE_AMR 02188 // use parent_ID of -1 to indicate a level 0 element 02189 if (elem->level() == 0) 02190 { 02191 conn.push_back(-1); 02192 conn.push_back(-1); 02193 } 02194 else 02195 { 02196 conn.push_back(elem->parent()->id()); 02197 conn.push_back(elem->parent()->which_child_am_i(elem)); 02198 } 02199 #else 02200 conn.push_back (-1); 02201 conn.push_back (-1); 02202 #endif 02203 02204 for (unsigned int n=0; n<elem->n_nodes(); n++) 02205 conn.push_back (elem->node(n)); 02206 02207 for (unsigned int n=0; n<elem->n_neighbors(); n++) 02208 { 02209 Elem *neigh = elem->neighbor(n); 02210 if (neigh) 02211 conn.push_back (neigh->id()); 02212 else 02213 conn.push_back (-1); 02214 } 02215 02216 elem->pack_indexing(std::back_inserter(conn)); 02217 } 02218 02219 02220 02221 Elem * Elem::PackedElem::unpack (MeshBase &mesh, Elem *parent) const 02222 { 02223 02224 Elem *elem = Elem::build(this->type(),parent).release(); 02225 libmesh_assert (elem); 02226 02227 #ifdef LIBMESH_ENABLE_AMR 02228 if (this->level() != 0) 02229 { 02230 libmesh_assert(parent); 02231 parent->add_child(elem, this->which_child_am_i()); 02232 libmesh_assert_equal_to (parent->type(), elem->type()); 02233 libmesh_assert_equal_to (parent->child(this->which_child_am_i()), elem); 02234 } 02235 #endif 02236 02237 // Assign the refinement flags and levels 02238 #ifdef LIBMESH_ENABLE_AMR 02239 elem->set_p_level(this->p_level()); 02240 elem->set_refinement_flag(this->refinement_flag()); 02241 elem->set_p_refinement_flag(this->p_refinement_flag()); 02242 libmesh_assert_equal_to (elem->level(), this->level()); 02243 02244 // If this element definitely should have children, assign 02245 // remote_elem for now; later unpacked elements may overwrite that. 02246 if (!elem->active()) 02247 for (unsigned int c=0; c != elem->n_children(); ++c) 02248 elem->add_child(const_cast<RemoteElem*>(remote_elem), c); 02249 #endif 02250 02251 // Assign the IDs 02252 elem->subdomain_id() = this->subdomain_id(); 02253 elem->processor_id() = this->processor_id(); 02254 elem->set_id() = this->id(); 02255 02256 // Assign the connectivity 02257 libmesh_assert_equal_to (elem->n_nodes(), this->n_nodes()); 02258 02259 for (unsigned int n=0; n<elem->n_nodes(); n++) 02260 elem->set_node(n) = mesh.node_ptr (this->node(n)); 02261 02262 // Assign the connectivity 02263 libmesh_assert_equal_to (elem->n_neighbors(), this->n_neighbors()); 02264 02265 for (unsigned int n=0; n<elem->n_neighbors(); n++) 02266 { 02267 dof_id_type neighbor_id = this->neighbor(n); 02268 02269 // We should only be unpacking elements sent by their owners, 02270 // and their owners should know all their neighbors 02271 libmesh_assert_not_equal_to (neighbor_id, remote_elem->id()); 02272 02273 if (neighbor_id == DofObject::invalid_id) 02274 continue; 02275 02276 Elem *neigh = mesh.query_elem(neighbor_id); 02277 if (!neigh) 02278 { 02279 elem->set_neighbor(n, const_cast<RemoteElem*>(remote_elem)); 02280 continue; 02281 } 02282 02283 // We never have neighbors more refined than us 02284 libmesh_assert_less_equal (neigh->level(), elem->level()); 02285 02286 // We never have subactive neighbors of non subactive elements 02287 libmesh_assert(!neigh->subactive() || elem->subactive()); 02288 02289 // If we have a neighbor less refined than us then it must not 02290 // have any more refined active descendants we could have 02291 // pointed to instead. 02292 libmesh_assert(neigh->level() == elem->level() || 02293 neigh->active()); 02294 02295 elem->set_neighbor(n, neigh); 02296 02297 // If neigh is at elem's level, then its family might have 02298 // remote_elem neighbor links which need to point to elem 02299 // instead, but if not, then we're done. 02300 if (neigh->level() != elem->level()) 02301 continue; 02302 02303 // What side of neigh is elem on? We can't use the usual Elem 02304 // method because we haven't finished restoring topology 02305 const AutoPtr<Elem> my_side = elem->side(n); 02306 unsigned int nn = 0; 02307 for (; nn != neigh->n_sides(); ++nn) 02308 { 02309 const AutoPtr<Elem> neigh_side = neigh->side(nn); 02310 if (*my_side == *neigh_side) 02311 break; 02312 } 02313 02314 // elem had better be on *some* side of neigh 02315 libmesh_assert_less (nn, neigh->n_sides()); 02316 02317 // Find any elements that ought to point to elem 02318 std::vector<const Elem*> neigh_family; 02319 #ifdef LIBMESH_ENABLE_AMR 02320 if (!neigh->subactive()) 02321 neigh->family_tree_by_side(neigh_family, nn); 02322 #else 02323 neigh_family.push_back(neigh); 02324 #endif 02325 02326 // And point them to elem 02327 for (unsigned int i = 0; i != neigh_family.size(); ++i) 02328 { 02329 Elem* neigh_family_member = const_cast<Elem*>(neigh_family[i]); 02330 02331 // The neighbor link ought to either be correct already or 02332 // ought to be to remote_elem 02333 libmesh_assert(neigh_family_member->neighbor(nn) == elem || 02334 neigh_family_member->neighbor(nn) == remote_elem); 02335 02336 neigh_family_member->set_neighbor(nn, elem); 02337 } 02338 } 02339 02340 elem->unpack_indexing(this->indices()); 02341 02342 return elem; 02343 } 02344 02345 02346 02347 unsigned int Elem::opposite_side(const unsigned int /*s*/) const 02348 { 02349 // If the subclass didn't rederive this, using it is an error 02350 libmesh_not_implemented(); 02351 } 02352 02353 02354 02355 unsigned int Elem::opposite_node(const unsigned int /*n*/, 02356 const unsigned int /*s*/) const 02357 { 02358 // If the subclass didn't rederive this, using it is an error 02359 libmesh_not_implemented(); 02360 } 02361 02362 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: