$extrastylesheet
#include <parmetis_partitioner.h>

Public Member Functions | |
| ParmetisPartitioner () | |
| virtual UniquePtr< Partitioner > | clone () const |
| void | partition (MeshBase &mesh, const unsigned int n) |
| void | partition (MeshBase &mesh) |
| void | repartition (MeshBase &mesh, const unsigned int n) |
| void | repartition (MeshBase &mesh) |
| virtual void | attach_weights (ErrorVector *) |
Static Public Member Functions | |
| static void | partition_unpartitioned_elements (MeshBase &mesh) |
| static void | partition_unpartitioned_elements (MeshBase &mesh, const unsigned int n) |
| 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 |
| vectormap< 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 |
The ParmetisPartitioner uses the Parmetis graph partitioner to partition the elements.
Definition at line 46 of file parmetis_partitioner.h.
| libMesh::ParmetisPartitioner::ParmetisPartitioner | ( | ) | [inline] |
| 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().
{
this->_do_repartition (mesh, n_sbdmns);
}
| 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 65 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::ParallelObject::comm(), libMesh::err, libMesh::Parallel::Communicator::get(), initialize(), libMesh::MIN_ELEM_PER_PROC, libMesh::ParallelObject::n_processors(), libMesh::Partitioner::partition(), libMesh::Partitioner::single_partition(), and libMesh::START_LOG().
Referenced by _do_partition().
{
libmesh_assert_greater (n_sbdmns, 0);
// Check for an easy return
if (n_sbdmns == 1)
{
this->single_partition(mesh);
return;
}
// This function must be run on all processors at once
libmesh_parallel_only(mesh.comm());
// What to do if the Parmetis library IS NOT present
#ifndef LIBMESH_HAVE_PARMETIS
libmesh_here();
libMesh::err << "ERROR: The library has been built without" << std::endl
<< "Parmetis support. Using a Metis" << std::endl
<< "partitioner instead!" << std::endl;
MetisPartitioner mp;
mp.partition (mesh, n_sbdmns);
// What to do if the Parmetis library IS present
#else
// Revert to METIS on one processor.
if (mesh.n_processors() == 1)
{
MetisPartitioner mp;
mp.partition (mesh, n_sbdmns);
return;
}
START_LOG("repartition()", "ParmetisPartitioner");
// Initialize the data structures required by ParMETIS
this->initialize (mesh, n_sbdmns);
// Make sure all processors have enough active local elements.
// Parmetis tends to crash when it's given only a couple elements
// per partition.
{
bool all_have_enough_elements = true;
for (processor_id_type pid=0; pid<_n_active_elem_on_proc.size(); pid++)
if (_n_active_elem_on_proc[pid] < MIN_ELEM_PER_PROC)
all_have_enough_elements = false;
// Parmetis will not work unless each processor has some
// elements. Specifically, it will abort when passed a NULL
// partition array on *any* of the processors.
if (!all_have_enough_elements)
{
// FIXME: revert to METIS, although this requires a serial mesh
MeshSerializer serialize(mesh);
STOP_LOG ("repartition()", "ParmetisPartitioner");
MetisPartitioner mp;
mp.partition (mesh, n_sbdmns);
return;
}
}
// build the graph corresponding to the mesh
this->build_graph (mesh);
// Partition the graph
std::vector<int> vsize(_vwgt.size(), 1);
float itr = 1000000.0;
MPI_Comm mpi_comm = mesh.comm().get();
// Call the ParMETIS adaptive repartitioning method. This respects the
// original partitioning when computing the new partitioning so as to
// minimize the required data redistribution.
Parmetis::ParMETIS_V3_AdaptiveRepart(_vtxdist.empty() ? NULL : &_vtxdist[0],
_xadj.empty() ? NULL : &_xadj[0],
_adjncy.empty() ? NULL : &_adjncy[0],
_vwgt.empty() ? NULL : &_vwgt[0],
vsize.empty() ? NULL : &vsize[0],
NULL,
&_wgtflag,
&_numflag,
&_ncon,
&_nparts,
_tpwgts.empty() ? NULL : &_tpwgts[0],
_ubvec.empty() ? NULL : &_ubvec[0],
&itr,
&_options[0],
&_edgecut,
_part.empty() ? NULL : &_part[0],
&mpi_comm);
// Assign the returned processor ids
this->assign_partitioning (mesh);
STOP_LOG ("repartition()", "ParmetisPartitioner");
#endif // #ifndef LIBMESH_HAVE_PARMETIS ... else ...
}
| void libMesh::ParmetisPartitioner::assign_partitioning | ( | MeshBase & | mesh | ) | [private] |
Assign the computed partitioning to the mesh.
Definition at line 513 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::ParallelObject::comm(), libMesh::vectormap< Key, Tp >::count(), libMesh::DofObject::id(), libMesh::libmesh_assert(), libMesh::MeshBase::n_active_local_elem(), libMesh::ParallelObject::n_processors(), libMesh::ParallelObject::processor_id(), libMesh::DofObject::processor_id(), and libMesh::Parallel::Communicator::send_receive().
Referenced by _do_repartition().
{
// This function must be run on all processors at once
libmesh_parallel_only(mesh.comm());
const dof_id_type
first_local_elem = _vtxdist[mesh.processor_id()];
std::vector<std::vector<dof_id_type> >
requested_ids(mesh.n_processors()),
requests_to_fill(mesh.n_processors());
MeshBase::element_iterator elem_it = mesh.active_elements_begin();
MeshBase::element_iterator elem_end = mesh.active_elements_end();
for (; elem_it != elem_end; ++elem_it)
{
Elem *elem = *elem_it;
// we need to get the index from the owning processor
// (note we cannot assign it now -- we are iterating
// over elements again and this will be bad!)
libmesh_assert_less (elem->processor_id(), requested_ids.size());
requested_ids[elem->processor_id()].push_back(elem->id());
}
// Trade with all processors (including self) to get their indices
for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
{
// Trade my requests with processor procup and procdown
const processor_id_type procup = (mesh.processor_id() + pid) %
mesh.n_processors();
const processor_id_type procdown = (mesh.n_processors() +
mesh.processor_id() - pid) %
mesh.n_processors();
mesh.comm().send_receive (procup, requested_ids[procup],
procdown, requests_to_fill[procdown]);
// we can overwrite these requested ids in-place.
for (std::size_t i=0; i<requests_to_fill[procdown].size(); i++)
{
const dof_id_type requested_elem_index =
requests_to_fill[procdown][i];
libmesh_assert(_global_index_by_pid_map.count(requested_elem_index));
const dof_id_type global_index_by_pid =
_global_index_by_pid_map[requested_elem_index];
const dof_id_type local_index =
global_index_by_pid - first_local_elem;
libmesh_assert_less (local_index, _part.size());
libmesh_assert_less (local_index, mesh.n_active_local_elem());
const unsigned int elem_procid =
static_cast<unsigned int>(_part[local_index]);
libmesh_assert_less (elem_procid, static_cast<unsigned int>(_nparts));
requests_to_fill[procdown][i] = elem_procid;
}
// Trade back
mesh.comm().send_receive (procdown, requests_to_fill[procdown],
procup, requested_ids[procup]);
}
// and finally assign the partitioning.
// note we are iterating in exactly the same order
// used to build up the request, so we can expect the
// required entries to be in the proper sequence.
elem_it = mesh.active_elements_begin();
elem_end = mesh.active_elements_end();
for (std::vector<unsigned int> counters(mesh.n_processors(), 0);
elem_it != elem_end; ++elem_it)
{
Elem *elem = *elem_it;
const processor_id_type current_pid = elem->processor_id();
libmesh_assert_less (counters[current_pid], requested_ids[current_pid].size());
const processor_id_type elem_procid =
requested_ids[current_pid][counters[current_pid]++];
libmesh_assert_less (elem_procid, static_cast<unsigned int>(_nparts));
elem->processor_id() = elem_procid;
}
}
| 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 131 of file partitioner.h.
{ libmesh_not_implemented(); }
| void libMesh::ParmetisPartitioner::build_graph | ( | const MeshBase & | mesh | ) | [private] |
Build the graph.
Definition at line 387 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::vectormap< Key, Tp >::count(), libMesh::DofObject::id(), libMesh::libmesh_assert(), libMesh::MeshBase::n_active_local_elem(), libMesh::Elem::n_neighbors(), libMesh::Elem::neighbor(), libMesh::ParallelObject::processor_id(), and libMesh::Elem::which_neighbor_am_i().
Referenced by _do_repartition().
{
// build the graph in distributed CSR format. Note that
// the edges in the graph will correspond to
// face neighbors
const dof_id_type n_active_local_elem = mesh.n_active_local_elem();
#ifdef LIBMESH_ENABLE_AMR
std::vector<const Elem*> neighbors_offspring;
#endif
std::vector<std::vector<dof_id_type> > graph(n_active_local_elem);
dof_id_type graph_size=0;
const dof_id_type first_local_elem = _vtxdist[mesh.processor_id()];
MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin();
const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end();
for (; elem_it != elem_end; ++elem_it)
{
const Elem* elem = *elem_it;
libmesh_assert (_global_index_by_pid_map.count(elem->id()));
const dof_id_type global_index_by_pid =
_global_index_by_pid_map[elem->id()];
const dof_id_type local_index =
global_index_by_pid - first_local_elem;
libmesh_assert_less (local_index, n_active_local_elem);
std::vector<dof_id_type> &graph_row = graph[local_index];
// Loop over the element's neighbors. An element
// adjacency corresponds to a face neighbor
for (unsigned int ms=0; ms<elem->n_neighbors(); ms++)
{
const Elem* neighbor = elem->neighbor(ms);
if (neighbor != NULL)
{
// If the neighbor is active treat it
// as a connection
if (neighbor->active())
{
libmesh_assert(_global_index_by_pid_map.count(neighbor->id()));
const dof_id_type neighbor_global_index_by_pid =
_global_index_by_pid_map[neighbor->id()];
graph_row.push_back(neighbor_global_index_by_pid);
graph_size++;
}
#ifdef LIBMESH_ENABLE_AMR
// Otherwise we need to find all of the
// neighbor's children that are connected to
// us and add them
else
{
// The side of the neighbor to which
// we are connected
const unsigned int ns =
neighbor->which_neighbor_am_i (elem);
libmesh_assert_less (ns, neighbor->n_neighbors());
// Get all the active children (& grandchildren, etc...)
// of the neighbor.
neighbor->active_family_tree (neighbors_offspring);
// Get all the neighbor's children that
// live on that side and are thus connected
// to us
for (unsigned int nc=0; nc<neighbors_offspring.size(); nc++)
{
const Elem* child =
neighbors_offspring[nc];
// This does not assume a level-1 mesh.
// Note that since children have sides numbered
// coincident with the parent then this is a sufficient test.
if (child->neighbor(ns) == elem)
{
libmesh_assert (child->active());
libmesh_assert (_global_index_by_pid_map.count(child->id()));
const dof_id_type child_global_index_by_pid =
_global_index_by_pid_map[child->id()];
graph_row.push_back(child_global_index_by_pid);
graph_size++;
}
}
}
#endif /* ifdef LIBMESH_ENABLE_AMR */
}
}
}
// Reserve space in the adjacency array
_xadj.clear();
_xadj.reserve (n_active_local_elem + 1);
_adjncy.clear();
_adjncy.reserve (graph_size);
for (std::size_t r=0; r<graph.size(); r++)
{
_xadj.push_back(_adjncy.size());
std::vector<dof_id_type> graph_row; // build this emtpy
graph_row.swap(graph[r]); // this will deallocate at the end of scope
_adjncy.insert(_adjncy.end(),
graph_row.begin(),
graph_row.end());
}
// The end of the adjacency array for the last elem
_xadj.push_back(_adjncy.size());
libmesh_assert_equal_to (_xadj.size(), n_active_local_elem+1);
libmesh_assert_equal_to (_adjncy.size(), graph_size);
}
| virtual UniquePtr<Partitioner> libMesh::ParmetisPartitioner::clone | ( | ) | const [inline, virtual] |
Creates a new partitioner of this type and returns it in an UniquePtr.
Implements libMesh::Partitioner.
Definition at line 59 of file parmetis_partitioner.h.
References ParmetisPartitioner().
{
return UniquePtr<Partitioner>(new ParmetisPartitioner());
}
| void libMesh::ParmetisPartitioner::initialize | ( | const MeshBase & | mesh, |
| const unsigned int | n_sbdmns | ||
| ) | [private] |
Initialize data structures.
Definition at line 179 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::ParallelObject::comm(), libMesh::vectormap< Key, Tp >::count(), end, libMesh::MeshCommunication::find_global_indices(), libMesh::DofObject::id(), libMesh::vectormap< Key, Tp >::insert(), libMesh::libmesh_assert(), std::min(), libMesh::MeshBase::n_active_local_elem(), libMesh::Elem::n_nodes(), libMesh::ParallelObject::n_processors(), and libMesh::ParallelObject::processor_id().
Referenced by _do_repartition().
{
const dof_id_type n_active_local_elem = mesh.n_active_local_elem();
// Set parameters.
_wgtflag = 2; // weights on vertices only
_ncon = 1; // one weight per vertex
_numflag = 0; // C-style 0-based numbering
_nparts = static_cast<int>(n_sbdmns); // number of subdomains to create
_edgecut = 0; // the numbers of edges cut by the
// partition
// Initialize data structures for ParMETIS
_vtxdist.resize (mesh.n_processors()+1); std::fill (_vtxdist.begin(), _vtxdist.end(), 0);
_tpwgts.resize (_nparts); std::fill (_tpwgts.begin(), _tpwgts.end(), 1./_nparts);
_ubvec.resize (_ncon); std::fill (_ubvec.begin(), _ubvec.end(), 1.05);
_part.resize (n_active_local_elem); std::fill (_part.begin(), _part.end(), 0);
_options.resize (5);
_vwgt.resize (n_active_local_elem);
// Set the options
_options[0] = 1; // don't use default options
_options[1] = 0; // default (level of timing)
_options[2] = 15; // random seed (default)
_options[3] = 2; // processor distribution and subdomain distribution are decoupled
// Find the number of active elements on each processor. We cannot use
// mesh.n_active_elem_on_proc(pid) since that only returns the number of
// elements assigned to pid which are currently stored on the calling
// processor. This will not in general be correct for parallel meshes
// when (pid!=mesh.processor_id()).
_n_active_elem_on_proc.resize(mesh.n_processors());
mesh.comm().allgather(n_active_local_elem, _n_active_elem_on_proc);
// count the total number of active elements in the mesh. Note we cannot
// use mesh.n_active_elem() in general since this only returns the number
// of active elements which are stored on the calling processor.
// We should not use n_active_elem for any allocation because that will
// be inheritly unscalable, but it can be useful for libmesh_assertions.
dof_id_type n_active_elem=0;
// Set up the vtxdist array. This will be the same on each processor.
// ***** Consult the Parmetis documentation. *****
libmesh_assert_equal_to (_vtxdist.size(),
cast_int<std::size_t>(mesh.n_processors()+1));
libmesh_assert_equal_to (_vtxdist[0], 0);
for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
{
_vtxdist[pid+1] = _vtxdist[pid] + _n_active_elem_on_proc[pid];
n_active_elem += _n_active_elem_on_proc[pid];
}
libmesh_assert_equal_to (_vtxdist.back(), static_cast<int>(n_active_elem));
// ParMetis expects the elements to be numbered in contiguous blocks
// by processor, i.e. [0, ne0), [ne0, ne0+ne1), ...
// Since we only partition active elements we should have no expectation
// that we currently have such a distribution. So we need to create it.
// Also, at the same time we are going to map all the active elements into a globally
// unique range [0,n_active_elem) which is *independent* of the current partitioning.
// This can be fed to ParMetis as the initial partitioning of the subdomains (decoupled
// from the partitioning of the objects themselves). This allows us to get the same
// resultant partitioning independed of the input partitioning.
MeshTools::BoundingBox bbox =
MeshTools::bounding_box(mesh);
_global_index_by_pid_map.clear();
// Maps active element ids into a contiguous range independent of partitioning.
// (only needs local scope)
vectormap<dof_id_type, dof_id_type> global_index_map;
{
std::vector<dof_id_type> global_index;
// create the mapping which is contiguous by processor
dof_id_type pid_offset=0;
for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
{
MeshBase::const_element_iterator it = mesh.active_pid_elements_begin(pid);
const MeshBase::const_element_iterator end = mesh.active_pid_elements_end(pid);
// note that we may not have all (or any!) the active elements which belong on this processor,
// but by calling this on all processors a unique range in [0,_n_active_elem_on_proc[pid])
// is constructed. Only the indices for the elements we pass in are returned in the array.
MeshCommunication().find_global_indices (mesh.comm(),
bbox, it, end,
global_index);
for (dof_id_type cnt=0; it != end; ++it)
{
const Elem *elem = *it;
libmesh_assert (!_global_index_by_pid_map.count(elem->id()));
libmesh_assert_less (cnt, global_index.size());
libmesh_assert_less (global_index[cnt], _n_active_elem_on_proc[pid]);
_global_index_by_pid_map.insert(std::make_pair(elem->id(), global_index[cnt++] + pid_offset));
}
pid_offset += _n_active_elem_on_proc[pid];
}
// create the unique mapping for all active elements independent of partitioning
{
MeshBase::const_element_iterator it = mesh.active_elements_begin();
const MeshBase::const_element_iterator end = mesh.active_elements_end();
// Calling this on all processors a unique range in [0,n_active_elem) is constructed.
// Only the indices for the elements we pass in are returned in the array.
MeshCommunication().find_global_indices (mesh.comm(),
bbox, it, end,
global_index);
for (dof_id_type cnt=0; it != end; ++it)
{
const Elem *elem = *it;
libmesh_assert (!global_index_map.count(elem->id()));
libmesh_assert_less (cnt, global_index.size());
libmesh_assert_less (global_index[cnt], n_active_elem);
global_index_map.insert(std::make_pair(elem->id(), global_index[cnt++]));
}
}
// really, shouldn't be close!
libmesh_assert_less_equal (global_index_map.size(), n_active_elem);
libmesh_assert_less_equal (_global_index_by_pid_map.size(), n_active_elem);
// At this point the two maps should be the same size. If they are not
// then the number of active elements is not the same as the sum over all
// processors of the number of active elements per processor, which means
// there must be some unpartitioned objects out there.
if (global_index_map.size() != _global_index_by_pid_map.size())
libmesh_error_msg("ERROR: ParmetisPartitioner cannot handle unpartitioned objects!");
}
// Finally, we need to initialize the vertex (partition) weights and the initial subdomain
// mapping. The subdomain mapping will be independent of the processor mapping, and is
// defined by a simple mapping of the global indices we just found.
{
std::vector<dof_id_type> subdomain_bounds(mesh.n_processors());
const dof_id_type first_local_elem = _vtxdist[mesh.processor_id()];
for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
{
dof_id_type tgt_subdomain_size = 0;
// watch out for the case that n_subdomains < n_processors
if (pid < static_cast<unsigned int>(_nparts))
{
tgt_subdomain_size = n_active_elem/std::min
(cast_int<int>(mesh.n_processors()), _nparts);
if (pid < n_active_elem%_nparts)
tgt_subdomain_size++;
}
if (pid == 0)
subdomain_bounds[0] = tgt_subdomain_size;
else
subdomain_bounds[pid] = subdomain_bounds[pid-1] + tgt_subdomain_size;
}
libmesh_assert_equal_to (subdomain_bounds.back(), n_active_elem);
MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin();
const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end();
for (; elem_it != elem_end; ++elem_it)
{
const Elem *elem = *elem_it;
libmesh_assert (_global_index_by_pid_map.count(elem->id()));
const dof_id_type global_index_by_pid =
_global_index_by_pid_map[elem->id()];
libmesh_assert_less (global_index_by_pid, n_active_elem);
const dof_id_type local_index =
global_index_by_pid - first_local_elem;
libmesh_assert_less (local_index, n_active_local_elem);
libmesh_assert_less (local_index, _vwgt.size());
// TODO:[BSK] maybe there is a better weight?
_vwgt[local_index] = elem->n_nodes();
// find the subdomain this element belongs in
libmesh_assert (global_index_map.count(elem->id()));
const dof_id_type global_index =
global_index_map[elem->id()];
libmesh_assert_less (global_index, subdomain_bounds.back());
const unsigned int subdomain_id =
std::distance(subdomain_bounds.begin(),
std::lower_bound(subdomain_bounds.begin(),
subdomain_bounds.end(),
global_index));
libmesh_assert_less (subdomain_id, static_cast<unsigned int>(_nparts));
libmesh_assert_less (local_index, _part.size());
_part[local_index] = subdomain_id;
}
}
}
| void libMesh::Partitioner::partition | ( | MeshBase & | mesh, |
| const unsigned int | n | ||
| ) | [inherited] |
Partition the MeshBase into n parts. 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 55 of file partitioner.C.
References libMesh::Partitioner::_do_partition(), libMesh::ParallelObject::comm(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_remote_elems(), mesh, 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::MetisPartitioner::_do_partition(), libMesh::SFCPartitioner::_do_partition(), _do_repartition(), and libMesh::Partitioner::partition().
{
libmesh_parallel_only(mesh.comm());
// BSK - temporary fix while redistribution is integrated 6/26/2008
// Uncomment this to not repartition in parallel
// if (!mesh.is_serial())
// return;
// we cannot partition into more pieces than we have
// active elements!
const unsigned int n_parts =
static_cast<unsigned int>
(std::min(mesh.n_active_elem(), static_cast<dof_id_type>(n)));
// Set the number of partitions in the mesh
mesh.set_n_partitions()=n_parts;
if (n_parts == 1)
{
this->single_partition (mesh);
return;
}
// First assign a temporary partitioning to any unpartitioned elements
Partitioner::partition_unpartitioned_elements(mesh, n_parts);
// Call the partitioning function
this->_do_partition(mesh,n_parts);
// Set the parent's processor ids
Partitioner::set_parent_processor_ids(mesh);
// Redistribute elements if necessary, before setting node processor
// ids, to make sure those will be set consistently
mesh.redistribute();
#ifdef DEBUG
MeshTools::libmesh_assert_valid_remote_elems(mesh);
// Messed up elem processor_id()s can leave us without the child
// elements we need to restrict vectors on a distributed mesh
MeshTools::libmesh_assert_valid_procids<Elem>(mesh);
#endif
// Set the node's processor ids
Partitioner::set_node_processor_ids(mesh);
#ifdef DEBUG
MeshTools::libmesh_assert_valid_procids<Elem>(mesh);
#endif
// Give derived Mesh classes a chance to update any cached data to
// reflect the new partitioning
mesh.update_post_partitioning();
}
| void libMesh::Partitioner::partition | ( | MeshBase & | mesh | ) | [inherited] |
Partition the MeshBase into mesh.n_processors() parts. 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::ParallelObject::n_processors(), and libMesh::Partitioner::partition().
| void libMesh::Partitioner::partition_unpartitioned_elements | ( | MeshBase & | mesh | ) | [static, inherited] |
This function
Definition at line 180 of file partitioner.C.
References libMesh::ParallelObject::n_processors().
Referenced by libMesh::Partitioner::partition(), and libMesh::Partitioner::repartition().
{
Partitioner::partition_unpartitioned_elements(mesh, mesh.n_processors());
}
| void libMesh::Partitioner::partition_unpartitioned_elements | ( | MeshBase & | mesh, |
| const unsigned int | n | ||
| ) | [static, inherited] |
Definition at line 187 of file partitioner.C.
References libMesh::MeshTools::bounding_box(), libMesh::ParallelObject::comm(), end, libMesh::MeshCommunication::find_global_indices(), libMesh::MeshTools::n_elem(), libMesh::ParallelObject::n_processors(), libMesh::DofObject::processor_id(), libMesh::MeshBase::unpartitioned_elements_begin(), and libMesh::MeshBase::unpartitioned_elements_end().
{
MeshBase::element_iterator it = mesh.unpartitioned_elements_begin();
const MeshBase::element_iterator end = mesh.unpartitioned_elements_end();
const dof_id_type n_unpartitioned_elements = MeshTools::n_elem (it, end);
// the unpartitioned elements must exist on all processors. If the range is empty on one
// it is empty on all, and we can quit right here.
if (!n_unpartitioned_elements) return;
// find the target subdomain sizes
std::vector<dof_id_type> subdomain_bounds(mesh.n_processors());
for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
{
dof_id_type tgt_subdomain_size = 0;
// watch out for the case that n_subdomains < n_processors
if (pid < n_subdomains)
{
tgt_subdomain_size = n_unpartitioned_elements/n_subdomains;
if (pid < n_unpartitioned_elements%n_subdomains)
tgt_subdomain_size++;
}
//libMesh::out << "pid, #= " << pid << ", " << tgt_subdomain_size << std::endl;
if (pid == 0)
subdomain_bounds[0] = tgt_subdomain_size;
else
subdomain_bounds[pid] = subdomain_bounds[pid-1] + tgt_subdomain_size;
}
libmesh_assert_equal_to (subdomain_bounds.back(), n_unpartitioned_elements);
// create the unique mapping for all unpartitioned elements independent of partitioning
// determine the global indexing for all the unpartitoned elements
std::vector<dof_id_type> global_indices;
// Calling this on all processors a unique range in [0,n_unpartitioned_elements) is constructed.
// Only the indices for the elements we pass in are returned in the array.
MeshCommunication().find_global_indices (mesh.comm(),
MeshTools::bounding_box(mesh), it, end,
global_indices);
for (dof_id_type cnt=0; it != end; ++it)
{
Elem *elem = *it;
libmesh_assert_less (cnt, global_indices.size());
const dof_id_type global_index =
global_indices[cnt++];
libmesh_assert_less (global_index, subdomain_bounds.back());
libmesh_assert_less (global_index, n_unpartitioned_elements);
const processor_id_type subdomain_id =
cast_int<processor_id_type>
(std::distance(subdomain_bounds.begin(),
std::upper_bound(subdomain_bounds.begin(),
subdomain_bounds.end(),
global_index)));
libmesh_assert_less (subdomain_id, n_subdomains);
elem->processor_id() = subdomain_id;
//libMesh::out << "assigning " << global_index << " to " << subdomain_id << std::endl;
}
}
| void libMesh::Partitioner::repartition | ( | MeshBase & | mesh, |
| const unsigned int | n | ||
| ) | [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(mesh,n)
Definition at line 122 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().
Referenced by libMesh::Partitioner::repartition().
{
// we cannot partition into more pieces than we have
// active elements!
const unsigned int n_parts =
static_cast<unsigned int>
(std::min(mesh.n_active_elem(), static_cast<dof_id_type>(n)));
// Set the number of partitions in the mesh
mesh.set_n_partitions()=n_parts;
if (n_parts == 1)
{
this->single_partition (mesh);
return;
}
// First assign a temporary partitioning to any unpartitioned elements
Partitioner::partition_unpartitioned_elements(mesh, n_parts);
// Call the partitioning function
this->_do_repartition(mesh,n_parts);
// Set the parent's processor ids
Partitioner::set_parent_processor_ids(mesh);
// Set the node's processor ids
Partitioner::set_node_processor_ids(mesh);
}
| void libMesh::Partitioner::repartition | ( | MeshBase & | mesh | ) | [inherited] |
Repartitions the MeshBase into mesh.n_processors() parts. This is required since some partitoning algorithms can repartition more efficiently than computing a new partitioning from scratch.
Definition at line 115 of file partitioner.C.
References libMesh::ParallelObject::n_processors(), and libMesh::Partitioner::repartition().
{
this->repartition(mesh,mesh.n_processors());
}
| 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 439 of file partitioner.C.
References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::ParallelObject::comm(), libMesh::Elem::get_node(), libMesh::DofObject::id(), libMesh::DofObject::invalid_processor_id, libMesh::DofObject::invalidate_processor_id(), libMesh::libmesh_assert(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), mesh, std::min(), libMesh::MeshTools::n_elem(), libMesh::Elem::n_nodes(), libMesh::ParallelObject::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::ParallelObject::processor_id(), libMesh::DofObject::processor_id(), libMesh::Parallel::Communicator::send_receive(), libMesh::START_LOG(), 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::MeshBase::partition(), libMesh::XdrIO::read(), libMesh::Partitioner::repartition(), and libMesh::BoundaryInfo::sync().
{
START_LOG("set_node_processor_ids()","Partitioner");
// This function must be run on all processors at once
libmesh_parallel_only(mesh.comm());
// If we have any unpartitioned elements at this
// stage there is a problem
libmesh_assert (MeshTools::n_elem(mesh.unpartitioned_elements_begin(),
mesh.unpartitioned_elements_end()) == 0);
// const dof_id_type orig_n_local_nodes = mesh.n_local_nodes();
// libMesh::err << "[" << mesh.processor_id() << "]: orig_n_local_nodes="
// << orig_n_local_nodes << std::endl;
// Build up request sets. Each node is currently owned by a processor because
// it is connected to an element owned by that processor. However, during the
// repartitioning phase that element may have been assigned a new processor id, but
// it is still resident on the original processor. We need to know where to look
// for new ids before assigning new ids, otherwise we may be asking the wrong processors
// for the wrong information.
//
// The only remaining issue is what to do with unpartitioned nodes. Since they are required
// to live on all processors we can simply rely on ourselves to number them properly.
std::vector<std::vector<dof_id_type> >
requested_node_ids(mesh.n_processors());
// Loop over all the nodes, count the ones on each processor. We can skip ourself
std::vector<dof_id_type> ghost_nodes_from_proc(mesh.n_processors(), 0);
MeshBase::node_iterator node_it = mesh.nodes_begin();
const MeshBase::node_iterator node_end = mesh.nodes_end();
for (; node_it != node_end; ++node_it)
{
Node *node = *node_it;
libmesh_assert(node);
const processor_id_type current_pid = node->processor_id();
if (current_pid != mesh.processor_id() &&
current_pid != DofObject::invalid_processor_id)
{
libmesh_assert_less (current_pid, ghost_nodes_from_proc.size());
ghost_nodes_from_proc[current_pid]++;
}
}
// We know how many objects live on each processor, so reserve()
// space for each.
for (processor_id_type pid=0; pid != mesh.n_processors(); ++pid)
requested_node_ids[pid].reserve(ghost_nodes_from_proc[pid]);
// We need to get the new pid for each node from the processor
// which *currently* owns the node. We can safely skip ourself
for (node_it = mesh.nodes_begin(); node_it != node_end; ++node_it)
{
Node *node = *node_it;
libmesh_assert(node);
const processor_id_type current_pid = node->processor_id();
if (current_pid != mesh.processor_id() &&
current_pid != DofObject::invalid_processor_id)
{
libmesh_assert_less (current_pid, requested_node_ids.size());
libmesh_assert_less (requested_node_ids[current_pid].size(),
ghost_nodes_from_proc[current_pid]);
requested_node_ids[current_pid].push_back(node->id());
}
// Unset any previously-set node processor ids
node->invalidate_processor_id();
}
// Loop over all the active elements
MeshBase::element_iterator elem_it = mesh.active_elements_begin();
const MeshBase::element_iterator elem_end = mesh.active_elements_end();
for ( ; elem_it != elem_end; ++elem_it)
{
Elem* elem = *elem_it;
libmesh_assert(elem);
libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id);
// For each node, set the processor ID to the min of
// its current value and this Element's processor id.
//
// TODO: we would probably get better parallel partitioning if
// we did something like "min for even numbered nodes, max for
// odd numbered". We'd need to be careful about how that would
// affect solution ordering for I/O, though.
for (unsigned int n=0; n<elem->n_nodes(); ++n)
elem->get_node(n)->processor_id() = std::min(elem->get_node(n)->processor_id(),
elem->processor_id());
}
// And loop over the subactive elements, but don't reassign
// nodes that are already active on another processor.
MeshBase::element_iterator sub_it = mesh.subactive_elements_begin();
const MeshBase::element_iterator sub_end = mesh.subactive_elements_end();
for ( ; sub_it != sub_end; ++sub_it)
{
Elem* elem = *sub_it;
libmesh_assert(elem);
libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id);
for (unsigned int n=0; n<elem->n_nodes(); ++n)
if (elem->get_node(n)->processor_id() == DofObject::invalid_processor_id)
elem->get_node(n)->processor_id() = elem->processor_id();
}
// Same for the inactive elements -- we will have already gotten most of these
// nodes, *except* for the case of a parent with a subset of children which are
// ghost elements. In that case some of the parent nodes will not have been
// properly handled yet
MeshBase::element_iterator not_it = mesh.not_active_elements_begin();
const MeshBase::element_iterator not_end = mesh.not_active_elements_end();
for ( ; not_it != not_end; ++not_it)
{
Elem* elem = *not_it;
libmesh_assert(elem);
libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id);
for (unsigned int n=0; n<elem->n_nodes(); ++n)
if (elem->get_node(n)->processor_id() == DofObject::invalid_processor_id)
elem->get_node(n)->processor_id() = elem->processor_id();
}
// We can't assert that all nodes are connected to elements, because
// a ParallelMesh with NodeConstraints might have pulled in some
// remote nodes solely for evaluating those constraints.
// MeshTools::libmesh_assert_connected_nodes(mesh);
// For such nodes, we'll do a sanity check later when making sure
// that we successfully reset their processor ids to something
// valid.
// Next set node ids from other processors, excluding self
for (processor_id_type p=1; p != mesh.n_processors(); ++p)
{
// Trade my requests with processor procup and procdown
processor_id_type procup = cast_int<processor_id_type>
((mesh.processor_id() + p) % mesh.n_processors());
processor_id_type procdown = cast_int<processor_id_type>
((mesh.n_processors() + mesh.processor_id() - p) %
mesh.n_processors());
std::vector<dof_id_type> request_to_fill;
mesh.comm().send_receive(procup, requested_node_ids[procup],
procdown, request_to_fill);
// Fill those requests in-place
for (std::size_t i=0; i != request_to_fill.size(); ++i)
{
Node *node = mesh.node_ptr(request_to_fill[i]);
libmesh_assert(node);
const processor_id_type new_pid = node->processor_id();
// We may have an invalid processor_id() on nodes that have been
// "detatched" from coarsened-away elements but that have not yet
// themselves been removed.
// libmesh_assert_not_equal_to (new_pid, DofObject::invalid_processor_id);
// libmesh_assert_less (new_pid, mesh.n_partitions()); // this is the correct test --
request_to_fill[i] = new_pid; // the number of partitions may
} // not equal the number of processors
// Trade back the results
std::vector<dof_id_type> filled_request;
mesh.comm().send_receive(procdown, request_to_fill,
procup, filled_request);
libmesh_assert_equal_to (filled_request.size(), requested_node_ids[procup].size());
// And copy the id changes we've now been informed of
for (std::size_t i=0; i != filled_request.size(); ++i)
{
Node *node = mesh.node_ptr(requested_node_ids[procup][i]);
libmesh_assert(node);
// this is the correct test -- the number of partitions may
// not equal the number of processors
// But: we may have an invalid processor_id() on nodes that
// have been "detatched" from coarsened-away elements but
// that have not yet themselves been removed.
// libmesh_assert_less (filled_request[i], mesh.n_partitions());
node->processor_id(cast_int<processor_id_type>(filled_request[i]));
}
}
#ifdef DEBUG
MeshTools::libmesh_assert_valid_procids<Node>(mesh);
#endif
STOP_LOG("set_node_processor_ids()","Partitioner");
}
| 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.
Definition at line 261 of file partitioner.C.
References libMesh::Elem::active_family_tree(), libMesh::Elem::child(), libMesh::Partitioner::communication_blocksize, end, libMesh::DofObject::id(), libMesh::DofObject::invalid_processor_id, libMesh::DofObject::invalidate_processor_id(), libMesh::Elem::is_remote(), libMesh::libmesh_assert(), mesh, std::min(), libMesh::Elem::n_children(), libMesh::MeshTools::n_elem(), libMesh::Elem::parent(), libMesh::processor_id(), libMesh::DofObject::processor_id(), libMesh::START_LOG(), and libMesh::Elem::total_family_tree().
Referenced by libMesh::Partitioner::partition(), and libMesh::Partitioner::repartition().
{
START_LOG("set_parent_processor_ids()","Partitioner");
#ifdef LIBMESH_ENABLE_AMR
// If the mesh is serial we have access to all the elements,
// in particular all the active ones. We can therefore set
// the parent processor ids indirecly through their children, and
// set the subactive processor ids while examining their active
// ancestors.
// By convention a parent is assigned to the minimum processor
// of all its children, and a subactive is assigned to the processor
// of its active ancestor.
if (mesh.is_serial())
{
// Loop over all the active elements in the mesh
MeshBase::element_iterator it = mesh.active_elements_begin();
const MeshBase::element_iterator end = mesh.active_elements_end();
for ( ; it!=end; ++it)
{
Elem *child = *it;
// First set descendents
std::vector<const Elem*> subactive_family;
child->total_family_tree(subactive_family);
for (unsigned int i = 0; i != subactive_family.size(); ++i)
const_cast<Elem*>(subactive_family[i])->processor_id() = child->processor_id();
// Then set ancestors
Elem *parent = child->parent();
while (parent)
{
// invalidate the parent id, otherwise the min below
// will not work if the current parent id is less
// than all the children!
parent->invalidate_processor_id();
for(unsigned int c=0; c<parent->n_children(); c++)
{
child = parent->child(c);
libmesh_assert(child);
libmesh_assert(!child->is_remote());
libmesh_assert_not_equal_to (child->processor_id(), DofObject::invalid_processor_id);
parent->processor_id() = std::min(parent->processor_id(),
child->processor_id());
}
parent = parent->parent();
}
}
}
// When the mesh is parallel we cannot guarantee that parents have access to
// all their children.
else
{
// Setting subactive processor ids is easy: we can guarantee
// that children have access to all their parents.
// Loop over all the active elements in the mesh
MeshBase::element_iterator it = mesh.active_elements_begin();
const MeshBase::element_iterator end = mesh.active_elements_end();
for ( ; it!=end; ++it)
{
Elem *child = *it;
std::vector<const Elem*> subactive_family;
child->total_family_tree(subactive_family);
for (unsigned int i = 0; i != subactive_family.size(); ++i)
const_cast<Elem*>(subactive_family[i])->processor_id() = child->processor_id();
}
// When the mesh is parallel we cannot guarantee that parents have access to
// all their children.
// We will use a brute-force approach here. Each processor finds its parent
// elements and sets the parent pid to the minimum of its
// semilocal descendants.
// A global reduction is then performed to make sure the true minimum is found.
// As noted, this is required because we cannot guarantee that a parent has
// access to all its children on any single processor.
libmesh_parallel_only(mesh.comm());
libmesh_assert(MeshTools::n_elem(mesh.unpartitioned_elements_begin(),
mesh.unpartitioned_elements_end()) == 0);
const dof_id_type max_elem_id = mesh.max_elem_id();
std::vector<processor_id_type>
parent_processor_ids (std::min(communication_blocksize,
max_elem_id));
for (dof_id_type blk=0, last_elem_id=0; last_elem_id<max_elem_id; blk++)
{
last_elem_id =
std::min(static_cast<dof_id_type>((blk+1)*communication_blocksize),
max_elem_id);
const dof_id_type first_elem_id = blk*communication_blocksize;
std::fill (parent_processor_ids.begin(),
parent_processor_ids.end(),
DofObject::invalid_processor_id);
// first build up local contributions to parent_processor_ids
MeshBase::element_iterator not_it = mesh.ancestor_elements_begin();
const MeshBase::element_iterator not_end = mesh.ancestor_elements_end();
bool have_parent_in_block = false;
for ( ; not_it != not_end; ++not_it)
{
Elem *parent = *not_it;
const dof_id_type parent_idx = parent->id();
libmesh_assert_less (parent_idx, max_elem_id);
if ((parent_idx >= first_elem_id) &&
(parent_idx < last_elem_id))
{
have_parent_in_block = true;
processor_id_type parent_pid = DofObject::invalid_processor_id;
std::vector<const Elem*> active_family;
parent->active_family_tree(active_family);
for (unsigned int i = 0; i != active_family.size(); ++i)
parent_pid = std::min (parent_pid, active_family[i]->processor_id());
const dof_id_type packed_idx = parent_idx - first_elem_id;
libmesh_assert_less (packed_idx, parent_processor_ids.size());
parent_processor_ids[packed_idx] = parent_pid;
}
}
// then find the global minimum
mesh.comm().min (parent_processor_ids);
// and assign the ids, if we have a parent in this block.
if (have_parent_in_block)
for (not_it = mesh.ancestor_elements_begin();
not_it != not_end; ++not_it)
{
Elem *parent = *not_it;
const dof_id_type parent_idx = parent->id();
if ((parent_idx >= first_elem_id) &&
(parent_idx < last_elem_id))
{
const dof_id_type packed_idx = parent_idx - first_elem_id;
libmesh_assert_less (packed_idx, parent_processor_ids.size());
const processor_id_type parent_pid =
parent_processor_ids[packed_idx];
libmesh_assert_not_equal_to (parent_pid, DofObject::invalid_processor_id);
parent->processor_id() = parent_pid;
}
}
}
}
#endif // LIBMESH_ENABLE_AMR
STOP_LOG("set_parent_processor_ids()","Partitioner");
}
| 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 157 of file partitioner.C.
References libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), libMesh::MeshBase::nodes_begin(), libMesh::MeshBase::nodes_end(), and libMesh::START_LOG().
Referenced by libMesh::LinearPartitioner::_do_partition(), libMesh::MetisPartitioner::_do_partition(), libMesh::SFCPartitioner::_do_partition(), libMesh::CentroidPartitioner::_do_partition(), _do_repartition(), libMesh::Partitioner::partition(), and libMesh::Partitioner::repartition().
{
START_LOG("single_partition()","Partitioner");
// Loop over all the elements and assign them to processor 0.
MeshBase::element_iterator elem_it = mesh.elements_begin();
const MeshBase::element_iterator elem_end = mesh.elements_end();
for ( ; elem_it != elem_end; ++elem_it)
(*elem_it)->processor_id() = 0;
// For a single partition, all the nodes are on processor 0
MeshBase::node_iterator node_it = mesh.nodes_begin();
const MeshBase::node_iterator node_end = mesh.nodes_end();
for ( ; node_it != node_end; ++node_it)
(*node_it)->processor_id() = 0;
STOP_LOG("single_partition()","Partitioner");
}
std::vector<int> libMesh::ParmetisPartitioner::_adjncy [private] |
Definition at line 121 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and build_graph().
int libMesh::ParmetisPartitioner::_edgecut [private] |
Definition at line 132 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
vectormap<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 113 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 108 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
int libMesh::ParmetisPartitioner::_ncon [private] |
Definition at line 129 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
int libMesh::ParmetisPartitioner::_nparts [private] |
Definition at line 131 of file parmetis_partitioner.h.
Referenced by _do_repartition(), assign_partitioning(), and initialize().
int libMesh::ParmetisPartitioner::_numflag [private] |
Definition at line 130 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<int> libMesh::ParmetisPartitioner::_options [private] |
Definition at line 125 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<int> libMesh::ParmetisPartitioner::_part [private] |
Definition at line 122 of file parmetis_partitioner.h.
Referenced by _do_repartition(), assign_partitioning(), and initialize().
std::vector<float> libMesh::ParmetisPartitioner::_tpwgts [private] |
Definition at line 123 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<float> libMesh::ParmetisPartitioner::_ubvec [private] |
Definition at line 124 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 119 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 126 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 168 of file partitioner.h.
Referenced by libMesh::MetisPartitioner::_do_partition(), and libMesh::MetisPartitioner::attach_weights().
int libMesh::ParmetisPartitioner::_wgtflag [private] |
Definition at line 128 of file parmetis_partitioner.h.
Referenced by _do_repartition(), and initialize().
std::vector<int> libMesh::ParmetisPartitioner::_xadj [private] |
Definition at line 120 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 163 of file partitioner.h.
Referenced by libMesh::Partitioner::set_parent_processor_ids().