]> www.ginac.de Git - ginac.git/commitdiff
Use neseted initializer lists to construct matrix objects.
authorRichard Kreckel <kreckel@ginac.de>
Sat, 28 Nov 2015 14:43:46 +0000 (15:43 +0100)
committerRichard Kreckel <kreckel@ginac.de>
Sun, 29 Nov 2015 10:32:23 +0000 (11:32 +0100)
Add constructor of initializer_list<initializer_list<ex>> to matrix.
Use this syntax where, previously, ctor from comma-separated list of
elements was used. Deprecate the ctor from comma-separated list.

Note: The output format '[[a,b],[c,d]]' and ginsh syntax are
unchanged because lists are printed '{a,b,c}' and a matrix is not a
list of lists.

check/exam_archive.cpp
check/exam_clifford.cpp
check/exam_indexed.cpp
check/exam_matrices [new file with mode: 0755]
check/exam_matrices.cpp
check/time_lw_O.cpp
doc/tutorial/ginac.texi
ginac/matrix.cpp
ginac/matrix.h

index 8a9bf5f91cdf7ae7ec12f31a8df7edcd01aa824c..f6ca95e933e73eb00cd04a22644bf0ee180f9923 100644 (file)
@@ -43,7 +43,7 @@ unsigned exam_archive()
          + lorentz_g(
              varidx(lst{x, -11*y, acos(2*x).series(x==3-5*I, 3)} * color_ONE()
                * metric_tensor(varidx(log(cos(128.0/(x*y))), 5), varidx(2, 5)), zeta(3)),
          + lorentz_g(
              varidx(lst{x, -11*y, acos(2*x).series(x==3-5*I, 3)} * color_ONE()
                * metric_tensor(varidx(log(cos(128.0/(x*y))), 5), varidx(2, 5)), zeta(3)),
-             varidx(diag_matrix(lst{-1, Euler, atan(x/y==-15*I/17)})
+             varidx(diag_matrix({-1, Euler, atan(x/y==-15*I/17)})
                * delta_tensor(idx(x, 2), idx(wild(7), 3)), zeta(3), true),
              true
            )
                * delta_tensor(idx(x, 2), idx(wild(7), 3)), zeta(3), true),
              true
            )
index 0430e2daccd4df8bedb1ec1a55297cf3d1b7408d..ec64ba9604711198c32db37b3aa04a40e7f0401a 100644 (file)
@@ -432,14 +432,13 @@ template <typename IDX> unsigned clifford_check6(const matrix &A)
        result += check_equal_lst(clifford_to_lst(lst_to_clifford(e, c), c, true), e);
 
 /* Moebius map (both forms) checks for symmetric metrics only */
        result += check_equal_lst(clifford_to_lst(lst_to_clifford(e, c), c, true), e);
 
 /* Moebius map (both forms) checks for symmetric metrics only */
-       matrix M1(2, 2),  M2(2, 2);
        c = clifford_unit(nu, A);
 
        e = clifford_moebius_map(0, dirac_ONE(), 
        c = clifford_unit(nu, A);
 
        e = clifford_moebius_map(0, dirac_ONE(), 
-                                                        dirac_ONE(), 0, lst{t, x, y, z}, A);
+                                dirac_ONE(), 0, lst{t, x, y, z}, A);
 /* this is just the inversion*/
 /* this is just the inversion*/
-       M1 = 0, dirac_ONE(),
-               dirac_ONE(), 0;
+       matrix M1 = {{0, dirac_ONE()},
+                    {dirac_ONE(), 0}};
        e1 = clifford_moebius_map(M1, lst{t, x, y, z}, A);
 /* the inversion again*/
        result += check_equal_lst(e, e1);
        e1 = clifford_moebius_map(M1, lst{t, x, y, z}, A);
 /* the inversion again*/
        result += check_equal_lst(e, e1);
@@ -448,10 +447,10 @@ template <typename IDX> unsigned clifford_check6(const matrix &A)
        result += check_equal_lst(e, e1);
 
        e = clifford_moebius_map(dirac_ONE(), lst_to_clifford(lst{1, 2, 3, 4}, nu, A),
        result += check_equal_lst(e, e1);
 
        e = clifford_moebius_map(dirac_ONE(), lst_to_clifford(lst{1, 2, 3, 4}, nu, A),
-                                                        0, dirac_ONE(), lst{t, x, y, z}, A);
+                                0, dirac_ONE(), lst{t, x, y, z}, A);
 /*this is just a shift*/
 /*this is just a shift*/
-       M2 = dirac_ONE(), lst_to_clifford(lst{1, 2, 3, 4}, c),
-               0, dirac_ONE();
+       matrix M2 = {{dirac_ONE(), lst_to_clifford(lst{1, 2, 3, 4}, c),},
+                    {0, dirac_ONE()}};
        e1 = clifford_moebius_map(M2, lst{t, x, y, z}, c);
 /* the same shift*/
        result += check_equal_lst(e, e1);
        e1 = clifford_moebius_map(M2, lst{t, x, y, z}, c);
 /* the same shift*/
        result += check_equal_lst(e, e1);
@@ -546,7 +545,7 @@ static unsigned clifford_check8()
        realsymbol a("a");
        varidx mu(symbol("mu", "\\mu"), 1);
 
        realsymbol a("a");
        varidx mu(symbol("mu", "\\mu"), 1);
 
-       ex e = clifford_unit(mu, diag_matrix(lst{-1})), e0 = e.subs(mu==0);
+       ex e = clifford_unit(mu, diag_matrix({-1})), e0 = e.subs(mu==0);
        result += ( exp(a*e0)*e0*e0 == -exp(e0*a) ) ? 0 : 1;
 
        return result;
        result += ( exp(a*e0)*e0*e0 == -exp(e0*a) ) ? 0 : 1;
 
        return result;
@@ -565,45 +564,45 @@ unsigned exam_clifford()
        result += clifford_check5(); cout << '.' << flush;
 
        // anticommuting, symmetric examples
        result += clifford_check5(); cout << '.' << flush;
 
        // anticommuting, symmetric examples
-       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix(lst{-1, 1, 1, 1})));
-       result += clifford_check6<idx>(ex_to<matrix>(diag_matrix(lst{-1, 1, 1, 1})));; cout << '.' << flush;
-       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix(lst{-1, -1, -1, -1})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix(lst{-1, -1, -1, -1})));; cout << '.' << flush;
-       result += clifford_check6<idx>(ex_to<matrix>(diag_matrix(lst{-1, 1, 1, -1})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix(lst{-1, 1, 1, -1})));; cout << '.' << flush;
-       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix(lst{-1, 0, 1, -1})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix(lst{-1, 0, 1, -1})));; cout << '.' << flush;
-       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix(lst{-3, 0, 2, -1})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix(lst{-3, 0, 2, -1})));; cout << '.' << flush;
+       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix({-1, 1, 1, 1})));
+       result += clifford_check6<idx>(ex_to<matrix>(diag_matrix({-1, 1, 1, 1})));; cout << '.' << flush;
+       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix({-1, -1, -1, -1})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix({-1, -1, -1, -1})));; cout << '.' << flush;
+       result += clifford_check6<idx>(ex_to<matrix>(diag_matrix({-1, 1, 1, -1})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix({-1, 1, 1, -1})));; cout << '.' << flush;
+       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix({-1, 0, 1, -1})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix({-1, 0, 1, -1})));; cout << '.' << flush;
+       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix({-3, 0, 2, -1})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix({-3, 0, 2, -1})));; cout << '.' << flush;
 
        realsymbol s("s"), t("t"); // symbolic entries in matric
 
        realsymbol s("s"), t("t"); // symbolic entries in matric
-       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix(lst{-1, 1, s, t})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix(lst{-1, 1, s, t})));; cout << '.' << flush;
+       result += clifford_check6<varidx>(ex_to<matrix>(diag_matrix({-1, 1, s, t})))+clifford_check6<idx>(ex_to<matrix>(diag_matrix({-1, 1, s, t})));; cout << '.' << flush;
 
        matrix A(4, 4);
 
        matrix A(4, 4);
-       A = 1,  0,  0,  0, // anticommuting, not symmetric, Tr=0
-           0, -1,  0,  0,
-           0,  0,  0, -1,
-           0,  0,  1,  0;
+       A = {{1,  0,  0,  0}, // anticommuting, not symmetric, Tr=0
+            {0, -1,  0,  0},
+            {0,  0,  0, -1},
+            {0,  0,  1,  0}};
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
-       A = 1,  0,  0,  0, // anticommuting, not symmetric, Tr=2
-           0,  1,  0,  0,
-           0,  0,  0, -1,
-           0,  0,  1,  0;
+       A = {{1,  0,  0,  0}, // anticommuting, not symmetric, Tr=2
+            {0,  1,  0,  0},
+            {0,  0,  0, -1},
+            {0,  0,  1,  0}};
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
-       A = 1,  0,  0,  0, // not anticommuting, symmetric, Tr=0
-           0, -1,  0,  0,
-           0,  0,  0, -1,
-           0,  0, -1,  0;
+       A = {{1,  0,  0,  0}, // not anticommuting, symmetric, Tr=0
+            {0, -1,  0,  0},
+            {0,  0,  0, -1},
+            {0,  0, -1,  0}};
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
-       A = 1,  0,  0,  0, // not anticommuting, symmetric, Tr=2
-           0,  1,  0,  0,
-           0,  0,  0, -1,
-           0,  0, -1,  0;
+       A = {{1,  0,  0,  0}, // not anticommuting, symmetric, Tr=2
+            {0,  1,  0,  0},
+            {0,  0,  0, -1},
+            {0,  0, -1,  0}};
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
-       A = 1,  1,  0,  0, // not anticommuting, not symmetric, Tr=4
-           0,  1,  1,  0,
-           0,  0,  1,  1,
-           0,  0,  0,  1;
+       A = {{1,  1,  0,  0}, // not anticommuting, not symmetric, Tr=4
+            {0,  1,  1,  0},
+            {0,  0,  1,  1},
+            {0,  0,  0,  1}};
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
        symbol dim("D");
        result += clifford_check6<varidx>(A)+clifford_check6<idx>(A);; cout << '.' << flush;
 
        symbol dim("D");
index cf58ac7fbe7ded16cb5f78ca97d92bad3f010ce1..efb12c6368fac5507627f82e91af84634285a6ed 100644 (file)
@@ -244,18 +244,16 @@ static unsigned edyn_check()
        symbol Bx("Bx"), By("By"), Bz("Bz");
 
        // Lorentz transformation matrix (boost along x axis)
        symbol Bx("Bx"), By("By"), Bz("Bz");
 
        // Lorentz transformation matrix (boost along x axis)
-       matrix L(4, 4);
-       L =       gamma, -beta*gamma, 0, 0,
-           -beta*gamma,       gamma, 0, 0,
-                     0,           0, 1, 0,
-                     0,           0, 0, 1;
+       matrix L = {{      gamma, -beta*gamma, 0, 0},
+                   {-beta*gamma,       gamma, 0, 0},
+                   {          0,           0, 1, 0},
+                   {          0,           0, 0, 1}};
 
        // Electromagnetic field tensor
 
        // Electromagnetic field tensor
-       matrix F(4, 4);
-       F =  0, -Ex, -Ey, -Ez,
-               Ex,   0, -Bz,  By,
-               Ey,  Bz,   0, -Bx,
-               Ez, -By,  Bx,   0;
+       matrix F = {{ 0, -Ex, -Ey, -Ez},
+                   {Ex,   0, -Bz,  By},
+                   {Ey,  Bz,   0, -Bx},
+                   {Ez, -By,  Bx,   0}};
 
        // Indices
        symbol s_mu("mu"), s_nu("nu"), s_rho("rho"), s_sigma("sigma");
 
        // Indices
        symbol s_mu("mu"), s_nu("nu"), s_rho("rho"), s_sigma("sigma");
diff --git a/check/exam_matrices b/check/exam_matrices
new file mode 100755 (executable)
index 0000000..92ffc81
Binary files /dev/null and b/check/exam_matrices differ
index b5b27ffc37a068719c598b696c2fa93eb34e9d99..6e06330af8159fc95ad77a2cdf4a7d81c4e2b0b2 100644 (file)
@@ -37,7 +37,7 @@ static unsigned matrix_determinants()
        symbol g("g"), h("h"), i("i");
        
        // check symbolic trivial matrix determinant
        symbol g("g"), h("h"), i("i");
        
        // check symbolic trivial matrix determinant
-       m1.set(0,0,a);
+       m1 = matrix{{a}};
        det = m1.determinant();
        if (det != a) {
                clog << "determinant of 1x1 matrix " << m1
        det = m1.determinant();
        if (det != a) {
                clog << "determinant of 1x1 matrix " << m1
@@ -46,8 +46,8 @@ static unsigned matrix_determinants()
        }
        
        // check generic dense symbolic 2x2 matrix determinant
        }
        
        // check generic dense symbolic 2x2 matrix determinant
-       m2.set(0,0,a).set(0,1,b);
-       m2.set(1,0,c).set(1,1,d);
+       m2 = matrix{{a, b},
+                   {c, d}};
        det = m2.determinant();
        if (det != (a*d-b*c)) {
                clog << "determinant of 2x2 matrix " << m2
        det = m2.determinant();
        if (det != (a*d-b*c)) {
                clog << "determinant of 2x2 matrix " << m2
@@ -56,9 +56,9 @@ static unsigned matrix_determinants()
        }
        
        // check generic dense symbolic 3x3 matrix determinant
        }
        
        // check generic dense symbolic 3x3 matrix determinant
-       m3.set(0,0,a).set(0,1,b).set(0,2,c);
-       m3.set(1,0,d).set(1,1,e).set(1,2,f);
-       m3.set(2,0,g).set(2,1,h).set(2,2,i);
+       m3 = matrix{{a, b, c},
+                   {d, e, f},
+                   {g, h, i}};
        det = m3.determinant();
        if (det != (a*e*i - a*f*h - d*b*i + d*c*h + g*b*f - g*c*e)) {
                clog << "determinant of 3x3 matrix " << m3
        det = m3.determinant();
        if (det != (a*e*i - a*f*h - d*b*i + d*c*h + g*b*f - g*c*e)) {
                clog << "determinant of 3x3 matrix " << m3
@@ -67,9 +67,9 @@ static unsigned matrix_determinants()
        }
        
        // check dense numeric 3x3 matrix determinant
        }
        
        // check dense numeric 3x3 matrix determinant
-       m3.set(0,0,numeric(0)).set(0,1,numeric(-1)).set(0,2,numeric(3));
-       m3.set(1,0,numeric(3)).set(1,1,numeric(-2)).set(1,2,numeric(2));
-       m3.set(2,0,numeric(3)).set(2,1,numeric(4)).set(2,2,numeric(-2));
+       m3 = matrix{{0, -1,  3},
+                   {3, -2,  2},
+                   {3,  4, -2}};
        det = m3.determinant();
        if (det != 42) {
                clog << "determinant of 3x3 matrix " << m3
        det = m3.determinant();
        if (det != 42) {
                clog << "determinant of 3x3 matrix " << m3
@@ -78,8 +78,8 @@ static unsigned matrix_determinants()
        }
        
        // check dense symbolic 2x2 matrix determinant
        }
        
        // check dense symbolic 2x2 matrix determinant
-       m2.set(0,0,a/(a-b)).set(0,1,1);
-       m2.set(1,0,b/(a-b)).set(1,1,1);
+       m2 = matrix{{a/(a-b), 1},
+                   {b/(a-b), 1}};
        det = m2.determinant();
        if (det != 1) {
                if (det.normal() == 1)  // only half wrong
        det = m2.determinant();
        if (det != 1) {
                if (det.normal() == 1)  // only half wrong
@@ -101,9 +101,9 @@ static unsigned matrix_determinants()
        }
        
        // check characteristic polynomial
        }
        
        // check characteristic polynomial
-       m3.set(0,0,a).set(0,1,-2).set(0,2,2);
-       m3.set(1,0,3).set(1,1,a-1).set(1,2,2);
-       m3.set(2,0,3).set(2,1,4).set(2,2,a-3);
+       m3 = matrix{{a, -2,   2},
+                   {3, a-1,  2},
+                   {3,  4,  a-3}};
        ex p = m3.charpoly(a);
        if (p != 0) {
                clog << "charpoly of 3x3 matrix " << m3
        ex p = m3.charpoly(a);
        if (p != 0) {
                clog << "charpoly of 3x3 matrix " << m3
@@ -135,10 +135,9 @@ static unsigned matrix_invert1()
 static unsigned matrix_invert2()
 {
        unsigned result = 0;
 static unsigned matrix_invert2()
 {
        unsigned result = 0;
-       matrix m(2,2);
        symbol a("a"), b("b"), c("c"), d("d");
        symbol a("a"), b("b"), c("c"), d("d");
-       m.set(0,0,a).set(0,1,b);
-       m.set(1,0,c).set(1,1,d);
+       matrix m = {{a, b},
+                   {c, d}};
        matrix m_i = m.inverse();
        ex det = m.determinant();
        
        matrix m_i = m.inverse();
        ex det = m.determinant();
        
@@ -157,13 +156,12 @@ static unsigned matrix_invert2()
 static unsigned matrix_invert3()
 {
        unsigned result = 0;
 static unsigned matrix_invert3()
 {
        unsigned result = 0;
-       matrix m(3,3);
        symbol a("a"), b("b"), c("c");
        symbol d("d"), e("e"), f("f");
        symbol g("g"), h("h"), i("i");
        symbol a("a"), b("b"), c("c");
        symbol d("d"), e("e"), f("f");
        symbol g("g"), h("h"), i("i");
-       m.set(0,0,a).set(0,1,b).set(0,2,c);
-       m.set(1,0,d).set(1,1,e).set(1,2,f);
-       m.set(2,0,g).set(2,1,h).set(2,2,i);
+       matrix m = {{a, b, c},
+                   {d, e, f},
+                   {g, h, i}};
        matrix m_i = m.inverse();
        ex det = m.determinant();
        
        matrix m_i = m.inverse();
        ex det = m.determinant();
        
@@ -187,36 +185,32 @@ static unsigned matrix_invert3()
 static unsigned matrix_solve2()
 {
        // check the solution of the multiple system A*X = B:
 static unsigned matrix_solve2()
 {
        // check the solution of the multiple system A*X = B:
-       //       [ 1  2 -1 ] [ x0 y0 ]   [ 4 0 ]
-       //       [ 1  4 -2 ]*[ x1 y1 ] = [ 7 0 ]
-       //       [ a -2  2 ] [ x2 y2 ]   [ a 4 ]
+       //       [ 1  2 -1 ] [ x0 y0 ]   [ 4 0 ]
+       //       [ 1  4 -2 ]*[ x1 y1 ] = [ 7 0 ]
+       //       [ a -2  2 ] [ x2 y2 ]   [ a 4 ]
        unsigned result = 0;
        symbol a("a");
        symbol x0("x0"), x1("x1"), x2("x2");
        symbol y0("y0"), y1("y1"), y2("y2");
        unsigned result = 0;
        symbol a("a");
        symbol x0("x0"), x1("x1"), x2("x2");
        symbol y0("y0"), y1("y1"), y2("y2");
-       matrix A(3,3);
-       A.set(0,0,1).set(0,1,2).set(0,2,-1);
-       A.set(1,0,1).set(1,1,4).set(1,2,-2);
-       A.set(2,0,a).set(2,1,-2).set(2,2,2);
-       matrix B(3,2);
-       B.set(0,0,4).set(1,0,7).set(2,0,a);
-       B.set(0,1,0).set(1,1,0).set(2,1,4);
-       matrix X(3,2);
-       X.set(0,0,x0).set(1,0,x1).set(2,0,x2);
-       X.set(0,1,y0).set(1,1,y1).set(2,1,y2);
-       matrix cmp(3,2);
-       cmp.set(0,0,1).set(1,0,3).set(2,0,3);
-       cmp.set(0,1,0).set(1,1,2).set(2,1,4);
+       matrix A = {{1,  2, -1},
+                   {1,  4, -2},
+                   {a, -2,  2}};
+       matrix B = {{4, 0},
+                   {7, 0},
+                   {a, 4}};
+       matrix X = {{x0 ,y0},
+                   {x1, y1},
+                   {x2, y2}};
+       matrix cmp = {{1, 0},
+                     {3, 2},
+                     {3, 4}};
        matrix sol(A.solve(X, B));
        matrix sol(A.solve(X, B));
-       for (unsigned ro=0; ro<3; ++ro)
-               for (unsigned co=0; co<2; ++co)
-                       if (cmp(ro,co) != sol(ro,co))
-                               result = 1;
-       if (result) {
+       if (cmp != sol) {
                clog << "Solving " << A << " * " << X << " == " << B << endl
                     << "erroneously returned " << sol << endl;
                clog << "Solving " << A << " * " << X << " == " << B << endl
                     << "erroneously returned " << sol << endl;
+               result = 1;
        }
        }
-       
+
        return result;
 }
 
        return result;
 }
 
@@ -224,16 +218,12 @@ static unsigned matrix_evalm()
 {
        unsigned result = 0;
 
 {
        unsigned result = 0;
 
-       matrix S(2, 2, lst{
-               1, 2,
-               3, 4
-       }), T(2, 2, lst{
-               1, 1,
-               2, -1
-       }), R(2, 2, lst{
-               27, 14,
-               36, 26
-       });
+       matrix S {{1, 2},
+                 {3, 4}};
+       matrix T {{1, 1},
+                 {2, -1}};
+       matrix R {{27, 14},
+                 {36, 26}};
 
        ex e = ((S + T) * (S + 2*T));
        ex f = e.evalm();
 
        ex e = ((S + T) * (S + 2*T));
        ex f = e.evalm();
@@ -258,18 +248,18 @@ static unsigned matrix_rank()
        }
 
        // a trivial rank one example
        }
 
        // a trivial rank one example
-       m = 1, 0, 0,
-           2, 0, 0,
-           3, 0, 0;
+       m = {{1, 0, 0},
+            {2, 0, 0},
+            {3, 0, 0}};
        if (m.rank() != 1) {
                clog << "The rank of " << m << " was not computed correctly." << endl;
                ++result;
        }
 
        // an example from Maple's help with rank two
        if (m.rank() != 1) {
                clog << "The rank of " << m << " was not computed correctly." << endl;
                ++result;
        }
 
        // an example from Maple's help with rank two
-       m = x,  1,  0,
-           0,  0,  1,
-          x*y, y,  1;
+       m = {{x,   1,  0},
+            {0,   0,  1},
+            {x*y, y,  1}};
        if (m.rank() != 2) {
                clog << "The rank of " << m << " was not computed correctly." << endl;
                ++result;
        if (m.rank() != 2) {
                clog << "The rank of " << m << " was not computed correctly." << endl;
                ++result;
@@ -288,10 +278,9 @@ static unsigned matrix_rank()
 static unsigned matrix_misc()
 {
        unsigned result = 0;
 static unsigned matrix_misc()
 {
        unsigned result = 0;
-       matrix m1(2,2);
        symbol a("a"), b("b"), c("c"), d("d"), e("e"), f("f");
        symbol a("a"), b("b"), c("c"), d("d"), e("e"), f("f");
-       m1.set(0,0,a).set(0,1,b);
-       m1.set(1,0,c).set(1,1,d);
+       matrix m1 = {{a, b},
+                    {c, d}};
        ex tr = trace(m1);
        
        // check a simple trace
        ex tr = trace(m1);
        
        // check a simple trace
@@ -308,10 +297,9 @@ static unsigned matrix_misc()
                         << " erroneously returned " << m2 << endl;
                ++result;
        }
                         << " erroneously returned " << m2 << endl;
                ++result;
        }
-       matrix m3(3,2);
-       m3.set(0,0,a).set(0,1,b);
-       m3.set(1,0,c).set(1,1,d);
-       m3.set(2,0,e).set(2,1,f);
+       matrix m3 = {{a, b},
+                    {c, d},
+                    {e, f}};
        if (transpose(transpose(m3)) != m3) {
                clog << "transposing 3x2 matrix " << m3 << " twice"
                     << " erroneously returned " << transpose(transpose(m3)) << endl;
        if (transpose(transpose(m3)) != m3) {
                clog << "transposing 3x2 matrix " << m3 << " twice"
                     << " erroneously returned " << transpose(transpose(m3)) << endl;
index 715d5019ef66114861287970a063126e26942d0c..b7e41c0ef5477bcc7754bd13743b78e2e1cd5d6c 100644 (file)
@@ -37,66 +37,63 @@ static const symbol c1("c1"), c2("c2"), c3("c3"), c4("c4"), c5("c5"), c6("c6");
 
 static const ex det1()
 {
 
 static const ex det1()
 {
-       matrix d1(15,15);
-       d1 = a6, a5, a4, a3, a2, a1, 0,  0,  0,  0,  0,  0,  0,  0,  0,
-            0,  0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,
-            0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,  0,
-            0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0,  0,
-            0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0,
-            0,  0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1,
-            0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0,  0,
-            0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0,
-            0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,  0,
-            0,  0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,
-            0,  0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1,
-            0,  0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1,
-            0,  0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,
-            0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,  0,
-            0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0;
+       matrix d1 = {{a6, a5, a4, a3, a2, a1, 0,  0,  0,  0,  0,  0,  0,  0,  0},
+                    {0,  0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0},
+                    {0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0,  0},
+                    {0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0},
+                    {0,  0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1},
+                    {0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0,  0},
+                    {0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0},
+                    {0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0},
+                    {0,  0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1},
+                    {0,  0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1},
+                    {0,  0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0},
+                    {0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0}};
 
        return d1.determinant();
 }
 
 static const ex det2()
 {
 
        return d1.determinant();
 }
 
 static const ex det2()
 {
-       matrix d2(15,15);
-       d2 = b6, b5, b4, b3, b2, b1, 0,  0,  0,  0,  0,  0,  0,  0,  0,
-            0,  0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,
-            0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,  0,
-            0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0,  0,
-            0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0,
-            0,  0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1,
-            0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0,  0,
-            0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0,
-            0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,  0,
-            0,  0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,
-            0,  0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1,
-            0,  0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1,
-            0,  0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,
-            0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,  0,
-            0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0;
+       matrix d2 = {{b6, b5, b4, b3, b2, b1, 0,  0,  0,  0,  0,  0,  0,  0,  0},
+                    {0,  0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0},
+                    {0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0,  0},
+                    {0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0},
+                    {0,  0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1},
+                    {0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0,  0},
+                    {0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0},
+                    {0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0},
+                    {0,  0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1},
+                    {0,  0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1},
+                    {0,  0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0},
+                    {0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0}};
 
        return d2.determinant();
 }
 
 static const ex det3()
 {
 
        return d2.determinant();
 }
 
 static const ex det3()
 {
-       matrix d3(15,15);
-       d3 = c6, c5, c4, c3, c2, c1, 0,  0,  0,  0,  0,  0,  0,  0,  0,
-            0,  0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,
-            0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,  0,
-            0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0,  0,
-            0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0,
-            0,  0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1,
-            0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0,  0,
-            0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0,
-            0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,  0,
-            0,  0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,
-            0,  0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1,
-            0,  0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1,
-            0,  0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,
-            0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,  0,
-            0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0;
+       matrix d3 = {{c6, c5, c4, c3, c2, c1, 0,  0,  0,  0,  0,  0,  0,  0,  0},
+                    {0,  0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0},
+                    {0,  c6, 0,  c5, c4, 0,  c3, c2, c1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0,  0},
+                    {0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1, 0},
+                    {0,  0,  0,  0,  0,  c6, 0,  0,  c5, c4, 0,  0,  c3, c2, c1},
+                    {0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0,  0},
+                    {0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1, 0},
+                    {0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  a6, 0,  a5, a4, 0,  a3, a2, a1, 0,  0,  0,  0,  0},
+                    {0,  0,  0,  0,  0,  a6, 0,  0,  a5, a4, 0,  0,  a3, a2, a1},
+                    {0,  0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1},
+                    {0,  0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0},
+                    {0,  b6, 0,  b5, b4, 0,  b3, b2, b1, 0,  0,  0,  0,  0,  0},
+                    {0,  0,  0,  0,  b6, 0,  0,  b5, b4, 0,  0,  b3, b2, b1, 0}};
 
        return d3.determinant();
 }
 
        return d3.determinant();
 }
index 7f6c9535f918704d96f76fee51285d6a9955111c..0d86eeee4d673b48d20c210989cb3a498165b802 100644 (file)
@@ -1940,9 +1940,17 @@ matrix::matrix(unsigned r, unsigned c);
 creates a matrix with @samp{r} rows and @samp{c} columns with all elements
 set to zero.
 
 creates a matrix with @samp{r} rows and @samp{c} columns with all elements
 set to zero.
 
-The fastest way to create a matrix with preinitialized elements is to assign
-a list of comma-separated expressions to an empty matrix (see below for an
-example). But you can also specify the elements as a (flat) list with
+The easiest way to create a matrix is using an initializer list of
+initializer lists, all of the same size:
+
+@example
+@{
+    matrix m = @{@{1, -a@},
+                @{a,  1@}@};
+@}
+@end example
+
+You can also specify the elements as a (flat) list with
 
 @example
 matrix::matrix(unsigned r, unsigned c, const lst & l);
 
 @example
 matrix::matrix(unsigned r, unsigned c, const lst & l);
@@ -1965,6 +1973,7 @@ matrices:
 @cindex @code{symbolic_matrix()}
 @example
 ex diag_matrix(const lst & l);
 @cindex @code{symbolic_matrix()}
 @example
 ex diag_matrix(const lst & l);
+ex diag_matrix(initializer_list<ex> l);
 ex unit_matrix(unsigned x);
 ex unit_matrix(unsigned r, unsigned c);
 ex symbolic_matrix(unsigned r, unsigned c, const string & base_name);
 ex unit_matrix(unsigned x);
 ex unit_matrix(unsigned r, unsigned c);
 ex symbolic_matrix(unsigned r, unsigned c, const string & base_name);
@@ -1972,7 +1981,7 @@ ex symbolic_matrix(unsigned r, unsigned c, const string & base_name,
                    const string & tex_base_name);
 @end example
 
                    const string & tex_base_name);
 @end example
 
-@code{diag_matrix()} constructs a diagonal matrix given the list of diagonal
+@code{diag_matrix()} constructs a square diagonal matrix given the diagonal
 elements. @code{unit_matrix()} creates an @samp{x} by @samp{x} (or @samp{r}
 by @samp{c}) unit matrix. And finally, @code{symbolic_matrix} constructs a
 matrix filled with newly generated symbols made of the specified base name
 elements. @code{unit_matrix()} creates an @samp{x} by @samp{x} (or @samp{r}
 by @samp{c}) unit matrix. And finally, @code{symbolic_matrix} constructs a
 matrix filled with newly generated symbols made of the specified base name
@@ -1999,10 +2008,9 @@ that specify which row and column to remove:
 
 @example
 @{
 
 @example
 @{
-    matrix m(3,3);
-    m = 11, 12, 13,
-        21, 22, 23,
-        31, 32, 33;
+    matrix m = @{@{11, 12, 13@},
+                @{21, 22, 23@},
+                @{31, 32, 33@}@};
     cout << reduced_matrix(m, 1, 1) << endl;
     // -> [[11,13],[31,33]]
     cout << sub_matrix(m, 1, 2, 1, 2) << endl;
     cout << reduced_matrix(m, 1, 1) << endl;
     // -> [[11,13],[31,33]]
     cout << sub_matrix(m, 1, 2, 1, 2) << endl;
@@ -2028,9 +2036,8 @@ Here are a couple of examples for constructing matrices:
 @{
     symbol a("a"), b("b");
 
 @{
     symbol a("a"), b("b");
 
-    matrix M(2, 2);
-    M = a, 0,
-        0, b;
+    matrix M = @{@{a, 0@},
+                @{0, b@}@};
     cout << M << endl;
      // -> [[a,0],[0,b]]
 
     cout << M << endl;
      // -> [[a,0],[0,b]]
 
@@ -2082,13 +2089,12 @@ and @math{C}:
 
 @example
 @{
 
 @example
 @{
-    matrix A(2, 2), B(2, 2), C(2, 2);
-    A =  1, 2,
-         3, 4;
-    B = -1, 0,
-         2, 1;
-    C =  8, 4,
-         2, 1;
+    matrix A = @{@{ 1, 2@},
+                @{ 3, 4@}@};
+    matrix B = @{@{-1, 0@},
+                @{ 2, 1@}@};
+    matrix C = @{@{ 8, 4@},
+                @{ 2, 1@}@};
 
     matrix result = A.mul(B).sub(C.mul_scalar(2));
     cout << result << endl;
 
     matrix result = A.mul(B).sub(C.mul_scalar(2));
     cout << result << endl;
@@ -2928,10 +2934,9 @@ and scalar products):
     symbol x("x"), y("y");
 
     // A is a 2x2 matrix, X is a 2x1 vector
     symbol x("x"), y("y");
 
     // A is a 2x2 matrix, X is a 2x1 vector
-    matrix A(2, 2), X(2, 1);
-    A = 1, 2,
-        3, 4;
-    X = x, y;
+    matrix A = @{@{1, 2@},
+                @{3, 4@}@};
+    matrix X = @{@{x, y@}@};
 
     cout << indexed(A, i, i) << endl;
      // -> 5
 
     cout << indexed(A, i, i) << endl;
      // -> 5
@@ -3413,7 +3418,7 @@ The previous code may be rewritten with the help of @code{lst_to_clifford()} as
     ...
     idx i(symbol("i"), 4);
     realsymbol s("s");
     ...
     idx i(symbol("i"), 4);
     realsymbol s("s");
-    ex M = diag_matrix(lst@{1, -1, 0, s@});
+    ex M = diag_matrix(@{1, -1, 0, s@});
     ex e0 = lst_to_clifford(lst@{1, 0, 0, 0@}, i, M);
     ex e1 = lst_to_clifford(lst@{0, 1, 0, 0@}, i, M);
     ex e2 = lst_to_clifford(lst@{0, 0, 1, 0@}, i, M);
     ex e0 = lst_to_clifford(lst@{1, 0, 0, 0@}, i, M);
     ex e1 = lst_to_clifford(lst@{0, 1, 0, 0@}, i, M);
     ex e2 = lst_to_clifford(lst@{0, 0, 1, 0@}, i, M);
index 2217bae85a3f4e64816d98bb4e70945959eba1c3..50f4631546a5795216287f27c9b9f272731fefa6 100644 (file)
@@ -73,20 +73,6 @@ matrix::matrix(unsigned r, unsigned c) : row(r), col(c), m(r*c, _ex0)
        setflag(status_flags::not_shareable);
 }
 
        setflag(status_flags::not_shareable);
 }
 
-// protected
-
-/** Ctor from representation, for internal use only. */
-matrix::matrix(unsigned r, unsigned c, const exvector & m2)
-  : row(r), col(c), m(m2)
-{
-       setflag(status_flags::not_shareable);
-}
-matrix::matrix(unsigned r, unsigned c, exvector && m2)
-  : row(r), col(c), m(std::move(m2))
-{
-       setflag(status_flags::not_shareable);
-}
-
 /** Construct matrix from (flat) list of elements. If the list has fewer
  *  elements than the matrix, the remaining matrix elements are set to zero.
  *  If the list has more elements than the matrix, the excessive elements are
 /** Construct matrix from (flat) list of elements. If the list has fewer
  *  elements than the matrix, the remaining matrix elements are set to zero.
  *  If the list has more elements than the matrix, the excessive elements are
@@ -107,6 +93,40 @@ matrix::matrix(unsigned r, unsigned c, const lst & l)
        }
 }
 
        }
 }
 
+/** Construct a matrix from an 2 dimensional initializer list.
+ *  Throws an exception if some row has a different length than all the others.
+ */
+matrix::matrix(std::initializer_list<std::initializer_list<ex>> l)
+  : row(l.size()), col(l.begin()->size())
+{
+       setflag(status_flags::not_shareable);
+
+       m.reserve(row*col);
+       for (const auto & r : l) {
+               unsigned c = 0;
+               for (const auto & e : r) {
+                       m.push_back(e);
+                       ++c;
+               }
+               if (c != col)
+                       throw std::invalid_argument("matrix::matrix{{}}: wrong dimension");
+       }
+}
+
+// protected
+
+/** Ctor from representation, for internal use only. */
+matrix::matrix(unsigned r, unsigned c, const exvector & m2)
+  : row(r), col(c), m(m2)
+{
+       setflag(status_flags::not_shareable);
+}
+matrix::matrix(unsigned r, unsigned c, exvector && m2)
+  : row(r), col(c), m(std::move(m2))
+{
+       setflag(status_flags::not_shareable);
+}
+
 //////////
 // archiving
 //////////
 //////////
 // archiving
 //////////
@@ -1585,6 +1605,23 @@ ex diag_matrix(const lst & l)
        return M;
 }
 
        return M;
 }
 
+ex diag_matrix(std::initializer_list<ex> l)
+{
+       size_t dim = l.size();
+
+       // Allocate and fill matrix
+       matrix &M = *new matrix(dim, dim);
+       M.setflag(status_flags::dynallocated);
+
+       unsigned i = 0;
+       for (auto & it : l) {
+               M(i, i) = it;
+               ++i;
+       }
+
+       return M;
+}
+
 ex unit_matrix(unsigned r, unsigned c)
 {
        matrix &Id = *new matrix(r, c);
 ex unit_matrix(unsigned r, unsigned c)
 {
        matrix &Id = *new matrix(r, c);
index 6081e900257131db359b0f40fc89e0b9a9e6b19e..988e257b0ddce7726e325038fd1ede00366917d9 100644 (file)
@@ -26,6 +26,7 @@
 #include "basic.h"
 #include "ex.h"
 #include "archive.h"
 #include "basic.h"
 #include "ex.h"
 #include "archive.h"
+#include "compiler.h"
 
 #include <string>
 #include <vector>
 
 #include <string>
 #include <vector>
@@ -99,15 +100,9 @@ class matrix : public basic
 public:
        matrix(unsigned r, unsigned c);
        matrix(unsigned r, unsigned c, const lst & l);
 public:
        matrix(unsigned r, unsigned c);
        matrix(unsigned r, unsigned c, const lst & l);
+       matrix(std::initializer_list<std::initializer_list<ex>> l);
 
 
-       // First step of initialization of matrix with a comma-separated sequence
-       // of expressions. Subsequent steps are handled by matrix_init<>::operator,().
-       matrix_init<ex, exvector::iterator> operator=(const ex & x)
-       {
-               m[0] = x;
-               return matrix_init<ex, exvector::iterator>(++m.begin());
-       }
-
+       matrix_init<ex, exvector::iterator> operator=(const ex & x) deprecated;
 protected:
        matrix(unsigned r, unsigned c, const exvector & m2);
        matrix(unsigned r, unsigned c, exvector && m2);
 protected:
        matrix(unsigned r, unsigned c, const exvector & m2);
        matrix(unsigned r, unsigned c, exvector && m2);
@@ -179,6 +174,13 @@ protected:
 };
 GINAC_DECLARE_UNARCHIVER(matrix); 
 
 };
 GINAC_DECLARE_UNARCHIVER(matrix); 
 
+// First step of initialization of matrix with a comma-separated sequence
+// of expressions. Subsequent steps are handled by matrix_init<>::operator,().
+inline matrix_init<ex, exvector::iterator> matrix::operator=(const ex & x)
+{
+       m[0] = x;
+       return matrix_init<ex, exvector::iterator>(++m.begin());
+}
 
 // wrapper functions around member functions
 
 
 // wrapper functions around member functions
 
@@ -225,6 +227,7 @@ extern ex lst_to_matrix(const lst & l);
 
 /** Convert list of diagonal elements to matrix. */
 extern ex diag_matrix(const lst & l);
 
 /** Convert list of diagonal elements to matrix. */
 extern ex diag_matrix(const lst & l);
+extern ex diag_matrix(std::initializer_list<ex> l);
 
 /** Create an r times c unit matrix. */
 extern ex unit_matrix(unsigned r, unsigned c);
 
 /** Create an r times c unit matrix. */
 extern ex unit_matrix(unsigned r, unsigned c);