libMesh::ParmetisPartitioner Class Reference
#include <parmetis_partitioner.h>

Public Member Functions | |
| ParmetisPartitioner () | |
| virtual AutoPtr< Partitioner > | clone () const |
| void | partition (MeshBase &mesh, const unsigned int n=libMesh::n_processors()) |
| void | repartition (MeshBase &mesh, const unsigned int n=libMesh::n_processors()) |
| virtual void | attach_weights (ErrorVector *) |
Static Public Member Functions | |
| static void | partition_unpartitioned_elements (MeshBase &mesh, const unsigned int n=libMesh::n_processors()) |
| static void | set_parent_processor_ids (MeshBase &mesh) |
| static void | set_node_processor_ids (MeshBase &mesh) |
Protected Member Functions | |
| virtual void | _do_repartition (MeshBase &mesh, const unsigned int n) |
| virtual void | _do_partition (MeshBase &mesh, const unsigned int n) |
| void | single_partition (MeshBase &mesh) |
Protected Attributes | |
| ErrorVector * | _weights |
Static Protected Attributes | |
| static const dof_id_type | communication_blocksize = 1000000 |
Private Member Functions | |
| void | initialize (const MeshBase &mesh, const unsigned int n_sbdmns) |
| void | build_graph (const MeshBase &mesh) |
| void | assign_partitioning (MeshBase &mesh) |
Private Attributes | |
| std::vector< dof_id_type > | _n_active_elem_on_proc |
| std::map< dof_id_type, dof_id_type > | _global_index_by_pid_map |
| std::vector< int > | _vtxdist |
| std::vector< int > | _xadj |
| std::vector< int > | _adjncy |
| std::vector< int > | _part |
| std::vector< float > | _tpwgts |
| std::vector< float > | _ubvec |
| std::vector< int > | _options |
| std::vector< int > | _vwgt |
| int | _wgtflag |
| int | _ncon |
| int | _numflag |
| int | _nparts |
| int | _edgecut |
Detailed Description
The ParmetisPartitioner uses the Parmetis graph partitioner to partition the elements.
Definition at line 44 of file parmetis_partitioner.h.
Constructor & Destructor Documentation
| libMesh::ParmetisPartitioner::ParmetisPartitioner | ( | ) | [inline] |
Member Function Documentation
| void libMesh::ParmetisPartitioner::_do_partition | ( | MeshBase & | mesh, | |
| const unsigned int | n | |||
| ) | [protected, virtual] |
Partition the MeshBase into n subdomains.
Implements libMesh::Partitioner.
Definition at line 57 of file parmetis_partitioner.C.
References _do_repartition().
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 }
| void libMesh::ParmetisPartitioner::_do_repartition | ( | MeshBase & | mesh, | |
| const unsigned int | n | |||
| ) | [protected, virtual] |
Parmetis can handle dynamically repartitioning a mesh such that the redistribution costs are minimized. This method takes a previously partitioned domain (which may have then been adaptively refined) and repartitions it.
Reimplemented from libMesh::Partitioner.
Definition at line 120 of file parmetis_partitioner.C.
References _adjncy, _edgecut, _n_active_elem_on_proc, _ncon, _nparts, _numflag, _options, _part, _tpwgts, _ubvec, _vtxdist, _vwgt, _wgtflag, _xadj, assign_partitioning(), build_graph(), libMesh::COMM_WORLD, libMesh::err, initialize(), libMesh::MIN_ELEM_PER_PROC, libMesh::n_processors(), libMesh::Partitioner::partition(), and libMesh::Partitioner::single_partition().
Referenced by _do_partition().
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 }
| void libMesh::ParmetisPartitioner::assign_partitioning | ( | MeshBase & | mesh | ) | [private] |
Assign the computed partitioning to the mesh.
Definition at line 569 of file parmetis_partitioner.C.
References _global_index_by_pid_map, _nparts, _part, _vtxdist, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::CommWorld, libMesh::DofObject::id(), libMesh::MeshBase::n_active_local_elem(), libMesh::n_processors(), libMesh::DofObject::processor_id(), libMesh::processor_id(), and libMesh::Parallel::Communicator::send_receive().
Referenced by _do_repartition().
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 }
| virtual void libMesh::Partitioner::attach_weights | ( | ErrorVector * | ) | [inline, virtual, inherited] |
Attach weights that can be used for partitioning. This ErrorVector should be _exactly_ the same on every processor and should have mesh->max_elem_id() entries.
Reimplemented in libMesh::MetisPartitioner.
Definition at line 118 of file partitioner.h.
| void libMesh::ParmetisPartitioner::build_graph | ( | const MeshBase & | mesh | ) | [private] |
Build the graph.
Definition at line 445 of file parmetis_partitioner.C.
References _adjncy, _global_index_by_pid_map, _vtxdist, _xadj, libMesh::Elem::active(), libMesh::Elem::active_family_tree(), libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::DofObject::id(), libMesh::MeshBase::n_active_local_elem(), libMesh::Elem::n_neighbors(), libMesh::Elem::neighbor(), libMesh::processor_id(), and libMesh::Elem::which_neighbor_am_i().
Referenced by _do_repartition().
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 }
| virtual AutoPtr<Partitioner> libMesh::ParmetisPartitioner::clone | ( | ) | const [inline, virtual] |
Creates a new partitioner of this type and returns it in an AutoPtr.
Implements libMesh::Partitioner.
Definition at line 57 of file parmetis_partitioner.h.
References ParmetisPartitioner().
00057 { 00058 AutoPtr<Partitioner> cloned_partitioner 00059 (new ParmetisPartitioner()); 00060 return cloned_partitioner; 00061 }
| void libMesh::ParmetisPartitioner::initialize | ( | const MeshBase & | mesh, | |
| const unsigned int | n_sbdmns | |||
| ) | [private] |
Initialize data structures.
Definition at line 234 of file parmetis_partitioner.C.
References _edgecut, _global_index_by_pid_map, _n_active_elem_on_proc, _ncon, _nparts, _numflag, _options, _part, _tpwgts, _ubvec, _vtxdist, _vwgt, _wgtflag, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::MeshBase::active_pid_elements_begin(), libMesh::MeshBase::active_pid_elements_end(), libMesh::Parallel::Communicator::allgather(), libMesh::MeshTools::bounding_box(), libMesh::CommWorld, end, libMesh::err, libMesh::DofObject::id(), std::min(), libMesh::MeshBase::n_active_local_elem(), libMesh::Elem::n_nodes(), libMesh::n_processors(), and libMesh::processor_id().
Referenced by _do_repartition().
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 }
| void libMesh::Partitioner::partition | ( | MeshBase & | mesh, | |
| const unsigned int | n = libMesh::n_processors() | |||
| ) | [inherited] |
Partition the MeshBase into n parts. If the user does not specify a number of pieces into which the mesh should be partitioned, then the default behavior of the partitioner is to partition according to the number of processors defined in libMesh::n_processors(). The partitioner currently does not modify the subdomain_id of each element. This number is reserved for things like material properties, etc.
Definition at line 48 of file partitioner.C.
References libMesh::Partitioner::_do_partition(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_remote_elems(), std::min(), libMesh::MeshBase::n_active_elem(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::MeshBase::redistribute(), libMesh::MeshBase::set_n_partitions(), libMesh::Partitioner::set_node_processor_ids(), libMesh::Partitioner::set_parent_processor_ids(), libMesh::Partitioner::single_partition(), and libMesh::MeshBase::update_post_partitioning().
Referenced by libMesh::SFCPartitioner::_do_partition(), libMesh::MetisPartitioner::_do_partition(), and _do_repartition().
00050 { 00051 parallel_only(); 00052 00053 // BSK - temporary fix while redistribution is integrated 6/26/2008 00054 // Uncomment this to not repartition in parallel 00055 // if (!mesh.is_serial()) 00056 // return; 00057 00058 // we cannot partition into more pieces than we have 00059 // active elements! 00060 const unsigned int n_parts = 00061 static_cast<unsigned int> 00062 (std::min(mesh.n_active_elem(), static_cast<dof_id_type>(n))); 00063 00064 // Set the number of partitions in the mesh 00065 mesh.set_n_partitions()=n_parts; 00066 00067 if (n_parts == 1) 00068 { 00069 this->single_partition (mesh); 00070 return; 00071 } 00072 00073 // First assign a temporary partitioning to any unpartitioned elements 00074 Partitioner::partition_unpartitioned_elements(mesh, n_parts); 00075 00076 // Call the partitioning function 00077 this->_do_partition(mesh,n_parts); 00078 00079 // Set the parent's processor ids 00080 Partitioner::set_parent_processor_ids(mesh); 00081 00082 // Redistribute elements if necessary, before setting node processor 00083 // ids, to make sure those will be set consistently 00084 mesh.redistribute(); 00085 00086 #ifdef DEBUG 00087 MeshTools::libmesh_assert_valid_remote_elems(mesh); 00088 00089 // Messed up elem processor_id()s can leave us without the child 00090 // elements we need to restrict vectors on a distributed mesh 00091 MeshTools::libmesh_assert_valid_procids<Elem>(mesh); 00092 #endif 00093 00094 // Set the node's processor ids 00095 Partitioner::set_node_processor_ids(mesh); 00096 00097 #ifdef DEBUG 00098 MeshTools::libmesh_assert_valid_procids<Elem>(mesh); 00099 #endif 00100 00101 // Give derived Mesh classes a chance to update any cached data to 00102 // reflect the new partitioning 00103 mesh.update_post_partitioning(); 00104 }
| void libMesh::Partitioner::partition_unpartitioned_elements | ( | MeshBase & | mesh, | |
| const unsigned int | n = libMesh::n_processors() | |||
| ) | [static, inherited] |
This function
Definition at line 168 of file partitioner.C.
References libMesh::MeshTools::bounding_box(), end, libMesh::MeshTools::n_elem(), libMesh::n_processors(), libMesh::DofObject::processor_id(), libMesh::MeshBase::unpartitioned_elements_begin(), and libMesh::MeshBase::unpartitioned_elements_end().
Referenced by libMesh::Partitioner::partition(), and libMesh::Partitioner::repartition().
00170 { 00171 MeshBase::element_iterator it = mesh.unpartitioned_elements_begin(); 00172 const MeshBase::element_iterator end = mesh.unpartitioned_elements_end(); 00173 00174 const dof_id_type n_unpartitioned_elements = MeshTools::n_elem (it, end); 00175 00176 // the unpartitioned elements must exist on all processors. If the range is empty on one 00177 // it is empty on all, and we can quit right here. 00178 if (!n_unpartitioned_elements) return; 00179 00180 // find the target subdomain sizes 00181 std::vector<dof_id_type> subdomain_bounds(libMesh::n_processors()); 00182 00183 for (processor_id_type pid=0; pid<libMesh::n_processors(); pid++) 00184 { 00185 dof_id_type tgt_subdomain_size = 0; 00186 00187 // watch out for the case that n_subdomains < n_processors 00188 if (pid < n_subdomains) 00189 { 00190 tgt_subdomain_size = n_unpartitioned_elements/n_subdomains; 00191 00192 if (pid < n_unpartitioned_elements%n_subdomains) 00193 tgt_subdomain_size++; 00194 00195 } 00196 00197 //libMesh::out << "pid, #= " << pid << ", " << tgt_subdomain_size << std::endl; 00198 if (pid == 0) 00199 subdomain_bounds[0] = tgt_subdomain_size; 00200 else 00201 subdomain_bounds[pid] = subdomain_bounds[pid-1] + tgt_subdomain_size; 00202 } 00203 00204 libmesh_assert_equal_to (subdomain_bounds.back(), n_unpartitioned_elements); 00205 00206 // create the unique mapping for all unpartitioned elements independent of partitioning 00207 // determine the global indexing for all the unpartitoned elements 00208 std::vector<dof_id_type> global_indices; 00209 00210 // Calling this on all processors a unique range in [0,n_unpartitioned_elements) is constructed. 00211 // Only the indices for the elements we pass in are returned in the array. 00212 MeshCommunication().find_global_indices (MeshTools::bounding_box(mesh), it, end, 00213 global_indices); 00214 00215 for (dof_id_type cnt=0; it != end; ++it) 00216 { 00217 Elem *elem = *it; 00218 00219 libmesh_assert_less (cnt, global_indices.size()); 00220 const dof_id_type global_index = 00221 global_indices[cnt++]; 00222 00223 libmesh_assert_less (global_index, subdomain_bounds.back()); 00224 libmesh_assert_less (global_index, n_unpartitioned_elements); 00225 00226 const processor_id_type subdomain_id = 00227 libmesh_cast_int<processor_id_type> 00228 (std::distance(subdomain_bounds.begin(), 00229 std::upper_bound(subdomain_bounds.begin(), 00230 subdomain_bounds.end(), 00231 global_index))); 00232 libmesh_assert_less (subdomain_id, n_subdomains); 00233 00234 elem->processor_id() = subdomain_id; 00235 //libMesh::out << "assigning " << global_index << " to " << subdomain_id << std::endl; 00236 } 00237 }
| void libMesh::Partitioner::repartition | ( | MeshBase & | mesh, | |
| const unsigned int | n = libMesh::n_processors() | |||
| ) | [inherited] |
Repartitions the MeshBase into n parts. This is required since some partitoning algorithms can repartition more efficiently than computing a new partitioning from scratch. The default behavior is to simply call this->partition(n)
Definition at line 110 of file partitioner.C.
References libMesh::Partitioner::_do_repartition(), std::min(), libMesh::MeshBase::n_active_elem(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::MeshBase::set_n_partitions(), libMesh::Partitioner::set_node_processor_ids(), libMesh::Partitioner::set_parent_processor_ids(), and libMesh::Partitioner::single_partition().
00112 { 00113 // we cannot partition into more pieces than we have 00114 // active elements! 00115 const unsigned int n_parts = 00116 static_cast<unsigned int> 00117 (std::min(mesh.n_active_elem(), static_cast<dof_id_type>(n))); 00118 00119 // Set the number of partitions in the mesh 00120 mesh.set_n_partitions()=n_parts; 00121 00122 if (n_parts == 1) 00123 { 00124 this->single_partition (mesh); 00125 return; 00126 } 00127 00128 // First assign a temporary partitioning to any unpartitioned elements 00129 Partitioner::partition_unpartitioned_elements(mesh, n_parts); 00130 00131 // Call the partitioning function 00132 this->_do_repartition(mesh,n_parts); 00133 00134 // Set the parent's processor ids 00135 Partitioner::set_parent_processor_ids(mesh); 00136 00137 // Set the node's processor ids 00138 Partitioner::set_node_processor_ids(mesh); 00139 }
| void libMesh::Partitioner::set_node_processor_ids | ( | MeshBase & | mesh | ) | [static, inherited] |
This function is called after partitioning to set the processor IDs for the nodes. By definition, a Node's processor ID is the minimum processor ID for all of the elements which share the node.
Definition at line 419 of file partitioner.C.
References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::CommWorld, libMesh::Elem::get_node(), libMesh::DofObject::id(), libMesh::DofObject::invalid_processor_id, libMesh::DofObject::invalidate_processor_id(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), std::min(), libMesh::MeshTools::n_elem(), libMesh::Elem::n_nodes(), libMesh::MeshBase::n_partitions(), libMesh::n_processors(), libMesh::MeshBase::node_ptr(), libMesh::MeshBase::nodes_begin(), libMesh::MeshBase::nodes_end(), libMesh::MeshBase::not_active_elements_begin(), libMesh::MeshBase::not_active_elements_end(), libMesh::processor_id(), libMesh::DofObject::processor_id(), libMesh::Parallel::Communicator::send_receive(), libMesh::MeshBase::subactive_elements_begin(), libMesh::MeshBase::subactive_elements_end(), libMesh::MeshBase::unpartitioned_elements_begin(), and libMesh::MeshBase::unpartitioned_elements_end().
Referenced by libMesh::UnstructuredMesh::all_first_order(), libMesh::Partitioner::partition(), libMesh::XdrIO::read(), libMesh::Partitioner::repartition(), and libMesh::BoundaryInfo::sync().
00420 { 00421 START_LOG("set_node_processor_ids()","Partitioner"); 00422 00423 // This function must be run on all processors at once 00424 parallel_only(); 00425 00426 // If we have any unpartitioned elements at this 00427 // stage there is a problem 00428 libmesh_assert (MeshTools::n_elem(mesh.unpartitioned_elements_begin(), 00429 mesh.unpartitioned_elements_end()) == 0); 00430 00431 00432 // const dof_id_type orig_n_local_nodes = mesh.n_local_nodes(); 00433 00434 // libMesh::err << "[" << libMesh::processor_id() << "]: orig_n_local_nodes=" 00435 // << orig_n_local_nodes << std::endl; 00436 00437 // Build up request sets. Each node is currently owned by a processor because 00438 // it is connected to an element owned by that processor. However, during the 00439 // repartitioning phase that element may have been assigned a new processor id, but 00440 // it is still resident on the original processor. We need to know where to look 00441 // for new ids before assigning new ids, otherwise we may be asking the wrong processors 00442 // for the wrong information. 00443 // 00444 // The only remaining issue is what to do with unpartitioned nodes. Since they are required 00445 // to live on all processors we can simply rely on ourselves to number them properly. 00446 std::vector<std::vector<dof_id_type> > 00447 requested_node_ids(libMesh::n_processors()); 00448 00449 // Loop over all the nodes, count the ones on each processor. We can skip ourself 00450 std::vector<dof_id_type> ghost_nodes_from_proc(libMesh::n_processors(), 0); 00451 00452 MeshBase::node_iterator node_it = mesh.nodes_begin(); 00453 const MeshBase::node_iterator node_end = mesh.nodes_end(); 00454 00455 for (; node_it != node_end; ++node_it) 00456 { 00457 Node *node = *node_it; 00458 libmesh_assert(node); 00459 const processor_id_type current_pid = node->processor_id(); 00460 if (current_pid != libMesh::processor_id() && 00461 current_pid != DofObject::invalid_processor_id) 00462 { 00463 libmesh_assert_less (current_pid, ghost_nodes_from_proc.size()); 00464 ghost_nodes_from_proc[current_pid]++; 00465 } 00466 } 00467 00468 // We know how many objects live on each processor, so reserve() 00469 // space for each. 00470 for (processor_id_type pid=0; pid != libMesh::n_processors(); ++pid) 00471 requested_node_ids[pid].reserve(ghost_nodes_from_proc[pid]); 00472 00473 // We need to get the new pid for each node from the processor 00474 // which *currently* owns the node. We can safely skip ourself 00475 for (node_it = mesh.nodes_begin(); node_it != node_end; ++node_it) 00476 { 00477 Node *node = *node_it; 00478 libmesh_assert(node); 00479 const processor_id_type current_pid = node->processor_id(); 00480 if (current_pid != libMesh::processor_id() && 00481 current_pid != DofObject::invalid_processor_id) 00482 { 00483 libmesh_assert_less (current_pid, requested_node_ids.size()); 00484 libmesh_assert_less (requested_node_ids[current_pid].size(), 00485 ghost_nodes_from_proc[current_pid]); 00486 requested_node_ids[current_pid].push_back(node->id()); 00487 } 00488 00489 // Unset any previously-set node processor ids 00490 node->invalidate_processor_id(); 00491 } 00492 00493 // Loop over all the active elements 00494 MeshBase::element_iterator elem_it = mesh.active_elements_begin(); 00495 const MeshBase::element_iterator elem_end = mesh.active_elements_end(); 00496 00497 for ( ; elem_it != elem_end; ++elem_it) 00498 { 00499 Elem* elem = *elem_it; 00500 libmesh_assert(elem); 00501 00502 libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id); 00503 00504 // For each node, set the processor ID to the min of 00505 // its current value and this Element's processor id. 00506 // 00507 // TODO: we would probably get better parallel partitioning if 00508 // we did something like "min for even numbered nodes, max for 00509 // odd numbered". We'd need to be careful about how that would 00510 // affect solution ordering for I/O, though. 00511 for (unsigned int n=0; n<elem->n_nodes(); ++n) 00512 elem->get_node(n)->processor_id() = std::min(elem->get_node(n)->processor_id(), 00513 elem->processor_id()); 00514 } 00515 00516 // And loop over the subactive elements, but don't reassign 00517 // nodes that are already active on another processor. 00518 MeshBase::element_iterator sub_it = mesh.subactive_elements_begin(); 00519 const MeshBase::element_iterator sub_end = mesh.subactive_elements_end(); 00520 00521 for ( ; sub_it != sub_end; ++sub_it) 00522 { 00523 Elem* elem = *sub_it; 00524 libmesh_assert(elem); 00525 00526 libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id); 00527 00528 for (unsigned int n=0; n<elem->n_nodes(); ++n) 00529 if (elem->get_node(n)->processor_id() == DofObject::invalid_processor_id) 00530 elem->get_node(n)->processor_id() = elem->processor_id(); 00531 } 00532 00533 // Same for the inactive elements -- we will have already gotten most of these 00534 // nodes, *except* for the case of a parent with a subset of children which are 00535 // ghost elements. In that case some of the parent nodes will not have been 00536 // properly handled yet 00537 MeshBase::element_iterator not_it = mesh.not_active_elements_begin(); 00538 const MeshBase::element_iterator not_end = mesh.not_active_elements_end(); 00539 00540 for ( ; not_it != not_end; ++not_it) 00541 { 00542 Elem* elem = *not_it; 00543 libmesh_assert(elem); 00544 00545 libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id); 00546 00547 for (unsigned int n=0; n<elem->n_nodes(); ++n) 00548 if (elem->get_node(n)->processor_id() == DofObject::invalid_processor_id) 00549 elem->get_node(n)->processor_id() = elem->processor_id(); 00550 } 00551 00552 // We can't assert that all nodes are connected to elements, because 00553 // a ParallelMesh with NodeConstraints might have pulled in some 00554 // remote nodes solely for evaluating those constraints. 00555 // MeshTools::libmesh_assert_connected_nodes(mesh); 00556 00557 // For such nodes, we'll do a sanity check later when making sure 00558 // that we successfully reset their processor ids to something 00559 // valid. 00560 00561 // Next set node ids from other processors, excluding self 00562 for (processor_id_type p=1; p != libMesh::n_processors(); ++p) 00563 { 00564 // Trade my requests with processor procup and procdown 00565 processor_id_type procup = (libMesh::processor_id() + p) % 00566 libMesh::n_processors(); 00567 processor_id_type procdown = (libMesh::n_processors() + 00568 libMesh::processor_id() - p) % 00569 libMesh::n_processors(); 00570 std::vector<dof_id_type> request_to_fill; 00571 CommWorld.send_receive(procup, requested_node_ids[procup], 00572 procdown, request_to_fill); 00573 00574 // Fill those requests in-place 00575 for (std::size_t i=0; i != request_to_fill.size(); ++i) 00576 { 00577 Node *node = mesh.node_ptr(request_to_fill[i]); 00578 libmesh_assert(node); 00579 const processor_id_type new_pid = node->processor_id(); 00580 libmesh_assert_not_equal_to (new_pid, DofObject::invalid_processor_id); 00581 libmesh_assert_less (new_pid, mesh.n_partitions()); // this is the correct test -- 00582 request_to_fill[i] = new_pid; // the number of partitions may 00583 } // not equal the number of processors 00584 00585 // Trade back the results 00586 std::vector<dof_id_type> filled_request; 00587 CommWorld.send_receive(procdown, request_to_fill, 00588 procup, filled_request); 00589 libmesh_assert_equal_to (filled_request.size(), requested_node_ids[procup].size()); 00590 00591 // And copy the id changes we've now been informed of 00592 for (std::size_t i=0; i != filled_request.size(); ++i) 00593 { 00594 Node *node = mesh.node_ptr(requested_node_ids[procup][i]); 00595 libmesh_assert(node); 00596 libmesh_assert_less (filled_request[i], mesh.n_partitions()); // this is the correct test -- 00597 node->processor_id(filled_request[i]); // the number of partitions may 00598 } // not equal the number of processors 00599 } 00600 00601 #ifdef DEBUG 00602 MeshTools::libmesh_assert_valid_procids<Node>(mesh); 00603 #endif 00604 00605 STOP_LOG("set_node_processor_ids()","Partitioner"); 00606 }
| static void libMesh::Partitioner::set_parent_processor_ids | ( | MeshBase & | mesh | ) | [static, inherited] |
This function is called after partitioning to set the processor IDs for the inactive parent elements. A Parent's processor ID is the same as its first child.
Referenced by libMesh::Partitioner::partition(), and libMesh::Partitioner::repartition().
| void libMesh::Partitioner::single_partition | ( | MeshBase & | mesh | ) | [protected, inherited] |
Trivially "partitions" the mesh for one processor. Simply loops through the elements and assigns all of them to processor 0. Is is provided as a separate function so that derived classes may use it without reimplementing it.
Definition at line 145 of file partitioner.C.
References libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), libMesh::MeshBase::nodes_begin(), and libMesh::MeshBase::nodes_end().
Referenced by libMesh::SFCPartitioner::_do_partition(), libMesh::MetisPartitioner::_do_partition(), libMesh::LinearPartitioner::_do_partition(), libMesh::CentroidPartitioner::_do_partition(), _do_repartition(), libMesh::Partitioner::partition(), and libMesh::Partitioner::repartition().
00146 { 00147 START_LOG("single_partition()","Partitioner"); 00148 00149 // Loop over all the elements and assign them to processor 0. 00150 MeshBase::element_iterator elem_it = mesh.elements_begin(); 00151 const MeshBase::element_iterator elem_end = mesh.elements_end(); 00152 00153 for ( ; elem_it != elem_end; ++elem_it) 00154 (*elem_it)->processor_id() = 0; 00155 00156 // For a single partition, all the nodes are on processor 0 00157 MeshBase::node_iterator node_it = mesh.nodes_begin(); 00158 const MeshBase::node_iterator node_end = mesh.nodes_end(); 00159 00160 for ( ; node_it != node_end; ++node_it) 00161 (*node_it)->processor_id() = 0; 00162 00163 STOP_LOG("single_partition()","Partitioner"); 00164 }
Member Data Documentation
std::vector<int> libMesh::ParmetisPartitioner::_adjncy [private] |
Definition at line 120 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and build_graph().
int libMesh::ParmetisPartitioner::_edgecut [private] |
Definition at line 131 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::map<dof_id_type, dof_id_type> libMesh::ParmetisPartitioner::_global_index_by_pid_map [private] |
Maps active element ids into a contiguous range, as needed by ParMETIS.
Definition at line 112 of file parmetis_partitioner.h.
Referenced by assign_partitioning(), build_graph(), and initialize().
std::vector<dof_id_type> libMesh::ParmetisPartitioner::_n_active_elem_on_proc [private] |
The number of active elements on each processor. Note that ParMETIS requires that each processor have some active elements, it will abort if any processor passes a NULL _part array.
Definition at line 107 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
int libMesh::ParmetisPartitioner::_ncon [private] |
Definition at line 128 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
int libMesh::ParmetisPartitioner::_nparts [private] |
Definition at line 130 of file parmetis_partitioner.h.
Referenced by _do_repartition(), assign_partitioning(), and initialize().
int libMesh::ParmetisPartitioner::_numflag [private] |
Definition at line 129 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<int> libMesh::ParmetisPartitioner::_options [private] |
Definition at line 124 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<int> libMesh::ParmetisPartitioner::_part [private] |
Definition at line 121 of file parmetis_partitioner.h.
Referenced by _do_repartition(), assign_partitioning(), and initialize().
std::vector<float> libMesh::ParmetisPartitioner::_tpwgts [private] |
Definition at line 122 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<float> libMesh::ParmetisPartitioner::_ubvec [private] |
Definition at line 123 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<int> libMesh::ParmetisPartitioner::_vtxdist [private] |
Data structures used by ParMETIS to describe the connectivity graph of the mesh. Consult the ParMETIS documentation.
Definition at line 118 of file parmetis_partitioner.h.
Referenced by _do_repartition(), assign_partitioning(), build_graph(), and initialize().
std::vector<int> libMesh::ParmetisPartitioner::_vwgt [private] |
Definition at line 125 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
ErrorVector* libMesh::Partitioner::_weights [protected, inherited] |
The weights that might be used for partitioning.
Definition at line 155 of file partitioner.h.
Referenced by libMesh::MetisPartitioner::_do_partition(), and libMesh::MetisPartitioner::attach_weights().
int libMesh::ParmetisPartitioner::_wgtflag [private] |
Definition at line 127 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<int> libMesh::ParmetisPartitioner::_xadj [private] |
Definition at line 119 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and build_graph().
const dof_id_type libMesh::Partitioner::communication_blocksize = 1000000 [static, protected, inherited] |
The blocksize to use when doing blocked parallel communication. This limits the maximum vector size which can be used in a single communication step.
Definition at line 150 of file partitioner.h.
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:32 UTC
Hosted By: