mesh_modification.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 <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC 00022 #include <cmath> // for std::acos() 00023 #include <algorithm> 00024 #include <limits> 00025 #include <map> 00026 00027 // Local includes 00028 #include "libmesh/boundary_info.h" 00029 #include "libmesh/face_tri3.h" 00030 #include "libmesh/face_tri6.h" 00031 #include "libmesh/libmesh_logging.h" 00032 #include "libmesh/location_maps.h" 00033 #include "libmesh/mesh_communication.h" 00034 #include "libmesh/mesh_modification.h" 00035 #include "libmesh/mesh_tools.h" 00036 #include "libmesh/parallel.h" 00037 #include "libmesh/remote_elem.h" 00038 #include "libmesh/string_to_enum.h" 00039 #include "libmesh/unstructured_mesh.h" 00040 00041 namespace libMesh 00042 { 00043 00044 00045 // ------------------------------------------------------------ 00046 // MeshTools::Modification functions for mesh modification 00047 void MeshTools::Modification::distort (MeshBase& mesh, 00048 const Real factor, 00049 const bool perturb_boundary) 00050 { 00051 libmesh_assert (mesh.n_nodes()); 00052 libmesh_assert (mesh.n_elem()); 00053 libmesh_assert ((factor >= 0.) && (factor <= 1.)); 00054 00055 START_LOG("distort()", "MeshTools::Modification"); 00056 00057 00058 00059 // First find nodes on the boundary and flag them 00060 // so that we don't move them 00061 // on_boundary holds false (not on boundary) and true (on boundary) 00062 std::vector<bool> on_boundary (mesh.n_nodes(), false); 00063 00064 if (!perturb_boundary) MeshTools::find_boundary_nodes (mesh, on_boundary); 00065 00066 // Now calculate the minimum distance to 00067 // neighboring nodes for each node. 00068 // hmin holds these distances. 00069 std::vector<float> hmin (mesh.n_nodes(), 00070 std::numeric_limits<float>::max()); 00071 00072 MeshBase::element_iterator el = mesh.active_elements_begin(); 00073 const MeshBase::element_iterator end = mesh.active_elements_end(); 00074 00075 for (; el!=end; ++el) 00076 for (unsigned int n=0; n<(*el)->n_nodes(); n++) 00077 hmin[(*el)->node(n)] = std::min(hmin[(*el)->node(n)], 00078 static_cast<float>((*el)->hmin())); 00079 00080 00081 // Now actually move the nodes 00082 { 00083 const unsigned int seed = 123456; 00084 00085 // seed the random number generator. 00086 // We'll loop from 1 to n_nodes on every processor, even those 00087 // that don't have a particular node, so that the pseudorandom 00088 // numbers will be the same everywhere. 00089 std::srand(seed); 00090 00091 // If the node is on the boundary or 00092 // the node is not used by any element (hmin[n]<1.e20) 00093 // then we should not move it. 00094 // [Note: Testing for (in)equality might be wrong 00095 // (different types, namely float and double)] 00096 for (unsigned int n=0; n<mesh.n_nodes(); n++) 00097 if (!on_boundary[n] && (hmin[n] < 1.e20) ) 00098 { 00099 // the direction, random but unit normalized 00100 00101 Point dir( static_cast<Real>(std::rand())/static_cast<Real>(RAND_MAX), 00102 (mesh.mesh_dimension() > 1) ? 00103 static_cast<Real>(std::rand())/static_cast<Real>(RAND_MAX) 00104 : 0., 00105 ((mesh.mesh_dimension() == 3) ? 00106 static_cast<Real>(std::rand())/static_cast<Real>(RAND_MAX) 00107 : 0.) 00108 ); 00109 00110 dir(0) = (dir(0)-.5)*2.; 00111 if (mesh.mesh_dimension() > 1) 00112 dir(1) = (dir(1)-.5)*2.; 00113 if (mesh.mesh_dimension() == 3) 00114 dir(2) = (dir(2)-.5)*2.; 00115 00116 dir = dir.unit(); 00117 00118 Node *node = mesh.node_ptr(n); 00119 if (!node) 00120 continue; 00121 00122 (*node)(0) += dir(0)*factor*hmin[n]; 00123 if (mesh.mesh_dimension() > 1) 00124 (*node)(1) += dir(1)*factor*hmin[n]; 00125 if (mesh.mesh_dimension() == 3) 00126 (*node)(2) += dir(2)*factor*hmin[n]; 00127 } 00128 } 00129 00130 00131 // All done 00132 STOP_LOG("distort()", "MeshTools::Modification"); 00133 } 00134 00135 00136 00137 void MeshTools::Modification::translate (MeshBase& mesh, 00138 const Real xt, 00139 const Real yt, 00140 const Real zt) 00141 { 00142 const Point p(xt, yt, zt); 00143 00144 const MeshBase::node_iterator nd_end = mesh.nodes_end(); 00145 00146 for (MeshBase::node_iterator nd = mesh.nodes_begin(); 00147 nd != nd_end; ++nd) 00148 **nd += p; 00149 } 00150 00151 00152 // void MeshTools::Modification::rotate2D (MeshBase& mesh, 00153 // const Real alpha) 00154 // { 00155 // libmesh_assert_not_equal_to (mesh.mesh_dimension(), 1); 00156 00157 // const Real pi = std::acos(-1); 00158 // const Real a = alpha/180.*pi; 00159 // for (unsigned int n=0; n<mesh.n_nodes(); n++) 00160 // { 00161 // const Point p = mesh.node(n); 00162 // const Real x = p(0); 00163 // const Real y = p(1); 00164 // const Real z = p(2); 00165 // mesh.node(n) = Point(std::cos(a)*x - std::sin(a)*y, 00166 // std::sin(a)*x + std::cos(a)*y, 00167 // z); 00168 // } 00169 00170 // } 00171 00172 00173 00174 void MeshTools::Modification::rotate (MeshBase& mesh, 00175 const Real phi, 00176 const Real theta, 00177 const Real psi) 00178 { 00179 libmesh_assert_not_equal_to (mesh.mesh_dimension(), 1); 00180 00181 const Real p = -phi/180.*libMesh::pi; 00182 const Real t = -theta/180.*libMesh::pi; 00183 const Real s = -psi/180.*libMesh::pi; 00184 const Real sp = std::sin(p), cp = std::cos(p); 00185 const Real st = std::sin(t), ct = std::cos(t); 00186 const Real ss = std::sin(s), cs = std::cos(s); 00187 00188 // We follow the convention described at http://mathworld.wolfram.com/EulerAngles.html 00189 // (equations 6-14 give the entries of the composite transformation matrix). 00190 // The rotations are performed sequentially about the z, x, and z axes, in that order. 00191 // A positive angle yields a counter-clockwise rotation about the axis in question. 00192 const MeshBase::node_iterator nd_end = mesh.nodes_end(); 00193 00194 for (MeshBase::node_iterator nd = mesh.nodes_begin(); 00195 nd != nd_end; ++nd) 00196 { 00197 const Point pt = **nd; 00198 const Real x = pt(0); 00199 const Real y = pt(1); 00200 const Real z = pt(2); 00201 **nd = Point(( cp*cs-sp*ct*ss)*x + ( sp*cs+cp*ct*ss)*y + (st*ss)*z, 00202 (-cp*ss-sp*ct*cs)*x + (-sp*ss+cp*ct*cs)*y + (st*cs)*z, 00203 ( sp*st)*x + (-cp*st)*y + (ct)*z ); 00204 } 00205 } 00206 00207 00208 void MeshTools::Modification::scale (MeshBase& mesh, 00209 const Real xs, 00210 const Real ys, 00211 const Real zs) 00212 { 00213 const Real x_scale = xs; 00214 Real y_scale = ys; 00215 Real z_scale = zs; 00216 00217 if (ys == 0.) 00218 { 00219 libmesh_assert_equal_to (zs, 0.); 00220 00221 y_scale = z_scale = x_scale; 00222 } 00223 00224 // Scale the x coordinate in all dimensions 00225 const MeshBase::node_iterator nd_end = mesh.nodes_end(); 00226 00227 for (MeshBase::node_iterator nd = mesh.nodes_begin(); 00228 nd != nd_end; ++nd) 00229 (**nd)(0) *= x_scale; 00230 00231 00232 // Only scale the y coordinate in 2 and 3D 00233 if (mesh.spatial_dimension() < 2) 00234 return; 00235 00236 for (MeshBase::node_iterator nd = mesh.nodes_begin(); 00237 nd != nd_end; ++nd) 00238 (**nd)(1) *= y_scale; 00239 00240 // Only scale the z coordinate in 3D 00241 if (mesh.spatial_dimension() < 3) 00242 return; 00243 00244 for (MeshBase::node_iterator nd = mesh.nodes_begin(); 00245 nd != nd_end; ++nd) 00246 (**nd)(2) *= z_scale; 00247 } 00248 00249 00250 00251 00252 // ------------------------------------------------------------ 00253 // UnstructuredMesh class member functions for mesh modification 00254 void UnstructuredMesh::all_first_order () 00255 { 00256 /* 00257 * when the mesh is not prepared, 00258 * at least renumber the nodes and 00259 * elements, so that the node ids 00260 * are correct 00261 */ 00262 if (!this->_is_prepared) 00263 this->renumber_nodes_and_elements (); 00264 00265 START_LOG("all_first_order()", "Mesh"); 00266 00270 std::vector<bool> node_touched_by_me(this->max_node_id(), false); 00271 00277 element_iterator endit = elements_end(); 00278 for (element_iterator it = elements_begin(); 00279 it != endit; ++it) 00280 { 00281 Elem* so_elem = *it; 00282 00283 libmesh_assert(so_elem); 00284 00285 /* 00286 * build the first-order equivalent, add to 00287 * the new_elements list. 00288 */ 00289 Elem* lo_elem = Elem::build 00290 (Elem::first_order_equivalent_type 00291 (so_elem->type()), so_elem->parent()).release(); 00292 00293 for (unsigned int s=0; s != so_elem->n_sides(); ++s) 00294 if (so_elem->neighbor(s) == remote_elem) 00295 lo_elem->set_neighbor(s, const_cast<RemoteElem*>(remote_elem)); 00296 00297 #ifdef LIBMESH_ENABLE_AMR 00298 /* 00299 * Reset the parent links of any child elements 00300 */ 00301 if (so_elem->has_children()) 00302 for (unsigned int c=0; c != so_elem->n_children(); ++c) 00303 { 00304 so_elem->child(c)->set_parent(lo_elem); 00305 lo_elem->add_child(so_elem->child(c), c); 00306 } 00307 00308 /* 00309 * Reset the child link of any parent element 00310 */ 00311 if (so_elem->parent()) 00312 { 00313 unsigned int c = 00314 so_elem->parent()->which_child_am_i(so_elem); 00315 lo_elem->parent()->replace_child(lo_elem, c); 00316 } 00317 00318 /* 00319 * Copy as much data to the new element as makes sense 00320 */ 00321 lo_elem->set_p_level(so_elem->p_level()); 00322 lo_elem->set_refinement_flag(so_elem->refinement_flag()); 00323 lo_elem->set_p_refinement_flag(so_elem->p_refinement_flag()); 00324 #endif 00325 00326 libmesh_assert_equal_to (lo_elem->n_vertices(), so_elem->n_vertices()); 00327 00328 /* 00329 * By definition the vertices of the linear and 00330 * second order element are identically numbered. 00331 * transfer these. 00332 */ 00333 for (unsigned int v=0; v < so_elem->n_vertices(); v++) 00334 { 00335 lo_elem->set_node(v) = so_elem->get_node(v); 00336 node_touched_by_me[lo_elem->node(v)] = true; 00337 } 00338 00345 libmesh_assert_equal_to (lo_elem->n_sides(), so_elem->n_sides()); 00346 00347 for (unsigned int s=0; s<so_elem->n_sides(); s++) 00348 { 00349 const std::vector<boundary_id_type> boundary_ids = 00350 this->boundary_info->raw_boundary_ids (so_elem, s); 00351 00352 this->boundary_info->add_side (lo_elem, s, boundary_ids); 00353 } 00354 00355 /* 00356 * The new first-order element is ready. 00357 * Inserting it into the mesh will replace and delete 00358 * the second-order element. 00359 */ 00360 lo_elem->set_id(so_elem->id()); 00361 lo_elem->processor_id() = so_elem->processor_id(); 00362 lo_elem->subdomain_id() = so_elem->subdomain_id(); 00363 this->insert_elem(lo_elem); 00364 } 00365 00366 const MeshBase::node_iterator nd_end = this->nodes_end(); 00367 MeshBase::node_iterator nd = this->nodes_begin(); 00368 while (nd != nd_end) 00369 { 00370 Node *the_node = *nd; 00371 ++nd; 00372 if (!node_touched_by_me[the_node->id()]) 00373 this->delete_node(the_node); 00374 } 00375 00376 STOP_LOG("all_first_order()", "Mesh"); 00377 00378 // On hanging nodes that used to also be second order nodes, we 00379 // might now have an invalid nodal processor_id() 00380 Partitioner::set_node_processor_ids(*this); 00381 00382 // delete or renumber nodes, etc 00383 this->prepare_for_use(/*skip_renumber =*/ false); 00384 } 00385 00386 00387 00388 void UnstructuredMesh::all_second_order (const bool full_ordered) 00389 { 00390 // This function must be run on all processors at once 00391 parallel_only(); 00392 00393 /* 00394 * when the mesh is not prepared, 00395 * at least renumber the nodes and 00396 * elements, so that the node ids 00397 * are correct 00398 */ 00399 if (!this->_is_prepared) 00400 this->renumber_nodes_and_elements (); 00401 00402 /* 00403 * If the mesh is empty 00404 * then we have nothing to do 00405 */ 00406 if (!this->n_elem()) 00407 return; 00408 00409 /* 00410 * If the mesh is already second order 00411 * then we have nothing to do. 00412 * We have to test for this in a round-about way to avoid 00413 * a bug on distributed parallel meshes with more processors 00414 * than elements. 00415 */ 00416 bool already_second_order = false; 00417 if (this->elements_begin() != this->elements_end() && 00418 (*(this->elements_begin()))->default_order() != FIRST) 00419 already_second_order = true; 00420 CommWorld.max(already_second_order); 00421 if (already_second_order) 00422 return; 00423 00424 START_LOG("all_second_order()", "Mesh"); 00425 00426 /* 00427 * this map helps in identifying second order 00428 * nodes. Namely, a second-order node: 00429 * - edge node 00430 * - face node 00431 * - bubble node 00432 * is uniquely defined through a set of adjacent 00433 * vertices. This set of adjacent vertices is 00434 * used to identify already added higher-order 00435 * nodes. We are safe to use node id's since we 00436 * make sure that these are correctly numbered. 00437 */ 00438 std::map<std::vector<dof_id_type>, Node*> adj_vertices_to_so_nodes; 00439 00440 /* 00441 * for speed-up of the \p add_point() method, we 00442 * can reserve memory. Guess the number of additional 00443 * nodes for different dimensions 00444 */ 00445 switch (this->mesh_dimension()) 00446 { 00447 case 1: 00448 /* 00449 * in 1D, there can only be order-increase from Edge2 00450 * to Edge3. Something like 1/2 of n_nodes() have 00451 * to be added 00452 */ 00453 this->reserve_nodes(static_cast<unsigned int> 00454 (1.5*static_cast<double>(this->n_nodes()))); 00455 break; 00456 00457 case 2: 00458 /* 00459 * in 2D, either refine from Tri3 to Tri6 (double the nodes) 00460 * or from Quad4 to Quad8 (again, double) or Quad9 (2.25 that much) 00461 */ 00462 this->reserve_nodes(static_cast<unsigned int> 00463 (2*static_cast<double>(this->n_nodes()))); 00464 break; 00465 00466 00467 case 3: 00468 /* 00469 * in 3D, either refine from Tet4 to Tet10 (factor = 2.5) up to 00470 * Hex8 to Hex27 (something > 3). Since in 3D there _are_ already 00471 * quite some nodes, and since we do not want to overburden the memory by 00472 * a too conservative guess, use the lower bound 00473 */ 00474 this->reserve_nodes(static_cast<unsigned int> 00475 (2.5*static_cast<double>(this->n_nodes()))); 00476 break; 00477 00478 default: 00479 // Hm? 00480 libmesh_error(); 00481 } 00482 00483 00484 00485 /* 00486 * form a vector that will hold the node id's of 00487 * the vertices that are adjacent to the son-th 00488 * second-order node. Pull this outside of the 00489 * loop so that silly compilers don't repeatedly 00490 * create and destroy the vector. 00491 */ 00492 std::vector<dof_id_type> adjacent_vertices_ids; 00493 00500 const_element_iterator endit = elements_end(); 00501 for (const_element_iterator it = elements_begin(); 00502 it != endit; ++it) 00503 { 00504 // the linear-order element 00505 const Elem* lo_elem = *it; 00506 00507 libmesh_assert(lo_elem); 00508 00509 // make sure it is linear order 00510 if (lo_elem->default_order() != FIRST) 00511 { 00512 libMesh::err << "ERROR: This is not a linear element: type=" 00513 << lo_elem->type() << std::endl; 00514 libmesh_error(); 00515 } 00516 00517 // this does _not_ work for refined elements 00518 libmesh_assert_equal_to (lo_elem->level (), 0); 00519 00520 /* 00521 * build the second-order equivalent, add to 00522 * the new_elements list. Note that this here 00523 * is the only point where \p full_ordered 00524 * is necessary. The remaining code works well 00525 * for either type of seconrd-order equivalent, e.g. 00526 * Hex20 or Hex27, as equivalents for Hex8 00527 */ 00528 Elem* so_elem = 00529 Elem::build (Elem::second_order_equivalent_type(lo_elem->type(), 00530 full_ordered) ).release(); 00531 00532 libmesh_assert_equal_to (lo_elem->n_vertices(), so_elem->n_vertices()); 00533 00534 00535 /* 00536 * By definition the vertices of the linear and 00537 * second order element are identically numbered. 00538 * transfer these. 00539 */ 00540 for (unsigned int v=0; v < lo_elem->n_vertices(); v++) 00541 so_elem->set_node(v) = lo_elem->get_node(v); 00542 00543 /* 00544 * Now handle the additional mid-side nodes. This 00545 * is simply handled through a map that remembers 00546 * the already-added nodes. This map maps the global 00547 * ids of the vertices (that uniquely define this 00548 * higher-order node) to the new node. 00549 * Notation: son = second-order node 00550 */ 00551 const unsigned int son_begin = so_elem->n_vertices(); 00552 const unsigned int son_end = so_elem->n_nodes(); 00553 00554 00555 for (unsigned int son=son_begin; son<son_end; son++) 00556 { 00557 const unsigned int n_adjacent_vertices = 00558 so_elem->n_second_order_adjacent_vertices(son); 00559 00560 adjacent_vertices_ids.resize(n_adjacent_vertices); 00561 00562 for (unsigned int v=0; v<n_adjacent_vertices; v++) 00563 adjacent_vertices_ids[v] = 00564 so_elem->node( so_elem->second_order_adjacent_vertex(son,v) ); 00565 00566 /* 00567 * \p adjacent_vertices_ids is now in order of the current 00568 * side. sort it, so that comparisons with the 00569 * \p adjacent_vertices_ids created through other elements' 00570 * sides can match 00571 */ 00572 std::sort(adjacent_vertices_ids.begin(), 00573 adjacent_vertices_ids.end()); 00574 00575 00576 // does this set of vertices already has a mid-node added? 00577 std::pair<std::map<std::vector<dof_id_type>, Node*>::iterator, 00578 std::map<std::vector<dof_id_type>, Node*>::iterator> 00579 pos = adj_vertices_to_so_nodes.equal_range (adjacent_vertices_ids); 00580 00581 // no, not added yet 00582 if (pos.first == pos.second) 00583 { 00584 /* 00585 * for this set of vertices, there is no 00586 * second_order node yet. Add it. 00587 * 00588 * compute the location of the new node as 00589 * the average over the adjacent vertices. 00590 */ 00591 Point new_location = this->point(adjacent_vertices_ids[0]); 00592 for (unsigned int v=1; v<n_adjacent_vertices; v++) 00593 new_location += this->point(adjacent_vertices_ids[v]); 00594 00595 new_location /= static_cast<Real>(n_adjacent_vertices); 00596 00597 /* Add the new point to the mesh, giving it a globally 00598 * well-defined processor id. 00599 */ 00600 Node* so_node = this->add_point 00601 (new_location, DofObject::invalid_id, 00602 this->node(adjacent_vertices_ids[0]).processor_id()); 00603 00604 /* 00605 * insert the new node with its defining vertex 00606 * set into the map, and relocate pos to this 00607 * new entry, so that the so_elem can use 00608 * \p pos for inserting the node 00609 */ 00610 adj_vertices_to_so_nodes.insert(pos.first, 00611 std::make_pair(adjacent_vertices_ids, 00612 so_node)); 00613 00614 so_elem->set_node(son) = so_node; 00615 } 00616 // yes, already added. 00617 else 00618 { 00619 libmesh_assert(pos.first->second); 00620 00621 so_elem->set_node(son) = pos.first->second; 00622 } 00623 } 00624 00625 00637 libmesh_assert_equal_to (lo_elem->n_sides(), so_elem->n_sides()); 00638 00639 for (unsigned int s=0; s<lo_elem->n_sides(); s++) 00640 { 00641 const std::vector<boundary_id_type> boundary_ids = 00642 this->boundary_info->raw_boundary_ids (lo_elem, s); 00643 00644 this->boundary_info->add_side (so_elem, s, boundary_ids); 00645 00646 if (lo_elem->neighbor(s) == remote_elem) 00647 so_elem->set_neighbor(s, const_cast<RemoteElem*>(remote_elem)); 00648 } 00649 00650 /* 00651 * The new second-order element is ready. 00652 * Inserting it into the mesh will replace and delete 00653 * the first-order element. 00654 */ 00655 so_elem->set_id(lo_elem->id()); 00656 so_elem->processor_id() = lo_elem->processor_id(); 00657 so_elem->subdomain_id() = lo_elem->subdomain_id(); 00658 this->insert_elem(so_elem); 00659 } 00660 00661 // we can clear the map 00662 adj_vertices_to_so_nodes.clear(); 00663 00664 00665 STOP_LOG("all_second_order()", "Mesh"); 00666 00667 // In a ParallelMesh our ghost node processor ids may be bad and 00668 // the ids of nodes touching remote elements may be inconsistent. 00669 // Fix them. 00670 if (!this->is_serial()) 00671 { 00672 LocationMap<Node> loc_map; 00673 MeshCommunication().make_nodes_parallel_consistent 00674 (*this, loc_map); 00675 } 00676 00677 // renumber nodes, elements etc 00678 this->prepare_for_use(/*skip_renumber =*/ false); 00679 } 00680 00681 00682 00683 00684 00685 00686 void MeshTools::Modification::all_tri (MeshBase& mesh) 00687 { 00688 // The number of elements in the original mesh before any additions 00689 // or deletions. 00690 const dof_id_type n_orig_elem = mesh.n_elem(); 00691 00692 // We store pointers to the newly created elements in a vector 00693 // until they are ready to be added to the mesh. This is because 00694 // adding new elements on the fly can cause reallocation and invalidation 00695 // of existing iterators. 00696 std::vector<Elem*> new_elements; 00697 new_elements.reserve (2*n_orig_elem); 00698 00699 // If the original mesh has boundary data, we carry that over 00700 // to the new mesh with triangular elements. 00701 const bool mesh_has_boundary_data = (mesh.boundary_info->n_boundary_ids() > 0); 00702 00703 // Temporary vectors to store the new boundary element pointers, side numbers, and boundary ids 00704 std::vector<Elem*> new_bndry_elements; 00705 std::vector<unsigned short int> new_bndry_sides; 00706 std::vector<boundary_id_type> new_bndry_ids; 00707 00708 // Iterate over the elements, splitting QUADS into 00709 // pairs of conforming triangles. 00710 // FIXME: This algorithm does not work on refined grids! 00711 { 00712 MeshBase::element_iterator el = mesh.elements_begin(); 00713 const MeshBase::element_iterator end = mesh.elements_end(); 00714 00715 for (; el!=end; ++el) 00716 { 00717 Elem* elem = *el; 00718 00719 const ElemType etype = elem->type(); 00720 00721 // all_tri currently only works on coarse meshes 00722 libmesh_assert (!elem->parent()); 00723 00724 // We split the quads using the shorter of the two diagonals 00725 // to maintain the best angle properties. 00726 bool edge_swap = false; 00727 00728 // True if we actually split the current element. 00729 bool split_elem = false; 00730 00731 // The two new triangular elements we will split the quad into. 00732 Elem* tri0 = NULL; 00733 Elem* tri1 = NULL; 00734 00735 00736 switch (etype) 00737 { 00738 case QUAD4: 00739 { 00740 split_elem = true; 00741 00742 tri0 = new Tri3; 00743 tri1 = new Tri3; 00744 00745 // Check for possible edge swap 00746 if ((elem->point(0) - elem->point(2)).size() < 00747 (elem->point(1) - elem->point(3)).size()) 00748 { 00749 tri0->set_node(0) = elem->get_node(0); 00750 tri0->set_node(1) = elem->get_node(1); 00751 tri0->set_node(2) = elem->get_node(2); 00752 00753 tri1->set_node(0) = elem->get_node(0); 00754 tri1->set_node(1) = elem->get_node(2); 00755 tri1->set_node(2) = elem->get_node(3); 00756 } 00757 00758 else 00759 { 00760 edge_swap=true; 00761 00762 tri0->set_node(0) = elem->get_node(0); 00763 tri0->set_node(1) = elem->get_node(1); 00764 tri0->set_node(2) = elem->get_node(3); 00765 00766 tri1->set_node(0) = elem->get_node(1); 00767 tri1->set_node(1) = elem->get_node(2); 00768 tri1->set_node(2) = elem->get_node(3); 00769 } 00770 00771 00772 break; 00773 } 00774 00775 case QUAD8: 00776 { 00777 split_elem = true; 00778 00779 tri0 = new Tri6; 00780 tri1 = new Tri6; 00781 00782 Node* new_node = mesh.add_point( (mesh.node(elem->node(0)) + 00783 mesh.node(elem->node(1)) + 00784 mesh.node(elem->node(2)) + 00785 mesh.node(elem->node(3)) / 4) 00786 ); 00787 00788 // Check for possible edge swap 00789 if ((elem->point(0) - elem->point(2)).size() < 00790 (elem->point(1) - elem->point(3)).size()) 00791 { 00792 tri0->set_node(0) = elem->get_node(0); 00793 tri0->set_node(1) = elem->get_node(1); 00794 tri0->set_node(2) = elem->get_node(2); 00795 tri0->set_node(3) = elem->get_node(4); 00796 tri0->set_node(4) = elem->get_node(5); 00797 tri0->set_node(5) = new_node; 00798 00799 tri1->set_node(0) = elem->get_node(0); 00800 tri1->set_node(1) = elem->get_node(2); 00801 tri1->set_node(2) = elem->get_node(3); 00802 tri1->set_node(3) = new_node; 00803 tri1->set_node(4) = elem->get_node(6); 00804 tri1->set_node(5) = elem->get_node(7); 00805 00806 } 00807 00808 else 00809 { 00810 edge_swap=true; 00811 00812 tri0->set_node(0) = elem->get_node(3); 00813 tri0->set_node(1) = elem->get_node(0); 00814 tri0->set_node(2) = elem->get_node(1); 00815 tri0->set_node(3) = elem->get_node(7); 00816 tri0->set_node(4) = elem->get_node(4); 00817 tri0->set_node(5) = new_node; 00818 00819 tri1->set_node(0) = elem->get_node(1); 00820 tri1->set_node(1) = elem->get_node(2); 00821 tri1->set_node(2) = elem->get_node(3); 00822 tri1->set_node(3) = elem->get_node(5); 00823 tri1->set_node(4) = elem->get_node(6); 00824 tri1->set_node(5) = new_node; 00825 } 00826 00827 break; 00828 } 00829 00830 case QUAD9: 00831 { 00832 split_elem = true; 00833 00834 tri0 = new Tri6; 00835 tri1 = new Tri6; 00836 00837 // Check for possible edge swap 00838 if ((elem->point(0) - elem->point(2)).size() < 00839 (elem->point(1) - elem->point(3)).size()) 00840 { 00841 tri0->set_node(0) = elem->get_node(0); 00842 tri0->set_node(1) = elem->get_node(1); 00843 tri0->set_node(2) = elem->get_node(2); 00844 tri0->set_node(3) = elem->get_node(4); 00845 tri0->set_node(4) = elem->get_node(5); 00846 tri0->set_node(5) = elem->get_node(8); 00847 00848 tri1->set_node(0) = elem->get_node(0); 00849 tri1->set_node(1) = elem->get_node(2); 00850 tri1->set_node(2) = elem->get_node(3); 00851 tri1->set_node(3) = elem->get_node(8); 00852 tri1->set_node(4) = elem->get_node(6); 00853 tri1->set_node(5) = elem->get_node(7); 00854 } 00855 00856 else 00857 { 00858 edge_swap=true; 00859 00860 tri0->set_node(0) = elem->get_node(0); 00861 tri0->set_node(1) = elem->get_node(1); 00862 tri0->set_node(2) = elem->get_node(3); 00863 tri0->set_node(3) = elem->get_node(4); 00864 tri0->set_node(4) = elem->get_node(8); 00865 tri0->set_node(5) = elem->get_node(7); 00866 00867 tri1->set_node(0) = elem->get_node(1); 00868 tri1->set_node(1) = elem->get_node(2); 00869 tri1->set_node(2) = elem->get_node(3); 00870 tri1->set_node(3) = elem->get_node(5); 00871 tri1->set_node(4) = elem->get_node(6); 00872 tri1->set_node(5) = elem->get_node(8); 00873 } 00874 00875 break; 00876 } 00877 // No need to split elements that are already triangles 00878 case TRI3: 00879 case TRI6: 00880 continue; 00881 // Try to ignore non-2D elements for now 00882 default: 00883 { 00884 libMesh::err << "Warning, encountered non-2D element " 00885 << Utility::enum_to_string<ElemType>(etype) 00886 << " in MeshTools::Modification::all_tri(), hope that's OK..." 00887 << std::endl; 00888 } 00889 } // end switch (etype) 00890 00891 00892 00893 if (split_elem) 00894 { 00895 // Be sure the correct ID's are also set for tri0 and 00896 // tri1. 00897 tri0->processor_id() = elem->processor_id(); 00898 tri0->subdomain_id() = elem->subdomain_id(); 00899 tri1->processor_id() = elem->processor_id(); 00900 tri1->subdomain_id() = elem->subdomain_id(); 00901 00902 if (mesh_has_boundary_data) 00903 { 00904 for (unsigned int sn=0; sn<elem->n_sides(); ++sn) 00905 { 00906 const std::vector<boundary_id_type>& bc_ids = mesh.boundary_info->boundary_ids(*el, sn); 00907 for (std::vector<boundary_id_type>::const_iterator id_it=bc_ids.begin(); id_it!=bc_ids.end(); ++id_it) 00908 { 00909 const boundary_id_type b_id = *id_it; 00910 00911 if (b_id != BoundaryInfo::invalid_id) 00912 { 00913 // Add the boundary ID to the list of new boundary ids 00914 new_bndry_ids.push_back(b_id); 00915 00916 // Convert the boundary side information of the old element to 00917 // boundary side information for the new element. 00918 if (!edge_swap) 00919 { 00920 switch (sn) 00921 { 00922 case 0: 00923 { 00924 // New boundary side is Tri 0, side 0 00925 new_bndry_elements.push_back(tri0); 00926 new_bndry_sides.push_back(0); 00927 break; 00928 } 00929 case 1: 00930 { 00931 // New boundary side is Tri 0, side 1 00932 new_bndry_elements.push_back(tri0); 00933 new_bndry_sides.push_back(1); 00934 break; 00935 } 00936 case 2: 00937 { 00938 // New boundary side is Tri 1, side 1 00939 new_bndry_elements.push_back(tri1); 00940 new_bndry_sides.push_back(1); 00941 break; 00942 } 00943 case 3: 00944 { 00945 // New boundary side is Tri 1, side 2 00946 new_bndry_elements.push_back(tri1); 00947 new_bndry_sides.push_back(2); 00948 break; 00949 } 00950 00951 default: 00952 { 00953 libMesh::err << "Quad4/8/9 cannot have more than 4 sides." << std::endl; 00954 libmesh_error(); 00955 } 00956 } 00957 } 00958 00959 else // edge_swap==true 00960 { 00961 switch (sn) 00962 { 00963 case 0: 00964 { 00965 // New boundary side is Tri 0, side 0 00966 new_bndry_elements.push_back(tri0); 00967 new_bndry_sides.push_back(0); 00968 break; 00969 } 00970 case 1: 00971 { 00972 // New boundary side is Tri 1, side 0 00973 new_bndry_elements.push_back(tri1); 00974 new_bndry_sides.push_back(0); 00975 break; 00976 } 00977 case 2: 00978 { 00979 // New boundary side is Tri 1, side 1 00980 new_bndry_elements.push_back(tri1); 00981 new_bndry_sides.push_back(1); 00982 break; 00983 } 00984 case 3: 00985 { 00986 // New boundary side is Tri 0, side 2 00987 new_bndry_elements.push_back(tri0); 00988 new_bndry_sides.push_back(2); 00989 break; 00990 } 00991 00992 default: 00993 { 00994 libMesh::err << "Quad4/8/9 cannot have more than 4 sides." << std::endl; 00995 libmesh_error(); 00996 } 00997 } 00998 } // end edge_swap==true 00999 } // end if (b_id != BoundaryInfo::invalid_id) 01000 } // end for loop over boundary IDs 01001 } // end for loop over sides 01002 01003 // Remove the original element from the BoundaryInfo structure. 01004 mesh.boundary_info->remove(elem); 01005 01006 } // end if (mesh_has_boundary_data) 01007 01008 01009 // On a distributed mesh, we need to preserve remote_elem 01010 // links, since prepare_for_use can't reconstruct them for 01011 // us. 01012 for (unsigned int sn=0; sn<elem->n_sides(); ++sn) 01013 { 01014 if (elem->neighbor(sn) == remote_elem) 01015 { 01016 // Create a remote_elem link on one of the new 01017 // elements corresponding to the link from the old 01018 // element. 01019 if (!edge_swap) 01020 { 01021 switch (sn) 01022 { 01023 case 0: 01024 { 01025 // New remote side is Tri 0, side 0 01026 tri0->set_neighbor(0, const_cast<RemoteElem*>(remote_elem)); 01027 break; 01028 } 01029 case 1: 01030 { 01031 // New remote side is Tri 0, side 1 01032 tri0->set_neighbor(1, const_cast<RemoteElem*>(remote_elem)); 01033 break; 01034 } 01035 case 2: 01036 { 01037 // New remote side is Tri 1, side 1 01038 tri1->set_neighbor(1, const_cast<RemoteElem*>(remote_elem)); 01039 break; 01040 } 01041 case 3: 01042 { 01043 // New remote side is Tri 1, side 2 01044 tri1->set_neighbor(2, const_cast<RemoteElem*>(remote_elem)); 01045 break; 01046 } 01047 01048 default: 01049 { 01050 libMesh::err << "Quad4/8/9 cannot have more than 4 sides." << std::endl; 01051 libmesh_error(); 01052 } 01053 } 01054 } 01055 01056 else // edge_swap==true 01057 { 01058 switch (sn) 01059 { 01060 case 0: 01061 { 01062 // New remote side is Tri 0, side 0 01063 tri0->set_neighbor(0, const_cast<RemoteElem*>(remote_elem)); 01064 break; 01065 } 01066 case 1: 01067 { 01068 // New remote side is Tri 1, side 0 01069 tri1->set_neighbor(0, const_cast<RemoteElem*>(remote_elem)); 01070 break; 01071 } 01072 case 2: 01073 { 01074 // New remote side is Tri 1, side 1 01075 tri1->set_neighbor(1, const_cast<RemoteElem*>(remote_elem)); 01076 break; 01077 } 01078 case 3: 01079 { 01080 // New remote side is Tri 0, side 2 01081 tri0->set_neighbor(2, const_cast<RemoteElem*>(remote_elem)); 01082 break; 01083 } 01084 01085 default: 01086 { 01087 libMesh::err << "Quad4/8/9 cannot have more than 4 sides." << std::endl; 01088 libmesh_error(); 01089 } 01090 } 01091 } // end edge_swap==true 01092 } // end if (elem->neighbor(sn) == remote_elem) 01093 } // end for loop over sides 01094 01095 // Determine new IDs for the split elements which will be 01096 // the same on all processors, therefore keeping the Mesh 01097 // in sync. Note: we offset the new IDs by n_orig_elem to 01098 // avoid overwriting any of the original IDs, this assumes 01099 // they were contiguously-numbered to begin with... 01100 tri0->set_id( n_orig_elem + 2*elem->id() + 0 ); 01101 tri1->set_id( n_orig_elem + 2*elem->id() + 1 ); 01102 01103 // Add the newly-created triangles to the temporary vector of new elements. 01104 new_elements.push_back(tri0); 01105 new_elements.push_back(tri1); 01106 01107 // Delete the original element 01108 mesh.delete_elem(elem); 01109 } // end if (split_elem) 01110 } // End for loop over elements 01111 } // end scope 01112 01113 01114 // Now, iterate over the new elements vector, and add them each to 01115 // the Mesh. 01116 { 01117 std::vector<Elem*>::iterator el = new_elements.begin(); 01118 const std::vector<Elem*>::iterator end = new_elements.end(); 01119 for (; el != end; ++el) 01120 mesh.add_elem(*el); 01121 } 01122 01123 if (mesh_has_boundary_data) 01124 { 01125 // By this time, we should have removed all of the original boundary sides 01126 // - except on a hybrid mesh, where we can't "start from a blank slate"! - RHS 01127 // libmesh_assert_equal_to (mesh.boundary_info->n_boundary_conds(), 0); 01128 01129 // Clear the boundary info, to be sure and start from a blank slate. 01130 // mesh.boundary_info->clear(); 01131 01132 // If the old mesh had boundary data, the new mesh better have some. 01133 libmesh_assert_greater (new_bndry_elements.size(), 0); 01134 01135 // We should also be sure that the lengths of the new boundary data vectors 01136 // are all the same. 01137 libmesh_assert_equal_to (new_bndry_elements.size(), new_bndry_sides.size()); 01138 libmesh_assert_equal_to (new_bndry_sides.size(), new_bndry_ids.size()); 01139 01140 // Add the new boundary info to the mesh 01141 for (unsigned int s=0; s<new_bndry_elements.size(); ++s) 01142 mesh.boundary_info->add_side(new_bndry_elements[s], 01143 new_bndry_sides[s], 01144 new_bndry_ids[s]); 01145 } 01146 01147 01148 // Prepare the newly created mesh for use. 01149 mesh.prepare_for_use(/*skip_renumber =*/ false); 01150 01151 // Let the new_elements and new_bndry_elements vectors go out of scope. 01152 } 01153 01154 01155 void MeshTools::Modification::smooth (MeshBase& mesh, 01156 const unsigned int n_iterations, 01157 const Real power) 01158 { 01162 libmesh_assert_equal_to (mesh.mesh_dimension(), 2); 01163 01164 /* 01165 * find the boundary nodes 01166 */ 01167 std::vector<bool> on_boundary; 01168 MeshTools::find_boundary_nodes(mesh, on_boundary); 01169 01170 for (unsigned int iter=0; iter<n_iterations; iter++) 01171 01172 { 01173 /* 01174 * loop over the mesh refinement level 01175 */ 01176 unsigned int n_levels = MeshTools::n_levels(mesh); 01177 for (unsigned int refinement_level=0; refinement_level != n_levels; 01178 refinement_level++) 01179 { 01180 // initialize the storage (have to do it on every level to get empty vectors 01181 std::vector<Point> new_positions; 01182 std::vector<Real> weight; 01183 new_positions.resize(mesh.n_nodes()); 01184 weight.resize(mesh.n_nodes()); 01185 01186 { 01187 /* 01188 * Loop over the elements to calculate new node positions 01189 */ 01190 MeshBase::element_iterator el = mesh.level_elements_begin(refinement_level); 01191 const MeshBase::element_iterator end = mesh.level_elements_end(refinement_level); 01192 01193 for (; el != end; ++el) 01194 { 01195 /* 01196 * Constant handle for the element 01197 */ 01198 const Elem* elem = *el; 01199 01200 /* 01201 * We relax all nodes on level 0 first 01202 * If the element is refined (level > 0), we interpolate the 01203 * parents nodes with help of the embedding matrix 01204 */ 01205 if (refinement_level == 0) 01206 { 01207 for (unsigned int s=0; s<elem->n_neighbors(); s++) 01208 { 01209 /* 01210 * Only operate on sides which are on the 01211 * boundary or for which the current element's 01212 * id is greater than its neighbor's. 01213 * Sides get only built once. 01214 */ 01215 if ((elem->neighbor(s) != NULL) && 01216 (elem->id() > elem->neighbor(s)->id()) ) 01217 { 01218 AutoPtr<Elem> side(elem->build_side(s)); 01219 01220 Node* node0 = side->get_node(0); 01221 Node* node1 = side->get_node(1); 01222 01223 Real node_weight = 1.; 01224 // calculate the weight of the nodes 01225 if (power > 0) 01226 { 01227 Point diff = (*node0)-(*node1); 01228 node_weight = std::pow( diff.size(), power ); 01229 } 01230 01231 const dof_id_type id0 = node0->id(), id1 = node1->id(); 01232 new_positions[id0].add_scaled( *node1, node_weight ); 01233 new_positions[id1].add_scaled( *node0, node_weight ); 01234 weight[id0] += node_weight; 01235 weight[id1] += node_weight; 01236 } 01237 } // element neighbor loop 01238 } 01239 #ifdef LIBMESH_ENABLE_AMR 01240 else // refinement_level > 0 01241 { 01242 /* 01243 * Find the positions of the hanging nodes of refined elements. 01244 * We do this by calculating their position based on the parent 01245 * (one level less refined) element, and the embedding matrix 01246 */ 01247 01248 const Elem* parent = elem->parent(); 01249 01250 /* 01251 * find out which child I am 01252 */ 01253 for (unsigned int c=0; c < parent->n_children(); c++) 01254 { 01255 if (parent->child(c) == elem) 01256 { 01257 /* 01258 *loop over the childs (that is, the current elements) nodes 01259 */ 01260 for (unsigned int nc=0; nc < elem->n_nodes(); nc++) 01261 { 01262 /* 01263 * the new position of the node 01264 */ 01265 Point point; 01266 for (unsigned int n=0; n<parent->n_nodes(); n++) 01267 { 01268 /* 01269 * The value from the embedding matrix 01270 */ 01271 const float em_val = parent->embedding_matrix(c,nc,n); 01272 01273 if (em_val != 0.) 01274 point.add_scaled (parent->point(n), em_val); 01275 } 01276 01277 const dof_id_type id = elem->get_node(nc)->id(); 01278 new_positions[id] = point; 01279 weight[id] = 1.; 01280 } 01281 01282 } // if parent->child == elem 01283 } // for parent->n_children 01284 } // if element refinement_level 01285 #endif // #ifdef LIBMESH_ENABLE_AMR 01286 01287 } // element loop 01288 01289 /* 01290 * finally reposition the vertex nodes 01291 */ 01292 for (unsigned int nid=0; nid<mesh.n_nodes(); ++nid) 01293 if (!on_boundary[nid] && weight[nid] > 0.) 01294 mesh.node(nid) = new_positions[nid]/weight[nid]; 01295 } 01296 01297 { 01298 /* 01299 * Now handle the additional second_order nodes by calculating 01300 * their position based on the vertex postitions 01301 * we do a second loop over the level elements 01302 */ 01303 MeshBase::element_iterator el = mesh.level_elements_begin(refinement_level); 01304 const MeshBase::element_iterator end = mesh.level_elements_end(refinement_level); 01305 01306 for (; el != end; ++el) 01307 { 01308 /* 01309 * Constant handle for the element 01310 */ 01311 const Elem* elem = *el; 01312 const unsigned int son_begin = elem->n_vertices(); 01313 const unsigned int son_end = elem->n_nodes(); 01314 for (unsigned int n=son_begin; n<son_end; n++) 01315 { 01316 const unsigned int n_adjacent_vertices = 01317 elem->n_second_order_adjacent_vertices(n); 01318 01319 Point point; 01320 for (unsigned int v=0; v<n_adjacent_vertices; v++) 01321 point.add(elem->point( elem->second_order_adjacent_vertex(n,v) )); 01322 01323 const dof_id_type id = elem->get_node(n)->id(); 01324 mesh.node(id) = point/n_adjacent_vertices; 01325 } 01326 } 01327 } 01328 01329 } // refinement_level loop 01330 01331 } // end iteration 01332 } 01333 01334 01335 01336 #ifdef LIBMESH_ENABLE_AMR 01337 void MeshTools::Modification::flatten(MeshBase& mesh) 01338 { 01339 // Algorithm: 01340 // .) For each active element in the mesh: construct a 01341 // copy which is the same in every way *except* it is 01342 // a level 0 element. Store the pointers to these in 01343 // a separate vector. Save any boundary information as well. 01344 // Delete the active element from the mesh. 01345 // .) Loop over all (remaining) elements in the mesh, delete them. 01346 // .) Add the level-0 copies back to the mesh 01347 01348 // Temporary storage for new element pointers 01349 std::vector<Elem*> new_elements; 01350 01351 // BoundaryInfo Storage for element ids, sides, and BC ids 01352 std::vector<Elem*> saved_boundary_elements; 01353 std::vector<boundary_id_type> saved_bc_ids; 01354 std::vector<unsigned short int> saved_bc_sides; 01355 01356 // Reserve a reasonable amt. of space for each 01357 new_elements.reserve(mesh.n_active_elem()); 01358 saved_boundary_elements.reserve(mesh.boundary_info->n_boundary_conds()); 01359 saved_bc_ids.reserve(mesh.boundary_info->n_boundary_conds()); 01360 saved_bc_sides.reserve(mesh.boundary_info->n_boundary_conds()); 01361 { 01362 MeshBase::element_iterator it = mesh.active_elements_begin(); 01363 const MeshBase::element_iterator end = mesh.active_elements_end(); 01364 01365 for (; it != end; ++it) 01366 { 01367 Elem* elem = *it; 01368 01369 // Make a new element of the same type 01370 Elem* copy = Elem::build(elem->type()).release(); 01371 01372 // Set node pointers (they still point to nodes in the original mesh) 01373 for(unsigned int n=0; n<elem->n_nodes(); n++) 01374 copy->set_node(n) = elem->get_node(n); 01375 01376 // Copy over ids 01377 copy->processor_id() = elem->processor_id(); 01378 copy->subdomain_id() = elem->subdomain_id(); 01379 01380 // Retain the original element's ID as well, otherwise ParallelMesh will 01381 // try to create one for you... 01382 copy->set_id( elem->id() ); 01383 01384 // This element could have boundary info or ParallelMesh 01385 // remote_elem links as well. We need to save the (elem, 01386 // side, bc_id) triples and those links 01387 for (unsigned int s=0; s<elem->n_sides(); s++) 01388 { 01389 if (elem->neighbor(s) == remote_elem) 01390 copy->set_neighbor(s, const_cast<RemoteElem*>(remote_elem)); 01391 01392 const std::vector<boundary_id_type>& bc_ids = mesh.boundary_info->boundary_ids(elem,s); 01393 for (std::vector<boundary_id_type>::const_iterator id_it=bc_ids.begin(); id_it!=bc_ids.end(); ++id_it) 01394 { 01395 const boundary_id_type bc_id = *id_it; 01396 01397 if (bc_id != BoundaryInfo::invalid_id) 01398 { 01399 saved_boundary_elements.push_back(copy); 01400 saved_bc_ids.push_back(bc_id); 01401 saved_bc_sides.push_back(s); 01402 } 01403 } 01404 } 01405 01406 01407 // We're done with this element 01408 mesh.delete_elem(elem); 01409 01410 // But save the copy 01411 new_elements.push_back(copy); 01412 } 01413 01414 // Make sure we saved the same number of boundary conditions 01415 // in each vector. 01416 libmesh_assert_equal_to (saved_boundary_elements.size(), saved_bc_ids.size()); 01417 libmesh_assert_equal_to (saved_bc_ids.size(), saved_bc_sides.size()); 01418 } 01419 01420 01421 // Loop again, delete any remaining elements 01422 { 01423 MeshBase::element_iterator it = mesh.elements_begin(); 01424 const MeshBase::element_iterator end = mesh.elements_end(); 01425 01426 for (; it != end; ++it) 01427 mesh.delete_elem( *it ); 01428 } 01429 01430 01431 // Add the copied (now level-0) elements back to the mesh 01432 { 01433 for (std::vector<Elem*>::iterator it = new_elements.begin(); 01434 it != new_elements.end(); 01435 ++it) 01436 { 01437 dof_id_type orig_id = (*it)->id(); 01438 01439 Elem* added_elem = mesh.add_elem(*it); 01440 01441 dof_id_type added_id = added_elem->id(); 01442 01443 // If the Elem, as it was re-added to the mesh, now has a 01444 // different ID (this is unlikely, so it's just an assert) 01445 // the boundary information will no longer be correct. 01446 libmesh_assert_equal_to (orig_id, added_id); 01447 } 01448 } 01449 01450 // Finally, also add back the saved boundary information 01451 for (unsigned int e=0; e<saved_boundary_elements.size(); ++e) 01452 mesh.boundary_info->add_side(saved_boundary_elements[e], 01453 saved_bc_sides[e], 01454 saved_bc_ids[e]); 01455 01456 // Trim unused and renumber nodes and elements 01457 mesh.prepare_for_use(/*skip_renumber =*/ false); 01458 } 01459 #endif // #ifdef LIBMESH_ENABLE_AMR 01460 01461 01462 01463 void MeshTools::Modification::change_boundary_id (MeshBase& mesh, 01464 const boundary_id_type old_id, 01465 const boundary_id_type new_id) 01466 { 01467 // Only level-0 elements store BCs. Loop over them. 01468 MeshBase::element_iterator el = mesh.level_elements_begin(0); 01469 const MeshBase::element_iterator end_el = mesh.level_elements_end(0); 01470 for (; el != end_el; ++el) 01471 { 01472 Elem *elem = *el; 01473 unsigned int n_sides = elem->n_sides(); 01474 for (unsigned int s=0; s != n_sides; ++s) 01475 { 01476 const std::vector<boundary_id_type>& old_ids = mesh.boundary_info->boundary_ids(elem, s); 01477 if (std::find(old_ids.begin(), old_ids.end(), old_id) != old_ids.end()) 01478 { 01479 std::vector<boundary_id_type> new_ids(old_ids); 01480 std::replace(new_ids.begin(), new_ids.end(), old_id, new_id); 01481 mesh.boundary_info->remove_side(elem, s); 01482 mesh.boundary_info->add_side(elem, s, new_ids); 01483 } 01484 } 01485 } 01486 } 01487 01488 01489 01490 void MeshTools::Modification::change_subdomain_id (MeshBase& mesh, 01491 const subdomain_id_type old_id, 01492 const subdomain_id_type new_id) 01493 { 01494 MeshBase::element_iterator el = mesh.elements_begin(); 01495 const MeshBase::element_iterator end_el = mesh.elements_end(); 01496 01497 for (; el != end_el; ++el) 01498 { 01499 Elem *elem = *el; 01500 01501 if (elem->subdomain_id() == old_id) 01502 elem->subdomain_id() = new_id; 01503 } 01504 } 01505 01506 01507 } // namespace libMesh 01508 01509 01510 01511
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: