]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Fix a bug in graph coloring and add a test for it.
authorturcksin <turcksin@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 1 Oct 2013 17:12:43 +0000 (17:12 +0000)
committerturcksin <turcksin@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 1 Oct 2013 17:12:43 +0000 (17:12 +0000)
git-svn-id: https://svn.dealii.org/trunk@31058 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/base/graph_coloring.h
tests/base/graph_coloring_04.cc [new file with mode: 0644]
tests/base/graph_coloring_04/cmp/generic [new file with mode: 0644]

index 2349fbe99435b89b3764c3b452ea925ed27092e5..fcd4d0b0c9c30f2916ac00f57a7419c63fb9e8dc 100644 (file)
 #include <boost/unordered_map.hpp>
 #include <boost/unordered_set.hpp>
 
+#include <set>
 #include <vector>
 
+
 DEAL_II_NAMESPACE_OPEN
 
 /// This namespace contains the functions necessary to color a graph.
@@ -120,7 +122,7 @@ namespace graph_coloring {
     // Number of zones composing the partitioning.
     const unsigned int partition_size(partition.size());
     std::vector<unsigned int> sorted_vertices(partition_size);
-    std::vector<unsigned int> degrees(partition_size);
+    std::vector<int> degrees(partition_size);
     std::vector<std::vector<types::global_dof_index> > conflict_indices(partition_size);
     std::vector<std::vector<unsigned int> > graph(partition_size);
 
@@ -132,7 +134,7 @@ namespace graph_coloring {
       std::sort(conflict_indices[i].begin(),conflict_indices[i].end());
     }
 
-    // Compute the degree of each vertex of the graph  using the
+    // Compute the degree of each vertex of the graph using the
     // intersection of the conflict indices.
     std::vector<types::global_dof_index> conflict_indices_intersection;
     std::vector<types::global_dof_index>::iterator intersection_it;
@@ -156,14 +158,14 @@ namespace graph_coloring {
       }
 
     // Sort the vertices by decreasing degree.
-    std::vector<unsigned int>::iterator degrees_it;
+    std::vector<int>::iterator degrees_it;
     for (unsigned int i=0; i<partition_size; ++i)
     {
       // Find the largest element.
       degrees_it = std::max_element(degrees.begin(),degrees.end());
       sorted_vertices[i] = degrees_it-degrees.begin();
-      // Zero the largest element.
-      *degrees_it = 0;
+      // Put the largest element to -1 so it cannot be chosen again.
+      *degrees_it = -1;
     }
 
     // Color the graph.
diff --git a/tests/base/graph_coloring_04.cc b/tests/base/graph_coloring_04.cc
new file mode 100644 (file)
index 0000000..a417c2c
--- /dev/null
@@ -0,0 +1,108 @@
+// ---------------------------------------------------------------------
+// $Id: graph_coloring_03.cc 30045 2013-07-18 19:18:00Z maier $
+//
+// Copyright (C) 2003 - 2013 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Check that graph coloring colors every cells.
+
+
+#include "../tests.h"
+#include <fstream>
+#include <vector>
+
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/graph_coloring.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_refinement.h>
+#include <deal.II/lac/vector.h>
+
+
+
+template <int dim>
+std::vector<types::global_dof_index> get_conflict_indices(
+    const typename DoFHandler<dim>::active_cell_iterator &it)
+{
+  std::vector<types::global_dof_index> local_dof_indices(it->get_fe().dofs_per_cell);
+  it->get_dof_indices(local_dof_indices);
+
+  return local_dof_indices;
+}
+
+
+template <int dim>
+void check ()
+{
+  Triangulation<dim> triangulation;
+  GridGenerator::hyper_shell (triangulation,
+      Point<dim>(),
+      1, 2,
+      (dim==3) ? 96 : 12,
+      true);
+
+  triangulation.refine_global (3);
+
+  for (unsigned int i=0; i<3; ++i)
+  {
+    Vector<float> estimated_error_per_cell (triangulation.n_active_cells());
+    for (unsigned int i=0; i< estimated_error_per_cell.size(); ++i)
+      estimated_error_per_cell(i) = i;
+    GridRefinement::
+      refine_and_coarsen_fixed_fraction (triangulation,
+          estimated_error_per_cell,
+          0.3, 0.1);
+    triangulation.execute_coarsening_and_refinement();
+  }
+  
+  FE_Q<dim> fe(1);
+  DoFHandler<dim> stokes_dof_handler (triangulation);
+  stokes_dof_handler.distribute_dofs (fe);
+
+  for (typename DoFHandler<dim>::active_cell_iterator cell=stokes_dof_handler.begin_active();
+      cell != stokes_dof_handler.end(); ++cell)
+    cell->clear_user_flag ();
+  std::vector<std::vector<typename DoFHandler<dim>::active_cell_iterator> > coloring 
+    = graph_coloring::make_graph_coloring(stokes_dof_handler.begin_active(),
+        stokes_dof_handler.end(),
+        std_cxx1x::function<std::vector<types::global_dof_index> (
+          const typename DoFHandler<dim>::active_cell_iterator &)>(&get_conflict_indices<dim>));
+
+  for (unsigned int c=0; c<coloring.size(); ++c)
+    for (unsigned int i=0; i<coloring[c].size(); ++i)
+      coloring[c][i]->set_user_flag();
+
+  for (typename DoFHandler<dim>::active_cell_iterator cell=stokes_dof_handler.begin_active();
+      cell != stokes_dof_handler.end(); ++cell)
+    Assert (cell->user_flag_set () == true, ExcInternalError());
+
+  deallog<<"OK"<<std::endl;
+}
+
+int main (int argc, char *argv[])
+{
+  std::ofstream logfile("graph_coloring_04/output");
+  deallog.attach(logfile);
+  deallog.depth_console (0);
+
+  check<2> ();
+
+  return 0;
+}
diff --git a/tests/base/graph_coloring_04/cmp/generic b/tests/base/graph_coloring_04/cmp/generic
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::OK

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.