$extrastylesheet
#include <metis_partitioner.h>

Public Member Functions | |
| MetisPartitioner () | |
| virtual UniquePtr< Partitioner > | clone () const |
| virtual void | attach_weights (ErrorVector *weights) |
| void | partition (MeshBase &mesh, const unsigned int n) |
| void | partition (MeshBase &mesh) |
| void | repartition (MeshBase &mesh, const unsigned int n) |
| void | repartition (MeshBase &mesh) |
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_partition (MeshBase &mesh, const unsigned int n) |
| void | single_partition (MeshBase &mesh) |
| virtual void | _do_repartition (MeshBase &mesh, const unsigned int n) |
Protected Attributes | |
| ErrorVector * | _weights |
Static Protected Attributes | |
| static const dof_id_type | communication_blocksize = 1000000 |
The MetisPartitioner uses the Metis graph partitioner to partition the elements.
Definition at line 42 of file metis_partitioner.h.
| libMesh::MetisPartitioner::MetisPartitioner | ( | ) | [inline] |
| void libMesh::MetisPartitioner::_do_partition | ( | MeshBase & | mesh, |
| const unsigned int | n | ||
| ) | [protected, virtual] |
Partition the MeshBase into n subdomains.
Implements libMesh::Partitioner.
Definition at line 58 of file metis_partitioner.C.
References libMesh::Partitioner::_weights, libMesh::Elem::active(), libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::Elem::active_family_tree(), libMesh::MeshTools::bounding_box(), libMesh::Parallel::Communicator::broadcast(), libMesh::ParallelObject::comm(), libMesh::vectormap< Key, Tp >::count(), end, libMesh::err, libMesh::MeshCommunication::find_global_indices(), libMesh::DofObject::id(), libMesh::vectormap< Key, Tp >::insert(), libMesh::MeshBase::is_serial(), libMesh::libmesh_assert(), std::max(), libMesh::MeshBase::n_active_elem(), libMesh::Elem::n_neighbors(), libMesh::Elem::n_nodes(), libMesh::Elem::neighbor(), libMesh::METIS_CSR_Graph::offsets, libMesh::Partitioner::partition(), libMesh::METIS_CSR_Graph::prep_n_nonzeros(), libMesh::METIS_CSR_Graph::prepare_for_use(), libMesh::ParallelObject::processor_id(), libMesh::DofObject::processor_id(), libMesh::Partitioner::single_partition(), libMesh::START_LOG(), libMesh::METIS_CSR_Graph::vals, and libMesh::Elem::which_neighbor_am_i().
{
libmesh_assert_greater (n_pieces, 0);
libmesh_assert (mesh.is_serial());
// Check for an easy return
if (n_pieces == 1)
{
this->single_partition (mesh);
return;
}
// What to do if the Metis library IS NOT present
#ifndef LIBMESH_HAVE_METIS
libmesh_here();
libMesh::err << "ERROR: The library has been built without" << std::endl
<< "Metis support. Using a space-filling curve" << std::endl
<< "partitioner instead!" << std::endl;
SFCPartitioner sfcp;
sfcp.partition (mesh, n_pieces);
// What to do if the Metis library IS present
#else
START_LOG("partition()", "MetisPartitioner");
const dof_id_type n_active_elem = mesh.n_active_elem();
// build the graph
// std::vector<int> options(5);
std::vector<int> vwgt(n_active_elem);
std::vector<int> part(n_active_elem);
int
n = static_cast<int>(n_active_elem), // number of "nodes" (elements)
// in the graph
// wgtflag = 2, // weights on vertices only,
// // none on edges
// numflag = 0, // C-style 0-based numbering
nparts = static_cast<int>(n_pieces), // number of subdomains to create
edgecut = 0; // the numbers of edges cut by the
// resulting partition
// Set the options
// options[0] = 0; // use default options
// Metis will only consider the active elements.
// We need to map the active element ids into a
// contiguous range. Further, we want the unique range indexing to be
// independednt of the element ordering, otherwise a circular dependency
// can result in which the partitioning depends on the ordering which
// depends on the partitioning...
vectormap<dof_id_type, dof_id_type> global_index_map;
global_index_map.reserve (n_active_elem);
{
std::vector<dof_id_type> global_index;
MeshBase::element_iterator it = mesh.active_elements_begin();
const MeshBase::element_iterator end = mesh.active_elements_end();
MeshCommunication().find_global_indices (mesh.comm(),
MeshTools::bounding_box(mesh),
it, end, global_index);
libmesh_assert_equal_to (global_index.size(), n_active_elem);
for (std::size_t cnt=0; it != end; ++it)
{
const Elem *elem = *it;
global_index_map.insert (std::make_pair(elem->id(), global_index[cnt++]));
}
libmesh_assert_equal_to (global_index_map.size(), n_active_elem);
}
// Invoke METIS, but only on processor 0.
// Then broadcast the resulting decomposition
if (mesh.processor_id() == 0)
{
METIS_CSR_Graph csr_graph;
csr_graph.offsets.resize(n_active_elem+1, 0);
// Local scope for these
{
// build the graph in CSR format. Note that
// the edges in the graph will correspond to
// face neighbors
#ifdef LIBMESH_ENABLE_AMR
std::vector<const Elem*> neighbors_offspring;
#endif
MeshBase::element_iterator elem_it = mesh.active_elements_begin();
const MeshBase::element_iterator elem_end = mesh.active_elements_end();
#ifndef NDEBUG
std::size_t graph_size=0;
#endif
// (1) first pass - get the row sizes for each element by counting the number
// of face neighbors. Also populate the vwght array if necessary
for (; elem_it != elem_end; ++elem_it)
{
const Elem* elem = *elem_it;
const dof_id_type elem_global_index =
global_index_map[elem->id()];
libmesh_assert_less (elem_global_index, vwgt.size());
// maybe there is a better weight?
// The weight is used to define what a balanced graph is
if(!_weights)
vwgt[elem_global_index] = elem->n_nodes();
else
vwgt[elem_global_index] = static_cast<int>((*_weights)[elem->id()]);
unsigned int num_neighbors = 0;
// 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())
num_neighbors++;
#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());
num_neighbors++;
}
}
}
#endif /* ifdef LIBMESH_ENABLE_AMR */
}
}
csr_graph.prep_n_nonzeros(elem_global_index, num_neighbors);
#ifndef NDEBUG
graph_size += num_neighbors;
#endif
}
csr_graph.prepare_for_use();
// (2) second pass - fill the compressed adjacency array
elem_it = mesh.active_elements_begin();
for (; elem_it != elem_end; ++elem_it)
{
const Elem* elem = *elem_it;
const dof_id_type elem_global_index =
global_index_map[elem->id()];
unsigned int connection=0;
// 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())
csr_graph(elem_global_index, connection++) = global_index_map[neighbor->id()];
#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());
csr_graph(elem_global_index, connection++) = global_index_map[child->id()];
}
}
}
#endif /* ifdef LIBMESH_ENABLE_AMR */
}
}
}
// We create a non-empty vals for a disconnected graph, to
// work around a segfault from METIS.
libmesh_assert_equal_to (csr_graph.vals.size(),
std::max(graph_size,std::size_t(1)));
} // done building the graph
int ncon = 1;
// Select which type of partitioning to create
// Use recursive if the number of partitions is less than or equal to 8
if (n_pieces <= 8)
Metis::METIS_PartGraphRecursive(&n, &ncon, &csr_graph.offsets[0], &csr_graph.vals[0], &vwgt[0], NULL,
NULL, &nparts, NULL, NULL, NULL,
&edgecut, &part[0]);
// Otherwise use kway
else
Metis::METIS_PartGraphKway(&n, &ncon, &csr_graph.offsets[0], &csr_graph.vals[0], &vwgt[0], NULL,
NULL, &nparts, NULL, NULL, NULL,
&edgecut, &part[0]);
} // end processor 0 part
// Broadcase the resutling partition
mesh.comm().broadcast(part);
// Assign the returned processor ids. The part array contains
// the processor id for each active element, but in terms of
// the contiguous indexing we defined above
{
MeshBase::element_iterator it = mesh.active_elements_begin();
const MeshBase::element_iterator end = mesh.active_elements_end();
for (; it!=end; ++it)
{
Elem* elem = *it;
libmesh_assert (global_index_map.count(elem->id()));
const dof_id_type elem_global_index =
global_index_map[elem->id()];
libmesh_assert_less (elem_global_index, part.size());
const processor_id_type elem_procid =
static_cast<processor_id_type>(part[elem_global_index]);
elem->processor_id() = elem_procid;
}
}
STOP_LOG("partition()", "MetisPartitioner");
#endif
}
| virtual void libMesh::Partitioner::_do_repartition | ( | MeshBase & | mesh, |
| const unsigned int | n | ||
| ) | [inline, protected, virtual, inherited] |
This is the actual re-partitioning method which can be overloaded in derived classes. Note that the default behavior is to simply call the partition function.
Reimplemented in libMesh::ParmetisPartitioner.
Definition at line 156 of file partitioner.h.
References libMesh::Partitioner::_do_partition().
Referenced by libMesh::Partitioner::repartition().
{ this->_do_partition (mesh, n); }
| virtual void libMesh::MetisPartitioner::attach_weights | ( | ErrorVector * | ) | [inline, virtual] |
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 from libMesh::Partitioner.
Definition at line 60 of file metis_partitioner.h.
References libMesh::Partitioner::_weights.
{ _weights = weights; }
| virtual UniquePtr<Partitioner> libMesh::MetisPartitioner::clone | ( | ) | const [inline, virtual] |
Creates a new partitioner of this type and returns it in an UniquePtr.
Implements libMesh::Partitioner.
Definition at line 55 of file metis_partitioner.h.
References MetisPartitioner().
{
return UniquePtr<Partitioner>(new MetisPartitioner());
}
| 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 _do_partition(), libMesh::SFCPartitioner::_do_partition(), libMesh::ParmetisPartitioner::_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(), _do_partition(), libMesh::SFCPartitioner::_do_partition(), libMesh::CentroidPartitioner::_do_partition(), libMesh::ParmetisPartitioner::_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");
}
ErrorVector* libMesh::Partitioner::_weights [protected, inherited] |
The weights that might be used for partitioning.
Definition at line 168 of file partitioner.h.
Referenced by _do_partition(), and attach_weights().
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().