diff --git a/dune/solvers/test/common.hh b/dune/solvers/test/common.hh
index 11bb9ea8d92c75fb31474a5f93558ae79c2b9eed..286185404d5c9fc491b83f9127e2c17ac89f4ace 100644
--- a/dune/solvers/test/common.hh
+++ b/dune/solvers/test/common.hh
@@ -76,6 +76,7 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix)
     typedef typename GridView::IndexSet IndexSet;
     typedef typename GridView::template Codim<0>::Iterator ElementIterator;
     typedef typename GridView::template Codim<0>::Entity Element;
+    typedef typename GridView::template Codim<0>::Entity::Geometry Geometry;
     typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache;
     typedef typename FiniteElementCache::FiniteElementType FiniteElement;
 
@@ -93,6 +94,7 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix)
     for (; it != end; ++it)
     {
         Element& e = *it;
+        const Geometry geometry = e.geometry();
         const FiniteElement& fe = cache.get(e.type());
 
         int localSize = fe.localBasis().size();
@@ -115,10 +117,10 @@ void assemblePQ1Stiffness(const GridView& gridView, Matrix& matrix)
             const LocalCoordinate& quadPos = quad[pt].position();
 
             // get transposed inverse of Jacobian of transformation
-            const JacobianInverseTransposed& invJacobian = e.geometry().jacobianInverseTransposed(quadPos);
+            const JacobianInverseTransposed invJacobian = geometry.jacobianInverseTransposed(quadPos);
 
             // get integration factor
-            const double integrationElement = e.geometry().integrationElement(quadPos);
+            const double integrationElement = geometry.integrationElement(quadPos);
 
             // evaluate gradients of shape functions
             fe.localBasis().evaluateJacobian(quadPos, referenceGradients);
@@ -156,6 +158,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix)
     typedef typename GridView::IndexSet IndexSet;
     typedef typename GridView::template Codim<0>::Iterator ElementIterator;
     typedef typename GridView::template Codim<0>::Entity Element;
+    typedef typename GridView::template Codim<0>::Entity::Geometry Geometry;
     typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache;
     typedef typename FiniteElementCache::FiniteElementType FiniteElement;
 
@@ -171,6 +174,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix)
     for (; it != end; ++it)
     {
         Element& e = *it;
+        const Geometry geometry = e.geometry();
         const FiniteElement& fe = cache.get(e.type());
 
         int localSize = fe.localBasis().size();
@@ -192,7 +196,7 @@ void assemblePQ1Mass(const GridView& gridView, Matrix& matrix)
             const LocalCoordinate& quadPos = quad[pt].position();
 
             // get integration factor
-            const double integrationElement = e.geometry().integrationElement(quadPos);
+            const double integrationElement = geometry.integrationElement(quadPos);
 
             // evaluate basis functions
             fe.localBasis().evaluateFunction(quadPos, values);
@@ -227,6 +231,7 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f)
     typedef typename GridView::IndexSet IndexSet;
     typedef typename GridView::template Codim<0>::Iterator ElementIterator;
     typedef typename GridView::template Codim<0>::Entity Element;
+    typedef typename GridView::template Codim<0>::Entity::Geometry Geometry;
     typedef typename Dune::PQkLocalFiniteElementCache<double, double, dim, 1> FiniteElementCache;
     typedef typename FiniteElementCache::FiniteElementType FiniteElement;
 
@@ -243,6 +248,7 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f)
     for (; it != end; ++it)
     {
         Element& e = *it;
+        const Geometry geometry = e.geometry();
         const FiniteElement& fe = cache.get(e.type());
 
         int localSize = fe.localBasis().size();
@@ -264,14 +270,14 @@ void assemblePQ1RHS(const GridView& gridView, Vector& r, const Function& f)
             const LocalCoordinate& quadPos = quad[pt].position();
 
             // get integration factor
-            const double integrationElement = e.geometry().integrationElement(quadPos);
+            const double integrationElement = geometry.integrationElement(quadPos);
 
             // evaluate basis functions
             fe.localBasis().evaluateFunction(quadPos, values);
 
             // evaluate function
             FunctionRangeType fAtPos;
-            f.evaluate(e.geometry().global(quadPos), fAtPos);
+            f.evaluate(geometry.global(quadPos), fAtPos);
 
             // add vector entries
             double z = quad[pt].weight() * integrationElement;