Reference documentation for deal.II version 9.1.0-pre
sparsity_tools.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/exceptions.h>
18 #include <deal.II/base/std_cxx14/memory.h>
19 
20 #include <deal.II/lac/exceptions.h>
21 #include <deal.II/lac/sparsity_pattern.h>
22 #include <deal.II/lac/sparsity_tools.h>
23 
24 #include <algorithm>
25 #include <functional>
26 #include <set>
27 
28 #ifdef DEAL_II_WITH_MPI
29 # include <deal.II/base/mpi.h>
30 # include <deal.II/base/utilities.h>
31 
32 # include <deal.II/lac/block_sparsity_pattern.h>
33 # include <deal.II/lac/dynamic_sparsity_pattern.h>
34 #endif
35 
36 #ifdef DEAL_II_WITH_METIS
37 extern "C"
38 {
39 # include <metis.h>
40 }
41 #endif
42 
43 #ifdef DEAL_II_TRILINOS_WITH_ZOLTAN
44 # include <zoltan_cpp.h>
45 #endif
46 
47 #include <string>
48 
49 
50 DEAL_II_NAMESPACE_OPEN
51 
52 namespace SparsityTools
53 {
54  namespace
55  {
56  void
57  partition_metis(const SparsityPattern & sparsity_pattern,
58  const std::vector<unsigned int> &cell_weights,
59  const unsigned int n_partitions,
60  std::vector<unsigned int> & partition_indices)
61  {
62  // Make sure that METIS is actually
63  // installed and detected
64 #ifndef DEAL_II_WITH_METIS
65  (void)sparsity_pattern;
66  (void)cell_weights;
67  (void)n_partitions;
68  (void)partition_indices;
70 #else
71 
72  // generate the data structures for
73  // METIS. Note that this is particularly
74  // simple, since METIS wants exactly our
75  // compressed row storage format. we only
76  // have to set up a few auxiliary arrays
77  idx_t n = static_cast<signed int>(sparsity_pattern.n_rows()),
78  ncon = 1, // number of balancing constraints (should be >0)
79  nparts =
80  static_cast<int>(n_partitions), // number of subdomains to create
81  dummy; // the numbers of edges cut by the
82  // resulting partition
83 
84  // We can not partition n items into more than n parts. METIS will
85  // generate non-sensical output (everything is owned by a single process)
86  // and complain with a message (but won't return an error code!):
87  // ***Cannot bisect a graph with 0 vertices!
88  // ***You are trying to partition a graph into too many parts!
89  nparts = std::min(n, nparts);
90 
91  // use default options for METIS
92  idx_t options[METIS_NOPTIONS];
93  METIS_SetDefaultOptions(options);
94 
95  // one more nuisance: we have to copy our own data to arrays that store
96  // signed integers :-(
97  std::vector<idx_t> int_rowstart(1);
98  int_rowstart.reserve(sparsity_pattern.n_rows() + 1);
99  std::vector<idx_t> int_colnums;
100  int_colnums.reserve(sparsity_pattern.n_nonzero_elements());
101  for (SparsityPattern::size_type row = 0; row < sparsity_pattern.n_rows();
102  ++row)
103  {
104  for (SparsityPattern::iterator col = sparsity_pattern.begin(row);
105  col < sparsity_pattern.end(row);
106  ++col)
107  int_colnums.push_back(col->column());
108  int_rowstart.push_back(int_colnums.size());
109  }
110 
111  std::vector<idx_t> int_partition_indices(sparsity_pattern.n_rows());
112 
113  // Setup cell weighting option
114  std::vector<idx_t> int_cell_weights;
115  if (cell_weights.size() > 0)
116  {
117  Assert(cell_weights.size() == sparsity_pattern.n_rows(),
118  ExcDimensionMismatch(cell_weights.size(),
119  sparsity_pattern.n_rows()));
120  int_cell_weights.resize(cell_weights.size());
121  std::copy(cell_weights.begin(),
122  cell_weights.end(),
123  int_cell_weights.begin());
124  }
125  // Set a pointer to the optional cell weighting information.
126  // METIS expects a null pointer if there are no weights to be considered.
127  idx_t *const p_int_cell_weights =
128  (cell_weights.size() > 0 ? int_cell_weights.data() : nullptr);
129 
130 
131  // Make use of METIS' error code.
132  int ierr;
133 
134  // Select which type of partitioning to create
135 
136  // Use recursive if the number of partitions is less than or equal to 8
137  if (nparts <= 8)
138  ierr = METIS_PartGraphRecursive(&n,
139  &ncon,
140  int_rowstart.data(),
141  int_colnums.data(),
142  p_int_cell_weights,
143  nullptr,
144  nullptr,
145  &nparts,
146  nullptr,
147  nullptr,
148  options,
149  &dummy,
150  int_partition_indices.data());
151 
152  // Otherwise use kway
153  else
154  ierr = METIS_PartGraphKway(&n,
155  &ncon,
156  int_rowstart.data(),
157  int_colnums.data(),
158  p_int_cell_weights,
159  nullptr,
160  nullptr,
161  &nparts,
162  nullptr,
163  nullptr,
164  options,
165  &dummy,
166  int_partition_indices.data());
167 
168  // If metis returns normally, an error code METIS_OK=1 is returned from
169  // the above functions (see metish.h)
170  AssertThrow(ierr == 1, ExcMETISError(ierr));
171 
172  // now copy back generated indices into the output array
173  std::copy(int_partition_indices.begin(),
174  int_partition_indices.end(),
175  partition_indices.begin());
176 #endif
177  }
178 
179 
180 // Query functions unused if zoltan is not installed
181 #ifdef DEAL_II_TRILINOS_WITH_ZOLTAN
182  // Query functions for partition_zoltan
183  int
184  get_number_of_objects(void *data, int *ierr)
185  {
186  SparsityPattern *graph = reinterpret_cast<SparsityPattern *>(data);
187 
188  *ierr = ZOLTAN_OK;
189 
190  return graph->n_rows();
191  }
192 
193 
194  void
195  get_object_list(void *data,
196  int /*sizeGID*/,
197  int /*sizeLID*/,
198  ZOLTAN_ID_PTR globalID,
199  ZOLTAN_ID_PTR localID,
200  int /*wgt_dim*/,
201  float * /*obj_wgts*/,
202  int *ierr)
203  {
204  SparsityPattern *graph = reinterpret_cast<SparsityPattern *>(data);
205  *ierr = ZOLTAN_OK;
206 
207  Assert(globalID != nullptr, ExcInternalError());
208  Assert(localID != nullptr, ExcInternalError());
209 
210  // set global degrees of freedom
211  auto n_dofs = graph->n_rows();
212 
213  for (unsigned int i = 0; i < n_dofs; i++)
214  {
215  globalID[i] = i;
216  localID[i] = i; // Same as global ids.
217  }
218  }
219 
220 
221  void
222  get_num_edges_list(void *data,
223  int /*sizeGID*/,
224  int /*sizeLID*/,
225  int num_obj,
226  ZOLTAN_ID_PTR globalID,
227  ZOLTAN_ID_PTR /*localID*/,
228  int *numEdges,
229  int *ierr)
230  {
231  SparsityPattern *graph = reinterpret_cast<SparsityPattern *>(data);
232 
233  *ierr = ZOLTAN_OK;
234 
235  Assert(numEdges != nullptr, ExcInternalError());
236 
237  for (int i = 0; i < num_obj; ++i)
238  {
239  if (graph->exists(i, i)) // Check if diagonal element is present
240  numEdges[i] = graph->row_length(globalID[i]) - 1;
241  else
242  numEdges[i] = graph->row_length(globalID[i]);
243  }
244  }
245 
246 
247 
248  void
249  get_edge_list(void *data,
250  int /*sizeGID*/,
251  int /*sizeLID*/,
252  int num_obj,
253  ZOLTAN_ID_PTR /*globalID*/,
254  ZOLTAN_ID_PTR /*localID*/,
255  int * /*num_edges*/,
256  ZOLTAN_ID_PTR nborGID,
257  int * nborProc,
258  int /*wgt_dim*/,
259  float * /*ewgts*/,
260  int *ierr)
261  {
262  SparsityPattern *graph = reinterpret_cast<SparsityPattern *>(data);
263  *ierr = ZOLTAN_OK;
264 
265  ZOLTAN_ID_PTR nextNborGID = nborGID;
266  int * nextNborProc = nborProc;
267 
268  // Loop through rows corresponding to indices in globalID implicitly
269  for (SparsityPattern::size_type i = 0;
270  i < static_cast<SparsityPattern::size_type>(num_obj);
271  ++i)
272  {
273  // Loop through each column to find neighbours
274  for (SparsityPattern::iterator col = graph->begin(i);
275  col < graph->end(i);
276  ++col)
277  // Ignore diagonal entries. Not needed for partitioning.
278  if (i != col->column())
279  {
280  Assert(nextNborGID != nullptr, ExcInternalError());
281  Assert(nextNborProc != nullptr, ExcInternalError());
282 
283  *nextNborGID++ = col->column();
284  *nextNborProc++ = 0; // All the vertices on processor 0
285  }
286  }
287  }
288 #endif
289 
290 
291  void
292  partition_zoltan(const SparsityPattern & sparsity_pattern,
293  const std::vector<unsigned int> &cell_weights,
294  const unsigned int n_partitions,
295  std::vector<unsigned int> & partition_indices)
296  {
297  // Make sure that ZOLTAN is actually
298  // installed and detected
299 #ifndef DEAL_II_TRILINOS_WITH_ZOLTAN
300  (void)sparsity_pattern;
301  (void)cell_weights;
302  (void)n_partitions;
303  (void)partition_indices;
305 #else
306 
307  Assert(
308  cell_weights.size() == 0,
309  ExcMessage(
310  "The cell weighting functionality for Zoltan has not yet been implemented."));
311  (void)cell_weights;
312 
313  // MPI environment must have been initialized by this point.
314  std::unique_ptr<Zoltan> zz =
315  std_cxx14::make_unique<Zoltan>(MPI_COMM_SELF);
316 
317  // General parameters
318  // DEBUG_LEVEL call must precede the call to LB_METHOD
319  zz->Set_Param("DEBUG_LEVEL", "0"); // set level of debug info
320  zz->Set_Param(
321  "LB_METHOD",
322  "GRAPH"); // graph based partition method (LB-load balancing)
323  zz->Set_Param("NUM_LOCAL_PARTS",
324  std::to_string(n_partitions)); // set number of partitions
325 
326  // The PHG partitioner is a hypergraph partitioner that Zoltan could use
327  // for graph partitioning.
328  // If number of vertices in hyperedge divided by total vertices in
329  // hypergraph exceeds PHG_EDGE_SIZE_THRESHOLD,
330  // then the hyperedge will be omitted as such (dense) edges will likely
331  // incur high communication costs regardless of the partition.
332  // PHG_EDGE_SIZE_THRESHOLD value is raised to 0.5 from the default
333  // value of 0.25 so that the PHG partitioner doesn't throw warning saying
334  // "PHG_EDGE_SIZE_THRESHOLD is low ..." after removing all dense edges.
335  // For instance, in two dimensions if the triangulation being partitioned
336  // is two quadrilaterals sharing an edge and if PHG_EDGE_SIZE_THRESHOLD
337  // value is set to 0.25, PHG will remove all the edges throwing the
338  // above warning.
339  zz->Set_Param("PHG_EDGE_SIZE_THRESHOLD", "0.5");
340 
341  // Need a non-const object equal to sparsity_pattern
342  SparsityPattern graph;
343  graph.copy_from(sparsity_pattern);
344 
345  // Set query functions
346  zz->Set_Num_Obj_Fn(get_number_of_objects, &graph);
347  zz->Set_Obj_List_Fn(get_object_list, &graph);
348  zz->Set_Num_Edges_Multi_Fn(get_num_edges_list, &graph);
349  zz->Set_Edge_List_Multi_Fn(get_edge_list, &graph);
350 
351  // Variables needed by partition function
352  int changes = 0;
353  int num_gid_entries = 1;
354  int num_lid_entries = 1;
355  int num_import = 0;
356  ZOLTAN_ID_PTR import_global_ids = nullptr;
357  ZOLTAN_ID_PTR import_local_ids = nullptr;
358  int * import_procs = nullptr;
359  int * import_to_part = nullptr;
360  int num_export = 0;
361  ZOLTAN_ID_PTR export_global_ids = nullptr;
362  ZOLTAN_ID_PTR export_local_ids = nullptr;
363  int * export_procs = nullptr;
364  int * export_to_part = nullptr;
365 
366  // call partitioner
367  const int rc = zz->LB_Partition(changes,
368  num_gid_entries,
369  num_lid_entries,
370  num_import,
371  import_global_ids,
372  import_local_ids,
373  import_procs,
374  import_to_part,
375  num_export,
376  export_global_ids,
377  export_local_ids,
378  export_procs,
379  export_to_part);
380  (void)rc;
381 
382  // check for error code in partitioner
383  Assert(rc == ZOLTAN_OK, ExcInternalError());
384 
385  // By default, all indices belong to part 0. After zoltan partition
386  // some are migrated to different part ID, which is stored in
387  // export_to_part array.
388  std::fill(partition_indices.begin(), partition_indices.end(), 0);
389 
390  // copy from export_to_part to partition_indices, whose part_ids != 0.
391  for (int i = 0; i < num_export; i++)
392  partition_indices[export_local_ids[i]] = export_to_part[i];
393 #endif
394  }
395  } // namespace
396 
397 
398  void
399  partition(const SparsityPattern & sparsity_pattern,
400  const unsigned int n_partitions,
401  std::vector<unsigned int> &partition_indices,
402  const Partitioner partitioner)
403  {
404  std::vector<unsigned int> cell_weights;
405 
406  // Call the other more general function
407  partition(sparsity_pattern,
408  cell_weights,
409  n_partitions,
410  partition_indices,
411  partitioner);
412  }
413 
414 
415  void
416  partition(const SparsityPattern & sparsity_pattern,
417  const std::vector<unsigned int> &cell_weights,
418  const unsigned int n_partitions,
419  std::vector<unsigned int> & partition_indices,
420  const Partitioner partitioner)
421  {
422  Assert(sparsity_pattern.n_rows() == sparsity_pattern.n_cols(),
423  ExcNotQuadratic());
424  Assert(sparsity_pattern.is_compressed(),
426 
427  Assert(n_partitions > 0, ExcInvalidNumberOfPartitions(n_partitions));
428  Assert(partition_indices.size() == sparsity_pattern.n_rows(),
429  ExcInvalidArraySize(partition_indices.size(),
430  sparsity_pattern.n_rows()));
431 
432  // check for an easy return
433  if (n_partitions == 1 || (sparsity_pattern.n_rows() == 1))
434  {
435  std::fill_n(partition_indices.begin(), partition_indices.size(), 0U);
436  return;
437  }
438 
439  if (partitioner == Partitioner::metis)
440  partition_metis(sparsity_pattern,
441  cell_weights,
442  n_partitions,
443  partition_indices);
444  else if (partitioner == Partitioner::zoltan)
445  partition_zoltan(sparsity_pattern,
446  cell_weights,
447  n_partitions,
448  partition_indices);
449  else
450  AssertThrow(false, ExcInternalError());
451  }
452 
453 
454  unsigned int
455  color_sparsity_pattern(const SparsityPattern & sparsity_pattern,
456  std::vector<unsigned int> &color_indices)
457  {
458  // Make sure that ZOLTAN is actually
459  // installed and detected
460 #ifndef DEAL_II_TRILINOS_WITH_ZOLTAN
461  (void)sparsity_pattern;
462  (void)color_indices;
464  return 0;
465 #else
466  // coloring algorithm is run in serial by each processor.
467  std::unique_ptr<Zoltan> zz = std_cxx14::make_unique<Zoltan>(MPI_COMM_SELF);
468 
469  // Coloring parameters
470  // DEBUG_LEVEL must precede all other calls
471  zz->Set_Param("DEBUG_LEVEL", "0"); // level of debug info
472  zz->Set_Param("COLORING_PROBLEM", "DISTANCE-1"); // Standard coloring
473  zz->Set_Param("NUM_GID_ENTRIES", "1"); // 1 entry represents global ID
474  zz->Set_Param("NUM_LID_ENTRIES", "1"); // 1 entry represents local ID
475  zz->Set_Param("OBJ_WEIGHT_DIM", "0"); // object weights not used
476  zz->Set_Param("RECOLORING_NUM_OF_ITERATIONS", "0");
477 
478  // Zoltan::Color function requires a non-const SparsityPattern object
479  SparsityPattern graph;
480  graph.copy_from(sparsity_pattern);
481 
482  // Set query functions required by coloring function
483  zz->Set_Num_Obj_Fn(get_number_of_objects, &graph);
484  zz->Set_Obj_List_Fn(get_object_list, &graph);
485  zz->Set_Num_Edges_Multi_Fn(get_num_edges_list, &graph);
486  zz->Set_Edge_List_Multi_Fn(get_edge_list, &graph);
487 
488  // Variables needed by coloring function
489  int num_gid_entries = 1;
490  const int num_objects = graph.n_rows();
491 
492  // Preallocate input variables. Element type fixed by ZOLTAN.
493  std::vector<ZOLTAN_ID_TYPE> global_ids(num_objects);
494  std::vector<int> color_exp(num_objects);
495 
496  // Set ids for which coloring needs to be done
497  for (int i = 0; i < num_objects; i++)
498  global_ids[i] = i;
499 
500  // Call ZOLTAN coloring algorithm
501  int rc = zz->Color(num_gid_entries,
502  num_objects,
503  global_ids.data(),
504  color_exp.data());
505 
506  (void)rc;
507  // Check for error code
508  Assert(rc == ZOLTAN_OK, ExcInternalError());
509 
510  // Allocate and assign color indices
511  color_indices.resize(num_objects);
512  Assert(color_exp.size() == color_indices.size(),
513  ExcDimensionMismatch(color_exp.size(), color_indices.size()));
514 
515  std::copy(color_exp.begin(), color_exp.end(), color_indices.begin());
516 
517  unsigned int n_colors =
518  *(std::max_element(color_indices.begin(), color_indices.end()));
519  return n_colors;
520 #endif
521  }
522 
523 
524  namespace internal
525  {
532  find_unnumbered_starting_index(
533  const DynamicSparsityPattern & sparsity,
534  const std::vector<DynamicSparsityPattern::size_type> &new_indices)
535  {
536  DynamicSparsityPattern::size_type starting_point =
538  DynamicSparsityPattern::size_type min_coordination = sparsity.n_rows();
539  for (DynamicSparsityPattern::size_type row = 0; row < sparsity.n_rows();
540  ++row)
541  // look over all as-yet unnumbered indices
542  if (new_indices[row] == numbers::invalid_size_type)
543  {
544  if (sparsity.row_length(row) < min_coordination)
545  {
546  min_coordination = sparsity.row_length(row);
547  starting_point = row;
548  }
549  }
550 
551  // now we still have to care for the case that no unnumbered dof has a
552  // coordination number less than sparsity.n_rows(). this rather exotic
553  // case only happens if we only have one cell, as far as I can see,
554  // but there may be others as well.
555  //
556  // if that should be the case, we can chose an arbitrary dof as
557  // starting point, e.g. the first unnumbered one
558  if (starting_point == numbers::invalid_size_type)
559  {
560  for (DynamicSparsityPattern::size_type i = 0; i < new_indices.size();
561  ++i)
562  if (new_indices[i] == numbers::invalid_size_type)
563  {
564  starting_point = i;
565  break;
566  }
567 
568  Assert(starting_point != numbers::invalid_size_type,
569  ExcInternalError());
570  }
571 
572  return starting_point;
573  }
574  } // namespace internal
575 
576 
577 
578  void
580  const DynamicSparsityPattern & sparsity,
581  std::vector<DynamicSparsityPattern::size_type> & new_indices,
582  const std::vector<DynamicSparsityPattern::size_type> &starting_indices)
583  {
584  Assert(sparsity.n_rows() == sparsity.n_cols(),
585  ExcDimensionMismatch(sparsity.n_rows(), sparsity.n_cols()));
586  Assert(sparsity.n_rows() == new_indices.size(),
587  ExcDimensionMismatch(sparsity.n_rows(), new_indices.size()));
588  Assert(starting_indices.size() <= sparsity.n_rows(),
589  ExcMessage(
590  "You can't specify more starting indices than there are rows"));
591  Assert(sparsity.row_index_set().size() == 0 ||
592  sparsity.row_index_set().size() == sparsity.n_rows(),
593  ExcMessage(
594  "Only valid for sparsity patterns which store all rows."));
595  for (SparsityPattern::size_type i = 0; i < starting_indices.size(); ++i)
596  Assert(starting_indices[i] < sparsity.n_rows(),
597  ExcMessage("Invalid starting index: All starting indices need "
598  "to be between zero and the number of rows in the "
599  "sparsity pattern."));
600 
601  // store the indices of the dofs renumbered in the last round. Default to
602  // starting points
603  std::vector<DynamicSparsityPattern::size_type> last_round_dofs(
604  starting_indices);
605 
606  // initialize the new_indices array with invalid values
607  std::fill(new_indices.begin(),
608  new_indices.end(),
610 
611  // if no starting indices were given: find dof with lowest coordination
612  // number
613  if (last_round_dofs.empty())
614  last_round_dofs.push_back(
615  internal::find_unnumbered_starting_index(sparsity, new_indices));
616 
617  // store next free dof index
618  DynamicSparsityPattern::size_type next_free_number = 0;
619 
620  // enumerate the first round dofs
621  for (DynamicSparsityPattern::size_type i = 0; i != last_round_dofs.size();
622  ++i)
623  new_indices[last_round_dofs[i]] = next_free_number++;
624 
625  // now do as many steps as needed to renumber all dofs
626  while (true)
627  {
628  // store the indices of the dofs to be renumbered in the next round
629  std::vector<DynamicSparsityPattern::size_type> next_round_dofs;
630 
631  // find all neighbors of the dofs numbered in the last round
633  i < last_round_dofs.size();
634  ++i)
636  sparsity.begin(last_round_dofs[i]);
637  j < sparsity.end(last_round_dofs[i]);
638  ++j)
639  next_round_dofs.push_back(j->column());
640 
641  // sort dof numbers
642  std::sort(next_round_dofs.begin(), next_round_dofs.end());
643 
644  // delete multiple entries
645  std::vector<DynamicSparsityPattern::size_type>::iterator end_sorted;
646  end_sorted =
647  std::unique(next_round_dofs.begin(), next_round_dofs.end());
648  next_round_dofs.erase(end_sorted, next_round_dofs.end());
649 
650  // eliminate dofs which are already numbered
651  for (int s = next_round_dofs.size() - 1; s >= 0; --s)
652  if (new_indices[next_round_dofs[s]] != numbers::invalid_size_type)
653  next_round_dofs.erase(next_round_dofs.begin() + s);
654 
655  // check whether there are any new dofs in the list. if there are
656  // none, then we have completely numbered the current component of the
657  // graph. check if there are as yet unnumbered components of the graph
658  // that we would then have to do next
659  if (next_round_dofs.empty())
660  {
661  if (std::find(new_indices.begin(),
662  new_indices.end(),
663  numbers::invalid_size_type) == new_indices.end())
664  // no unnumbered indices, so we can leave now
665  break;
666 
667  // otherwise find a valid starting point for the next component of
668  // the graph and continue with numbering that one. we only do so
669  // if no starting indices were provided by the user (see the
670  // documentation of this function) so produce an error if we got
671  // here and starting indices were given
672  Assert(starting_indices.empty(),
673  ExcMessage("The input graph appears to have more than one "
674  "component, but as stated in the documentation "
675  "we only want to reorder such graphs if no "
676  "starting indices are given. The function was "
677  "called with starting indices, however."))
678 
679  next_round_dofs.push_back(
680  internal::find_unnumbered_starting_index(sparsity,
681  new_indices));
682  }
683 
684 
685 
686  // store for each coordination number the dofs with these coordination
687  // number
688  std::multimap<DynamicSparsityPattern::size_type, int>
689  dofs_by_coordination;
690 
691  // find coordination number for each of these dofs
692  for (std::vector<DynamicSparsityPattern::size_type>::iterator s =
693  next_round_dofs.begin();
694  s != next_round_dofs.end();
695  ++s)
696  {
697  const DynamicSparsityPattern::size_type coordination =
698  sparsity.row_length(*s);
699 
700  // insert this dof at its coordination number
701  const std::pair<const DynamicSparsityPattern::size_type, int>
702  new_entry(coordination, *s);
703  dofs_by_coordination.insert(new_entry);
704  }
705 
706  // assign new DoF numbers to the elements of the present front:
707  std::multimap<DynamicSparsityPattern::size_type, int>::iterator i;
708  for (i = dofs_by_coordination.begin(); i != dofs_by_coordination.end();
709  ++i)
710  new_indices[i->second] = next_free_number++;
711 
712  // after that: copy this round's dofs for the next round
713  last_round_dofs = next_round_dofs;
714  }
715 
716  // test for all indices numbered. this mostly tests whether the
717  // front-marching-algorithm (which Cuthill-McKee actually is) has reached
718  // all points.
719  Assert((std::find(new_indices.begin(),
720  new_indices.end(),
721  numbers::invalid_size_type) == new_indices.end()) &&
722  (next_free_number == sparsity.n_rows()),
723  ExcInternalError());
724  }
725 
726 
727 
728  namespace internal
729  {
730  void
732  const DynamicSparsityPattern & connectivity,
733  std::vector<DynamicSparsityPattern::size_type> &renumbering)
734  {
735  AssertDimension(connectivity.n_rows(), connectivity.n_cols());
736  AssertDimension(connectivity.n_rows(), renumbering.size());
737  Assert(connectivity.row_index_set().size() == 0 ||
738  connectivity.row_index_set().size() == connectivity.n_rows(),
739  ExcMessage(
740  "Only valid for sparsity patterns which store all rows."));
741 
742  std::vector<types::global_dof_index> touched_nodes(
743  connectivity.n_rows(), numbers::invalid_dof_index);
744  std::vector<unsigned int> row_lengths(connectivity.n_rows());
745  std::set<types::global_dof_index> current_neighbors;
746  std::vector<std::vector<types::global_dof_index>> groups;
747 
748  // First collect the number of neighbors for each node. We use this
749  // field to find next nodes with the minimum number of non-touched
750  // neighbors in the field n_remaining_neighbors, so we will count down
751  // on this field. We also cache the row lengths because we need this
752  // data frequently and getting it from the sparsity pattern is more
753  // expensive.
754  for (types::global_dof_index row = 0; row < connectivity.n_rows(); ++row)
755  {
756  row_lengths[row] = connectivity.row_length(row);
757  Assert(row_lengths[row] > 0, ExcInternalError());
758  }
759  std::vector<unsigned int> n_remaining_neighbors(row_lengths);
760 
761  // This outer loop is typically traversed only once, unless the global
762  // graph is not connected
763  while (true)
764  {
765  // Find cell with the minimal number of neighbors (typically a
766  // corner node when based on FEM meshes). If no cell is left, we are
767  // done. Together with the outer while loop, this loop can possibly
768  // be of quadratic complexity in the number of disconnected
769  // partitions, i.e. up to connectivity.n_rows() in the worst case,
770  // but that is not the usual use case of this loop and thus not
771  // optimized for.
772  std::pair<types::global_dof_index, types::global_dof_index>
773  min_neighbors(numbers::invalid_dof_index,
775  for (types::global_dof_index i = 0; i < touched_nodes.size(); ++i)
776  if (touched_nodes[i] == numbers::invalid_dof_index)
777  if (row_lengths[i] < min_neighbors.second)
778  {
779  min_neighbors = std::make_pair(i, n_remaining_neighbors[i]);
780  if (n_remaining_neighbors[i] <= 1)
781  break;
782  }
783  if (min_neighbors.first == numbers::invalid_dof_index)
784  break;
785 
786  Assert(min_neighbors.second > 0, ExcInternalError());
787 
788  current_neighbors.clear();
789  current_neighbors.insert(min_neighbors.first);
790  while (!current_neighbors.empty())
791  {
792  // Find node with minimum number of untouched neighbors among the
793  // next set of possible neighbors
794  min_neighbors = std::make_pair(numbers::invalid_dof_index,
796  for (std::set<types::global_dof_index>::iterator it =
797  current_neighbors.begin();
798  it != current_neighbors.end();
799  ++it)
800  {
801  Assert(touched_nodes[*it] == numbers::invalid_dof_index,
802  ExcInternalError());
803  if (n_remaining_neighbors[*it] < min_neighbors.second)
804  min_neighbors =
805  std::make_pair(*it, n_remaining_neighbors[*it]);
806  }
807 
808  // Among the set of nodes with the minimal number of neighbors,
809  // choose the one with the largest number of touched neighbors,
810  // i.e., the one with the largest row length
811  const types::global_dof_index best_row_length =
812  min_neighbors.second;
813  for (std::set<types::global_dof_index>::iterator it =
814  current_neighbors.begin();
815  it != current_neighbors.end();
816  ++it)
817  if (n_remaining_neighbors[*it] == best_row_length)
818  if (row_lengths[*it] > min_neighbors.second)
819  min_neighbors = std::make_pair(*it, row_lengths[*it]);
820 
821  // Add the pivot and all direct neighbors of the pivot node not
822  // yet touched to the list of new entries.
823  groups.emplace_back();
824  std::vector<types::global_dof_index> &next_group = groups.back();
825 
826  next_group.push_back(min_neighbors.first);
827  touched_nodes[min_neighbors.first] = groups.size() - 1;
829  connectivity.begin(min_neighbors.first);
830  it != connectivity.end(min_neighbors.first);
831  ++it)
832  if (touched_nodes[it->column()] == numbers::invalid_dof_index)
833  {
834  next_group.push_back(it->column());
835  touched_nodes[it->column()] = groups.size() - 1;
836  }
837 
838  // Add all neighbors of the current list not yet touched to the
839  // set of possible next pivots. The added node is no longer a
840  // valid neighbor (here we assume symmetry of the
841  // connectivity). Delete the entries of the current list from
842  // the set of possible next pivots.
843  for (unsigned int i = 0; i < next_group.size(); ++i)
844  {
846  connectivity.begin(next_group[i]);
847  it != connectivity.end(next_group[i]);
848  ++it)
849  {
850  if (touched_nodes[it->column()] ==
852  current_neighbors.insert(it->column());
853  n_remaining_neighbors[it->column()]--;
854  }
855  current_neighbors.erase(next_group[i]);
856  }
857  }
858  }
859 
860  // Sanity check: for all nodes, there should not be any neighbors left
861  for (types::global_dof_index row = 0; row < connectivity.n_rows(); ++row)
862  Assert(n_remaining_neighbors[row] == 0, ExcInternalError());
863 
864  // If the number of groups is smaller than the number of nodes, we
865  // continue by recursively calling this method
866  if (groups.size() < connectivity.n_rows())
867  {
868  // Form the connectivity of the groups
869  DynamicSparsityPattern connectivity_next(groups.size(),
870  groups.size());
871  for (types::global_dof_index i = 0; i < groups.size(); ++i)
872  for (types::global_dof_index col = 0; col < groups[i].size(); ++col)
874  connectivity.begin(groups[i][col]);
875  it != connectivity.end(groups[i][col]);
876  ++it)
877  connectivity_next.add(i, touched_nodes[it->column()]);
878 
879  // Recursively call the reordering
880  std::vector<types::global_dof_index> renumbering_next(groups.size());
881  reorder_hierarchical(connectivity_next, renumbering_next);
882 
883  // Renumber the indices group by group according to the incoming
884  // ordering for the groups
885  for (types::global_dof_index i = 0, count = 0; i < groups.size(); ++i)
886  for (types::global_dof_index col = 0;
887  col < groups[renumbering_next[i]].size();
888  ++col, ++count)
889  renumbering[count] = groups[renumbering_next[i]][col];
890  }
891  else
892  {
893  // All groups should have size one and no more recursion is possible,
894  // so use the numbering of the groups
895  for (types::global_dof_index i = 0, count = 0; i < groups.size(); ++i)
896  for (types::global_dof_index col = 0; col < groups[i].size();
897  ++col, ++count)
898  renumbering[count] = groups[i][col];
899  }
900  }
901  } // namespace internal
902 
903  void
905  const DynamicSparsityPattern & connectivity,
906  std::vector<DynamicSparsityPattern::size_type> &renumbering)
907  {
908  // the internal renumbering keeps the numbering the wrong way around (but
909  // we cannot invert the numbering inside that method because it is used
910  // recursively), so invert it here
911  internal::reorder_hierarchical(connectivity, renumbering);
912  renumbering = Utilities::invert_permutation(renumbering);
913  }
914 
915 
916 
917 #ifdef DEAL_II_WITH_MPI
918  void
921  const std::vector<DynamicSparsityPattern::size_type> &rows_per_cpu,
922  const MPI_Comm & mpi_comm,
923  const IndexSet & myrange)
924  {
925  const unsigned int myid = Utilities::MPI::this_mpi_process(mpi_comm);
926  std::vector<DynamicSparsityPattern::size_type> start_index(
927  rows_per_cpu.size() + 1);
928  start_index[0] = 0;
929  for (DynamicSparsityPattern::size_type i = 0; i < rows_per_cpu.size(); ++i)
930  start_index[i + 1] = start_index[i] + rows_per_cpu[i];
931 
932  using map_vec_t = std::map<DynamicSparsityPattern::size_type,
933  std::vector<DynamicSparsityPattern::size_type>>;
934 
935  map_vec_t send_data;
936 
937  {
938  unsigned int dest_cpu = 0;
939 
940  DynamicSparsityPattern::size_type n_local_rel_rows = myrange.n_elements();
941  for (DynamicSparsityPattern::size_type row_idx = 0;
942  row_idx < n_local_rel_rows;
943  ++row_idx)
944  {
945  DynamicSparsityPattern::size_type row =
946  myrange.nth_index_in_set(row_idx);
947 
948  // calculate destination CPU
949  while (row >= start_index[dest_cpu + 1])
950  ++dest_cpu;
951 
952  // skip myself
953  if (dest_cpu == myid)
954  {
955  row_idx += rows_per_cpu[myid] - 1;
956  continue;
957  }
958 
959  DynamicSparsityPattern::size_type rlen = dsp.row_length(row);
960 
961  // skip empty lines
962  if (!rlen)
963  continue;
964 
965  // save entries
966  std::vector<DynamicSparsityPattern::size_type> &dst =
967  send_data[dest_cpu];
968 
969  dst.push_back(rlen); // number of entries
970  dst.push_back(row); // row index
971  for (DynamicSparsityPattern::size_type c = 0; c < rlen; ++c)
972  {
973  // columns
974  DynamicSparsityPattern::size_type column =
975  dsp.column_number(row, c);
976  dst.push_back(column);
977  }
978  }
979  }
980 
981  unsigned int num_receive = 0;
982  {
983  std::vector<unsigned int> send_to;
984  send_to.reserve(send_data.size());
985  for (map_vec_t::iterator it = send_data.begin(); it != send_data.end();
986  ++it)
987  send_to.push_back(it->first);
988 
989  num_receive =
991  send_to)
992  .size();
993  }
994 
995  std::vector<MPI_Request> requests(send_data.size());
996 
997 
998  // send data
999  {
1000  unsigned int idx = 0;
1001  for (map_vec_t::iterator it = send_data.begin(); it != send_data.end();
1002  ++it, ++idx)
1003  {
1004  const int ierr = MPI_Isend(&(it->second[0]),
1005  it->second.size(),
1006  DEAL_II_DOF_INDEX_MPI_TYPE,
1007  it->first,
1008  124,
1009  mpi_comm,
1010  &requests[idx]);
1011  AssertThrowMPI(ierr);
1012  }
1013  }
1014 
1015  {
1016  // receive
1017  std::vector<DynamicSparsityPattern::size_type> recv_buf;
1018  for (unsigned int index = 0; index < num_receive; ++index)
1019  {
1020  MPI_Status status;
1021  int len;
1022  int ierr = MPI_Probe(MPI_ANY_SOURCE, MPI_ANY_TAG, mpi_comm, &status);
1023  AssertThrowMPI(ierr);
1024  Assert(status.MPI_TAG == 124, ExcInternalError());
1025 
1026  ierr = MPI_Get_count(&status, DEAL_II_DOF_INDEX_MPI_TYPE, &len);
1027  AssertThrowMPI(ierr);
1028  recv_buf.resize(len);
1029  ierr = MPI_Recv(recv_buf.data(),
1030  len,
1031  DEAL_II_DOF_INDEX_MPI_TYPE,
1032  status.MPI_SOURCE,
1033  status.MPI_TAG,
1034  mpi_comm,
1035  &status);
1036  AssertThrowMPI(ierr);
1037 
1038  std::vector<DynamicSparsityPattern::size_type>::const_iterator ptr =
1039  recv_buf.begin();
1040  std::vector<DynamicSparsityPattern::size_type>::const_iterator end =
1041  recv_buf.end();
1042  while (ptr != end)
1043  {
1044  DynamicSparsityPattern::size_type num = *(ptr++);
1045  Assert(ptr != end, ExcInternalError());
1046  DynamicSparsityPattern::size_type row = *(ptr++);
1047  for (unsigned int c = 0; c < num; ++c)
1048  {
1049  Assert(ptr != end, ExcInternalError());
1050  dsp.add(row, *ptr);
1051  ++ptr;
1052  }
1053  }
1054  Assert(ptr == end, ExcInternalError());
1055  }
1056  }
1057 
1058  // complete all sends, so that we can safely destroy the buffers.
1059  if (requests.size())
1060  {
1061  const int ierr =
1062  MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
1063  AssertThrowMPI(ierr);
1064  }
1065  }
1066 
1067  void
1069  const std::vector<IndexSet> &owned_set_per_cpu,
1070  const MPI_Comm & mpi_comm,
1071  const IndexSet & myrange)
1072  {
1073  const unsigned int myid = Utilities::MPI::this_mpi_process(mpi_comm);
1074 
1075  using map_vec_t =
1077  std::vector<BlockDynamicSparsityPattern::size_type>>;
1078  map_vec_t send_data;
1079 
1080  {
1081  unsigned int dest_cpu = 0;
1082 
1083  BlockDynamicSparsityPattern::size_type n_local_rel_rows =
1084  myrange.n_elements();
1085  for (BlockDynamicSparsityPattern::size_type row_idx = 0;
1086  row_idx < n_local_rel_rows;
1087  ++row_idx)
1088  {
1089  BlockDynamicSparsityPattern::size_type row =
1090  myrange.nth_index_in_set(row_idx);
1091 
1092  // calculate destination CPU, note that we start the search
1093  // at last destination cpu, because even if the owned ranges
1094  // are not contiguous, they hopefully consist of large blocks
1095  while (!owned_set_per_cpu[dest_cpu].is_element(row))
1096  {
1097  ++dest_cpu;
1098  if (dest_cpu == owned_set_per_cpu.size()) // wrap around
1099  dest_cpu = 0;
1100  }
1101 
1102  // skip myself
1103  if (dest_cpu == myid)
1104  continue;
1105 
1106  BlockDynamicSparsityPattern::size_type rlen = dsp.row_length(row);
1107 
1108  // skip empty lines
1109  if (!rlen)
1110  continue;
1111 
1112  // save entries
1113  std::vector<BlockDynamicSparsityPattern::size_type> &dst =
1114  send_data[dest_cpu];
1115 
1116  dst.push_back(rlen); // number of entries
1117  dst.push_back(row); // row index
1118  for (BlockDynamicSparsityPattern::size_type c = 0; c < rlen; ++c)
1119  {
1120  // columns
1121  BlockDynamicSparsityPattern::size_type column =
1122  dsp.column_number(row, c);
1123  dst.push_back(column);
1124  }
1125  }
1126  }
1127 
1128  unsigned int num_receive = 0;
1129  {
1130  std::vector<unsigned int> send_to;
1131  send_to.reserve(send_data.size());
1132  for (map_vec_t::iterator it = send_data.begin(); it != send_data.end();
1133  ++it)
1134  send_to.push_back(it->first);
1135 
1136  num_receive =
1138  send_to)
1139  .size();
1140  }
1141 
1142  std::vector<MPI_Request> requests(send_data.size());
1143 
1144 
1145  // send data
1146  {
1147  unsigned int idx = 0;
1148  for (map_vec_t::iterator it = send_data.begin(); it != send_data.end();
1149  ++it, ++idx)
1150  {
1151  const int ierr = MPI_Isend(&(it->second[0]),
1152  it->second.size(),
1153  DEAL_II_DOF_INDEX_MPI_TYPE,
1154  it->first,
1155  124,
1156  mpi_comm,
1157  &requests[idx]);
1158  AssertThrowMPI(ierr);
1159  }
1160  }
1161 
1162  {
1163  // receive
1164  std::vector<BlockDynamicSparsityPattern::size_type> recv_buf;
1165  for (unsigned int index = 0; index < num_receive; ++index)
1166  {
1167  MPI_Status status;
1168  int len;
1169  int ierr = MPI_Probe(MPI_ANY_SOURCE, MPI_ANY_TAG, mpi_comm, &status);
1170  AssertThrowMPI(ierr);
1171  Assert(status.MPI_TAG == 124, ExcInternalError());
1172 
1173  ierr = MPI_Get_count(&status, DEAL_II_DOF_INDEX_MPI_TYPE, &len);
1174  AssertThrowMPI(ierr);
1175  recv_buf.resize(len);
1176  ierr = MPI_Recv(recv_buf.data(),
1177  len,
1178  DEAL_II_DOF_INDEX_MPI_TYPE,
1179  status.MPI_SOURCE,
1180  status.MPI_TAG,
1181  mpi_comm,
1182  &status);
1183  AssertThrowMPI(ierr);
1184 
1185  std::vector<BlockDynamicSparsityPattern::size_type>::const_iterator
1186  ptr = recv_buf.begin();
1187  std::vector<BlockDynamicSparsityPattern::size_type>::const_iterator
1188  end = recv_buf.end();
1189  while (ptr != end)
1190  {
1191  BlockDynamicSparsityPattern::size_type num = *(ptr++);
1192  Assert(ptr != end, ExcInternalError());
1193  BlockDynamicSparsityPattern::size_type row = *(ptr++);
1194  for (unsigned int c = 0; c < num; ++c)
1195  {
1196  Assert(ptr != end, ExcInternalError());
1197  dsp.add(row, *ptr);
1198  ++ptr;
1199  }
1200  }
1201  Assert(ptr == end, ExcInternalError());
1202  }
1203  }
1204 
1205  // complete all sends, so that we can safely destroy the buffers.
1206  if (requests.size())
1207  {
1208  const int ierr =
1209  MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
1210  AssertThrowMPI(ierr);
1211  }
1212  }
1213 #endif
1214 } // namespace SparsityTools
1215 
1216 DEAL_II_NAMESPACE_CLOSE
size_type row_length(const size_type row) const
const types::global_dof_index invalid_size_type
Definition: types.h:182
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
void add(const size_type i, const size_type j)
static::ExceptionBase & ExcMETISNotInstalled()
iterator begin() const
size_type n_elements() const
Definition: index_set.h:1743
void reorder_Cuthill_McKee(const DynamicSparsityPattern &sparsity, std::vector< DynamicSparsityPattern::size_type > &new_indices, const std::vector< DynamicSparsityPattern::size_type > &starting_indices=std::vector< DynamicSparsityPattern::size_type >())
#define AssertThrow(cond, exc)
Definition: exceptions.h:1329
size_type size() const
Definition: index_set.h:1611
unsigned long long int global_dof_index
Definition: types.h:72
static::ExceptionBase & ExcNotCompressed()
static::ExceptionBase & ExcMessage(std::string arg1)
unsigned int row_length(const size_type row) const
static::ExceptionBase & ExcInvalidNumberOfPartitions(int arg1)
void add(const size_type i, const size_type j)
const IndexSet & row_index_set() const
void distribute_sparsity_pattern(DynamicSparsityPattern &dsp, const std::vector< DynamicSparsityPattern::size_type > &rows_per_cpu, const MPI_Comm &mpi_comm, const IndexSet &myrange)
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices, const Partitioner partitioner=Partitioner::metis)
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcInvalidArraySize(int arg1, int arg2)
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
size_type n_cols() const
unsigned int row_length(const size_type row) const
iterator end() const
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
Definition: utilities.cc:503
types::global_dof_index size_type
size_type column_number(const size_type row, const size_type index) const
static::ExceptionBase & ExcNotQuadratic()
void reorder_hierarchical(const DynamicSparsityPattern &sparsity, std::vector< DynamicSparsityPattern::size_type > &new_indices)
types::global_dof_index size_type
static::ExceptionBase & ExcZOLTANNotInstalled()
unsigned int color_sparsity_pattern(const SparsityPattern &sparsity_pattern, std::vector< unsigned int > &color_indices)
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1443
size_type column_number(const size_type row, const unsigned int index) const
size_type n_rows() const
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:80
bool exists(const size_type i, const size_type j) const
const types::global_dof_index invalid_dof_index
Definition: types.h:188
size_type nth_index_in_set(const unsigned int local_index) const
Definition: index_set.h:1793
std::vector< unsigned int > compute_point_to_point_communication_pattern(const MPI_Comm &mpi_comm, const std::vector< unsigned int > &destinations)
Definition: mpi.cc:184
bool is_compressed() const
std::size_t n_nonzero_elements() const
static::ExceptionBase & ExcMETISError(int arg1)
static::ExceptionBase & ExcInternalError()