diff --git a/src/bisection-example-new.cc b/src/bisection-example-new.cc
index f5bc2f02fc7fbc4d0f8f39d07034cf4144135ee7..f075e31abfbb01a5fe4515c4945118cd11b40e92 100644
--- a/src/bisection-example-new.cc
+++ b/src/bisection-example-new.cc
@@ -71,14 +71,11 @@ class SampleFunctional {
       A_.mv(x, tmp3);
       double const rest_b = (b_ - tmp3) * descDir;
 
-      Dune::BlockVector<SmallVector> xx; // FIXME: we actually want x here
-      Dune::BlockVector<SmallVector> dd; // FIXME: we actually want descDir here
-
       typedef MyNonlinearity<dimension, SampleFunction> MyNonlinearityType;
       MyNonlinearityType phi;
       typedef DirectionalConvexFunction<MyNonlinearityType>
       MyDirectionalConvexFunctionType;
-      MyDirectionalConvexFunctionType rest(rest_A, rest_b, phi, xx, dd);
+      MyDirectionalConvexFunctionType rest(rest_A, rest_b, phi, x, descDir);
 
       // Experiment a bit
       Interval<double> D;
diff --git a/src/mynonlinearity.cc b/src/mynonlinearity.cc
index db756cbb49fdbb9fd80f47eb2adfb621ed2a9883..4efdd603c12ac0d2f3d0f748745dbab141268a77 100644
--- a/src/mynonlinearity.cc
+++ b/src/mynonlinearity.cc
@@ -21,19 +21,19 @@ class MyNonlinearity {
   typedef Dune::FieldVector<double, dimension> SmallVector;
   typedef Dune::FieldMatrix<double, dimension, dimension> SmallMatrix;
 
-  typedef Dune::BlockVector<SmallVector> VectorType;
-  typedef Dune::BlockVector<SmallMatrix> MatrixType;
+  typedef SmallVector VectorType;
+  typedef SmallMatrix MatrixType;
 
   void directionalSubDiff(VectorType u, VectorType v, Interval<double>& D) {
-    // if (u == SmallVector(0.0)) {
-    //   D[0] = D[1] = func_.rightDifferential(0);
-    // } else if (u * v > 0) {
-    //   D[1] = (v * u) * func_.rightDifferential(u.two_norm())/u.two_norm();
-    //   D[0] = (v * u) * func_.leftDifferential(u.two_norm())/u.two_norm();
-    // } else {
-    //   D[1] = (v * u) * func_.leftDifferential(u.two_norm())/u.two_norm();
-    //   D[0] = (v * u) * func_.rightDifferential(u.two_norm())/u.two_norm();
-    // }
+    if (u == SmallVector(0.0)) {
+      D[0] = D[1] = func_.rightDifferential(0);
+    } else if (u * v > 0) {
+      D[1] = (v * u) * func_.rightDifferential(u.two_norm()) / u.two_norm();
+      D[0] = (v * u) * func_.leftDifferential(u.two_norm()) / u.two_norm();
+    } else {
+      D[1] = (v * u) * func_.leftDifferential(u.two_norm()) / u.two_norm();
+      D[0] = (v * u) * func_.rightDifferential(u.two_norm()) / u.two_norm();
+    }
   }
 
   void directionalDomain(const VectorType&, const VectorType&,