diff --git a/dune/common/fmatrix.hh b/dune/common/fmatrix.hh
index 2c7e2112dff0ff32294b73b1a5d7991cb4794702..95190a0b55ee8d4f7413bd0c6356c6e530e1a3d2 100644
--- a/dune/common/fmatrix.hh
+++ b/dune/common/fmatrix.hh
@@ -7,6 +7,7 @@
 #include <cmath>
 #include <cstddef>
 #include <iostream>
+#include <initializer_list>
 
 #include <dune/common/exceptions.hh>
 #include <dune/common/fvector.hh>
@@ -96,6 +97,15 @@ namespace Dune
       DenseMatrixAssigner< FieldMatrix< K, ROWS, COLS >, Other >::apply( *this, other );
     }
 
+    FieldMatrix (std::initializer_list<std::initializer_list<K> > const &ll)
+    {
+      assert(ll.size() == rows);
+      size_t i = 0;
+      for (typename std::initializer_list<std::initializer_list<K> >::
+             const_iterator lit = ll.begin(); lit != ll.end(); ++lit)
+        _data[i++] = *lit;
+    }
+
     //===== assignment
     using Base::operator=;
 
diff --git a/dune/common/fvector.hh b/dune/common/fvector.hh
index 76c77bc2c2739774e631958dd4ddc0ca964d2f41..710de61fd4053f9adbcd695f3f274125e025c968 100644
--- a/dune/common/fvector.hh
+++ b/dune/common/fvector.hh
@@ -10,6 +10,7 @@
 #include <complex>
 #include <cstring>
 #include <utility>
+#include<initializer_list>
 
 #include <dune/common/std/constexpr.hh>
 
@@ -129,6 +130,15 @@ namespace Dune {
     FieldVector (const FieldVector & x) : _data(x._data)
     {}
 
+    FieldVector (std::initializer_list<K> const &l)
+    {
+      assert(l.size() == dimension);
+      size_t i = 0;
+      for (typename std::initializer_list<K>::const_iterator it = l.begin();
+           it != l.end(); ++it)
+        _data[i++] = *it;
+    }
+
     /**
      * \brief Copy constructor from a second vector of possibly different type
      *
diff --git a/dune/common/test/fmatrixtest.cc b/dune/common/test/fmatrixtest.cc
index 2fb4389e0e20d0e803f0cacfdf7d3937c3a7ee28..48e659cde3cc8c4ecc04992c599aab5614baad39 100644
--- a/dune/common/test/fmatrixtest.cc
+++ b/dune/common/test/fmatrixtest.cc
@@ -553,11 +553,24 @@ void test_interface()
   checkMatrixInterface< FMatrix, Traits >( m );
 }
 
+void test_initialisation()
+{
+  Dune::FieldMatrix<int, 2, 2> const A = {
+    { 1, 2 },
+    { 3, 4 }
+  };
+  assert(A[0][0] == 1);
+  assert(A[0][1] == 2);
+  assert(A[1][0] == 3);
+  assert(A[1][1] == 4);
+}
+
 int main()
 {
   try {
     test_nan();
     test_infinity_norms();
+    test_initialisation();
 
     // test 1 x 1 matrices
     test_interface<float, 1, 1>();
diff --git a/dune/common/test/fvectortest.cc b/dune/common/test/fvectortest.cc
index 81ac1913f5f6b386758204f7f43c566bfafdc03a..b1dd3d0f5f7748b17f8396a0eb00cf3cc20997a8 100644
--- a/dune/common/test/fvectortest.cc
+++ b/dune/common/test/fvectortest.cc
@@ -342,6 +342,15 @@ test_infinity_norms()
   assert(std::abs(v.infinity_norm_real()-14.0) < 1e-10); // max(7,14)
 }
 
+void
+test_initialisation()
+{
+  Dune::FieldVector<int, 2> const b = { 1, 2 };
+
+  assert(b[0] == 1);
+  assert(b[1] == 2);
+}
+
 int main()
 {
   try {
@@ -351,6 +360,7 @@ int main()
 
     test_nan();
     test_infinity_norms();
+    test_initialisation();
   } catch (Dune::Exception& e) {
     std::cerr << e << std::endl;
     return 1;