diff --git a/configure.ac b/configure.ac
index 6558696ba25ed72f1789f5149b0ae022427c4610..189e8506c4fac8f1a45c41149884d15202af9491 100644
--- a/configure.ac
+++ b/configure.ac
@@ -484,7 +484,41 @@ AC_CHECK_LIB(pthread, posix_fallocate,
 	     AC_DEFINE([HAVE_POSIX_FALLOCATE], [1], [The posix library implements file allocation functions.]),
 	     AC_MSG_WARN(POSIX implementation does not have file allocation functions.))
 
-# Check for PARMETIS.
+# Check for METIS.
+have_metis="no"
+AC_ARG_WITH([metis],
+    [AS_HELP_STRING([--with-metis=PATH],
+       [root directory where METIS is installed @<:@yes/no@:>@]
+    )],
+    [with_metis="$withval"],
+    [with_metis="no"]
+)
+
+METIS_LIBS=""
+if test "x$with_metis" != "xno"; then
+
+# Check if we have METIS.
+   if test "x$with_metis" != "xyes" -a "x$with_metis" != "x"; then
+      METIS_LIBS="-L$with_metis/lib -lmetis"
+      METIS_INCS="-I$with_metis/include"
+   else
+      METIS_LIBS="-lmetis"
+      METIS_INCS=""
+   fi
+   AC_CHECK_LIB([metis],[METIS_PartGraphKway], [have_metis="yes"],
+                [have_metis="no"], $METIS_LIBS)
+   if test "$have_metis" == "yes"; then
+      AC_DEFINE([HAVE_METIS],1,[The METIS library is present.])
+   else
+      AC_MSG_ERROR("Failed to find a METIS library")
+   fi
+fi
+
+AC_SUBST([METIS_LIBS])
+AC_SUBST([METIS_INCS])
+AM_CONDITIONAL([HAVEMETIS],[test -n "$METIS_LIBS"])
+
+# Check for ParMETIS note we can have both as ParMETIS uses METIS.
 have_parmetis="no"
 AC_ARG_WITH([parmetis],
     [AS_HELP_STRING([--with-parmetis=PATH],
@@ -494,19 +528,6 @@ AC_ARG_WITH([parmetis],
     [with_parmetis="no"]
 )
 
-# Don't proceed with old --with-metis option.
-used_metis="no"
-AC_ARG_WITH([metis],
-    [AS_HELP_STRING([--with-metis],
-       [no longer supported, use --with-parmetis]
-    )],
-    [used_metis="yes"],
-    [:]
-)
-if test "x$used_metis" != "xno"; then
-   AC_MSG_ERROR([Using unsupported --with-metis option, use --with-parmetis])
-fi
-
 if test "x$with_parmetis" != "xno"; then
 
 # Check if we have ParMETIS.
@@ -531,7 +552,7 @@ if test "x$with_parmetis" != "xno"; then
          PARMETIS_INCS=""
       fi
       AC_CHECK_LIB([parmetis],[ParMETIS_V3_RefineKway], [have_parmetis="yes"],
-                   [have_parmetis="no"], $PARMETIS_LIBS)
+                   [have_parmetis="no"], [$METIS_LIBS $PARMETIS_LIBS])
 
    fi
    if test "$have_parmetis" == "yes"; then
@@ -1369,7 +1390,7 @@ AC_MSG_RESULT([
    MPI enabled        : $enable_mpi
    HDF5 enabled       : $with_hdf5
     - parallel        : $have_parallel_hdf5
-   ParMetis enabled   : $have_parmetis
+   METIS/ParMETIS     : $have_metis / $have_parmetis
    FFTW3 enabled      : $have_fftw
    GSL enabled        : $have_gsl
    libNUMA enabled    : $have_numa
diff --git a/examples/Makefile.am b/examples/Makefile.am
index b9c601e039d3f7566e56b0220935e63e42994b28..b4be7751dbcb06f88095482753d0404a1826c7d6 100644
--- a/examples/Makefile.am
+++ b/examples/Makefile.am
@@ -27,8 +27,8 @@ AM_LDFLAGS = $(HDF5_LDFLAGS)
 EXTRA_LIBS = $(HDF5_LIBS) $(FFTW_LIBS) $(PROFILER_LIBS) $(TCMALLOC_LIBS) $(JEMALLOC_LIBS) $(TBBMALLOC_LIBS) $(GRACKLE_LIBS) $(GSL_LIBS)
 
 # MPI libraries.
-MPI_LIBS = $(PARMETIS_LIBS) $(MPI_THREAD_LIBS)
-MPI_FLAGS = -DWITH_MPI $(PARMETIS_INCS)
+MPI_LIBS = $(PARMETIS_LIBS) $(METIS_LIBS) $(MPI_THREAD_LIBS)
+MPI_FLAGS = -DWITH_MPI $(PARMETIS_INCS) $(METIS_INCS)
 
 # Programs.
 bin_PROGRAMS = swift
diff --git a/src/Makefile.am b/src/Makefile.am
index 899d12440396f16e9ee106963938c47b0100bb34..903757e1611f598f088b4bd1f2aa8fc1de54feb8 100644
--- a/src/Makefile.am
+++ b/src/Makefile.am
@@ -28,8 +28,8 @@ GIT_CMD = @GIT_CMD@
 EXTRA_LIBS = $(HDF5_LIBS) $(FFTW_LIBS) $(PROFILER_LIBS) $(TCMALLOC_LIBS) $(JEMALLOC_LIBS) $(TBBMALLOC_LIBS) $(GRACKLE_LIB) $(GSL_LIBS)
 
 # MPI libraries.
-MPI_LIBS = $(PARMETIS_LIBS) $(MPI_THREAD_LIBS)
-MPI_FLAGS = -DWITH_MPI $(PARMETIS_INCS)
+MPI_LIBS = $(PARMETIS_LIBS) $(METIS_LIBS) $(MPI_THREAD_LIBS)
+MPI_FLAGS = -DWITH_MPI $(PARMETIS_INCS) $(METIS_INCS)
 
 # Build the libswiftsim library
 lib_LTLIBRARIES = libswiftsim.la
diff --git a/src/common_io.c b/src/common_io.c
index 3fce8c89f26ab5c4752dff222a0af00837832635..9240146a2c122b8f33997311d84729e1c890e892 100644
--- a/src/common_io.c
+++ b/src/common_io.c
@@ -355,6 +355,10 @@ void io_write_code_description(hid_t h_file) {
 #endif
 #ifdef WITH_MPI
   io_write_attribute_s(h_grpcode, "MPI library", mpi_version());
+#ifdef HAVE_METIS
+  io_write_attribute_s(h_grpcode, "METIS library version",
+                       metis_version());
+#endif
 #ifdef HAVE_PARMETIS
   io_write_attribute_s(h_grpcode, "ParMETIS library version",
                        parmetis_version());
diff --git a/src/debug.c b/src/debug.c
index 3af4c82724d794616155786482bdf20a7aa602b2..f692364bd7d5fd500340b599118df251907d499f 100644
--- a/src/debug.c
+++ b/src/debug.c
@@ -408,16 +408,16 @@ void dumpCells(const char *prefix, int super, int active, int mpiactive,
   fclose(file);
 }
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
 
 /**
- * @brief Dump the ParMETIS graph in standard METIS format, simple format and
- * weights only, to a file.
+ * @brief Dump a graph in METIS standard format, simple format and weights
+ * only, to a file.
  *
  * The standard format output can be read into the METIS and some ParMETIS
- * command-line tools and tests. The simple format is just the cell
- * connectivity (this should not change between calls).  The weights format is
- * the standard one, minus the cell connectivity.
+ * command-line tools. The simple format is just the cell connectivity (this
+ * should not change between calls).  The weights format is the standard one,
+ * minus the cell connectivity.
  *
  * The output filenames are generated from the prefix and the sequence number
  * of calls. So the first is called {prefix}_std_001.dat,
@@ -433,10 +433,9 @@ void dumpCells(const char *prefix, int super, int active, int mpiactive,
  * @param vertexsizes size of vertices
  * @param edgeweights weights of edges
  */
-void dumpParMETISGraph(const char *prefix, idx_t nvertices,
-                       idx_t nvertexweights, idx_t *cellconruns, idx_t *cellcon,
-                       idx_t *vertexweights, idx_t *vertexsizes,
-                       idx_t *edgeweights) {
+void dumpMETISGraph(const char *prefix, idx_t nvertices, idx_t nvertexweights,
+                    idx_t *cellconruns, idx_t *cellcon, idx_t *vertexweights,
+                    idx_t *vertexsizes, idx_t *edgeweights) {
   FILE *stdfile = NULL;
   FILE *simplefile = NULL;
   FILE *weightfile = NULL;
@@ -547,7 +546,7 @@ void dumpParMETISGraph(const char *prefix, idx_t nvertices,
   }
 }
 
-#endif /* HAVE_PARMETIS */
+#endif /* HAVE_METIS || HAVE_PARMETIS */
 
 #ifdef HAVE_MPI
 /**
diff --git a/src/debug.h b/src/debug.h
index 0f978722e157fcefcc585f8882ecbbce5cd3920f..6d0ec665355192d2569aedd36ccd59dce52a56c9 100644
--- a/src/debug.h
+++ b/src/debug.h
@@ -39,10 +39,10 @@ int checkCellhdxmax(const struct cell *c, int *depth);
 void dumpCells(const char *prefix, int super, int active, int mpiactive,
                int pactive, struct space *s, int rank, int step);
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
-#include "parmetis.h"
-void dumpParMETISGraph(const char *prefix, idx_t nvtxs, idx_t ncon, idx_t *xadj,
-                       idx_t *adjncy, idx_t *vwgt, idx_t *vsize, idx_t *adjwgt);
+#if defined(WITH_MPI) && (defined(HAVE_METIS) ||defined(HAVE_PARMETIS))
+#include "metis.h"
+void dumpMETISGraph(const char *prefix, idx_t nvtxs, idx_t ncon, idx_t *xadj,
+                    idx_t *adjncy, idx_t *vwgt, idx_t *vsize, idx_t *adjwgt);
 #endif
 
 #ifdef HAVE_MPI
diff --git a/src/engine.c b/src/engine.c
index 2bcd18bd5e628574e3fdc6227066498d0652fae9..d2d105bed9c347e7210ef16e16dbfc137a012f6c 100644
--- a/src/engine.c
+++ b/src/engine.c
@@ -1116,7 +1116,7 @@ void engine_redistribute(struct engine *e) {
  */
 void engine_repartition(struct engine *e) {
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_PARMETIS) || defined(HAVE_METIS))
 
   ticks tic = getticks();
 
@@ -1133,7 +1133,7 @@ void engine_repartition(struct engine *e) {
   /* Clear the repartition flag. */
   e->forcerepart = 0;
 
-  /* Nothing to do if only using a single node. Also avoids ParMETIS
+  /* Nothing to do if only using a single node. Also avoids METIS
    * bug that doesn't handle this case well. */
   if (e->nr_nodes == 1) return;
 
@@ -1174,7 +1174,7 @@ void engine_repartition(struct engine *e) {
             clocks_getunit());
 #else
   if (e->reparttype->type != REPART_NONE)
-    error("SWIFT was not compiled with MPI and ParMETIS support.");
+    error("SWIFT was not compiled with MPI and METIS or ParMETIS support.");
 
   /* Clear the repartition flag. */
   e->forcerepart = 0;
diff --git a/src/partition.c b/src/partition.c
index f39bcb1e2baf7916df1e369a6681435560fe9d41..f1bb85a32c1336ab2e63037df3169aef705f4ae2 100644
--- a/src/partition.c
+++ b/src/partition.c
@@ -24,7 +24,7 @@
  *  a grid of cells into geometrically connected regions and distributing
  *  these around a number of MPI nodes.
  *
- *  Currently supported partitioning types: grid, vectorise and ParMETIS.
+ *  Currently supported partitioning types: grid, vectorise and METIS/ParMETIS.
  */
 
 /* Config parameters. */
@@ -40,10 +40,13 @@
 /* MPI headers. */
 #ifdef WITH_MPI
 #include <mpi.h>
-/* ParMETIS headers only used when MPI is also available. */
+/* METIS/ParMETIS headers only used when MPI is also available. */
 #ifdef HAVE_PARMETIS
 #include <parmetis.h>
 #endif
+#ifdef HAVE_METIS
+#include <metis.h>
+#endif
 #endif
 
 /* Local headers. */
@@ -55,22 +58,22 @@
 #include "space.h"
 #include "tools.h"
 
-/* Maximum weight used for ParMETIS. */
-#define parmetis_maxweight 10000.0f
+/* Maximum weight used for METIS. */
+#define metis_maxweight 10000.0f
 
 /* Simple descriptions of initial partition types for reports. */
 const char *initial_partition_name[] = {
     "axis aligned grids of cells", "vectorized point associated cells",
-    "memory balanced, using ParMETIS particle weighted cells",
-    "similar sized regions, using ParMETIS unweighted cells"};
+    "memory balanced, using METIS particle weighted cells",
+    "similar sized regions, using METIS unweighted cells"};
 
 /* Simple descriptions of repartition types for reports. */
 const char *repartition_name[] = {
     "no",
-    "ParMETIS edge and vertex task cost weights",
-    "ParMETIS task cost edge weights",
-    "ParMETIS task cost vertex weights",
-    "ParMETIS vertex task costs and edge delta timebin weights"
+    "METIS edge and vertex task cost weights",
+    "METIS task cost edge weights",
+    "METIS task cost vertex weights",
+    "METIS vertex task costs and edge delta timebin weights"
 };
 
 /* Local functions, if needed. */
@@ -154,19 +157,19 @@ static void split_vector(struct space *s, int nregions, int *samplecells) {
 }
 #endif
 
-/* ParMETIS support
- * ================
+/* METIS/ParMETIS support (optional)
+ * =================================
  *
- * ParMETIS partitions using a multi-level k-way scheme. We support using this
- * in a unweighted scheme, which works well and seems to be guaranteed, and a
- * weighted by the number of particles scheme. Note ParMETIS is optional.
+ * METIS/ParMETIS partitions using a multi-level k-way scheme. We support
+ * using this in a unweighted scheme, which works well and seems to be
+ * guaranteed, and a weighted by the number of particles scheme.
  *
  * Repartitioning is based on ParMETIS and uses weights determined from the
  * estimated costs that a cells tasks will take or the relative time bins of
  * the cells next updates.
  */
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
 /**
  * @brief Fill the adjncy array defining the graph of cells in a space.
  *
@@ -238,7 +241,7 @@ static void graph_init(struct space *s, idx_t *adjncy, idx_t *xadj) {
 }
 #endif
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
 /**
  * @brief Accumulate the counts of particles per cell.
  *
@@ -270,9 +273,9 @@ static void accumulate_counts(struct space *s, double *counts) {
 }
 #endif
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
 /**
- * @brief Apply ParMETIS cell-list partitioning to a cell structure.
+ * @brief Apply METIS cell-list partitioning to a cell structure.
  *
  * @param s the space containing the cells to split into regions.
  * @param nregions number of regions.
@@ -283,11 +286,11 @@ static void split_metis(struct space *s, int nregions, int *celllist) {
   for (int i = 0; i < s->nr_cells; i++) s->cells_top[i].nodeID = celllist[i];
 
   /* To check or visualise the partition dump all the cells. */
-  /* dumpCellRanks("parmetis_partition", s->cells_top, s->nr_cells);*/
+  /* dumpCellRanks("metis_partition", s->cells_top, s->nr_cells);*/
 }
 #endif
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
 
 /* qsort support. */
 struct indexval {
@@ -413,16 +416,19 @@ void permute_regions(int *newlist, int *oldlist, int nregions, int ncells,
  *        in CSR format, so same as adjncy array. Need to be in the range of
  *        idx_t.
  * @param refine whether to refine an existing partition, or create a new one.
+ * @param adaptive whether to use an adaptive reparitition of an existing
+ *        partition or simple refinement. Adaptive repartition is controlled
+ *        by the itr parameter.
  * @param itr the ratio of inter-process communication time to data
  *            redistribution time. Used to weight repartitioning edge cuts
- *            when refine is true.
+ *            when refine and adaptive are true.
  * @param celllist on exit this contains the ids of the selected regions,
  *        sizeof number of cells. If refine is 1, then this should contain
  *        the old partition on entry.
  */
 static void pick_parmetis(int nodeID, struct space *s, int nregions,
-                          double *vertexw, double *edgew, int refine, float itr,
-                          int *celllist) {
+                          double *vertexw, double *edgew, int refine,
+                          int adaptive, float itr, int *celllist) {
   int res;
   MPI_Comm comm;
   MPI_Comm_dup(MPI_COMM_WORLD, &comm);
@@ -726,13 +732,25 @@ static void pick_parmetis(int nodeID, struct space *s, int nregions,
      * present on their expected ranks. */
     options[3] = PARMETIS_PSR_UNCOUPLED;
 
-    /* Balance between cuts and movement. */
-    real_t itr_real_t = itr;
-    if (ParMETIS_V3_AdaptiveRepart(vtxdist, xadj, adjncy, weights_v, NULL,
-                                   weights_e, &wgtflag, &numflag, &ncon,
-                                   &nparts, tpwgts, ubvec, &itr_real_t, options,
-                                   &edgecut, regionid, &comm) != METIS_OK)
-      error("Call to ParMETIS_V3_AdaptiveRepart failed.");
+    /* Choice is whether to use an adaptive repartition or a simple
+     * refinement. */
+    if (adaptive) {
+
+        /* Balance between cuts and movement. */
+        real_t itr_real_t = itr;
+        if (ParMETIS_V3_AdaptiveRepart(vtxdist, xadj, adjncy, weights_v, NULL,
+                                       weights_e, &wgtflag, &numflag, &ncon,
+                                       &nparts, tpwgts, ubvec, &itr_real_t, options,
+                                       &edgecut, regionid, &comm) != METIS_OK)
+            error("Call to ParMETIS_V3_AdaptiveRepart failed.");
+    } else {
+
+        if (ParMETIS_V3_RefineKway(vtxdist, xadj, adjncy, weights_v, weights_e,
+                                   &wgtflag, &numflag, &ncon, &nparts, tpwgts,
+                                   ubvec, options, &edgecut, regionid,
+                                   &comm) != METIS_OK)
+            error("Call to ParMETIS_V3_RefineKway failed.");
+    }
   } else {
 
     /* Create a new partition. */
@@ -856,16 +874,13 @@ static void pick_parmetis(int nodeID, struct space *s, int nregions,
 }
 #endif
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && defined(HAVE_METIS) && !defined(HAVE_PARMETIS)
 /**
  * @brief Partition the given space into a number of connected regions.
  *
  * Split the space using METIS to derive a partitions using the given edge and
  * vertex weights. If no weights are given then an unweighted partition is
- * performed. Tests show that METIS gives a better global solution than
- * ParMETIS, at least with the options set in this routine (which cannot be
- * repeated using ParMETIS, the options array is not used in the same way),
- * so this should be used when creating initial partitions.
+ * performed.
  *
  * @param nodeID the rank of our node.
  * @param s the space of cells to partition.
@@ -1020,7 +1035,7 @@ static void pick_metis(int nodeID, struct space *s, int nregions,
 }
 #endif
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
 /**
  * @brief Repartition the cells amongst the nodes using weights of
  *        various kinds.
@@ -1035,10 +1050,10 @@ static void pick_metis(int nodeID, struct space *s, int nregions,
  * @param tasks the completed tasks from the last engine step for our node.
  * @param nr_tasks the number of tasks.
  */
-static void repart_edge_parmetis(int vweights, int eweights, int timebins,
-                                 struct repartition *repartition, int nodeID,
-                                 int nr_nodes, struct space *s,
-                                 struct task *tasks, int nr_tasks) {
+static void repart_edge_metis(int vweights, int eweights, int timebins,
+                              struct repartition *repartition, int nodeID,
+                              int nr_nodes, struct space *s,
+                              struct task *tasks, int nr_tasks) {
 
   /* Create weight arrays using task ticks for vertices and edges (edges
    * assume the same graph structure as used in the part_ calls). */
@@ -1196,6 +1211,7 @@ static void repart_edge_parmetis(int vweights, int eweights, int timebins,
   }
 
   /* Allocate cell list for the partition. If not already done. */
+#ifdef HAVE_PARMETIS
   int refine = 1;
   if (repartition->ncelllist != s->nr_cells) {
     refine = 0;
@@ -1206,9 +1222,10 @@ static void repart_edge_parmetis(int vweights, int eweights, int timebins,
       error("Failed to allocate celllist");
     repartition->ncelllist = s->nr_cells;
   }
+#endif
 
   /* We need to rescale the weights into the range of an integer for
-   * ParMETIS (that is the range of idx_t). Also we would like the range of
+   * METIS (that is the range of idx_t). Also we would like the range of
    * vertex and edges weights to be similar so they balance. */
   double wminv = 0.0;
   double wmaxv = 0.0;
@@ -1246,7 +1263,7 @@ static void repart_edge_parmetis(int vweights, int eweights, int timebins,
         }
         wmine = wminv;
         wmaxe = wmaxv;
-        
+
       } else {
         double wscale = 1.0;
         if ((wmaxv - wminv) > 0.0) {
@@ -1260,10 +1277,10 @@ static void repart_edge_parmetis(int vweights, int eweights, int timebins,
       }
     }
 
-    /* Scale to the ParMETIS range. */
+    /* Scale to the METIS range. */
     double wscale = 1.0;
     if ((wmaxv - wminv) > 0.0) {
-      wscale = (parmetis_maxweight - 1.0) / (wmaxv - wminv);
+      wscale = (metis_maxweight - 1.0) / (wmaxv - wminv);
     }
     for (int k = 0; k < nr_cells; k++) {
       weights_v[k] = (weights_v[k] - wminv) * wscale + 1.0;
@@ -1273,10 +1290,10 @@ static void repart_edge_parmetis(int vweights, int eweights, int timebins,
 
   if (eweights) {
 
-    /* Scale to the ParMETIS range. */
+    /* Scale to the METIS range. */
     double wscale = 1.0;
     if ((wmaxe - wmine) > 0.0) {
-      wscale = (parmetis_maxweight - 1.0) / (wmaxe - wmine);
+      wscale = (metis_maxweight - 1.0) / (wmaxe - wmine);
     }
     for (int k = 0; k < 26 * nr_cells; k++) {
       weights_e[k] = (weights_e[k] - wmine) * wscale + 1.0;
@@ -1284,14 +1301,13 @@ static void repart_edge_parmetis(int vweights, int eweights, int timebins,
   }
 
   /* And repartition/ partition, using both weights or not as requested. */
-  if (refine) {
-    pick_parmetis(nodeID, s, nr_nodes, weights_v, weights_e, refine,
-                  repartition->itr, repartition->celllist);
-  } else {
-
-    /* Not refining, so use METIS. */
-    pick_metis(nodeID, s, nr_nodes, weights_v, weights_e, repartition->celllist);
-  }
+#ifdef HAVE_PARMETIS
+  pick_parmetis(nodeID, s, nr_nodes, weights_v, weights_e, refine,
+                repartition->adaptive, repartition->itr,
+                repartition->celllist);
+#else
+  pick_metis(nodeID, s, nr_nodes, weights_v, weights_e, repartition->celllist);
+#endif
 
   /* Check that all cells have good values. All nodes have same copy, so just
    * check on one. */
@@ -1351,22 +1367,22 @@ void partition_repartition(struct repartition *reparttype, int nodeID,
                            int nr_nodes, struct space *s, struct task *tasks,
                            int nr_tasks) {
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
 
-  if (reparttype->type == REPART_PARMETIS_VERTEX_EDGE_COSTS) {
-      repart_edge_parmetis(1, 1, 0, reparttype, nodeID, nr_nodes, s, tasks,
-                           nr_tasks);
+  if (reparttype->type == REPART_METIS_VERTEX_EDGE_COSTS) {
+      repart_edge_metis(1, 1, 0, reparttype, nodeID, nr_nodes, s, tasks,
+                        nr_tasks);
 
-  } else if (reparttype->type == REPART_PARMETIS_EDGE_COSTS) {
-      repart_edge_parmetis(0, 1, 0, reparttype, nodeID, nr_nodes, s, tasks,
-                           nr_tasks);
+  } else if (reparttype->type == REPART_METIS_EDGE_COSTS) {
+      repart_edge_metis(0, 1, 0, reparttype, nodeID, nr_nodes, s, tasks,
+                        nr_tasks);
 
-  } else if (reparttype->type == REPART_PARMETIS_VERTEX_COSTS) {
-      repart_edge_parmetis(1, 0, 0, reparttype, nodeID, nr_nodes, s, tasks,
+  } else if (reparttype->type == REPART_METIS_VERTEX_COSTS) {
+      repart_edge_metis(1, 0, 0, reparttype, nodeID, nr_nodes, s, tasks,
                            nr_tasks);
 
-  } else if (reparttype->type == REPART_PARMETIS_VERTEX_COSTS_TIMEBINS) {
-      repart_edge_parmetis(1, 1, 1, reparttype, nodeID, nr_nodes, s, tasks,
+  } else if (reparttype->type == REPART_METIS_VERTEX_COSTS_TIMEBINS) {
+      repart_edge_metis(1, 1, 1, reparttype, nodeID, nr_nodes, s, tasks,
                            nr_tasks);
 
   } else if (reparttype->type == REPART_NONE) {
@@ -1376,7 +1392,7 @@ void partition_repartition(struct repartition *reparttype, int nodeID,
     error("Impossible repartition type");
   }
 #else
-  error("SWIFT was not compiled with ParMETIS support.");
+  error("SWIFT was not compiled with METIS or ParMETIS support.");
 #endif
 }
 
@@ -1433,18 +1449,18 @@ void partition_initial_partition(struct partition *initial_partition,
       return;
     }
 
-  } else if (initial_partition->type == INITPART_PARMETIS_WEIGHT ||
-             initial_partition->type == INITPART_PARMETIS_NOWEIGHT) {
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
-    /* Simple k-way partition selected by ParMETIS using cell particle counts
-     * as weights or not. Should be best when starting with a inhomogeneous
-     * dist.
+  } else if (initial_partition->type == INITPART_METIS_WEIGHT ||
+             initial_partition->type == INITPART_METIS_NOWEIGHT) {
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
+    /* Simple k-way partition selected by METIS using cell particle
+     * counts as weights or not. Should be best when starting with a
+     * inhomogeneous dist.
      */
 
     /* Space for particles per cell counts, which will be used as weights or
      * not. */
     double *weights = NULL;
-    if (initial_partition->type == INITPART_PARMETIS_WEIGHT) {
+    if (initial_partition->type == INITPART_METIS_WEIGHT) {
       if ((weights = (double *)malloc(sizeof(double) * s->nr_cells)) == NULL)
         error("Failed to allocate weights buffer.");
       bzero(weights, sizeof(double) * s->nr_cells);
@@ -1453,8 +1469,8 @@ void partition_initial_partition(struct partition *initial_partition,
       accumulate_counts(s, weights);
 
       /* Get all the counts from all the nodes. */
-      if (MPI_Allreduce(MPI_IN_PLACE, weights, s->nr_cells, MPI_DOUBLE, MPI_SUM,
-                        MPI_COMM_WORLD) != MPI_SUCCESS)
+      if (MPI_Allreduce(MPI_IN_PLACE, weights, s->nr_cells, MPI_DOUBLE,
+                        MPI_SUM, MPI_COMM_WORLD) != MPI_SUCCESS)
         error("Failed to allreduce particle cell weights.");
     }
 
@@ -1462,7 +1478,11 @@ void partition_initial_partition(struct partition *initial_partition,
     int *celllist = NULL;
     if ((celllist = (int *)malloc(sizeof(int) * s->nr_cells)) == NULL)
       error("Failed to allocate celllist");
+#ifdef HAVE_PARMETIS
+    pick_parmetis(nodeID, s, nr_nodes, weights, NULL, 0, 0, 0.0f, celllist);
+#else
     pick_metis(nodeID, s, nr_nodes, weights, NULL, celllist);
+#endif
 
     /* And apply to our cells */
     split_metis(s, nr_nodes, celllist);
@@ -1479,7 +1499,7 @@ void partition_initial_partition(struct partition *initial_partition,
     if (weights != NULL) free(weights);
     free(celllist);
 #else
-    error("SWIFT was not compiled with METIS support");
+    error("SWIFT was not compiled with METIS or ParMETIS support");
 #endif
 
   } else if (initial_partition->type == INITPART_VECTORIZE) {
@@ -1525,7 +1545,7 @@ void partition_init(struct partition *partition,
 #ifdef WITH_MPI
 
 /* Defaults make use of METIS if available */
-#ifdef HAVE_PARMETIS
+#if defined(HAVE_METIS) || defined(HAVE_PARMETIS)
   const char *default_repart = "costs/costs";
   const char *default_part = "memory";
 #else
@@ -1551,12 +1571,12 @@ void partition_init(struct partition *partition,
     case 'v':
       partition->type = INITPART_VECTORIZE;
       break;
-#ifdef HAVE_PARMETIS
+#if defined(HAVE_METIS) || defined(HAVE_PARMETIS)
     case 'r':
-      partition->type = INITPART_PARMETIS_NOWEIGHT;
+      partition->type = INITPART_METIS_NOWEIGHT;
       break;
     case 'm':
-      partition->type = INITPART_PARMETIS_WEIGHT;
+      partition->type = INITPART_METIS_WEIGHT;
       break;
     default:
       message("Invalid choice of initial partition type '%s'.", part_type);
@@ -1568,7 +1588,7 @@ void partition_init(struct partition *partition,
       message("Invalid choice of initial partition type '%s'.", part_type);
       error(
           "Permitted values are: 'grid' or 'vectorized' when compiled "
-          "without ParMETIS.");
+          "without METIS or ParMETIS.");
 #endif
   }
 
@@ -1585,18 +1605,18 @@ void partition_init(struct partition *partition,
   if (strcmp("none/none", part_type) == 0) {
     repartition->type = REPART_NONE;
 
-#ifdef HAVE_PARMETIS
+#if defined(HAVE_METIS) || defined(HAVE_PARMETIS)
   } else if (strcmp("costs/costs", part_type) == 0) {
-    repartition->type = REPART_PARMETIS_VERTEX_EDGE_COSTS;
+    repartition->type = REPART_METIS_VERTEX_EDGE_COSTS;
 
   } else if (strcmp("none/costs", part_type) == 0) {
-    repartition->type = REPART_PARMETIS_EDGE_COSTS;
+    repartition->type = REPART_METIS_EDGE_COSTS;
 
   } else if (strcmp("costs/none", part_type) == 0) {
-    repartition->type = REPART_PARMETIS_VERTEX_COSTS;
+    repartition->type = REPART_METIS_VERTEX_COSTS;
 
   } else if (strcmp("costs/time", part_type) == 0) {
-    repartition->type = REPART_PARMETIS_VERTEX_COSTS_TIMEBINS;
+    repartition->type = REPART_METIS_VERTEX_COSTS_TIMEBINS;
 
   } else {
     message("Invalid choice of re-partition type '%s'.", part_type);
@@ -1606,7 +1626,8 @@ void partition_init(struct partition *partition,
 #else
   } else {
     message("Invalid choice of re-partition type '%s'.", part_type);
-    error("Permitted values are: 'none/none' when compiled without ParMETIS.");
+    error("Permitted values are: 'none/none' when compiled without "
+          "METIS or ParMETIS.");
 #endif
   }
 
@@ -1630,6 +1651,10 @@ void partition_init(struct partition *partition,
         "Invalid DomainDecomposition:minfrac, must be greater than 0 and less "
         "than equal to 1");
 
+  /* Use adaptive or simple refinement when repartitioning. */
+  repartition->adaptive = 
+      parser_get_opt_param_int(params, "DomainDecomposition:adaptive", 1   );
+
   /* Ratio of interprocess communication time to data redistribution time. */
   repartition->itr =
       parser_get_opt_param_float(params, "DomainDecomposition:itr", 0.1f);
diff --git a/src/partition.h b/src/partition.h
index d8c68b9b62c9e86c892557234ce5d351c3e101bb..782d4f411be05c99177e30b780c9a3b8d48cc10f 100644
--- a/src/partition.h
+++ b/src/partition.h
@@ -27,8 +27,8 @@
 enum partition_type {
   INITPART_GRID = 0,
   INITPART_VECTORIZE,
-  INITPART_PARMETIS_WEIGHT,
-  INITPART_PARMETIS_NOWEIGHT
+  INITPART_METIS_WEIGHT,
+  INITPART_METIS_NOWEIGHT
 };
 
 /* Simple descriptions of types for reports. */
@@ -43,10 +43,10 @@ struct partition {
 /* Repartition type to use. */
 enum repartition_type {
   REPART_NONE = 0,
-  REPART_PARMETIS_VERTEX_EDGE_COSTS,
-  REPART_PARMETIS_EDGE_COSTS,
-  REPART_PARMETIS_VERTEX_COSTS,
-  REPART_PARMETIS_VERTEX_COSTS_TIMEBINS
+  REPART_METIS_VERTEX_EDGE_COSTS,
+  REPART_METIS_EDGE_COSTS,
+  REPART_METIS_VERTEX_COSTS,
+  REPART_METIS_VERTEX_COSTS_TIMEBINS
 };
 
 /* Repartition preferences. */
@@ -55,6 +55,7 @@ struct repartition {
   float trigger;
   float minfrac;
   float itr;
+  int adaptive;
 
   /* The partition as a cell-list. */
   int ncelllist;
diff --git a/src/scheduler.c b/src/scheduler.c
index 38313251b89f9f7dd6a6e8dd5d93a7eba3894f5a..3e20b669ac1c1f9e2a8b43d4e10eb61123396508 100644
--- a/src/scheduler.c
+++ b/src/scheduler.c
@@ -1190,7 +1190,7 @@ void scheduler_reweight(struct scheduler *s, int verbose) {
       if (t->unlock_tasks[j]->weight > t->weight)
         t->weight = t->unlock_tasks[j]->weight;
     float cost = 0.f;
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_PARMETIS) || defined(HAVE_METIS))
     int partcost = 1;
 #endif
 
@@ -1279,7 +1279,7 @@ void scheduler_reweight(struct scheduler *s, int verbose) {
         cost = wscale * count_i + wscale * gcount_i;
         break;
       case task_type_send:
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_PARMETIS) || defined(HAVE_METIS))
         partcost = 0;
 #endif
         if (count_i < 1e5)
@@ -1288,7 +1288,7 @@ void scheduler_reweight(struct scheduler *s, int verbose) {
           cost = 2e9;
         break;
       case task_type_recv:
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_PARMETIS) || defined(HAVE_METIS))
         partcost = 0;
 #endif
         if (count_i < 1e5)
@@ -1301,7 +1301,7 @@ void scheduler_reweight(struct scheduler *s, int verbose) {
         break;
     }
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_PARMETIS) || defined(HAVE_METIS))
     if (partcost) t->cost = cost;
 #endif
     t->weight += cost;
diff --git a/src/space.c b/src/space.c
index 25450ba41ebb45bb639a23003c4cf26916977b86..92bac8ebe8c11b93782996ca601f2bfed589247b 100644
--- a/src/space.c
+++ b/src/space.c
@@ -439,8 +439,8 @@ void space_regrid(struct space *s, int verbose) {
         /* Failed, try another technique that requires no settings. */
         message("Failed to get a new partition, trying less optimal method");
         struct partition initial_partition;
-#ifdef HAVE_PARMETIS
-        initial_partition.type = INITPART_PARMETIS_NOWEIGHT;
+#if defined(HAVE_PARMETIS) || defined(HAVE_METIS)
+        initial_partition.type = INITPART_METIS_NOWEIGHT;
 #else
         initial_partition.type = INITPART_VECTORIZE;
 #endif
diff --git a/src/task.h b/src/task.h
index eb352f272e2f3599b57d8c8b4559fbb6236442e4..788a1968a35d3b7e8d77e8fa690331d9082c3e1a 100644
--- a/src/task.h
+++ b/src/task.h
@@ -138,7 +138,7 @@ struct task {
   /*! Weight of the task */
   float weight;
 
-#if defined(WITH_MPI) && defined(HAVE_PARMETIS)
+#if defined(WITH_MPI) && (defined(HAVE_METIS) || defined(HAVE_PARMETIS))
   /*! Individual cost estimate for this task. */
   float cost;
 #endif
diff --git a/src/version.c b/src/version.c
index 5d7a15f149fc3594c22b37239d74c1b3a3df2df1..6fe8c38fc22f3c06fae42adbba83a65aff208bb9 100644
--- a/src/version.c
+++ b/src/version.c
@@ -24,6 +24,9 @@
 /* MPI headers. */
 #ifdef WITH_MPI
 #include <mpi.h>
+#ifdef HAVE_METIS
+#include <metis.h>
+#endif
 #ifdef HAVE_PARMETIS
 #include <parmetis.h>
 #endif
@@ -303,6 +306,23 @@ const char *hdf5_version(void) {
   return version;
 }
 
+/**
+ * @brief return the METIS version used when SWIFT was built.
+ *
+ * @result description of the METIS version.
+ */
+const char *metis_version(void) {
+
+  static char version[256] = {0};
+#if defined(WITH_MPI) && defined(HAVE_METIS)
+  sprintf(version, "%i.%i.%i", METIS_VER_MAJOR, METIS_VER_MINOR,
+          METIS_VER_SUBMINOR);
+#else
+  sprintf(version, "Unknown version");
+#endif
+  return version;
+}
+
 /**
  * @brief return the ParMETIS version used when SWIFT was built.
  *
@@ -421,6 +441,9 @@ void greetings(void) {
 #endif
 #ifdef WITH_MPI
   printf(" MPI library: %s\n", mpi_version());
+#ifdef HAVE_METIS
+  printf(" METIS library version: %s\n", metis_version());
+#endif
 #ifdef HAVE_PARMETIS
   printf(" ParMETIS library version: %s\n", parmetis_version());
 #endif
diff --git a/src/version.h b/src/version.h
index 0d8e77cd9d07dc5e3a2f7d44e27d7bf9da6f7872..75a371bd9e47b19c1887b556accd486da50d9cea 100644
--- a/src/version.h
+++ b/src/version.h
@@ -31,6 +31,7 @@ const char* compilation_cflags(void);
 const char* compiler_name(void);
 const char* compiler_version(void);
 const char* mpi_version(void);
+const char* metis_version(void);
 const char* parmetis_version(void);
 const char* hdf5_version(void);
 const char* fftw3_version(void);