libMesh::ParmetisPartitioner Class Reference

Partitioner which provides an interface to ParMETIS. More...

#include <parmetis_partitioner.h>

Inheritance diagram for libMesh::ParmetisPartitioner:

Public Member Functions

 ParmetisPartitioner ()
 
 ~ParmetisPartitioner ()
 
virtual UniquePtr< Partitionerclone () const libmesh_override
 
virtual void partition (MeshBase &mesh, const unsigned int n)
 
virtual void partition (MeshBase &mesh)
 
virtual void partition_range (MeshBase &, MeshBase::element_iterator, MeshBase::element_iterator, const unsigned int)
 
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) libmesh_override
 
virtual void _do_partition (MeshBase &mesh, const unsigned int n) libmesh_override
 
void single_partition (MeshBase &mesh)
 
void single_partition_range (MeshBase::element_iterator it, MeshBase::element_iterator end)
 

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
 
UniquePtr< ParmetisHelper_pmetis
 

Detailed Description

Partitioner which provides an interface to ParMETIS.

The ParmetisPartitioner uses the Parmetis graph partitioner to partition the elements.

Author
Benjamin S. Kirk
Date
2003

Definition at line 46 of file parmetis_partitioner.h.

Constructor & Destructor Documentation

libMesh::ParmetisPartitioner::ParmetisPartitioner ( )

Constructor.

Definition at line 64 of file parmetis_partitioner.C.

Referenced by clone().

66  : _pmetis(new ParmetisHelper)
67 #endif
68 {}
69 
70 
71 
73 {
74 }
UniquePtr< ParmetisHelper > _pmetis
libMesh::ParmetisPartitioner::~ParmetisPartitioner ( )

Destructor.

Member Function Documentation

void libMesh::ParmetisPartitioner::_do_partition ( MeshBase mesh,
const unsigned int  n 
)
protectedvirtual

Partition the MeshBase into n subdomains.

Implements libMesh::Partitioner.

Definition at line 78 of file parmetis_partitioner.C.

Referenced by clone().

80 {
81  this->_do_repartition (mesh, n_sbdmns);
82 }
virtual void _do_repartition(MeshBase &mesh, const unsigned int n) libmesh_override
MeshBase & mesh
void libMesh::ParmetisPartitioner::_do_repartition ( MeshBase mesh,
const unsigned int  n 
)
protectedvirtual

Parmetis can handle dynamically repartitioning a mesh such that the redistribution costs are minimized. This method takes a previously partitioned mesh (which may have then been adaptively refined) and repartitions it.

Reimplemented from libMesh::Partitioner.

Definition at line 86 of file parmetis_partitioner.C.

References libMesh::ParallelObject::comm(), libMesh::err, libMesh::Parallel::Communicator::get(), libmesh_nullptr, libMesh::ParallelObject::n_processors(), and libMesh::Partitioner::partition().

Referenced by clone().

88 {
89  libmesh_assert_greater (n_sbdmns, 0);
90 
91  // Check for an easy return
92  if (n_sbdmns == 1)
93  {
94  this->single_partition(mesh);
95  return;
96  }
97 
98  // This function must be run on all processors at once
99  libmesh_parallel_only(mesh.comm());
100 
101  // What to do if the Parmetis library IS NOT present
102 #ifndef LIBMESH_HAVE_PARMETIS
103 
104  libmesh_here();
105  libMesh::err << "ERROR: The library has been built without" << std::endl
106  << "Parmetis support. Using a Metis" << std::endl
107  << "partitioner instead!" << std::endl;
108 
109  MetisPartitioner mp;
110 
111  mp.partition (mesh, n_sbdmns);
112 
113  // What to do if the Parmetis library IS present
114 #else
115 
116  // Revert to METIS on one processor.
117  if (mesh.n_processors() == 1)
118  {
119  MetisPartitioner mp;
120  mp.partition (mesh, n_sbdmns);
121  return;
122  }
123 
124  LOG_SCOPE("repartition()", "ParmetisPartitioner");
125 
126  // Initialize the data structures required by ParMETIS
127  this->initialize (mesh, n_sbdmns);
128 
129  // Make sure all processors have enough active local elements.
130  // Parmetis tends to crash when it's given only a couple elements
131  // per partition.
132  {
133  bool all_have_enough_elements = true;
134  for (std::size_t pid=0; pid<_n_active_elem_on_proc.size(); pid++)
136  all_have_enough_elements = false;
137 
138  // Parmetis will not work unless each processor has some
139  // elements. Specifically, it will abort when passed a NULL
140  // partition array on *any* of the processors.
141  if (!all_have_enough_elements)
142  {
143  // FIXME: revert to METIS, although this requires a serial mesh
144  MeshSerializer serialize(mesh);
145  MetisPartitioner mp;
146  mp.partition (mesh, n_sbdmns);
147  return;
148  }
149  }
150 
151  // build the graph corresponding to the mesh
152  this->build_graph (mesh);
153 
154 
155  // Partition the graph
156  std::vector<Parmetis::idx_t> vsize(_pmetis->vwgt.size(), 1);
157  Parmetis::real_t itr = 1000000.0;
158  MPI_Comm mpi_comm = mesh.comm().get();
159 
160  // Call the ParMETIS adaptive repartitioning method. This respects the
161  // original partitioning when computing the new partitioning so as to
162  // minimize the required data redistribution.
163  Parmetis::ParMETIS_V3_AdaptiveRepart(_pmetis->vtxdist.empty() ? libmesh_nullptr : &_pmetis->vtxdist[0],
164  _pmetis->xadj.empty() ? libmesh_nullptr : &_pmetis->xadj[0],
165  _pmetis->adjncy.empty() ? libmesh_nullptr : &_pmetis->adjncy[0],
166  _pmetis->vwgt.empty() ? libmesh_nullptr : &_pmetis->vwgt[0],
167  vsize.empty() ? libmesh_nullptr : &vsize[0],
169  &_pmetis->wgtflag,
170  &_pmetis->numflag,
171  &_pmetis->ncon,
172  &_pmetis->nparts,
173  _pmetis->tpwgts.empty() ? libmesh_nullptr : &_pmetis->tpwgts[0],
174  _pmetis->ubvec.empty() ? libmesh_nullptr : &_pmetis->ubvec[0],
175  &itr,
176  &_pmetis->options[0],
177  &_pmetis->edgecut,
178  _pmetis->part.empty() ? libmesh_nullptr : &_pmetis->part[0],
179  &mpi_comm);
180 
181  // Assign the returned processor ids
182  this->assign_partitioning (mesh);
183 
184 #endif // #ifndef LIBMESH_HAVE_PARMETIS ... else ...
185 
186 }
void single_partition(MeshBase &mesh)
Definition: partitioner.C:151
void initialize(const MeshBase &mesh, const unsigned int n_sbdmns)
MeshBase & mesh
const class libmesh_nullptr_t libmesh_nullptr
std::vector< dof_id_type > _n_active_elem_on_proc
UniquePtr< ParmetisHelper > _pmetis
void assign_partitioning(MeshBase &mesh)
void build_graph(const MeshBase &mesh)
OStreamProxy err(std::cerr)
const unsigned int MIN_ELEM_PER_PROC
void libMesh::ParmetisPartitioner::assign_partitioning ( MeshBase mesh)
private

Assign the computed partitioning to the mesh.

Definition at line 631 of file parmetis_partitioner.C.

References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::ParallelObject::comm(), 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 clone().

632 {
633  LOG_SCOPE("assign_partitioning()", "ParmetisPartitioner");
634 
635  // This function must be run on all processors at once
636  libmesh_parallel_only(mesh.comm());
637 
638  const dof_id_type
639  first_local_elem = _pmetis->vtxdist[mesh.processor_id()];
640 
641 #ifndef NDEBUG
642  const dof_id_type n_active_local_elem = mesh.n_active_local_elem();
643 #endif
644 
645  std::vector<std::vector<dof_id_type> >
646  requested_ids(mesh.n_processors()),
647  requests_to_fill(mesh.n_processors());
648 
649  MeshBase::element_iterator elem_it = mesh.active_elements_begin();
650  MeshBase::element_iterator elem_end = mesh.active_elements_end();
651 
652  for (; elem_it != elem_end; ++elem_it)
653  {
654  Elem * elem = *elem_it;
655 
656  // we need to get the index from the owning processor
657  // (note we cannot assign it now -- we are iterating
658  // over elements again and this will be bad!)
659  libmesh_assert_less (elem->processor_id(), requested_ids.size());
660  requested_ids[elem->processor_id()].push_back(elem->id());
661  }
662 
663  // Trade with all processors (including self) to get their indices
664  for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
665  {
666  // Trade my requests with processor procup and procdown
667  const processor_id_type procup = (mesh.processor_id() + pid) % mesh.n_processors();
668  const processor_id_type procdown = (mesh.n_processors() +
669  mesh.processor_id() - pid) % mesh.n_processors();
670 
671  mesh.comm().send_receive (procup, requested_ids[procup],
672  procdown, requests_to_fill[procdown]);
673 
674  // we can overwrite these requested ids in-place.
675  for (std::size_t i=0; i<requests_to_fill[procdown].size(); i++)
676  {
677  const dof_id_type requested_elem_index =
678  requests_to_fill[procdown][i];
679 
680  libmesh_assert(_global_index_by_pid_map.count(requested_elem_index));
681 
682  const dof_id_type global_index_by_pid =
683  _global_index_by_pid_map[requested_elem_index];
684 
685  const dof_id_type local_index =
686  global_index_by_pid - first_local_elem;
687 
688  libmesh_assert_less (local_index, _pmetis->part.size());
689  libmesh_assert_less (local_index, n_active_local_elem);
690 
691  const unsigned int elem_procid =
692  static_cast<unsigned int>(_pmetis->part[local_index]);
693 
694  libmesh_assert_less (elem_procid, static_cast<unsigned int>(_pmetis->nparts));
695 
696  requests_to_fill[procdown][i] = elem_procid;
697  }
698 
699  // Trade back
700  mesh.comm().send_receive (procdown, requests_to_fill[procdown],
701  procup, requested_ids[procup]);
702  }
703 
704  // and finally assign the partitioning.
705  // note we are iterating in exactly the same order
706  // used to build up the request, so we can expect the
707  // required entries to be in the proper sequence.
708  elem_it = mesh.active_elements_begin();
709  elem_end = mesh.active_elements_end();
710 
711  for (std::vector<unsigned int> counters(mesh.n_processors(), 0);
712  elem_it != elem_end; ++elem_it)
713  {
714  Elem * elem = *elem_it;
715 
716  const processor_id_type current_pid = elem->processor_id();
717 
718  libmesh_assert_less (counters[current_pid], requested_ids[current_pid].size());
719 
720  const processor_id_type elem_procid =
721  requested_ids[current_pid][counters[current_pid]++];
722 
723  libmesh_assert_less (elem_procid, static_cast<unsigned int>(_pmetis->nparts));
724  elem->processor_id() = elem_procid;
725  }
726 }
difference_type count(const key_type &key) const
Definition: vectormap.h:209
MeshBase & mesh
uint8_t processor_id_type
Definition: id_types.h:99
vectormap< dof_id_type, dof_id_type > _global_index_by_pid_map
libmesh_assert(j)
UniquePtr< ParmetisHelper > _pmetis
uint8_t dof_id_type
Definition: id_types.h:64
virtual void libMesh::Partitioner::attach_weights ( ErrorVector )
inlinevirtualinherited

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 169 of file partitioner.h.

References libMesh::Partitioner::_do_partition(), end, libMesh::Partitioner::single_partition(), and libMesh::Partitioner::single_partition_range().

169 { libmesh_not_implemented(); }
void libMesh::ParmetisPartitioner::build_graph ( const MeshBase mesh)
private

Build the graph.

Definition at line 407 of file parmetis_partitioner.C.

References libMesh::Elem::active(), libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::Elem::active_family_tree(), libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::Elem::dim(), libMesh::Elem::find_interior_neighbors(), libMesh::DofObject::id(), libMesh::Elem::interior_parent(), libMesh::libmesh_assert(), libmesh_nullptr, libMesh::MeshBase::n_active_local_elem(), libMesh::Elem::n_neighbors(), libMesh::Elem::neighbor_ptr(), libMesh::ParallelObject::processor_id(), and libMesh::Elem::which_neighbor_am_i().

Referenced by clone().

408 {
409  LOG_SCOPE("build_graph()", "ParmetisPartitioner");
410 
411  // build the graph in distributed CSR format. Note that
412  // the edges in the graph will correspond to
413  // face neighbors
414  const dof_id_type n_active_local_elem = mesh.n_active_local_elem();
415 
416  // If we have boundary elements in this mesh, we want to account for
417  // the connectivity between them and interior elements. We can find
418  // interior elements from boundary elements, but we need to build up
419  // a lookup map to do the reverse.
420 
421  typedef LIBMESH_BEST_UNORDERED_MULTIMAP<const Elem *, const Elem *>
422  map_type;
423  map_type interior_to_boundary_map;
424 
425  {
426  MeshBase::const_element_iterator elem_it = mesh.active_elements_begin();
427  const MeshBase::const_element_iterator elem_end = mesh.active_elements_end();
428 
429  for (; elem_it != elem_end; ++elem_it)
430  {
431  const Elem * elem = *elem_it;
432 
433  // If we don't have an interior_parent then there's nothing to look us
434  // up.
435  if ((elem->dim() >= LIBMESH_DIM) ||
436  !elem->interior_parent())
437  continue;
438 
439  // get all relevant interior elements
440  std::set<const Elem *> neighbor_set;
441  elem->find_interior_neighbors(neighbor_set);
442 
443  std::set<const Elem *>::iterator n_it = neighbor_set.begin();
444  for (; n_it != neighbor_set.end(); ++n_it)
445  {
446  // FIXME - non-const versions of the Elem set methods
447  // would be nice
448  Elem * neighbor = const_cast<Elem *>(*n_it);
449 
450 #if defined(LIBMESH_HAVE_UNORDERED_MULTIMAP) || \
451  defined(LIBMESH_HAVE_TR1_UNORDERED_MAP) || \
452  defined(LIBMESH_HAVE_HASH_MAP) || \
453  defined(LIBMESH_HAVE_EXT_HASH_MAP)
454  interior_to_boundary_map.insert
455  (std::make_pair(neighbor, elem));
456 #else
457  interior_to_boundary_map.insert
458  (interior_to_boundary_map.begin(),
459  std::make_pair(neighbor, elem));
460 #endif
461  }
462  }
463  }
464 
465 #ifdef LIBMESH_ENABLE_AMR
466  std::vector<const Elem *> neighbors_offspring;
467 #endif
468 
469  std::vector<std::vector<dof_id_type> > graph(n_active_local_elem);
470  dof_id_type graph_size=0;
471 
472  const dof_id_type first_local_elem = _pmetis->vtxdist[mesh.processor_id()];
473 
474  MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin();
475  const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end();
476 
477  for (; elem_it != elem_end; ++elem_it)
478  {
479  const Elem * elem = *elem_it;
480 
482  const dof_id_type global_index_by_pid =
483  _global_index_by_pid_map[elem->id()];
484 
485  const dof_id_type local_index =
486  global_index_by_pid - first_local_elem;
487  libmesh_assert_less (local_index, n_active_local_elem);
488 
489  std::vector<dof_id_type> & graph_row = graph[local_index];
490 
491  // Loop over the element's neighbors. An element
492  // adjacency corresponds to a face neighbor
493  for (unsigned int ms=0; ms<elem->n_neighbors(); ms++)
494  {
495  const Elem * neighbor = elem->neighbor_ptr(ms);
496 
497  if (neighbor != libmesh_nullptr)
498  {
499  // If the neighbor is active treat it
500  // as a connection
501  if (neighbor->active())
502  {
504  const dof_id_type neighbor_global_index_by_pid =
505  _global_index_by_pid_map[neighbor->id()];
506 
507  graph_row.push_back(neighbor_global_index_by_pid);
508  graph_size++;
509  }
510 
511 #ifdef LIBMESH_ENABLE_AMR
512 
513  // Otherwise we need to find all of the
514  // neighbor's children that are connected to
515  // us and add them
516  else
517  {
518  // The side of the neighbor to which
519  // we are connected
520  const unsigned int ns =
521  neighbor->which_neighbor_am_i (elem);
522  libmesh_assert_less (ns, neighbor->n_neighbors());
523 
524  // Get all the active children (& grandchildren, etc...)
525  // of the neighbor
526 
527  // FIXME - this is the wrong thing, since we
528  // should be getting the active family tree on
529  // our side only. But adding too many graph
530  // links may cause hanging nodes to tend to be
531  // on partition interiors, which would reduce
532  // communication overhead for constraint
533  // equations, so we'll leave it.
534 
535  neighbor->active_family_tree (neighbors_offspring);
536 
537  // Get all the neighbor's children that
538  // live on that side and are thus connected
539  // to us
540  for (std::size_t nc=0; nc<neighbors_offspring.size(); nc++)
541  {
542  const Elem * child =
543  neighbors_offspring[nc];
544 
545  // This does not assume a level-1 mesh.
546  // Note that since children have sides numbered
547  // coincident with the parent then this is a sufficient test.
548  if (child->neighbor_ptr(ns) == elem)
549  {
550  libmesh_assert (child->active());
552  const dof_id_type child_global_index_by_pid =
553  _global_index_by_pid_map[child->id()];
554 
555  graph_row.push_back(child_global_index_by_pid);
556  graph_size++;
557  }
558  }
559  }
560 
561 #endif /* ifdef LIBMESH_ENABLE_AMR */
562 
563 
564  }
565  }
566 
567  if ((elem->dim() < LIBMESH_DIM) &&
568  elem->interior_parent())
569  {
570  // get all relevant interior elements
571  std::set<const Elem *> neighbor_set;
572  elem->find_interior_neighbors(neighbor_set);
573 
574  std::set<const Elem *>::iterator n_it = neighbor_set.begin();
575  for (; n_it != neighbor_set.end(); ++n_it)
576  {
577  // FIXME - non-const versions of the Elem set methods
578  // would be nice
579  Elem * neighbor = const_cast<Elem *>(*n_it);
580 
581  const dof_id_type neighbor_global_index_by_pid =
582  _global_index_by_pid_map[neighbor->id()];
583 
584  graph_row.push_back(neighbor_global_index_by_pid);
585  graph_size++;
586  }
587  }
588 
589  // Check for any boundary neighbors
590  typedef map_type::iterator map_it_type;
591  std::pair<map_it_type, map_it_type>
592  bounds = interior_to_boundary_map.equal_range(elem);
593 
594  for (map_it_type it = bounds.first; it != bounds.second; ++it)
595  {
596  const Elem * neighbor = it->second;
597 
598  const dof_id_type neighbor_global_index_by_pid =
599  _global_index_by_pid_map[neighbor->id()];
600 
601  graph_row.push_back(neighbor_global_index_by_pid);
602  graph_size++;
603  }
604  }
605 
606  // Reserve space in the adjacency array
607  _pmetis->xadj.clear();
608  _pmetis->xadj.reserve (n_active_local_elem + 1);
609  _pmetis->adjncy.clear();
610  _pmetis->adjncy.reserve (graph_size);
611 
612  for (std::size_t r=0; r<graph.size(); r++)
613  {
614  _pmetis->xadj.push_back(_pmetis->adjncy.size());
615  std::vector<dof_id_type> graph_row; // build this emtpy
616  graph_row.swap(graph[r]); // this will deallocate at the end of scope
617  _pmetis->adjncy.insert(_pmetis->adjncy.end(),
618  graph_row.begin(),
619  graph_row.end());
620  }
621 
622  // The end of the adjacency array for the last elem
623  _pmetis->xadj.push_back(_pmetis->adjncy.size());
624 
625  libmesh_assert_equal_to (_pmetis->xadj.size(), n_active_local_elem+1);
626  libmesh_assert_equal_to (_pmetis->adjncy.size(), graph_size);
627 }
difference_type count(const key_type &key) const
Definition: vectormap.h:209
MeshBase & mesh
const class libmesh_nullptr_t libmesh_nullptr
vectormap< dof_id_type, dof_id_type > _global_index_by_pid_map
libmesh_assert(j)
UniquePtr< ParmetisHelper > _pmetis
uint8_t dof_id_type
Definition: id_types.h:64
virtual UniquePtr<Partitioner> libMesh::ParmetisPartitioner::clone ( ) const
inlinevirtual

Creates a new partitioner of this type and returns it in a UniquePtr.

Implements libMesh::Partitioner.

Definition at line 64 of file parmetis_partitioner.h.

References _do_partition(), _do_repartition(), assign_partitioning(), build_graph(), initialize(), mesh, and ParmetisPartitioner().

65  {
66  return UniquePtr<Partitioner>(new ParmetisPartitioner());
67  }
void libMesh::ParmetisPartitioner::initialize ( const MeshBase mesh,
const unsigned int  n_sbdmns 
)
private

Initialize data structures.

Definition at line 193 of file parmetis_partitioner.C.

References 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 clone().

195 {
196  LOG_SCOPE("initialize()", "ParmetisPartitioner");
197 
198  const dof_id_type n_active_local_elem = mesh.n_active_local_elem();
199 
200  // Set parameters.
201  _pmetis->wgtflag = 2; // weights on vertices only
202  _pmetis->ncon = 1; // one weight per vertex
203  _pmetis->numflag = 0; // C-style 0-based numbering
204  _pmetis->nparts = static_cast<Parmetis::idx_t>(n_sbdmns); // number of subdomains to create
205  _pmetis->edgecut = 0; // the numbers of edges cut by the
206  // partition
207 
208  // Initialize data structures for ParMETIS
209  _pmetis->vtxdist.assign (mesh.n_processors()+1, 0);
210  _pmetis->tpwgts.assign (_pmetis->nparts, 1./_pmetis->nparts);
211  _pmetis->ubvec.assign (_pmetis->ncon, 1.05);
212  _pmetis->part.assign (n_active_local_elem, 0);
213  _pmetis->options.resize (5);
214  _pmetis->vwgt.resize (n_active_local_elem);
215 
216  // Set the options
217  _pmetis->options[0] = 1; // don't use default options
218  _pmetis->options[1] = 0; // default (level of timing)
219  _pmetis->options[2] = 15; // random seed (default)
220  _pmetis->options[3] = 2; // processor distribution and subdomain distribution are decoupled
221 
222  // Find the number of active elements on each processor. We cannot use
223  // mesh.n_active_elem_on_proc(pid) since that only returns the number of
224  // elements assigned to pid which are currently stored on the calling
225  // processor. This will not in general be correct for parallel meshes
226  // when (pid!=mesh.processor_id()).
227  _n_active_elem_on_proc.resize(mesh.n_processors());
228  mesh.comm().allgather(n_active_local_elem, _n_active_elem_on_proc);
229 
230  // count the total number of active elements in the mesh. Note we cannot
231  // use mesh.n_active_elem() in general since this only returns the number
232  // of active elements which are stored on the calling processor.
233  // We should not use n_active_elem for any allocation because that will
234  // be inheritly unscalable, but it can be useful for libmesh_assertions.
235  dof_id_type n_active_elem=0;
236 
237  // Set up the vtxdist array. This will be the same on each processor.
238  // ***** Consult the Parmetis documentation. *****
239  libmesh_assert_equal_to (_pmetis->vtxdist.size(),
240  cast_int<std::size_t>(mesh.n_processors()+1));
241  libmesh_assert_equal_to (_pmetis->vtxdist[0], 0);
242 
243  for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
244  {
245  _pmetis->vtxdist[pid+1] = _pmetis->vtxdist[pid] + _n_active_elem_on_proc[pid];
246  n_active_elem += _n_active_elem_on_proc[pid];
247  }
248  libmesh_assert_equal_to (_pmetis->vtxdist.back(), static_cast<Parmetis::idx_t>(n_active_elem));
249 
250  // ParMetis expects the elements to be numbered in contiguous blocks
251  // by processor, i.e. [0, ne0), [ne0, ne0+ne1), ...
252  // Since we only partition active elements we should have no expectation
253  // that we currently have such a distribution. So we need to create it.
254  // Also, at the same time we are going to map all the active elements into a globally
255  // unique range [0,n_active_elem) which is *independent* of the current partitioning.
256  // This can be fed to ParMetis as the initial partitioning of the subdomains (decoupled
257  // from the partitioning of the objects themselves). This allows us to get the same
258  // resultant partitioning independed of the input partitioning.
259  BoundingBox bbox =
261 
262  _global_index_by_pid_map.clear();
263 
264  // Maps active element ids into a contiguous range independent of partitioning.
265  // (only needs local scope)
266  vectormap<dof_id_type, dof_id_type> global_index_map;
267 
268  {
269  std::vector<dof_id_type> global_index;
270 
271  // create the mapping which is contiguous by processor
272  dof_id_type pid_offset=0;
273  for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
274  {
275  MeshBase::const_element_iterator it = mesh.active_pid_elements_begin(pid);
276  const MeshBase::const_element_iterator end = mesh.active_pid_elements_end(pid);
277 
278  // note that we may not have all (or any!) the active elements which belong on this processor,
279  // but by calling this on all processors a unique range in [0,_n_active_elem_on_proc[pid])
280  // is constructed. Only the indices for the elements we pass in are returned in the array.
281  MeshCommunication().find_global_indices (mesh.comm(),
282  bbox, it, end,
283  global_index);
284 
285  for (dof_id_type cnt=0; it != end; ++it)
286  {
287  const Elem * elem = *it;
288  // vectormap::count forces a sort, which is too expensive
289  // in a loop
290  // libmesh_assert (!_global_index_by_pid_map.count(elem->id()));
291  libmesh_assert_less (cnt, global_index.size());
292  libmesh_assert_less (global_index[cnt], _n_active_elem_on_proc[pid]);
293 
294  _global_index_by_pid_map.insert(std::make_pair(elem->id(), global_index[cnt++] + pid_offset));
295  }
296 
297  pid_offset += _n_active_elem_on_proc[pid];
298  }
299 
300  // create the unique mapping for all active elements independent of partitioning
301  {
302  MeshBase::const_element_iterator it = mesh.active_elements_begin();
303  const MeshBase::const_element_iterator end = mesh.active_elements_end();
304 
305  // Calling this on all processors a unique range in [0,n_active_elem) is constructed.
306  // Only the indices for the elements we pass in are returned in the array.
307  MeshCommunication().find_global_indices (mesh.comm(),
308  bbox, it, end,
309  global_index);
310 
311  for (dof_id_type cnt=0; it != end; ++it)
312  {
313  const Elem * elem = *it;
314  // vectormap::count forces a sort, which is too expensive
315  // in a loop
316  // libmesh_assert (!global_index_map.count(elem->id()));
317  libmesh_assert_less (cnt, global_index.size());
318  libmesh_assert_less (global_index[cnt], n_active_elem);
319 
320  global_index_map.insert(std::make_pair(elem->id(), global_index[cnt++]));
321  }
322  }
323  // really, shouldn't be close!
324  libmesh_assert_less_equal (global_index_map.size(), n_active_elem);
325  libmesh_assert_less_equal (_global_index_by_pid_map.size(), n_active_elem);
326 
327  // At this point the two maps should be the same size. If they are not
328  // then the number of active elements is not the same as the sum over all
329  // processors of the number of active elements per processor, which means
330  // there must be some unpartitioned objects out there.
331  if (global_index_map.size() != _global_index_by_pid_map.size())
332  libmesh_error_msg("ERROR: ParmetisPartitioner cannot handle unpartitioned objects!");
333  }
334 
335  // Finally, we need to initialize the vertex (partition) weights and the initial subdomain
336  // mapping. The subdomain mapping will be independent of the processor mapping, and is
337  // defined by a simple mapping of the global indices we just found.
338  {
339  std::vector<dof_id_type> subdomain_bounds(mesh.n_processors());
340 
341  const dof_id_type first_local_elem = _pmetis->vtxdist[mesh.processor_id()];
342 
343  for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
344  {
345  dof_id_type tgt_subdomain_size = 0;
346 
347  // watch out for the case that n_subdomains < n_processors
348  if (pid < static_cast<unsigned int>(_pmetis->nparts))
349  {
350  tgt_subdomain_size = n_active_elem/std::min
351  (cast_int<Parmetis::idx_t>(mesh.n_processors()), _pmetis->nparts);
352 
353  if (pid < n_active_elem%_pmetis->nparts)
354  tgt_subdomain_size++;
355  }
356  if (pid == 0)
357  subdomain_bounds[0] = tgt_subdomain_size;
358  else
359  subdomain_bounds[pid] = subdomain_bounds[pid-1] + tgt_subdomain_size;
360  }
361 
362  libmesh_assert_equal_to (subdomain_bounds.back(), n_active_elem);
363 
364  MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin();
365  const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end();
366 
367  for (; elem_it != elem_end; ++elem_it)
368  {
369  const Elem * elem = *elem_it;
370 
372  const dof_id_type global_index_by_pid =
373  _global_index_by_pid_map[elem->id()];
374  libmesh_assert_less (global_index_by_pid, n_active_elem);
375 
376  const dof_id_type local_index =
377  global_index_by_pid - first_local_elem;
378 
379  libmesh_assert_less (local_index, n_active_local_elem);
380  libmesh_assert_less (local_index, _pmetis->vwgt.size());
381 
382  // TODO:[BSK] maybe there is a better weight?
383  _pmetis->vwgt[local_index] = elem->n_nodes();
384 
385  // find the subdomain this element belongs in
386  libmesh_assert (global_index_map.count(elem->id()));
387  const dof_id_type global_index =
388  global_index_map[elem->id()];
389 
390  libmesh_assert_less (global_index, subdomain_bounds.back());
391 
392  const unsigned int subdomain_id =
393  std::distance(subdomain_bounds.begin(),
394  std::lower_bound(subdomain_bounds.begin(),
395  subdomain_bounds.end(),
396  global_index));
397  libmesh_assert_less (subdomain_id, static_cast<unsigned int>(_pmetis->nparts));
398  libmesh_assert_less (local_index, _pmetis->part.size());
399 
400  _pmetis->part[local_index] = subdomain_id;
401  }
402  }
403 }
difference_type count(const key_type &key) const
Definition: vectormap.h:209
MeshBase & mesh
uint8_t processor_id_type
Definition: id_types.h:99
std::vector< dof_id_type > _n_active_elem_on_proc
IterBase * end
vectormap< dof_id_type, dof_id_type > _global_index_by_pid_map
libmesh_assert(j)
BoundingBox bounding_box(const MeshBase &mesh)
Definition: mesh_tools.C:358
UniquePtr< ParmetisHelper > _pmetis
void insert(const value_type &x)
Definition: vectormap.h:116
long double min(long double a, double b)
uint8_t dof_id_type
Definition: id_types.h:64
void libMesh::Partitioner::partition ( MeshBase mesh,
const unsigned int  n 
)
virtualinherited

Partitions the MeshBase into n parts by setting processor_id() on Nodes and Elems.

NOTE: If you are implementing a new type of Partitioner, you most likely do not want to override the partition() function, see instead the protected virtual _do_partition() method below. The partition() function is responsible for doing a lot of libmesh-internals-specific setup and finalization before and after the _do_partition() function is called. The only responsibility of the _do_partition() function, on the other hand, is to set the processor IDs of the elements according to a specific partitioning algorithm. See, e.g. MetisPartitioner for an example.

Definition at line 49 of file partitioner.C.

References libMesh::Partitioner::_do_partition(), libMesh::ParallelObject::comm(), 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_repartition(), libMesh::Partitioner::partition(), and libMesh::Partitioner::~Partitioner().

51 {
52  libmesh_parallel_only(mesh.comm());
53 
54  // BSK - temporary fix while redistribution is integrated 6/26/2008
55  // Uncomment this to not repartition in parallel
56  // if (!mesh.is_serial())
57  // return;
58 
59  // we cannot partition into more pieces than we have
60  // active elements!
61  const unsigned int n_parts =
62  static_cast<unsigned int>
63  (std::min(mesh.n_active_elem(), static_cast<dof_id_type>(n)));
64 
65  // Set the number of partitions in the mesh
66  mesh.set_n_partitions()=n_parts;
67 
68  if (n_parts == 1)
69  {
70  this->single_partition (mesh);
71  return;
72  }
73 
74  // First assign a temporary partitioning to any unpartitioned elements
76 
77  // Call the partitioning function
78  this->_do_partition(mesh,n_parts);
79 
80  // Set the parent's processor ids
82 
83  // Redistribute elements if necessary, before setting node processor
84  // ids, to make sure those will be set consistently
85  mesh.redistribute();
86 
87 #ifdef DEBUG
89 
90  // Messed up elem processor_id()s can leave us without the child
91  // elements we need to restrict vectors on a distributed mesh
92  MeshTools::libmesh_assert_valid_procids<Elem>(mesh);
93 #endif
94 
95  // Set the node's processor ids
97 
98 #ifdef DEBUG
99  MeshTools::libmesh_assert_valid_procids<Elem>(mesh);
100 #endif
101 
102  // Give derived Mesh classes a chance to update any cached data to
103  // reflect the new partitioning
104  mesh.update_post_partitioning();
105 }
void single_partition(MeshBase &mesh)
Definition: partitioner.C:151
void libmesh_assert_valid_remote_elems(const MeshBase &mesh)
Definition: mesh_tools.C:1040
static void set_node_processor_ids(MeshBase &mesh)
Definition: partitioner.C:431
MeshBase & mesh
virtual void _do_partition(MeshBase &mesh, const unsigned int n)=0
static void partition_unpartitioned_elements(MeshBase &mesh)
Definition: partitioner.C:175
static void set_parent_processor_ids(MeshBase &mesh)
Definition: partitioner.C:256
long double min(long double a, double b)
uint8_t dof_id_type
Definition: id_types.h:64
void libMesh::Partitioner::partition ( MeshBase mesh)
virtualinherited

Partitions the MeshBase into mesh.n_processors() by setting processor_id() on Nodes and Elems.

NOTE: If you are implementing a new type of Partitioner, you most likely do not want to override the partition() function, see instead the protected virtual _do_partition() method below. The partition() function is responsible for doing a lot of libmesh-internals-specific setup and finalization before and after the _do_partition() function is called. The only responsibility of the _do_partition() function, on the other hand, is to set the processor IDs of the elements according to a specific partitioning algorithm. See, e.g. MetisPartitioner for an example.

Definition at line 42 of file partitioner.C.

References libMesh::ParallelObject::n_processors(), and libMesh::Partitioner::partition().

43 {
44  this->partition(mesh,mesh.n_processors());
45 }
MeshBase & mesh
virtual void partition(MeshBase &mesh, const unsigned int n)
Definition: partitioner.C:49
virtual void libMesh::Partitioner::partition_range ( MeshBase ,
MeshBase::element_iterator  ,
MeshBase::element_iterator  ,
const unsigned int   
)
inlinevirtualinherited

Partitions elements in the range (it, end) into n parts. The mesh from which the iterators are created must also be passed in, since it is a parallel object and has other useful information in it.

Although partition_range() is part of the public Partitioner interface, it should not generally be called by applications. Its main purpose is to support the SubdomainPartitioner, which uses it internally to individually partition ranges of elements before combining them into the final partitioning. Most of the time, the protected _do_partition() function is implemented in terms of partition_range() by passing a range which includes all the elements of the Mesh.

Reimplemented in libMesh::CentroidPartitioner, libMesh::MappedSubdomainPartitioner, libMesh::SFCPartitioner, libMesh::LinearPartitioner, and libMesh::MetisPartitioner.

Definition at line 119 of file partitioner.h.

References libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::Partitioner::repartition(), libMesh::Partitioner::set_node_processor_ids(), and libMesh::Partitioner::set_parent_processor_ids().

123  { libmesh_not_implemented(); }
void libMesh::Partitioner::partition_unpartitioned_elements ( MeshBase mesh)
staticinherited

These functions assign processor IDs to newly-created elements (in parallel) which are currently assigned to processor 0.

Definition at line 175 of file partitioner.C.

References libMesh::ParallelObject::n_processors().

Referenced by libMesh::Partitioner::partition(), libMesh::Partitioner::partition_range(), and libMesh::Partitioner::repartition().

176 {
178 }
MeshBase & mesh
static void partition_unpartitioned_elements(MeshBase &mesh)
Definition: partitioner.C:175
void libMesh::Partitioner::partition_unpartitioned_elements ( MeshBase mesh,
const unsigned int  n 
)
staticinherited

Definition at line 182 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().

184 {
185  MeshBase::element_iterator it = mesh.unpartitioned_elements_begin();
186  const MeshBase::element_iterator end = mesh.unpartitioned_elements_end();
187 
188  const dof_id_type n_unpartitioned_elements = MeshTools::n_elem (it, end);
189 
190  // the unpartitioned elements must exist on all processors. If the range is empty on one
191  // it is empty on all, and we can quit right here.
192  if (!n_unpartitioned_elements) return;
193 
194  // find the target subdomain sizes
195  std::vector<dof_id_type> subdomain_bounds(mesh.n_processors());
196 
197  for (processor_id_type pid=0; pid<mesh.n_processors(); pid++)
198  {
199  dof_id_type tgt_subdomain_size = 0;
200 
201  // watch out for the case that n_subdomains < n_processors
202  if (pid < n_subdomains)
203  {
204  tgt_subdomain_size = n_unpartitioned_elements/n_subdomains;
205 
206  if (pid < n_unpartitioned_elements%n_subdomains)
207  tgt_subdomain_size++;
208 
209  }
210 
211  //libMesh::out << "pid, #= " << pid << ", " << tgt_subdomain_size << std::endl;
212  if (pid == 0)
213  subdomain_bounds[0] = tgt_subdomain_size;
214  else
215  subdomain_bounds[pid] = subdomain_bounds[pid-1] + tgt_subdomain_size;
216  }
217 
218  libmesh_assert_equal_to (subdomain_bounds.back(), n_unpartitioned_elements);
219 
220  // create the unique mapping for all unpartitioned elements independent of partitioning
221  // determine the global indexing for all the unpartitoned elements
222  std::vector<dof_id_type> global_indices;
223 
224  // Calling this on all processors a unique range in [0,n_unpartitioned_elements) is constructed.
225  // Only the indices for the elements we pass in are returned in the array.
226  MeshCommunication().find_global_indices (mesh.comm(),
228  global_indices);
229 
230  for (dof_id_type cnt=0; it != end; ++it)
231  {
232  Elem * elem = *it;
233 
234  libmesh_assert_less (cnt, global_indices.size());
235  const dof_id_type global_index =
236  global_indices[cnt++];
237 
238  libmesh_assert_less (global_index, subdomain_bounds.back());
239  libmesh_assert_less (global_index, n_unpartitioned_elements);
240 
241  const processor_id_type subdomain_id =
242  cast_int<processor_id_type>
243  (std::distance(subdomain_bounds.begin(),
244  std::upper_bound(subdomain_bounds.begin(),
245  subdomain_bounds.end(),
246  global_index)));
247  libmesh_assert_less (subdomain_id, n_subdomains);
248 
249  elem->processor_id() = subdomain_id;
250  //libMesh::out << "assigning " << global_index << " to " << subdomain_id << std::endl;
251  }
252 }
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Definition: mesh_tools.C:631
MeshBase & mesh
uint8_t processor_id_type
Definition: id_types.h:99
IterBase * end
BoundingBox bounding_box(const MeshBase &mesh)
Definition: mesh_tools.C:358
uint8_t dof_id_type
Definition: id_types.h:64
void libMesh::Partitioner::repartition ( MeshBase mesh,
const unsigned int  n 
)
inherited

Repartitions the MeshBase into n parts. (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 116 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::partition_range(), and libMesh::Partitioner::repartition().

118 {
119  // we cannot partition into more pieces than we have
120  // active elements!
121  const unsigned int n_parts =
122  static_cast<unsigned int>
123  (std::min(mesh.n_active_elem(), static_cast<dof_id_type>(n)));
124 
125  // Set the number of partitions in the mesh
126  mesh.set_n_partitions()=n_parts;
127 
128  if (n_parts == 1)
129  {
130  this->single_partition (mesh);
131  return;
132  }
133 
134  // First assign a temporary partitioning to any unpartitioned elements
136 
137  // Call the partitioning function
138  this->_do_repartition(mesh,n_parts);
139 
140  // Set the parent's processor ids
142 
143  // Set the node's processor ids
145 }
void single_partition(MeshBase &mesh)
Definition: partitioner.C:151
static void set_node_processor_ids(MeshBase &mesh)
Definition: partitioner.C:431
MeshBase & mesh
virtual void _do_repartition(MeshBase &mesh, const unsigned int n)
Definition: partitioner.h:201
static void partition_unpartitioned_elements(MeshBase &mesh)
Definition: partitioner.C:175
static void set_parent_processor_ids(MeshBase &mesh)
Definition: partitioner.C:256
long double min(long double a, double b)
uint8_t dof_id_type
Definition: id_types.h:64
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 109 of file partitioner.C.

References libMesh::ParallelObject::n_processors(), and libMesh::Partitioner::repartition().

110 {
111  this->repartition(mesh,mesh.n_processors());
112 }
MeshBase & mesh
void repartition(MeshBase &mesh, const unsigned int n)
Definition: partitioner.C:116
void libMesh::Partitioner::set_node_processor_ids ( MeshBase mesh)
staticinherited

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 431 of file partitioner.C.

References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::ParallelObject::comm(), libMesh::DofObject::id(), libMesh::DofObject::invalid_processor_id, libMesh::DofObject::invalidate_processor_id(), libMesh::libmesh_assert(), mesh, std::min(), libMesh::MeshTools::n_elem(), libMesh::Elem::n_nodes(), libMesh::ParallelObject::n_processors(), libMesh::Elem::node_ptr(), libMesh::MeshBase::node_ref(), 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::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::Partitioner::partition_range(), libMesh::XdrIO::read(), libMesh::Partitioner::repartition(), and libMesh::BoundaryInfo::sync().

432 {
433  LOG_SCOPE("set_node_processor_ids()","Partitioner");
434 
435  // This function must be run on all processors at once
436  libmesh_parallel_only(mesh.comm());
437 
438  // If we have any unpartitioned elements at this
439  // stage there is a problem
440  libmesh_assert (MeshTools::n_elem(mesh.unpartitioned_elements_begin(),
441  mesh.unpartitioned_elements_end()) == 0);
442 
443 
444  // const dof_id_type orig_n_local_nodes = mesh.n_local_nodes();
445 
446  // libMesh::err << "[" << mesh.processor_id() << "]: orig_n_local_nodes="
447  // << orig_n_local_nodes << std::endl;
448 
449  // Build up request sets. Each node is currently owned by a processor because
450  // it is connected to an element owned by that processor. However, during the
451  // repartitioning phase that element may have been assigned a new processor id, but
452  // it is still resident on the original processor. We need to know where to look
453  // for new ids before assigning new ids, otherwise we may be asking the wrong processors
454  // for the wrong information.
455  //
456  // The only remaining issue is what to do with unpartitioned nodes. Since they are required
457  // to live on all processors we can simply rely on ourselves to number them properly.
458  std::vector<std::vector<dof_id_type> >
459  requested_node_ids(mesh.n_processors());
460 
461  // Loop over all the nodes, count the ones on each processor. We can skip ourself
462  std::vector<dof_id_type> ghost_nodes_from_proc(mesh.n_processors(), 0);
463 
464  MeshBase::node_iterator node_it = mesh.nodes_begin();
465  const MeshBase::node_iterator node_end = mesh.nodes_end();
466 
467  for (; node_it != node_end; ++node_it)
468  {
469  Node * node = *node_it;
470  libmesh_assert(node);
471  const processor_id_type current_pid = node->processor_id();
472  if (current_pid != mesh.processor_id() &&
473  current_pid != DofObject::invalid_processor_id)
474  {
475  libmesh_assert_less (current_pid, ghost_nodes_from_proc.size());
476  ghost_nodes_from_proc[current_pid]++;
477  }
478  }
479 
480  // We know how many objects live on each processor, so reserve()
481  // space for each.
482  for (processor_id_type pid=0; pid != mesh.n_processors(); ++pid)
483  requested_node_ids[pid].reserve(ghost_nodes_from_proc[pid]);
484 
485  // We need to get the new pid for each node from the processor
486  // which *currently* owns the node. We can safely skip ourself
487  for (node_it = mesh.nodes_begin(); node_it != node_end; ++node_it)
488  {
489  Node * node = *node_it;
490  libmesh_assert(node);
491  const processor_id_type current_pid = node->processor_id();
492  if (current_pid != mesh.processor_id() &&
493  current_pid != DofObject::invalid_processor_id)
494  {
495  libmesh_assert_less (current_pid, requested_node_ids.size());
496  libmesh_assert_less (requested_node_ids[current_pid].size(),
497  ghost_nodes_from_proc[current_pid]);
498  requested_node_ids[current_pid].push_back(node->id());
499  }
500 
501  // Unset any previously-set node processor ids
502  node->invalidate_processor_id();
503  }
504 
505  // Loop over all the active elements
506  MeshBase::element_iterator elem_it = mesh.active_elements_begin();
507  const MeshBase::element_iterator elem_end = mesh.active_elements_end();
508 
509  for ( ; elem_it != elem_end; ++elem_it)
510  {
511  Elem * elem = *elem_it;
512  libmesh_assert(elem);
513 
514  libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id);
515 
516  // For each node, set the processor ID to the min of
517  // its current value and this Element's processor id.
518  //
519  // TODO: we would probably get better parallel partitioning if
520  // we did something like "min for even numbered nodes, max for
521  // odd numbered". We'd need to be careful about how that would
522  // affect solution ordering for I/O, though.
523  for (unsigned int n=0; n<elem->n_nodes(); ++n)
524  elem->node_ptr(n)->processor_id() = std::min(elem->node_ptr(n)->processor_id(),
525  elem->processor_id());
526  }
527 
528  // And loop over the subactive elements, but don't reassign
529  // nodes that are already active on another processor.
530  MeshBase::element_iterator sub_it = mesh.subactive_elements_begin();
531  const MeshBase::element_iterator sub_end = mesh.subactive_elements_end();
532 
533  for ( ; sub_it != sub_end; ++sub_it)
534  {
535  Elem * elem = *sub_it;
536  libmesh_assert(elem);
537 
538  libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id);
539 
540  for (unsigned int n=0; n<elem->n_nodes(); ++n)
541  if (elem->node_ptr(n)->processor_id() == DofObject::invalid_processor_id)
542  elem->node_ptr(n)->processor_id() = elem->processor_id();
543  }
544 
545  // Same for the inactive elements -- we will have already gotten most of these
546  // nodes, *except* for the case of a parent with a subset of children which are
547  // ghost elements. In that case some of the parent nodes will not have been
548  // properly handled yet
549  MeshBase::element_iterator not_it = mesh.not_active_elements_begin();
550  const MeshBase::element_iterator not_end = mesh.not_active_elements_end();
551 
552  for ( ; not_it != not_end; ++not_it)
553  {
554  Elem * elem = *not_it;
555  libmesh_assert(elem);
556 
557  libmesh_assert_not_equal_to (elem->processor_id(), DofObject::invalid_processor_id);
558 
559  for (unsigned int n=0; n<elem->n_nodes(); ++n)
560  if (elem->node_ptr(n)->processor_id() == DofObject::invalid_processor_id)
561  elem->node_ptr(n)->processor_id() = elem->processor_id();
562  }
563 
564  // We can't assert that all nodes are connected to elements, because
565  // a DistributedMesh with NodeConstraints might have pulled in some
566  // remote nodes solely for evaluating those constraints.
567  // MeshTools::libmesh_assert_connected_nodes(mesh);
568 
569  // For such nodes, we'll do a sanity check later when making sure
570  // that we successfully reset their processor ids to something
571  // valid.
572 
573  // Next set node ids from other processors, excluding self
574  for (processor_id_type p=1; p != mesh.n_processors(); ++p)
575  {
576  // Trade my requests with processor procup and procdown
577  processor_id_type procup = cast_int<processor_id_type>
578  ((mesh.processor_id() + p) % mesh.n_processors());
579  processor_id_type procdown = cast_int<processor_id_type>
580  ((mesh.n_processors() + mesh.processor_id() - p) %
581  mesh.n_processors());
582  std::vector<dof_id_type> request_to_fill;
583  mesh.comm().send_receive(procup, requested_node_ids[procup],
584  procdown, request_to_fill);
585 
586  // Fill those requests in-place
587  for (std::size_t i=0; i != request_to_fill.size(); ++i)
588  {
589  Node & node = mesh.node_ref(request_to_fill[i]);
590  const processor_id_type new_pid = node.processor_id();
591 
592  // We may have an invalid processor_id() on nodes that have been
593  // "detatched" from coarsened-away elements but that have not yet
594  // themselves been removed.
595  // libmesh_assert_not_equal_to (new_pid, DofObject::invalid_processor_id);
596  // libmesh_assert_less (new_pid, mesh.n_partitions()); // this is the correct test --
597  request_to_fill[i] = new_pid; // the number of partitions may
598  } // not equal the number of processors
599 
600  // Trade back the results
601  std::vector<dof_id_type> filled_request;
602  mesh.comm().send_receive(procdown, request_to_fill,
603  procup, filled_request);
604  libmesh_assert_equal_to (filled_request.size(), requested_node_ids[procup].size());
605 
606  // And copy the id changes we've now been informed of
607  for (std::size_t i=0; i != filled_request.size(); ++i)
608  {
609  Node & node = mesh.node_ref(requested_node_ids[procup][i]);
610 
611  // this is the correct test -- the number of partitions may
612  // not equal the number of processors
613 
614  // But: we may have an invalid processor_id() on nodes that
615  // have been "detatched" from coarsened-away elements but
616  // that have not yet themselves been removed.
617  // libmesh_assert_less (filled_request[i], mesh.n_partitions());
618 
619  node.processor_id(cast_int<processor_id_type>(filled_request[i]));
620  }
621  }
622 
623 #ifdef DEBUG
624  MeshTools::libmesh_assert_valid_procids<Node>(mesh);
625 #endif
626 }
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Definition: mesh_tools.C:631
MeshBase & mesh
uint8_t processor_id_type
Definition: id_types.h:99
libmesh_assert(j)
static const processor_id_type invalid_processor_id
Definition: dof_object.h:346
long double min(long double a, double b)
void libMesh::Partitioner::set_parent_processor_ids ( MeshBase mesh)
staticinherited

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 256 of file partitioner.C.

References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::Elem::active_family_tree(), libMesh::MeshBase::ancestor_elements_begin(), libMesh::MeshBase::ancestor_elements_end(), libMesh::Elem::child_ptr(), libMesh::ParallelObject::comm(), libMesh::Partitioner::communication_blocksize, end, libMesh::DofObject::id(), libMesh::DofObject::invalid_processor_id, libMesh::DofObject::invalidate_processor_id(), libMesh::Elem::is_remote(), libMesh::MeshBase::is_serial(), libMesh::libmesh_assert(), libMesh::libmesh_ignore(), libMesh::MeshBase::max_elem_id(), std::min(), libMesh::Parallel::Communicator::min(), libMesh::Elem::n_children(), libMesh::MeshTools::n_elem(), libMesh::Elem::parent(), libMesh::processor_id(), libMesh::DofObject::processor_id(), libMesh::Elem::total_family_tree(), libMesh::MeshBase::unpartitioned_elements_begin(), and libMesh::MeshBase::unpartitioned_elements_end().

Referenced by libMesh::Partitioner::partition(), libMesh::Partitioner::partition_range(), and libMesh::Partitioner::repartition().

257 {
258  // Ignore the parameter when !LIBMESH_ENABLE_AMR
260 
261  LOG_SCOPE("set_parent_processor_ids()", "Partitioner");
262 
263 #ifdef LIBMESH_ENABLE_AMR
264 
265  // If the mesh is serial we have access to all the elements,
266  // in particular all the active ones. We can therefore set
267  // the parent processor ids indirecly through their children, and
268  // set the subactive processor ids while examining their active
269  // ancestors.
270  // By convention a parent is assigned to the minimum processor
271  // of all its children, and a subactive is assigned to the processor
272  // of its active ancestor.
273  if (mesh.is_serial())
274  {
275  // Loop over all the active elements in the mesh
276  MeshBase::element_iterator it = mesh.active_elements_begin();
277  const MeshBase::element_iterator end = mesh.active_elements_end();
278 
279  for ( ; it!=end; ++it)
280  {
281  Elem * child = *it;
282 
283  // First set descendents
284 
285  std::vector<const Elem *> subactive_family;
286  child->total_family_tree(subactive_family);
287  for (std::size_t i = 0; i != subactive_family.size(); ++i)
288  const_cast<Elem *>(subactive_family[i])->processor_id() = child->processor_id();
289 
290  // Then set ancestors
291 
292  Elem * parent = child->parent();
293 
294  while (parent)
295  {
296  // invalidate the parent id, otherwise the min below
297  // will not work if the current parent id is less
298  // than all the children!
299  parent->invalidate_processor_id();
300 
301  for (unsigned int c=0; c<parent->n_children(); c++)
302  {
303  child = parent->child_ptr(c);
304  libmesh_assert(child);
305  libmesh_assert(!child->is_remote());
306  libmesh_assert_not_equal_to (child->processor_id(), DofObject::invalid_processor_id);
307  parent->processor_id() = std::min(parent->processor_id(),
308  child->processor_id());
309  }
310  parent = parent->parent();
311  }
312  }
313  }
314 
315  // When the mesh is parallel we cannot guarantee that parents have access to
316  // all their children.
317  else
318  {
319  // Setting subactive processor ids is easy: we can guarantee
320  // that children have access to all their parents.
321 
322  // Loop over all the active elements in the mesh
323  MeshBase::element_iterator it = mesh.active_elements_begin();
324  const MeshBase::element_iterator end = mesh.active_elements_end();
325 
326  for ( ; it!=end; ++it)
327  {
328  Elem * child = *it;
329 
330  std::vector<const Elem *> subactive_family;
331  child->total_family_tree(subactive_family);
332  for (std::size_t i = 0; i != subactive_family.size(); ++i)
333  const_cast<Elem *>(subactive_family[i])->processor_id() = child->processor_id();
334  }
335 
336  // When the mesh is parallel we cannot guarantee that parents have access to
337  // all their children.
338 
339  // We will use a brute-force approach here. Each processor finds its parent
340  // elements and sets the parent pid to the minimum of its
341  // semilocal descendants.
342  // A global reduction is then performed to make sure the true minimum is found.
343  // As noted, this is required because we cannot guarantee that a parent has
344  // access to all its children on any single processor.
345  libmesh_parallel_only(mesh.comm());
346  libmesh_assert(MeshTools::n_elem(mesh.unpartitioned_elements_begin(),
347  mesh.unpartitioned_elements_end()) == 0);
348 
349  const dof_id_type max_elem_id = mesh.max_elem_id();
350 
351  std::vector<processor_id_type>
352  parent_processor_ids (std::min(communication_blocksize,
353  max_elem_id));
354 
355  for (dof_id_type blk=0, last_elem_id=0; last_elem_id<max_elem_id; blk++)
356  {
357  last_elem_id =
358  std::min(static_cast<dof_id_type>((blk+1)*communication_blocksize),
359  max_elem_id);
360  const dof_id_type first_elem_id = blk*communication_blocksize;
361 
362  std::fill (parent_processor_ids.begin(),
363  parent_processor_ids.end(),
365 
366  // first build up local contributions to parent_processor_ids
367  MeshBase::element_iterator not_it = mesh.ancestor_elements_begin();
368  const MeshBase::element_iterator not_end = mesh.ancestor_elements_end();
369 
370  bool have_parent_in_block = false;
371 
372  for ( ; not_it != not_end; ++not_it)
373  {
374  Elem * parent = *not_it;
375 
376  const dof_id_type parent_idx = parent->id();
377  libmesh_assert_less (parent_idx, max_elem_id);
378 
379  if ((parent_idx >= first_elem_id) &&
380  (parent_idx < last_elem_id))
381  {
382  have_parent_in_block = true;
384 
385  std::vector<const Elem *> active_family;
386  parent->active_family_tree(active_family);
387  for (std::size_t i = 0; i != active_family.size(); ++i)
388  parent_pid = std::min (parent_pid, active_family[i]->processor_id());
389 
390  const dof_id_type packed_idx = parent_idx - first_elem_id;
391  libmesh_assert_less (packed_idx, parent_processor_ids.size());
392 
393  parent_processor_ids[packed_idx] = parent_pid;
394  }
395  }
396 
397  // then find the global minimum
398  mesh.comm().min (parent_processor_ids);
399 
400  // and assign the ids, if we have a parent in this block.
401  if (have_parent_in_block)
402  for (not_it = mesh.ancestor_elements_begin();
403  not_it != not_end; ++not_it)
404  {
405  Elem * parent = *not_it;
406 
407  const dof_id_type parent_idx = parent->id();
408 
409  if ((parent_idx >= first_elem_id) &&
410  (parent_idx < last_elem_id))
411  {
412  const dof_id_type packed_idx = parent_idx - first_elem_id;
413  libmesh_assert_less (packed_idx, parent_processor_ids.size());
414 
415  const processor_id_type parent_pid =
416  parent_processor_ids[packed_idx];
417 
418  libmesh_assert_not_equal_to (parent_pid, DofObject::invalid_processor_id);
419 
420  parent->processor_id() = parent_pid;
421  }
422  }
423  }
424  }
425 
426 #endif // LIBMESH_ENABLE_AMR
427 }
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Definition: mesh_tools.C:631
MeshBase & mesh
uint8_t processor_id_type
Definition: id_types.h:99
IterBase * end
libmesh_assert(j)
static const processor_id_type invalid_processor_id
Definition: dof_object.h:346
void libmesh_ignore(const T &)
static const dof_id_type communication_blocksize
Definition: partitioner.h:208
processor_id_type processor_id()
Definition: libmesh_base.h:96
long double min(long double a, double b)
uint8_t dof_id_type
Definition: id_types.h:64
void libMesh::Partitioner::single_partition ( MeshBase mesh)
protectedinherited

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 151 of file partitioner.C.

References libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), and libMesh::Partitioner::single_partition_range().

Referenced by libMesh::SubdomainPartitioner::_do_partition(), libMesh::Partitioner::attach_weights(), libMesh::Partitioner::partition(), and libMesh::Partitioner::repartition().

152 {
153  this->single_partition_range(mesh.elements_begin(),
154  mesh.elements_end());
155 }
MeshBase & mesh
void single_partition_range(MeshBase::element_iterator it, MeshBase::element_iterator end)
Definition: partitioner.C:159
void libMesh::Partitioner::single_partition_range ( MeshBase::element_iterator  it,
MeshBase::element_iterator  end 
)
protectedinherited

Slightly generalized version of single_partition which acts on a range of elements defined by the pair of iterators (it, end).

Definition at line 159 of file partitioner.C.

References end, libMesh::Elem::n_nodes(), libMesh::Elem::node_ptr(), and libMesh::DofObject::processor_id().

Referenced by libMesh::Partitioner::attach_weights(), libMesh::LinearPartitioner::partition_range(), libMesh::MappedSubdomainPartitioner::partition_range(), libMesh::CentroidPartitioner::partition_range(), and libMesh::Partitioner::single_partition().

161 {
162  LOG_SCOPE("single_partition_range()", "Partitioner");
163 
164  for ( ; it != end; ++it)
165  {
166  Elem * elem = *it;
167  elem->processor_id() = 0;
168 
169  // Assign all this element's nodes to processor 0 as well.
170  for (unsigned int n=0; n<elem->n_nodes(); ++n)
171  elem->node_ptr(n)->processor_id() = 0;
172  }
173 }
IterBase * end

Member Data Documentation

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 118 of file parmetis_partitioner.h.

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 113 of file parmetis_partitioner.h.

UniquePtr<ParmetisHelper> libMesh::ParmetisPartitioner::_pmetis
private

Pointer to the Parmetis-specific data structures. Lets us avoid including parmetis.h here.

Definition at line 124 of file parmetis_partitioner.h.

ErrorVector* libMesh::Partitioner::_weights
protectedinherited

The weights that might be used for partitioning.

Definition at line 213 of file partitioner.h.

Referenced by libMesh::MetisPartitioner::attach_weights().

const dof_id_type libMesh::Partitioner::communication_blocksize = 1000000
staticprotectedinherited

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 208 of file partitioner.h.

Referenced by libMesh::Partitioner::set_parent_processor_ids().


The documentation for this class was generated from the following files: