libMesh
parmetis_partitioner.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // Local includes
21 #include "libmesh/parmetis_partitioner.h"
22 
23 // libMesh includes
24 #include "libmesh/libmesh_config.h"
25 #include "libmesh/mesh_base.h"
26 #include "libmesh/mesh_serializer.h"
27 #include "libmesh/mesh_tools.h"
28 #include "libmesh/mesh_communication.h"
29 #include "libmesh/metis_partitioner.h"
30 #include "libmesh/parallel_only.h"
31 #include "libmesh/libmesh_logging.h"
32 #include "libmesh/elem.h"
33 #include "libmesh/parmetis_helper.h"
34 
35 // TIMPI includes
36 #include "timpi/communicator.h" // also includes mpi.h
37 #include "timpi/parallel_implementation.h" // for min()
38 
39 // Include the ParMETIS header file.
40 #ifdef LIBMESH_HAVE_PARMETIS
41 
42 // Before we include a header wrapped in a namespace, we'd better make
43 // sure none of its dependencies end up in that namespace
44 #include <mpi.h>
45 
46 namespace Parmetis {
47 extern "C" {
48 # include "libmesh/ignore_warnings.h"
49 # include "parmetis.h"
50 # include "libmesh/restore_warnings.h"
51 }
52 }
53 
54 #endif
55 
56 
57 // C++ includes
58 #include <unordered_map>
59 
60 
61 namespace libMesh
62 {
63 
64 // Minimum elements on each processor required for us to choose
65 // Parmetis over Metis.
66 #ifdef LIBMESH_HAVE_PARMETIS
67 const unsigned int MIN_ELEM_PER_PROC = 4;
68 #endif
69 
70 // ------------------------------------------------------------
71 // ParmetisPartitioner implementation
73 #ifdef LIBMESH_HAVE_PARMETIS
74  : _pmetis(libmesh_make_unique<ParmetisHelper>())
75 #endif
76 {}
77 
78 
79 
80 ParmetisPartitioner::ParmetisPartitioner (const ParmetisPartitioner & other)
81  : Partitioner(other)
82 #ifdef LIBMESH_HAVE_PARMETIS
83  , _pmetis(libmesh_make_unique<ParmetisHelper>(*(other._pmetis)))
84 #endif
85 {
86 }
87 
88 
89 
90 ParmetisPartitioner::~ParmetisPartitioner() = default;
91 
92 
93 
94 void ParmetisPartitioner::_do_partition (MeshBase & mesh,
95  const unsigned int n_sbdmns)
96 {
97  this->_do_repartition (mesh, n_sbdmns);
98 }
99 
100 
101 
102 void ParmetisPartitioner::_do_repartition (MeshBase & mesh,
103  const unsigned int n_sbdmns)
104 {
105  // This function must be run on all processors at once
106  libmesh_parallel_only(mesh.comm());
107 
108  // Check for easy returns
109  if (!mesh.n_elem())
110  return;
111 
112  if (n_sbdmns == 1)
113  {
114  this->single_partition(mesh);
115  return;
116  }
117 
118  libmesh_assert_greater (n_sbdmns, 0);
119 
120  // What to do if the Parmetis library IS NOT present
121 #ifndef LIBMESH_HAVE_PARMETIS
122 
123  libmesh_do_once(
124  libMesh::out << "ERROR: The library has been built without" << std::endl
125  << "Parmetis support. Using a Metis" << std::endl
126  << "partitioner instead!" << std::endl;);
127 
128  MetisPartitioner mp;
129 
130  // Don't just call partition() here; that would end up calling
131  // post-element-partitioning work redundantly (and at the moment
132  // incorrectly)
133  mp.partition_range (mesh, mesh.active_elements_begin(),
134  mesh.active_elements_end(), n_sbdmns);
135 
136  // What to do if the Parmetis library IS present
137 #else
138 
139  // Revert to METIS on one processor.
140  if (mesh.n_processors() == 1)
141  {
142  // Make sure the mesh knows it's serial
143  mesh.allgather();
144 
145  MetisPartitioner mp;
146  // Don't just call partition() here; that would end up calling
147  // post-element-partitioning work redundantly (and at the moment
148  // incorrectly)
149  mp.partition_range (mesh, mesh.active_elements_begin(),
150  mesh.active_elements_end(), n_sbdmns);
151  return;
152  }
153 
154  LOG_SCOPE("repartition()", "ParmetisPartitioner");
155 
156  // Initialize the data structures required by ParMETIS
157  this->initialize (mesh, n_sbdmns);
158 
159  // build the graph corresponding to the mesh
160  this->build_graph (mesh);
161 
162  // Make sure all processors have enough active local elements and
163  // enough connectivity among them.
164  // Parmetis tends to die when it's given only a couple elements
165  // per partition or when it can't reach elements from each other.
166  {
167  bool ready_for_parmetis = true;
168  for (const auto & nelem : _n_active_elem_on_proc)
169  if (nelem < MIN_ELEM_PER_PROC)
170  ready_for_parmetis = false;
171 
172  std::size_t my_adjacency = _pmetis->adjncy.size();
173  mesh.comm().min(my_adjacency);
174  if (!my_adjacency)
175  ready_for_parmetis = false;
176 
177  // Parmetis will not work unless each processor has some
178  // elements. Specifically, it will abort when passed a nullptr
179  // partition or adjacency array on *any* of the processors.
180  if (!ready_for_parmetis)
181  {
182  // FIXME: revert to METIS, although this requires a serial mesh
183  MeshSerializer serialize(mesh);
184  MetisPartitioner mp;
185  mp.partition (mesh, n_sbdmns);
186  return;
187  }
188  }
189 
190 
191  // Partition the graph
192  std::vector<Parmetis::idx_t> vsize(_pmetis->vwgt.size(), 1);
193  Parmetis::real_t itr = 1000000.0;
194  MPI_Comm mpi_comm = mesh.comm().get();
195 
196  // Call the ParMETIS adaptive repartitioning method. This respects the
197  // original partitioning when computing the new partitioning so as to
198  // minimize the required data redistribution.
199  Parmetis::ParMETIS_V3_AdaptiveRepart(_pmetis->vtxdist.empty() ? nullptr : _pmetis->vtxdist.data(),
200  _pmetis->xadj.empty() ? nullptr : _pmetis->xadj.data(),
201  _pmetis->adjncy.empty() ? nullptr : _pmetis->adjncy.data(),
202  _pmetis->vwgt.empty() ? nullptr : _pmetis->vwgt.data(),
203  vsize.empty() ? nullptr : vsize.data(),
204  nullptr,
205  &_pmetis->wgtflag,
206  &_pmetis->numflag,
207  &_pmetis->ncon,
208  &_pmetis->nparts,
209  _pmetis->tpwgts.empty() ? nullptr : _pmetis->tpwgts.data(),
210  _pmetis->ubvec.empty() ? nullptr : _pmetis->ubvec.data(),
211  &itr,
212  _pmetis->options.data(),
213  &_pmetis->edgecut,
214  _pmetis->part.empty() ? nullptr : reinterpret_cast<Parmetis::idx_t *>(_pmetis->part.data()),
215  &mpi_comm);
216 
217  // Assign the returned processor ids
218  this->assign_partitioning (mesh, _pmetis->part);
219 
220 #endif // #ifndef LIBMESH_HAVE_PARMETIS ... else ...
221 
222 }
223 
224 
225 
226 // Only need to compile these methods if ParMETIS is present
227 #ifdef LIBMESH_HAVE_PARMETIS
228 
230  const unsigned int n_sbdmns)
231 {
232  LOG_SCOPE("initialize()", "ParmetisPartitioner");
233 
234  const dof_id_type n_active_local_elem = mesh.n_active_local_elem();
235  // Set parameters.
236  _pmetis->wgtflag = 2; // weights on vertices only
237  _pmetis->ncon = 1; // one weight per vertex
238  _pmetis->numflag = 0; // C-style 0-based numbering
239  _pmetis->nparts = static_cast<Parmetis::idx_t>(n_sbdmns); // number of subdomains to create
240  _pmetis->edgecut = 0; // the numbers of edges cut by the
241  // partition
242 
243  // Initialize data structures for ParMETIS
244  _pmetis->vtxdist.assign (mesh.n_processors()+1, 0);
245  _pmetis->tpwgts.assign (_pmetis->nparts, 1./_pmetis->nparts);
246  _pmetis->ubvec.assign (_pmetis->ncon, 1.05);
247  _pmetis->part.assign (n_active_local_elem, 0);
248  _pmetis->options.resize (5);
249  _pmetis->vwgt.resize (n_active_local_elem);
250 
251  // Set the options
252  _pmetis->options[0] = 1; // don't use default options
253  _pmetis->options[1] = 0; // default (level of timing)
254  _pmetis->options[2] = 15; // random seed (default)
255  _pmetis->options[3] = 2; // processor distribution and subdomain distribution are decoupled
256 
257  // ParMetis expects the elements to be numbered in contiguous blocks
258  // by processor, i.e. [0, ne0), [ne0, ne0+ne1), ...
259  // Since we only partition active elements we should have no expectation
260  // that we currently have such a distribution. So we need to create it.
261  // Also, at the same time we are going to map all the active elements into a globally
262  // unique range [0,n_active_elem) which is *independent* of the current partitioning.
263  // This can be fed to ParMetis as the initial partitioning of the subdomains (decoupled
264  // from the partitioning of the objects themselves). This allows us to get the same
265  // resultant partitioning independent of the input partitioning.
266  libMesh::BoundingBox bbox =
268 
269  _find_global_index_by_pid_map(mesh);
270 
271 
272  // count the total number of active elements in the mesh. Note we cannot
273  // use mesh.n_active_elem() in general since this only returns the number
274  // of active elements which are stored on the calling processor.
275  // We should not use n_active_elem for any allocation because that will
276  // be inherently unscalable, but it can be useful for libmesh_assertions.
277  dof_id_type n_active_elem=0;
278 
279  // Set up the vtxdist array. This will be the same on each processor.
280  // ***** Consult the Parmetis documentation. *****
281  libmesh_assert_equal_to (_pmetis->vtxdist.size(),
282  cast_int<std::size_t>(mesh.n_processors()+1));
283  libmesh_assert_equal_to (_pmetis->vtxdist[0], 0);
284 
285  for (auto pid : IntRange<processor_id_type>(0, mesh.n_processors()))
286  {
287  _pmetis->vtxdist[pid+1] = _pmetis->vtxdist[pid] + _n_active_elem_on_proc[pid];
288  n_active_elem += _n_active_elem_on_proc[pid];
289  }
290  libmesh_assert_equal_to (_pmetis->vtxdist.back(), static_cast<Parmetis::idx_t>(n_active_elem));
291 
292 
293  // Maps active element ids into a contiguous range independent of partitioning.
294  // (only needs local scope)
295  std::unordered_map<dof_id_type, dof_id_type> global_index_map;
296 
297  {
298  std::vector<dof_id_type> global_index;
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.
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 (auto pid : IntRange<processor_id_type>(0, mesh.n_processors()))
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  for (const auto & elem : mesh.active_local_element_ptr_range())
365  {
366  libmesh_assert (_global_index_by_pid_map.count(elem->id()));
367  const dof_id_type global_index_by_pid =
368  _global_index_by_pid_map[elem->id()];
369  libmesh_assert_less (global_index_by_pid, n_active_elem);
370 
371  const dof_id_type local_index =
372  global_index_by_pid - first_local_elem;
373 
374  libmesh_assert_less (local_index, n_active_local_elem);
375  libmesh_assert_less (local_index, _pmetis->vwgt.size());
376 
377  // TODO:[BSK] maybe there is a better weight?
378  _pmetis->vwgt[local_index] = elem->n_nodes();
379 
380  // find the subdomain this element belongs in
381  libmesh_assert (global_index_map.count(elem->id()));
382  const dof_id_type global_index =
383  global_index_map[elem->id()];
384 
385  libmesh_assert_less (global_index, subdomain_bounds.back());
386 
387  const unsigned int subdomain_id =
388  cast_int<unsigned int>
389  (std::distance(subdomain_bounds.begin(),
390  std::lower_bound(subdomain_bounds.begin(),
391  subdomain_bounds.end(),
392  global_index)));
393  libmesh_assert_less (subdomain_id, _pmetis->nparts);
394  libmesh_assert_less (local_index, _pmetis->part.size());
395 
396  _pmetis->part[local_index] = subdomain_id;
397  }
398  }
399 }
400 
401 
402 
403 void ParmetisPartitioner::build_graph (const MeshBase & mesh)
404 {
405  LOG_SCOPE("build_graph()", "ParmetisPartitioner");
406 
407  // build the graph in distributed CSR format. Note that
408  // the edges in the graph will correspond to
409  // face neighbors
410  const dof_id_type n_active_local_elem = mesh.n_active_local_elem();
411 
412  Partitioner::build_graph(mesh);
413 
414  dof_id_type graph_size=0;
415 
416  for (auto & row: _dual_graph)
417  graph_size += cast_int<dof_id_type>(row.size());
418 
419  // Reserve space in the adjacency array
420  _pmetis->xadj.clear();
421  _pmetis->xadj.reserve (n_active_local_elem + 1);
422  _pmetis->adjncy.clear();
423  _pmetis->adjncy.reserve (graph_size);
424 
425  for (auto & graph_row : _dual_graph)
426  {
427  _pmetis->xadj.push_back(cast_int<int>(_pmetis->adjncy.size()));
428  _pmetis->adjncy.insert(_pmetis->adjncy.end(),
429  graph_row.begin(),
430  graph_row.end());
431  }
432 
433  // The end of the adjacency array for the last elem
434  _pmetis->xadj.push_back(cast_int<int>(_pmetis->adjncy.size()));
435 
436  libmesh_assert_equal_to (_pmetis->xadj.size(), n_active_local_elem+1);
437  libmesh_assert_equal_to (_pmetis->adjncy.size(), graph_size);
438 }
439 
440 #endif // #ifdef LIBMESH_HAVE_PARMETIS
441 
442 } // namespace libMesh
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
Parmetis
Definition: parmetis_helper.h:35
libMesh::BoundingBox
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
end
IterBase * end
Also have a polymorphic pointer to the end object, this prevents iterating past the end.
Definition: variant_filter_iterator.h:343
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::Partitioner::partition
virtual void partition(MeshBase &mesh, const unsigned int n)
Partitions the MeshBase into n parts by setting processor_id() on Nodes and Elems.
Definition: partitioner.C:64
libMesh::MeshCommunication::find_global_indices
void find_global_indices(const Parallel::Communicator &communicator, const libMesh::BoundingBox &, const ForwardIterator &, const ForwardIterator &, std::vector< dof_id_type > &) const
This method determines a globally unique, partition-agnostic index for each object in the input range...
Definition: mesh_communication_global_indices.C:710
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::MeshTools::create_bounding_box
libMesh::BoundingBox create_bounding_box(const MeshBase &mesh)
The same functionality as the deprecated MeshTools::bounding_box().
Definition: mesh_tools.C:389
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::MetisPartitioner::partition_range
virtual void partition_range(MeshBase &mesh, MeshBase::element_iterator it, MeshBase::element_iterator end, const unsigned int n) override
Called by the SubdomainPartitioner to partition elements in the range (it, end).
Definition: metis_partitioner.C:58
libMesh::MeshSerializer
Temporarily serialize a DistributedMesh for output; a distributed mesh is allgathered by the MeshSeri...
Definition: mesh_serializer.h:42
libMesh::MetisPartitioner
The MetisPartitioner uses the Metis graph partitioner to partition the elements.
Definition: metis_partitioner.h:38
parallel_implementation.h
initialize
void initialize(EquationSystems &es, const std::string &system_name)
Definition: main.C:52
libMesh::MeshBase::const_element_iterator
The definition of the const_element_iterator struct.
Definition: mesh_base.h:1891
libMesh::ParmetisPartitioner::ParmetisPartitioner
ParmetisPartitioner()
Default and copy ctors.
libMesh::MIN_ELEM_PER_PROC
const unsigned int MIN_ELEM_PER_PROC
Definition: parmetis_partitioner.C:67
distance
Real distance(const Point &p)
Definition: subdomains_ex3.C:50
libMesh::DofObject::id
dof_id_type id() const
Definition: dof_object.h:767
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::MeshCommunication
This is the MeshCommunication class.
Definition: mesh_communication.h:50
communicator.h
libMesh::out
OStreamProxy out