]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Actually change the code that handles Metis functions.
authorToby D. Young <tyoung@ippt.pan.pl>
Mon, 13 Aug 2012 13:45:23 +0000 (13:45 +0000)
committerToby D. Young <tyoung@ippt.pan.pl>
Mon, 13 Aug 2012 13:45:23 +0000 (13:45 +0000)
git-svn-id: https://svn.dealii.org/trunk@25934 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/source/lac/sparsity_tools.cc

index baae3e719f108f5be8e51cad40dd14d06511e06c..a4751585b27769131f57120cf0a0d632da250db0 100644 (file)
@@ -62,6 +62,13 @@ used to store boundary indicators internally.
 
 
 <ol>
+<li>
+Changed: Support for the METIS 4.x has been replaced with support for
+METIS 5.x. Use <code>--with-metis=path/to/metis</code> to configure
+with METIS 5.x.
+<br>
+(Stefano Zampini, Toby D. Young 2012/08/13)
+
 <li>
 Changed: numerics/vectors.h is now called numerics/vector_tools.h and
 numerics/matrices.h is now called numerics/matrix_tools.h The old files are
index 22432f1757c7dc8ddb1e12c5db038e641740ce26..b49dfeefcae33f8584b0e27b3608528275f9df69 100644 (file)
 #endif
 
 #ifdef DEAL_II_USE_METIS
-// This is sorta stupid. what we really would like to do here is this:
-//   extern "C" {
-//   #  include <metis.h>
-//   }
-// The problem with this is that (i) metis.h declares a couple of
-// functions with different exception specifications than are
-// declared in standard header files (in particular the drand48
-// function every one who's worked with Metis seems to know about);
-// and (ii) metis.h may include <mpi.h> which if we run a C++
-// compiler may include mpicxx.h -- but that leads to trouble
-// because we have wrapped an
-//   extern "C" {...}
-// around the whole block that now also includes a C++ header :-(
-//
-// The correct solution, of course, would be if Metis would produce
-// a new dot release after 10 years that would remove the unnecessary
-// function prototypes, and that would wrap the contents into
-//   #ifdef __C_PLUSPLUS__
-//   extern "C" {
-//   #endif
-// But that appears to not be happening, and so we have to do the
-// following hack where we just forward declare everything ourselves :-(
-
 extern "C"
 {
-// the following is from Metis-4.0/Lib/struct.h:
-  typedef int idxtype;
-// and this is from proto.h
-  void
-  METIS_PartGraphKway(int *, idxtype *, idxtype *,
-                      idxtype *, idxtype *, int *,
-                      int *, int *, int *, int *, idxtype *);
-  void
-  METIS_PartGraphRecursive(int *, idxtype *, idxtype *,
-                           idxtype *, idxtype *, int *,
-                           int *, int *, int *, int *, idxtype *);
+#include <metis.h>
 }
 #endif
 
@@ -103,28 +70,33 @@ namespace SparsityTools
                                      // simple, since METIS wants exactly our
                                      // compressed row storage format. we only
                                      // have to set up a few auxiliary arrays
-    int
-      n = static_cast<signed int>(sparsity_pattern.n_rows()),
-      wgtflag = 0,                          // no weights on nodes or edges
-      numflag = 0,                          // C-style 0-based numbering
+    idx_t
+      n       = static_cast<signed int>(sparsity_pattern.n_rows()),
+      ncon    = 1,                              // number of balancing constraints (should be >0)
+      wgtflag = 0,                              // no weights on nodes or edges
+      numflag = 0,                              // C-style 0-based numbering
       nparts  = static_cast<int>(n_partitions), // number of subdomains to create
-      dummy;                                // the numbers of edges cut by the
-                                     //   resulting partition
+      dummy;                                    // the numbers of edges cut by the
+                                               // resulting partition
 
                                      // use default options for METIS
-    int options[5] = { 0,0,0,0,0 };
+    idx_t options[METIS_NOPTIONS];
+    METIS_SetDefaultOptions (options);
 
                                      // one more nuisance: we have to copy our
                                      // own data to arrays that store signed
                                      // integers :-(
-    std::vector<idxtype> int_rowstart (sparsity_pattern.get_rowstart_indices(),
-                                       sparsity_pattern.get_rowstart_indices() +
-                                       sparsity_pattern.n_rows()+1);
-    std::vector<idxtype> int_colnums (sparsity_pattern.get_column_numbers(),
-                                      sparsity_pattern.get_column_numbers()+
-                                      int_rowstart[sparsity_pattern.n_rows()]);
+    std::vector<idx_t> int_rowstart (sparsity_pattern.get_rowstart_indices(),
+                                    sparsity_pattern.get_rowstart_indices() +
+                                    sparsity_pattern.n_rows()+1);
+    std::vector<idx_t> int_colnums (sparsity_pattern.get_column_numbers(),
+                                   sparsity_pattern.get_column_numbers()+
+                                   int_rowstart[sparsity_pattern.n_rows()]);
+
+    std::vector<idx_t> int_partition_indices (sparsity_pattern.n_rows());
 
-    std::vector<idxtype> int_partition_indices (sparsity_pattern.n_rows());
+                                     // Make use of METIS' error code.
+    int ierr;
 
                                      // Select which type of partitioning to
                                      // create
@@ -132,17 +104,23 @@ namespace SparsityTools
                                      // Use recursive if the number of
                                      // partitions is less than or equal to 8
     if (n_partitions <= 8)
-      METIS_PartGraphRecursive(&n, &int_rowstart[0], &int_colnums[0],
-                               NULL, NULL,
-                               &wgtflag, &numflag, &nparts, &options[0],
-                               &dummy, &int_partition_indices[0]);
+      ierr = METIS_PartGraphRecursive(&n, &ncon, &int_rowstart[0], &int_colnums[0], 
+                                     NULL, NULL, NULL,
+                                     &nparts,NULL,NULL,&options[0],
+                                     &dummy,&int_partition_indices[0]);
 
-                                     // Otherwise  use kway
+                                     // Otherwise use kway
     else
-      METIS_PartGraphKway(&n, &int_rowstart[0], &int_colnums[0],
-                          NULL, NULL,
-                          &wgtflag, &numflag, &nparts, &options[0],
-                          &dummy, &int_partition_indices[0]);
+      ierr = METIS_PartGraphKway(&n, &ncon, &int_rowstart[0], &int_colnums[0],
+                                NULL, NULL, NULL,
+                                &nparts,NULL,NULL,&options[0],
+                                &dummy,&int_partition_indices[0]);
+
+                                     // If metis returns normally, an
+                                     // error code METIS_OK=1 is
+                                     // returned from the above
+                                     // functions (see metish.h)
+    AssertThrow (ierr == 1, ExcMETISError (ierr));
 
                                      // now copy back generated indices into the
                                      // output array

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.