parmetis_partitioner.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 00022 // Local Includes ----------------------------------- 00023 #include "libmesh/libmesh_config.h" 00024 #include "libmesh/mesh_base.h" 00025 #include "libmesh/parallel.h" // also includes mpi.h 00026 #include "libmesh/mesh_serializer.h" 00027 #include "libmesh/mesh_tools.h" 00028 #include "libmesh/mesh_communication.h" 00029 #include "libmesh/parmetis_partitioner.h" 00030 #include "libmesh/metis_partitioner.h" 00031 #include "libmesh/libmesh_logging.h" 00032 #include "libmesh/elem.h" 00033 00034 #ifdef LIBMESH_HAVE_PARMETIS 00035 00036 // Include the ParMETIS header files 00037 namespace Parmetis { 00038 extern "C" { 00039 # include "libmesh/ignore_warnings.h" 00040 # include "parmetis.h" 00041 # include "libmesh/restore_warnings.h" 00042 } 00043 } 00044 00045 #endif // #ifdef LIBMESH_HAVE_PARMETIS ... else ... 00046 00047 namespace libMesh 00048 { 00049 00050 // Minimum elements on each processor required for us to choose 00051 // Parmetis over Metis. 00052 const unsigned int MIN_ELEM_PER_PROC = 4; 00053 00054 00055 // ------------------------------------------------------------ 00056 // ParmetisPartitioner implementation 00057 void ParmetisPartitioner::_do_partition (MeshBase& mesh, 00058 const unsigned int n_sbdmns) 00059 { 00060 this->_do_repartition (mesh, n_sbdmns); 00061 00062 // libmesh_assert_greater (n_sbdmns, 0); 00063 00064 // // Check for an easy return 00065 // if (n_sbdmns == 1) 00066 // { 00067 // this->single_partition (mesh); 00068 // return; 00069 // } 00070 00071 // // This function must be run on all processors at once 00072 // parallel_only(); 00073 00074 // // What to do if the Parmetis library IS NOT present 00075 // #ifndef LIBMESH_HAVE_PARMETIS 00076 00077 // libmesh_here(); 00078 // libMesh::err << "ERROR: The library has been built without" << std::endl 00079 // << "Parmetis support. Using a Metis" << std::endl 00080 // << "partitioner instead!" << std::endl; 00081 00082 // MetisPartitioner mp; 00083 00084 // mp.partition (mesh, n_sbdmns); 00085 00086 // // What to do if the Parmetis library IS present 00087 // #else 00088 00089 // START_LOG("partition()", "ParmetisPartitioner"); 00090 00091 // // Initialize the data structures required by ParMETIS 00092 // this->initialize (mesh, n_sbdmns); 00093 00094 // // build the graph corresponding to the mesh 00095 // this->build_graph (mesh); 00096 00097 00098 // // Partition the graph 00099 // MPI_Comm mpi_comm = libMesh::COMM_WORLD; 00100 00101 // // Call the ParMETIS k-way partitioning algorithm. 00102 // Parmetis::ParMETIS_V3_PartKway(&_vtxdist[0], &_xadj[0], &_adjncy[0], &_vwgt[0], NULL, 00103 // &_wgtflag, &_numflag, &_ncon, &_nparts, &_tpwgts[0], 00104 // &_ubvec[0], &_options[0], &_edgecut, 00105 // &_part[0], 00106 // &mpi_comm); 00107 00108 // // Assign the returned processor ids 00109 // this->assign_partitioning (mesh); 00110 00111 00112 // STOP_LOG ("partition()", "ParmetisPartitioner"); 00113 00114 // #endif // #ifndef LIBMESH_HAVE_PARMETIS ... else ... 00115 00116 } 00117 00118 00119 00120 void ParmetisPartitioner::_do_repartition (MeshBase& mesh, 00121 const unsigned int n_sbdmns) 00122 { 00123 libmesh_assert_greater (n_sbdmns, 0); 00124 00125 // Check for an easy return 00126 if (n_sbdmns == 1) 00127 { 00128 this->single_partition(mesh); 00129 return; 00130 } 00131 00132 // This function must be run on all processors at once 00133 parallel_only(); 00134 00135 // What to do if the Parmetis library IS NOT present 00136 #ifndef LIBMESH_HAVE_PARMETIS 00137 00138 libmesh_here(); 00139 libMesh::err << "ERROR: The library has been built without" << std::endl 00140 << "Parmetis support. Using a Metis" << std::endl 00141 << "partitioner instead!" << std::endl; 00142 00143 MetisPartitioner mp; 00144 00145 mp.partition (mesh, n_sbdmns); 00146 00147 // What to do if the Parmetis library IS present 00148 #else 00149 00150 // Revert to METIS on one processor. 00151 if (libMesh::n_processors() == 1) 00152 { 00153 MetisPartitioner mp; 00154 mp.partition (mesh, n_sbdmns); 00155 return; 00156 } 00157 00158 START_LOG("repartition()", "ParmetisPartitioner"); 00159 00160 // Initialize the data structures required by ParMETIS 00161 this->initialize (mesh, n_sbdmns); 00162 00163 // Make sure all processors have enough active local elements. 00164 // Parmetis tends to crash when it's given only a couple elements 00165 // per partition. 00166 { 00167 bool all_have_enough_elements = true; 00168 for (processor_id_type pid=0; pid<_n_active_elem_on_proc.size(); pid++) 00169 if (_n_active_elem_on_proc[pid] < MIN_ELEM_PER_PROC) 00170 all_have_enough_elements = false; 00171 00172 // Parmetis will not work unless each processor has some 00173 // elements. Specifically, it will abort when passed a NULL 00174 // partition array on *any* of the processors. 00175 if (!all_have_enough_elements) 00176 { 00177 // FIXME: revert to METIS, although this requires a serial mesh 00178 MeshSerializer serialize(mesh); 00179 00180 STOP_LOG ("repartition()", "ParmetisPartitioner"); 00181 00182 MetisPartitioner mp; 00183 mp.partition (mesh, n_sbdmns); 00184 00185 return; 00186 } 00187 } 00188 00189 // build the graph corresponding to the mesh 00190 this->build_graph (mesh); 00191 00192 00193 // Partition the graph 00194 std::vector<int> vsize(_vwgt.size(), 1); 00195 float itr = 1000000.0; 00196 MPI_Comm mpi_comm = libMesh::COMM_WORLD; 00197 00198 // Call the ParMETIS adaptive repartitioning method. This respects the 00199 // original partitioning when computing the new partitioning so as to 00200 // minimize the required data redistribution. 00201 Parmetis::ParMETIS_V3_AdaptiveRepart(_vtxdist.empty() ? NULL : &_vtxdist[0], 00202 _xadj.empty() ? NULL : &_xadj[0], 00203 _adjncy.empty() ? NULL : &_adjncy[0], 00204 _vwgt.empty() ? NULL : &_vwgt[0], 00205 vsize.empty() ? NULL : &vsize[0], 00206 NULL, 00207 &_wgtflag, 00208 &_numflag, 00209 &_ncon, 00210 &_nparts, 00211 _tpwgts.empty() ? NULL : &_tpwgts[0], 00212 _ubvec.empty() ? NULL : &_ubvec[0], 00213 &itr, 00214 &_options[0], 00215 &_edgecut, 00216 _part.empty() ? NULL : &_part[0], 00217 &mpi_comm); 00218 00219 // Assign the returned processor ids 00220 this->assign_partitioning (mesh); 00221 00222 00223 STOP_LOG ("repartition()", "ParmetisPartitioner"); 00224 00225 #endif // #ifndef LIBMESH_HAVE_PARMETIS ... else ... 00226 00227 } 00228 00229 00230 00231 // Only need to compile these methods if ParMETIS is present 00232 #ifdef LIBMESH_HAVE_PARMETIS 00233 00234 void ParmetisPartitioner::initialize (const MeshBase& mesh, 00235 const unsigned int n_sbdmns) 00236 { 00237 const dof_id_type n_active_local_elem = mesh.n_active_local_elem(); 00238 00239 // Set parameters. 00240 _wgtflag = 2; // weights on vertices only 00241 _ncon = 1; // one weight per vertex 00242 _numflag = 0; // C-style 0-based numbering 00243 _nparts = static_cast<int>(n_sbdmns); // number of subdomains to create 00244 _edgecut = 0; // the numbers of edges cut by the 00245 // partition 00246 00247 // Initialize data structures for ParMETIS 00248 _vtxdist.resize (libMesh::n_processors()+1); std::fill (_vtxdist.begin(), _vtxdist.end(), 0); 00249 _tpwgts.resize (_nparts); std::fill (_tpwgts.begin(), _tpwgts.end(), 1./_nparts); 00250 _ubvec.resize (_ncon); std::fill (_ubvec.begin(), _ubvec.end(), 1.05); 00251 _part.resize (n_active_local_elem); std::fill (_part.begin(), _part.end(), 0); 00252 _options.resize (5); 00253 _vwgt.resize (n_active_local_elem); 00254 00255 // Set the options 00256 _options[0] = 1; // don't use default options 00257 _options[1] = 0; // default (level of timing) 00258 _options[2] = 15; // random seed (default) 00259 _options[3] = 2; // processor distribution and subdomain distribution are decoupled 00260 00261 // Find the number of active elements on each processor. We cannot use 00262 // mesh.n_active_elem_on_proc(pid) since that only returns the number of 00263 // elements assigned to pid which are currently stored on the calling 00264 // processor. This will not in general be correct for parallel meshes 00265 // when (pid!=libMesh::processor_id()). 00266 _n_active_elem_on_proc.resize(libMesh::n_processors()); 00267 CommWorld.allgather(n_active_local_elem, _n_active_elem_on_proc); 00268 00269 // count the total number of active elements in the mesh. Note we cannot 00270 // use mesh.n_active_elem() in general since this only returns the number 00271 // of active elements which are stored on the calling processor. 00272 // We should not use n_active_elem for any allocation because that will 00273 // be inheritly unscalable, but it can be useful for libmesh_assertions. 00274 dof_id_type n_active_elem=0; 00275 00276 // Set up the vtxdist array. This will be the same on each processor. 00277 // ***** Consult the Parmetis documentation. ***** 00278 libmesh_assert_equal_to (_vtxdist.size(), 00279 libmesh_cast_int<std::size_t>(libMesh::n_processors()+1)); 00280 libmesh_assert_equal_to (_vtxdist[0], 0); 00281 00282 for (processor_id_type pid=0; pid<libMesh::n_processors(); pid++) 00283 { 00284 _vtxdist[pid+1] = _vtxdist[pid] + _n_active_elem_on_proc[pid]; 00285 n_active_elem += _n_active_elem_on_proc[pid]; 00286 } 00287 libmesh_assert_equal_to (_vtxdist.back(), static_cast<int>(n_active_elem)); 00288 00289 // ParMetis expects the elements to be numbered in contiguous blocks 00290 // by processor, i.e. [0, ne0), [ne0, ne0+ne1), ... 00291 // Since we only partition active elements we should have no expectation 00292 // that we currently have such a distribution. So we need to create it. 00293 // Also, at the same time we are going to map all the active elements into a globally 00294 // unique range [0,n_active_elem) which is *independent* of the current partitioning. 00295 // This can be fed to ParMetis as the initial partitioning of the subdomains (decoupled 00296 // from the partitioning of the objects themselves). This allows us to get the same 00297 // resultant partitioning independed of the input partitioning. 00298 MeshTools::BoundingBox bbox = 00299 MeshTools::bounding_box(mesh); 00300 00301 _global_index_by_pid_map.clear(); 00302 00303 // Maps active element ids into a contiguous range independent of partitioning. 00304 // (only needs local scope) 00305 std::map<dof_id_type, dof_id_type> global_index_map; 00306 00307 { 00308 std::vector<dof_id_type> global_index; 00309 00310 // create the mapping which is contiguous by processor 00311 dof_id_type pid_offset=0; 00312 for (processor_id_type pid=0; pid<libMesh::n_processors(); pid++) 00313 { 00314 MeshBase::const_element_iterator it = mesh.active_pid_elements_begin(pid); 00315 const MeshBase::const_element_iterator end = mesh.active_pid_elements_end(pid); 00316 00317 // note that we may not have all (or any!) the active elements which belong on this processor, 00318 // but by calling this on all processors a unique range in [0,_n_active_elem_on_proc[pid]) 00319 // is constructed. Only the indices for the elements we pass in are returned in the array. 00320 MeshCommunication().find_global_indices (bbox, it, end, 00321 global_index); 00322 00323 for (dof_id_type cnt=0; it != end; ++it) 00324 { 00325 const Elem *elem = *it; 00326 libmesh_assert (!_global_index_by_pid_map.count(elem->id())); 00327 libmesh_assert_less (cnt, global_index.size()); 00328 libmesh_assert_less (global_index[cnt], _n_active_elem_on_proc[pid]); 00329 00330 _global_index_by_pid_map[elem->id()] = global_index[cnt++] + pid_offset; 00331 } 00332 00333 pid_offset += _n_active_elem_on_proc[pid]; 00334 } 00335 00336 // create the unique mapping for all active elements independent of partitioning 00337 { 00338 MeshBase::const_element_iterator it = mesh.active_elements_begin(); 00339 const MeshBase::const_element_iterator end = mesh.active_elements_end(); 00340 00341 // Calling this on all processors a unique range in [0,n_active_elem) is constructed. 00342 // Only the indices for the elements we pass in are returned in the array. 00343 MeshCommunication().find_global_indices (bbox, it, end, 00344 global_index); 00345 00346 for (dof_id_type cnt=0; it != end; ++it) 00347 { 00348 const Elem *elem = *it; 00349 libmesh_assert (!global_index_map.count(elem->id())); 00350 libmesh_assert_less (cnt, global_index.size()); 00351 libmesh_assert_less (global_index[cnt], n_active_elem); 00352 00353 global_index_map[elem->id()] = global_index[cnt++]; 00354 } 00355 } 00356 // really, shouldn't be close! 00357 libmesh_assert_less_equal (global_index_map.size(), n_active_elem); 00358 libmesh_assert_less_equal (_global_index_by_pid_map.size(), n_active_elem); 00359 00360 // At this point the two maps should be the same size. If they are not 00361 // then the number of active elements is not the same as the sum over all 00362 // processors of the number of active elements per processor, which means 00363 // there must be some unpartitioned objects out there. 00364 if (global_index_map.size() != _global_index_by_pid_map.size()) 00365 { 00366 libMesh::err << "ERROR: ParmetisPartitioner cannot handle unpartitioned objects!" 00367 << std::endl; 00368 libmesh_error(); 00369 } 00370 } 00371 00372 // Finally, we need to initialize the vertex (partition) weights and the initial subdomain 00373 // mapping. The subdomain mapping will be independent of the processor mapping, and is 00374 // defined by a simple mapping of the global indices we just found. 00375 { 00376 std::vector<dof_id_type> subdomain_bounds(libMesh::n_processors()); 00377 00378 const dof_id_type first_local_elem = _vtxdist[libMesh::processor_id()]; 00379 00380 for (processor_id_type pid=0; pid<libMesh::n_processors(); pid++) 00381 { 00382 dof_id_type tgt_subdomain_size = 0; 00383 00384 // watch out for the case that n_subdomains < n_processors 00385 if (pid < static_cast<unsigned int>(_nparts)) 00386 { 00387 tgt_subdomain_size = n_active_elem/std::min 00388 (libmesh_cast_int<int>(libMesh::n_processors()), 00389 _nparts); 00390 00391 if (pid < n_active_elem%_nparts) 00392 tgt_subdomain_size++; 00393 } 00394 if (pid == 0) 00395 subdomain_bounds[0] = tgt_subdomain_size; 00396 else 00397 subdomain_bounds[pid] = subdomain_bounds[pid-1] + tgt_subdomain_size; 00398 } 00399 00400 libmesh_assert_equal_to (subdomain_bounds.back(), n_active_elem); 00401 00402 MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin(); 00403 const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end(); 00404 00405 for (; elem_it != elem_end; ++elem_it) 00406 { 00407 const Elem *elem = *elem_it; 00408 00409 libmesh_assert (_global_index_by_pid_map.count(elem->id())); 00410 const dof_id_type global_index_by_pid = 00411 _global_index_by_pid_map[elem->id()]; 00412 libmesh_assert_less (global_index_by_pid, n_active_elem); 00413 00414 const dof_id_type local_index = 00415 global_index_by_pid - first_local_elem; 00416 00417 libmesh_assert_less (local_index, n_active_local_elem); 00418 libmesh_assert_less (local_index, _vwgt.size()); 00419 00420 // TODO:[BSK] maybe there is a better weight? 00421 _vwgt[local_index] = elem->n_nodes(); 00422 00423 // find the subdomain this element belongs in 00424 libmesh_assert (global_index_map.count(elem->id())); 00425 const dof_id_type global_index = 00426 global_index_map[elem->id()]; 00427 00428 libmesh_assert_less (global_index, subdomain_bounds.back()); 00429 00430 const unsigned int subdomain_id = 00431 std::distance(subdomain_bounds.begin(), 00432 std::lower_bound(subdomain_bounds.begin(), 00433 subdomain_bounds.end(), 00434 global_index)); 00435 libmesh_assert_less (subdomain_id, static_cast<unsigned int>(_nparts)); 00436 libmesh_assert_less (local_index, _part.size()); 00437 00438 _part[local_index] = subdomain_id; 00439 } 00440 } 00441 } 00442 00443 00444 00445 void ParmetisPartitioner::build_graph (const MeshBase& mesh) 00446 { 00447 // build the graph in distributed CSR format. Note that 00448 // the edges in the graph will correspond to 00449 // face neighbors 00450 const dof_id_type n_active_local_elem = mesh.n_active_local_elem(); 00451 00452 std::vector<const Elem*> neighbors_offspring; 00453 00454 std::vector<std::vector<dof_id_type> > graph(n_active_local_elem); 00455 dof_id_type graph_size=0; 00456 00457 const dof_id_type first_local_elem = _vtxdist[libMesh::processor_id()]; 00458 00459 MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin(); 00460 const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end(); 00461 00462 for (; elem_it != elem_end; ++elem_it) 00463 { 00464 const Elem* elem = *elem_it; 00465 00466 libmesh_assert (_global_index_by_pid_map.count(elem->id())); 00467 const dof_id_type global_index_by_pid = 00468 _global_index_by_pid_map[elem->id()]; 00469 00470 const dof_id_type local_index = 00471 global_index_by_pid - first_local_elem; 00472 libmesh_assert_less (local_index, n_active_local_elem); 00473 00474 std::vector<dof_id_type> &graph_row = graph[local_index]; 00475 00476 // Loop over the element's neighbors. An element 00477 // adjacency corresponds to a face neighbor 00478 for (unsigned int ms=0; ms<elem->n_neighbors(); ms++) 00479 { 00480 const Elem* neighbor = elem->neighbor(ms); 00481 00482 if (neighbor != NULL) 00483 { 00484 // If the neighbor is active treat it 00485 // as a connection 00486 if (neighbor->active()) 00487 { 00488 libmesh_assert(_global_index_by_pid_map.count(neighbor->id())); 00489 const dof_id_type neighbor_global_index_by_pid = 00490 _global_index_by_pid_map[neighbor->id()]; 00491 00492 graph_row.push_back(neighbor_global_index_by_pid); 00493 graph_size++; 00494 } 00495 00496 #ifdef LIBMESH_ENABLE_AMR 00497 00498 // Otherwise we need to find all of the 00499 // neighbor's children that are connected to 00500 // us and add them 00501 else 00502 { 00503 // The side of the neighbor to which 00504 // we are connected 00505 const unsigned int ns = 00506 neighbor->which_neighbor_am_i (elem); 00507 libmesh_assert_less (ns, neighbor->n_neighbors()); 00508 00509 // Get all the active children (& grandchildren, etc...) 00510 // of the neighbor. 00511 neighbor->active_family_tree (neighbors_offspring); 00512 00513 // Get all the neighbor's children that 00514 // live on that side and are thus connected 00515 // to us 00516 for (unsigned int nc=0; nc<neighbors_offspring.size(); nc++) 00517 { 00518 const Elem* child = 00519 neighbors_offspring[nc]; 00520 00521 // This does not assume a level-1 mesh. 00522 // Note that since children have sides numbered 00523 // coincident with the parent then this is a sufficient test. 00524 if (child->neighbor(ns) == elem) 00525 { 00526 libmesh_assert (child->active()); 00527 libmesh_assert (_global_index_by_pid_map.count(child->id())); 00528 const dof_id_type child_global_index_by_pid = 00529 _global_index_by_pid_map[child->id()]; 00530 00531 graph_row.push_back(child_global_index_by_pid); 00532 graph_size++; 00533 } 00534 } 00535 } 00536 00537 #endif /* ifdef LIBMESH_ENABLE_AMR */ 00538 00539 00540 } 00541 } 00542 } 00543 00544 // Reserve space in the adjacency array 00545 _xadj.clear(); 00546 _xadj.reserve (n_active_local_elem + 1); 00547 _adjncy.clear(); 00548 _adjncy.reserve (graph_size); 00549 00550 for (std::size_t r=0; r<graph.size(); r++) 00551 { 00552 _xadj.push_back(_adjncy.size()); 00553 std::vector<dof_id_type> graph_row; // build this emtpy 00554 graph_row.swap(graph[r]); // this will deallocate at the end of scope 00555 _adjncy.insert(_adjncy.end(), 00556 graph_row.begin(), 00557 graph_row.end()); 00558 } 00559 00560 // The end of the adjacency array for the last elem 00561 _xadj.push_back(_adjncy.size()); 00562 00563 libmesh_assert_equal_to (_xadj.size(), n_active_local_elem+1); 00564 libmesh_assert_equal_to (_adjncy.size(), graph_size); 00565 } 00566 00567 00568 00569 void ParmetisPartitioner::assign_partitioning (MeshBase& mesh) 00570 { 00571 // This function must be run on all processors at once 00572 parallel_only(); 00573 00574 const dof_id_type 00575 first_local_elem = _vtxdist[libMesh::processor_id()]; 00576 00577 std::vector<std::vector<dof_id_type> > 00578 requested_ids(libMesh::n_processors()), 00579 requests_to_fill(libMesh::n_processors()); 00580 00581 MeshBase::element_iterator elem_it = mesh.active_elements_begin(); 00582 MeshBase::element_iterator elem_end = mesh.active_elements_end(); 00583 00584 for (; elem_it != elem_end; ++elem_it) 00585 { 00586 Elem *elem = *elem_it; 00587 00588 // we need to get the index from the owning processor 00589 // (note we cannot assign it now -- we are iterating 00590 // over elements again and this will be bad!) 00591 libmesh_assert_less (elem->processor_id(), requested_ids.size()); 00592 requested_ids[elem->processor_id()].push_back(elem->id()); 00593 } 00594 00595 // Trade with all processors (including self) to get their indices 00596 for (processor_id_type pid=0; pid<libMesh::n_processors(); pid++) 00597 { 00598 // Trade my requests with processor procup and procdown 00599 const processor_id_type procup = (libMesh::processor_id() + pid) % 00600 libMesh::n_processors(); 00601 const processor_id_type procdown = (libMesh::n_processors() + 00602 libMesh::processor_id() - pid) % 00603 libMesh::n_processors(); 00604 00605 CommWorld.send_receive (procup, requested_ids[procup], 00606 procdown, requests_to_fill[procdown]); 00607 00608 // we can overwrite these requested ids in-place. 00609 for (std::size_t i=0; i<requests_to_fill[procdown].size(); i++) 00610 { 00611 const dof_id_type requested_elem_index = 00612 requests_to_fill[procdown][i]; 00613 00614 libmesh_assert(_global_index_by_pid_map.count(requested_elem_index)); 00615 00616 const dof_id_type global_index_by_pid = 00617 _global_index_by_pid_map[requested_elem_index]; 00618 00619 const dof_id_type local_index = 00620 global_index_by_pid - first_local_elem; 00621 00622 libmesh_assert_less (local_index, _part.size()); 00623 libmesh_assert_less (local_index, mesh.n_active_local_elem()); 00624 00625 const unsigned int elem_procid = 00626 static_cast<unsigned int>(_part[local_index]); 00627 00628 libmesh_assert_less (elem_procid, static_cast<unsigned int>(_nparts)); 00629 00630 requests_to_fill[procdown][i] = elem_procid; 00631 } 00632 00633 // Trade back 00634 CommWorld.send_receive (procdown, requests_to_fill[procdown], 00635 procup, requested_ids[procup]); 00636 } 00637 00638 // and finally assign the partitioning. 00639 // note we are iterating in exactly the same order 00640 // used to build up the request, so we can expect the 00641 // required entries to be in the proper sequence. 00642 elem_it = mesh.active_elements_begin(); 00643 elem_end = mesh.active_elements_end(); 00644 00645 for (std::vector<unsigned int> counters(libMesh::n_processors(), 0); 00646 elem_it != elem_end; ++elem_it) 00647 { 00648 Elem *elem = *elem_it; 00649 00650 const processor_id_type current_pid = elem->processor_id(); 00651 00652 libmesh_assert_less (counters[current_pid], requested_ids[current_pid].size()); 00653 00654 const processor_id_type elem_procid = 00655 requested_ids[current_pid][counters[current_pid]++]; 00656 00657 libmesh_assert_less (elem_procid, static_cast<unsigned int>(_nparts)); 00658 elem->processor_id() = elem_procid; 00659 } 00660 } 00661 00662 #endif // #ifdef LIBMESH_HAVE_PARMETIS 00663 00664 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: