]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Allow Cuthill_McKee() to also use starting indices in parallel.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 12 Jun 2015 02:27:24 +0000 (21:27 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 12 Jun 2015 02:31:12 +0000 (21:31 -0500)
This was previously not possible and likely yielded at best results that were not intended.

doc/news/changes.h
source/dofs/dof_renumbering.cc
tests/mpi/renumber_cuthill_mckee_02.cc [new file with mode: 0644]
tests/mpi/renumber_cuthill_mckee_02.mpirun=4.output [new file with mode: 0644]

index cdc18334f91773b05f6885814eea41c6d0a3d353..e2f743a95d0bda5f3884db07ee5f0b0b4b778fc2 100644 (file)
@@ -512,6 +512,12 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+  <li> Improved: DoFRenumbering::Cuthill_McKee() can now also
+  use starting indices for parallel triangulations.
+  <br>
+  (Wolfgang Bangerth, 2015/06/11)
+  </li>
+
   <li> Improved: VectorTools::interpolate now works with FE_Nothing.
   <br>
   (Angel Rodriguez, 2015/06/03)
index d0cd8f6d4f722bb00db2e389bf796f1b57973349..24ad40588d33edb27ffe00ea0534fdd85cf5ae17 100644 (file)
@@ -431,9 +431,22 @@ namespace DoFRenumbering
                                  true);
           }
 
+        // translate starting indices from global to local indices
+        std::vector<types::global_dof_index> local_starting_indices (starting_indices.size());
+        for (unsigned int i=0; i<starting_indices.size(); ++i)
+          {
+            Assert (locally_owned.is_element (starting_indices[i]),
+                    ExcMessage ("You specified global degree of freedom "
+                                + Utilities::int_to_string(starting_indices[i]) +
+                                " as a starting index, but this index is not among the "
+                                "locally owned ones on this processor."));
+            local_starting_indices[i] = locally_owned.index_within_set(starting_indices[i]);
+          }
+
+        // then do the renumbering on the locally owned portion
         AssertDimension(new_indices.size(), locally_owned.n_elements());
         SparsityTools::reorder_Cuthill_McKee (sparsity, new_indices,
-                                              starting_indices);
+                                              local_starting_indices);
         if (reversed_numbering)
           new_indices = Utilities::reverse_permutation (new_indices);
 
diff --git a/tests/mpi/renumber_cuthill_mckee_02.cc b/tests/mpi/renumber_cuthill_mckee_02.cc
new file mode 100644 (file)
index 0000000..1d0ea46
--- /dev/null
@@ -0,0 +1,122 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2009 - 2014 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Test that DofRenumbering::Cuthill_McKee works in parallel also when
+// a set of starting indices is given.
+
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/mpi.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_renumbering.h>
+#include <deal.II/fe/fe_q.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  unsigned int nprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  parallel::distributed::Triangulation<dim> tr (MPI_COMM_WORLD);
+
+  GridGenerator::hyper_cube (tr, -1.0, 1.0);
+  tr.refine_global (8-2*dim);
+
+  for (typename Triangulation<dim>::active_cell_iterator
+       cell = tr.begin_active();
+       cell != tr.end(); ++cell)
+    if (!cell->is_ghost() && !cell->is_artificial())
+      if (cell->center().norm() < 0.3)
+        {
+          cell->set_refine_flag();
+        }
+
+  tr.execute_coarsening_and_refinement ();
+
+  DoFHandler<dim> dofh (tr);
+
+  static const FE_Q<dim> fe (1);
+  dofh.distribute_dofs (fe);
+  std::vector<types::global_dof_index> renumbering(dofh.locally_owned_dofs().n_elements());
+  std::vector<types::global_dof_index> starting_indices;
+  starting_indices.push_back (dofh.locally_owned_dofs().nth_index_in_set(0));
+  DoFRenumbering::compute_Cuthill_McKee (renumbering, dofh, false, false, starting_indices);
+
+  // send everything to processor 0 for output
+  std::vector<types::global_dof_index> complete_renumbering(dofh.n_dofs());
+  std::copy(renumbering.begin(), renumbering.end(), complete_renumbering.begin());
+  unsigned int offset = renumbering.size();
+  for (unsigned int i=1; i<nprocs; ++i)
+    {
+      if (myid == i)
+        MPI_Send (&renumbering[0], renumbering.size(),
+                  Utilities::MPI::internal::mpi_type_id(&complete_renumbering[0]), 0, i,
+                  MPI_COMM_WORLD);
+      else if (myid == 0)
+        MPI_Recv (&complete_renumbering[offset],
+                  dofh.locally_owned_dofs_per_processor()[i].n_elements(),
+                  Utilities::MPI::internal::mpi_type_id(&complete_renumbering[0]), i, i,
+                  MPI_COMM_WORLD, MPI_STATUSES_IGNORE);
+      offset += dofh.locally_owned_dofs_per_processor()[i].n_elements();
+    }
+
+  if (myid == 0)
+    {
+      AssertDimension(offset, complete_renumbering.size());
+      for (unsigned int i=0; i<complete_renumbering.size(); ++i)
+        deallog << complete_renumbering[i] << std::endl;
+    }
+}
+
+
+int main (int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+
+
+  deallog.push (Utilities::int_to_string (myid));
+
+  if (myid == 0)
+    {
+      std::ofstream logfile ("output");
+      deallog.attach (logfile);
+      deallog.depth_console (0);
+      deallog.threshold_double (1.e-10);
+
+      deallog.push ("2d");
+      test<2>();
+      deallog.pop();
+      deallog.push ("3d");
+      test<3>();
+      deallog.pop();
+    }
+  else
+    {
+      test<2>();
+      test<3>();
+    }
+}
diff --git a/tests/mpi/renumber_cuthill_mckee_02.mpirun=4.output b/tests/mpi/renumber_cuthill_mckee_02.mpirun=4.output
new file mode 100644 (file)
index 0000000..8d70654
--- /dev/null
@@ -0,0 +1,471 @@
+
+DEAL:0:2d::0
+DEAL:0:2d::1
+DEAL:0:2d::2
+DEAL:0:2d::3
+DEAL:0:2d::4
+DEAL:0:2d::6
+DEAL:0:2d::5
+DEAL:0:2d::7
+DEAL:0:2d::8
+DEAL:0:2d::9
+DEAL:0:2d::11
+DEAL:0:2d::16
+DEAL:0:2d::18
+DEAL:0:2d::12
+DEAL:0:2d::19
+DEAL:0:2d::10
+DEAL:0:2d::13
+DEAL:0:2d::14
+DEAL:0:2d::17
+DEAL:0:2d::20
+DEAL:0:2d::21
+DEAL:0:2d::15
+DEAL:0:2d::22
+DEAL:0:2d::23
+DEAL:0:2d::24
+DEAL:0:2d::25
+DEAL:0:2d::27
+DEAL:0:2d::36
+DEAL:0:2d::38
+DEAL:0:2d::28
+DEAL:0:2d::39
+DEAL:0:2d::49
+DEAL:0:2d::53
+DEAL:0:2d::66
+DEAL:0:2d::68
+DEAL:0:2d::54
+DEAL:0:2d::69
+DEAL:0:2d::29
+DEAL:0:2d::40
+DEAL:0:2d::30
+DEAL:0:2d::41
+DEAL:0:2d::55
+DEAL:0:2d::70
+DEAL:0:2d::56
+DEAL:0:2d::71
+DEAL:0:2d::26
+DEAL:0:2d::31
+DEAL:0:2d::32
+DEAL:0:2d::37
+DEAL:0:2d::42
+DEAL:0:2d::43
+DEAL:0:2d::33
+DEAL:0:2d::34
+DEAL:0:2d::44
+DEAL:0:2d::45
+DEAL:0:2d::50
+DEAL:0:2d::57
+DEAL:0:2d::58
+DEAL:0:2d::67
+DEAL:0:2d::72
+DEAL:0:2d::73
+DEAL:0:2d::59
+DEAL:0:2d::60
+DEAL:0:2d::74
+DEAL:0:2d::75
+DEAL:0:2d::35
+DEAL:0:2d::46
+DEAL:0:2d::47
+DEAL:0:2d::48
+DEAL:0:2d::61
+DEAL:0:2d::76
+DEAL:0:2d::64
+DEAL:0:2d::80
+DEAL:0:2d::62
+DEAL:0:2d::65
+DEAL:0:2d::77
+DEAL:0:2d::81
+DEAL:0:2d::51
+DEAL:0:2d::52
+DEAL:0:2d::63
+DEAL:0:2d::82
+DEAL:0:2d::83
+DEAL:0:2d::84
+DEAL:0:2d::78
+DEAL:0:2d::85
+DEAL:0:2d::87
+DEAL:0:2d::91
+DEAL:0:2d::88
+DEAL:0:2d::79
+DEAL:0:2d::86
+DEAL:0:2d::92
+DEAL:0:2d::89
+DEAL:0:2d::90
+DEAL:0:2d::93
+DEAL:0:2d::95
+DEAL:0:2d::96
+DEAL:0:2d::94
+DEAL:0:2d::97
+DEAL:0:2d::98
+DEAL:0:2d::99
+DEAL:0:2d::100
+DEAL:0:2d::101
+DEAL:0:2d::103
+DEAL:0:2d::102
+DEAL:0:2d::104
+DEAL:0:2d::106
+DEAL:0:2d::108
+DEAL:0:2d::105
+DEAL:0:2d::109
+DEAL:0:2d::107
+DEAL:0:2d::110
+DEAL:0:2d::113
+DEAL:0:2d::115
+DEAL:0:2d::111
+DEAL:0:2d::112
+DEAL:0:2d::116
+DEAL:0:2d::117
+DEAL:0:2d::114
+DEAL:0:2d::118
+DEAL:0:2d::122
+DEAL:0:2d::124
+DEAL:0:2d::119
+DEAL:0:2d::125
+DEAL:0:2d::133
+DEAL:0:2d::134
+DEAL:0:2d::146
+DEAL:0:2d::148
+DEAL:0:2d::135
+DEAL:0:2d::149
+DEAL:0:2d::120
+DEAL:0:2d::126
+DEAL:0:2d::121
+DEAL:0:2d::127
+DEAL:0:2d::136
+DEAL:0:2d::150
+DEAL:0:2d::137
+DEAL:0:2d::151
+DEAL:0:2d::123
+DEAL:0:2d::128
+DEAL:0:2d::138
+DEAL:0:2d::145
+DEAL:0:2d::129
+DEAL:0:2d::130
+DEAL:0:2d::139
+DEAL:0:2d::140
+DEAL:0:2d::165
+DEAL:0:2d::158
+DEAL:0:2d::159
+DEAL:0:2d::174
+DEAL:0:2d::167
+DEAL:0:2d::168
+DEAL:0:2d::131
+DEAL:0:2d::132
+DEAL:0:2d::141
+DEAL:0:2d::142
+DEAL:0:2d::143
+DEAL:0:2d::152
+DEAL:0:2d::144
+DEAL:0:2d::153
+DEAL:0:2d::160
+DEAL:0:2d::161
+DEAL:0:2d::169
+DEAL:0:2d::170
+DEAL:0:2d::162
+DEAL:0:2d::154
+DEAL:0:2d::171
+DEAL:0:2d::166
+DEAL:0:2d::147
+DEAL:0:2d::155
+DEAL:0:2d::163
+DEAL:0:2d::172
+DEAL:0:2d::175
+DEAL:0:2d::156
+DEAL:0:2d::164
+DEAL:0:2d::157
+DEAL:0:2d::176
+DEAL:0:2d::178
+DEAL:0:2d::181
+DEAL:0:2d::182
+DEAL:0:2d::179
+DEAL:0:2d::177
+DEAL:0:2d::173
+DEAL:0:2d::180
+DEAL:0:2d::183
+DEAL:0:2d::184
+DEAL:0:2d::187
+DEAL:0:2d::185
+DEAL:0:2d::186
+DEAL:0:2d::189
+DEAL:0:2d::192
+DEAL:0:2d::199
+DEAL:0:2d::194
+DEAL:0:2d::201
+DEAL:0:2d::188
+DEAL:0:2d::190
+DEAL:0:2d::191
+DEAL:0:2d::193
+DEAL:0:2d::195
+DEAL:0:2d::196
+DEAL:0:2d::197
+DEAL:0:2d::202
+DEAL:0:2d::198
+DEAL:0:2d::203
+DEAL:0:2d::208
+DEAL:0:2d::220
+DEAL:0:2d::210
+DEAL:0:2d::231
+DEAL:0:2d::211
+DEAL:0:2d::221
+DEAL:0:2d::212
+DEAL:0:2d::222
+DEAL:0:2d::251
+DEAL:0:2d::244
+DEAL:0:2d::260
+DEAL:0:2d::253
+DEAL:0:2d::245
+DEAL:0:2d::254
+DEAL:0:2d::200
+DEAL:0:2d::204
+DEAL:0:2d::205
+DEAL:0:2d::209
+DEAL:0:2d::213
+DEAL:0:2d::214
+DEAL:0:2d::206
+DEAL:0:2d::207
+DEAL:0:2d::215
+DEAL:0:2d::216
+DEAL:0:2d::219
+DEAL:0:2d::223
+DEAL:0:2d::224
+DEAL:0:2d::232
+DEAL:0:2d::234
+DEAL:0:2d::235
+DEAL:0:2d::225
+DEAL:0:2d::226
+DEAL:0:2d::236
+DEAL:0:2d::237
+DEAL:0:2d::217
+DEAL:0:2d::227
+DEAL:0:2d::218
+DEAL:0:2d::228
+DEAL:0:2d::246
+DEAL:0:2d::255
+DEAL:0:2d::247
+DEAL:0:2d::256
+DEAL:0:2d::229
+DEAL:0:2d::230
+DEAL:0:2d::238
+DEAL:0:2d::239
+DEAL:0:2d::248
+DEAL:0:2d::257
+DEAL:0:2d::240
+DEAL:0:2d::252
+DEAL:0:2d::233
+DEAL:0:2d::241
+DEAL:0:2d::258
+DEAL:0:2d::249
+DEAL:0:2d::261
+DEAL:0:2d::264
+DEAL:0:2d::268
+DEAL:0:2d::267
+DEAL:0:2d::265
+DEAL:0:2d::242
+DEAL:0:2d::250
+DEAL:0:2d::262
+DEAL:0:2d::243
+DEAL:0:2d::263
+DEAL:0:2d::266
+DEAL:0:2d::259
+DEAL:0:2d::269
+DEAL:0:2d::271
+DEAL:0:2d::278
+DEAL:0:2d::277
+DEAL:0:2d::274
+DEAL:0:2d::282
+DEAL:0:2d::283
+DEAL:0:2d::279
+DEAL:0:2d::284
+DEAL:0:2d::291
+DEAL:0:2d::295
+DEAL:0:2d::285
+DEAL:0:2d::286
+DEAL:0:2d::296
+DEAL:0:2d::297
+DEAL:0:2d::292
+DEAL:0:2d::301
+DEAL:0:2d::298
+DEAL:0:2d::303
+DEAL:0:2d::311
+DEAL:0:2d::323
+DEAL:0:2d::313
+DEAL:0:2d::324
+DEAL:0:2d::299
+DEAL:0:2d::304
+DEAL:0:2d::300
+DEAL:0:2d::305
+DEAL:0:2d::314
+DEAL:0:2d::325
+DEAL:0:2d::315
+DEAL:0:2d::326
+DEAL:0:2d::302
+DEAL:0:2d::306
+DEAL:0:2d::312
+DEAL:0:2d::316
+DEAL:0:2d::307
+DEAL:0:2d::308
+DEAL:0:2d::317
+DEAL:0:2d::318
+DEAL:0:2d::327
+DEAL:0:2d::331
+DEAL:0:2d::337
+DEAL:0:2d::339
+DEAL:0:2d::332
+DEAL:0:2d::333
+DEAL:0:2d::340
+DEAL:0:2d::341
+DEAL:0:2d::309
+DEAL:0:2d::310
+DEAL:0:2d::319
+DEAL:0:2d::320
+DEAL:0:2d::321
+DEAL:0:2d::328
+DEAL:0:2d::322
+DEAL:0:2d::329
+DEAL:0:2d::334
+DEAL:0:2d::335
+DEAL:0:2d::342
+DEAL:0:2d::343
+DEAL:0:2d::336
+DEAL:0:2d::330
+DEAL:0:2d::344
+DEAL:0:2d::338
+DEAL:0:2d::289
+DEAL:0:2d::280
+DEAL:0:2d::293
+DEAL:0:2d::287
+DEAL:0:2d::272
+DEAL:0:2d::270
+DEAL:0:2d::275
+DEAL:0:2d::294
+DEAL:0:2d::288
+DEAL:0:2d::290
+DEAL:0:2d::276
+DEAL:0:2d::273
+DEAL:0:2d::281
+DEAL:0:3d::0
+DEAL:0:3d::1
+DEAL:0:3d::2
+DEAL:0:3d::4
+DEAL:0:3d::3
+DEAL:0:3d::5
+DEAL:0:3d::6
+DEAL:0:3d::7
+DEAL:0:3d::11
+DEAL:0:3d::20
+DEAL:0:3d::21
+DEAL:0:3d::26
+DEAL:0:3d::8
+DEAL:0:3d::12
+DEAL:0:3d::13
+DEAL:0:3d::22
+DEAL:0:3d::14
+DEAL:0:3d::23
+DEAL:0:3d::9
+DEAL:0:3d::15
+DEAL:0:3d::16
+DEAL:0:3d::24
+DEAL:0:3d::17
+DEAL:0:3d::25
+DEAL:0:3d::10
+DEAL:0:3d::18
+DEAL:0:3d::19
+DEAL:0:3d::27
+DEAL:0:3d::31
+DEAL:0:3d::32
+DEAL:0:3d::35
+DEAL:0:3d::36
+DEAL:0:3d::40
+DEAL:0:3d::41
+DEAL:0:3d::44
+DEAL:0:3d::28
+DEAL:0:3d::33
+DEAL:0:3d::37
+DEAL:0:3d::42
+DEAL:0:3d::29
+DEAL:0:3d::34
+DEAL:0:3d::38
+DEAL:0:3d::43
+DEAL:0:3d::30
+DEAL:0:3d::39
+DEAL:0:3d::45
+DEAL:0:3d::47
+DEAL:0:3d::48
+DEAL:0:3d::51
+DEAL:0:3d::55
+DEAL:0:3d::61
+DEAL:0:3d::46
+DEAL:0:3d::49
+DEAL:0:3d::50
+DEAL:0:3d::52
+DEAL:0:3d::56
+DEAL:0:3d::62
+DEAL:0:3d::53
+DEAL:0:3d::57
+DEAL:0:3d::58
+DEAL:0:3d::54
+DEAL:0:3d::59
+DEAL:0:3d::60
+DEAL:0:3d::63
+DEAL:0:3d::67
+DEAL:0:3d::69
+DEAL:0:3d::73
+DEAL:0:3d::64
+DEAL:0:3d::68
+DEAL:0:3d::70
+DEAL:0:3d::74
+DEAL:0:3d::65
+DEAL:0:3d::71
+DEAL:0:3d::66
+DEAL:0:3d::72
+DEAL:0:3d::75
+DEAL:0:3d::77
+DEAL:0:3d::78
+DEAL:0:3d::81
+DEAL:0:3d::85
+DEAL:0:3d::91
+DEAL:0:3d::83
+DEAL:0:3d::86
+DEAL:0:3d::87
+DEAL:0:3d::76
+DEAL:0:3d::79
+DEAL:0:3d::80
+DEAL:0:3d::82
+DEAL:0:3d::88
+DEAL:0:3d::92
+DEAL:0:3d::84
+DEAL:0:3d::89
+DEAL:0:3d::90
+DEAL:0:3d::93
+DEAL:0:3d::97
+DEAL:0:3d::99
+DEAL:0:3d::103
+DEAL:0:3d::94
+DEAL:0:3d::100
+DEAL:0:3d::95
+DEAL:0:3d::98
+DEAL:0:3d::101
+DEAL:0:3d::104
+DEAL:0:3d::96
+DEAL:0:3d::102
+DEAL:0:3d::105
+DEAL:0:3d::109
+DEAL:0:3d::113
+DEAL:0:3d::106
+DEAL:0:3d::110
+DEAL:0:3d::114
+DEAL:0:3d::107
+DEAL:0:3d::111
+DEAL:0:3d::115
+DEAL:0:3d::108
+DEAL:0:3d::112
+DEAL:0:3d::116
+DEAL:0:3d::117
+DEAL:0:3d::121
+DEAL:0:3d::118
+DEAL:0:3d::122
+DEAL:0:3d::119
+DEAL:0:3d::123
+DEAL:0:3d::120
+DEAL:0:3d::124

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.