]> www.ginac.de Git - ginac.git/blobdiff - ginac/matrix.cpp
- Complete revamp of methods in class matrix. Some redundant (and poor)
[ginac.git] / ginac / matrix.cpp
index feae32fee687341c7a92a5d8ea4e6693dc61dc88..6240ce20e9e36a14557e38a9918841ecf63923c7 100644 (file)
@@ -108,7 +108,7 @@ matrix::matrix(unsigned r, unsigned c)
     m.resize(r*c, _ex0());
 }
 
- // protected
+// protected
 
 /** Ctor from representation, for internal use only. */
 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
@@ -152,7 +152,7 @@ void matrix::archive(archive_node &n) const
     exvector::const_iterator i = m.begin(), iend = m.end();
     while (i != iend) {
         n.add_ex("m", *i);
-        i++;
+        ++i;
     }
 }
 
@@ -174,15 +174,13 @@ void matrix::print(std::ostream & os, unsigned upper_precedence) const
     os << "[[ ";
     for (unsigned r=0; r<row-1; ++r) {
         os << "[[";
-        for (unsigned c=0; c<col-1; ++c) {
+        for (unsigned c=0; c<col-1; ++c)
             os << m[r*col+c] << ",";
-        }
         os << m[col*(r+1)-1] << "]], ";
     }
     os << "[[";
-    for (unsigned c=0; c<col-1; ++c) {
+    for (unsigned c=0; c<col-1; ++c)
         os << m[(row-1)*col+c] << ",";
-    }
     os << m[row*col-1] << "]] ]]";
 }
 
@@ -192,15 +190,13 @@ void matrix::printraw(std::ostream & os) const
     os << "matrix(" << row << "," << col <<",";
     for (unsigned r=0; r<row-1; ++r) {
         os << "(";
-        for (unsigned c=0; c<col-1; ++c) {
+        for (unsigned c=0; c<col-1; ++c)
             os << m[r*col+c] << ",";
-        }
         os << m[col*(r-1)-1] << "),";
     }
     os << "(";
-    for (unsigned c=0; c<col-1; ++c) {
+    for (unsigned c=0; c<col-1; ++c)
         os << m[(row-1)*col+c] << ",";
-    }
     os << m[row*col-1] << "))";
 }
 
@@ -219,6 +215,9 @@ ex matrix::op(int i) const
 /** returns matrix entry at position (i/col, i%col). */
 ex & matrix::let_op(int i)
 {
+    GINAC_ASSERT(i>=0);
+    GINAC_ASSERT(i<nops());
+    
     return m[i];
 }
 
@@ -226,9 +225,9 @@ ex & matrix::let_op(int i)
 ex matrix::expand(unsigned options) const
 {
     exvector tmp(row*col);
-    for (unsigned i=0; i<row*col; ++i) {
-        tmp[i]=m[i].expand(options);
-    }
+    for (unsigned i=0; i<row*col; ++i)
+        tmp[i] = m[i].expand(options);
+    
     return matrix(row, col, tmp);
 }
 
@@ -242,9 +241,9 @@ bool matrix::has(const ex & other) const
     if (is_equal(*other.bp)) return true;
     
     // search all the elements
-    for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r) {
+    for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r)
         if ((*r).has(other)) return true;
-    }
+    
     return false;
 }
 
@@ -263,12 +262,10 @@ ex matrix::eval(int level) const
     
     // eval() entry by entry
     exvector m2(row*col);
-    --level;    
-    for (unsigned r=0; r<row; ++r) {
-        for (unsigned c=0; c<col; ++c) {
+    --level;
+    for (unsigned r=0; r<row; ++r)
+        for (unsigned c=0; c<col; ++c)
             m2[r*col+c] = m[r*col+c].eval(level);
-        }
-    }
     
     return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
                                                status_flags::evaluated );
@@ -291,11 +288,10 @@ ex matrix::evalf(int level) const
     // evalf() entry by entry
     exvector m2(row*col);
     --level;
-    for (unsigned r=0; r<row; ++r) {
-        for (unsigned c=0; c<col; ++c) {
+    for (unsigned r=0; r<row; ++r)
+        for (unsigned c=0; c<col; ++c)
             m2[r*col+c] = m[r*col+c].evalf(level);
-        }
-    }
+    
     return matrix(row, col, m2);
 }
 
@@ -343,11 +339,9 @@ matrix matrix::add(const matrix & other) const
     exvector sum(this->m);
     exvector::iterator i;
     exvector::const_iterator ci;
-    for (i=sum.begin(), ci=other.m.begin();
-         i!=sum.end();
-         ++i, ++ci) {
+    for (i=sum.begin(), ci=other.m.begin(); i!=sum.end(); ++i, ++ci)
         (*i) += (*ci);
-    }
+    
     return matrix(row,col,sum);
 }
 
@@ -363,11 +357,9 @@ matrix matrix::sub(const matrix & other) const
     exvector dif(this->m);
     exvector::iterator i;
     exvector::const_iterator ci;
-    for (i=dif.begin(), ci=other.m.begin();
-         i!=dif.end();
-         ++i, ++ci) {
+    for (i=dif.begin(), ci=other.m.begin(); i!=dif.end(); ++i, ++ci)
         (*i) -= (*ci);
-    }
+    
     return matrix(row,col,dif);
 }
 
@@ -377,17 +369,17 @@ matrix matrix::sub(const matrix & other) const
  *  @exception logic_error (incompatible matrices) */
 matrix matrix::mul(const matrix & other) const
 {
-    if (col != other.row)
+    if (this->cols() != other.rows())
         throw (std::logic_error("matrix::mul(): incompatible matrices"));
     
-    exvector prod(row*other.col);
+    exvector prod(this->rows()*other.cols());
     
-    for (unsigned r1=0; r1<rows(); ++r1) {
-        for (unsigned c=0; c<cols(); ++c) {
+    for (unsigned r1=0; r1<this->rows(); ++r1) {
+        for (unsigned c=0; c<this->cols(); ++c) {
             if (m[r1*col+c].is_zero())
                 continue;
-            for (unsigned r2=0; r2<other.col; ++r2)
-                prod[r1*other.col+r2] += m[r1*col+c] * other.m[c*other.col+r2];
+            for (unsigned r2=0; r2<other.cols(); ++r2)
+                prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]).expand();
         }
     }
     return matrix(row, other.col, prod);
@@ -426,19 +418,19 @@ matrix & matrix::set(unsigned ro, unsigned co, ex value)
  *  represents the transposed. */
 matrix matrix::transpose(void) const
 {
-    exvector trans(col*row);
+    exvector trans(this->cols()*this->rows());
     
-    for (unsigned r=0; r<col; ++r)
-        for (unsigned c=0; c<row; ++c)
-            trans[r*row+c] = m[c*col+r];
+    for (unsigned r=0; r<this->cols(); ++r)
+        for (unsigned c=0; c<this->rows(); ++c)
+            trans[r*this->rows()+c] = m[c*this->cols()+r];
     
-    return matrix(col,row,trans);
+    return matrix(this->cols(),this->rows(),trans);
 }
 
 
 /** Determinant of square matrix.  This routine doesn't actually calculate the
  *  determinant, it only implements some heuristics about which algorithm to
- *  call.  If all the elements of the matrix are elements of an integral domain
+ *  run.  If all the elements of the matrix are elements of an integral domain
  *  the determinant is also in that integral domain and the result is expanded
  *  only.  If one or more elements are from a quotient field the determinant is
  *  usually also in that quotient field and the result is normalized before it
@@ -446,85 +438,127 @@ matrix matrix::transpose(void) const
  *  [[a/(a-b),1],[b/(a-b),1]] is returned as unity.  (In this respect, it
  *  behaves like MapleV and unlike Mathematica.)
  *
+ *  @param     algo allows to chose an algorithm
  *  @return    the determinant as a new expression
- *  @exception logic_error (matrix not square) */
-ex matrix::determinant(void) const
+ *  @exception logic_error (matrix not square)
+ *  @see       determinant_algo */
+ex matrix::determinant(unsigned algo) const
 {
     if (row!=col)
         throw (std::logic_error("matrix::determinant(): matrix not square"));
     GINAC_ASSERT(row*col==m.capacity());
-    if (this->row==1)  // continuation would be pointless
-        return m[0];
     
-    // Gather some information about the matrix:
+    // Gather some statistical information about this matrix:
     bool numeric_flag = true;
     bool normal_flag = false;
-    unsigned sparse_count = 0;  // count non-zero elements
+    unsigned sparse_count = 0;  // counts non-zero elements
     for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r) {
-        if (!(*r).is_zero())
+        lst srl;  // symbol replacement list
+        ex rtest = (*r).to_rational(srl);
+        if (!rtest.is_zero())
             ++sparse_count;
-        if (!(*r).info(info_flags::numeric))
+        if (!rtest.info(info_flags::numeric))
             numeric_flag = false;
-        if ((*r).info(info_flags::rational_function) &&
-            !(*r).info(info_flags::crational_polynomial))
+        if (!rtest.info(info_flags::crational_polynomial) &&
+             rtest.info(info_flags::rational_function))
             normal_flag = true;
     }
     
-    // Purely numeric matrix handled by Gauss elimination
-    if (numeric_flag) {
-        ex det = 1;
-        matrix tmp(*this);
-        int sign = tmp.gauss_elimination();
-        for (int d=0; d<row; ++d)
-            det *= tmp.m[d*col+d];
-        return (sign*det);
+    // Here is the heuristics in case this routine has to decide:
+    if (algo == determinant_algo::automatic) {
+        // Minor expansion is generally a good guess:
+        algo = determinant_algo::laplace;
+        // Does anybody know when a matrix is really sparse?
+        // Maybe <~row/2.236 nonzero elements average in a row?
+        if (row>3 && 5*sparse_count<=row*col)
+            algo = determinant_algo::bareiss;
+        // Purely numeric matrix can be handled by Gauss elimination.
+        // This overrides any prior decisions.
+        if (numeric_flag)
+            algo = determinant_algo::gauss;
     }
     
-    // Does anybody know when a matrix is really sparse?
-    // Maybe <~row/2.2 nonzero elements average in a row?
-    if (5*sparse_count<=row*col) {
-        // copy *this:
-        matrix tmp(*this);
-        int sign;
-        sign = tmp.fraction_free_elimination(true);
+    // Trap the trivial case here, since some algorithms don't like it
+    if (this->row==1) {
+        // for consistency with non-trivial determinants...
         if (normal_flag)
-            return (sign*tmp.m[row*col-1]).normal();
+            return m[0].normal();
         else
-            return (sign*tmp.m[row*col-1]).expand();
+            return m[0].expand();
     }
     
-    // Now come the minor expansion schemes.  We always develop such that the
-    // smallest minors (i.e, the trivial 1x1 ones) are on the rightmost column.
-    // For this to be efficient it turns out that the emptiest columns (i.e.
-    // the ones with most zeros) should be the ones on the right hand side.
-    // Therefore we presort the columns of the matrix:
-    typedef std::pair<unsigned,unsigned> uintpair;  // # of zeros, column
-    std::vector<uintpair> c_zeros;  // number of zeros in column
-    for (unsigned c=0; c<col; ++c) {
-        unsigned acc = 0;
-        for (unsigned r=0; r<row; ++r)
-            if (m[r*col+c].is_zero())
-                ++acc;
-        c_zeros.push_back(uintpair(acc,c));
-    }
-    sort(c_zeros.begin(),c_zeros.end());
-    // unfortunately std::vector<uintpair> can't be used for permutation_sign.
-    std::vector<unsigned> pre_sort;
-    for (std::vector<uintpair>::iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
-        pre_sort.push_back(i->second);
-    int sign = permutation_sign(pre_sort);
-    exvector result(row*col);  // represents sorted matrix
-    unsigned c = 0;
-    for (std::vector<unsigned>::iterator i=pre_sort.begin();
-         i!=pre_sort.end();
-         ++i,++c) {
-        for (unsigned r=0; r<row; ++r)
-            result[r*col+c] = m[r*col+(*i)];
+    // Compute the determinant
+    switch(algo) {
+        case determinant_algo::gauss: {
+            ex det = 1;
+            matrix tmp(*this);
+            int sign = tmp.gauss_elimination(true);
+            for (unsigned d=0; d<row; ++d)
+                det *= tmp.m[d*col+d];
+            if (normal_flag)
+                return (sign*det).normal();
+            else
+                return (sign*det).normal().expand();
+        }
+        case determinant_algo::bareiss: {
+            matrix tmp(*this);
+            int sign;
+            sign = tmp.fraction_free_elimination(true);
+            if (normal_flag)
+                return (sign*tmp.m[row*col-1]).normal();
+            else
+                return (sign*tmp.m[row*col-1]).expand();
+        }
+        case determinant_algo::divfree: {
+            matrix tmp(*this);
+            int sign;
+            sign = tmp.division_free_elimination(true);
+            if (sign==0)
+                return _ex0();
+            ex det = tmp.m[row*col-1];
+            // factor out accumulated bogus slag
+            for (unsigned d=0; d<row-2; ++d)
+                for (unsigned j=0; j<row-d-2; ++j)
+                    det = (det/tmp.m[d*col+d]).normal();
+            return (sign*det);
+        }
+        case determinant_algo::laplace:
+        default: {
+            // This is the minor expansion scheme.  We always develop such
+            // that the smallest minors (i.e, the trivial 1x1 ones) are on the
+            // rightmost column.  For this to be efficient it turns out that
+            // the emptiest columns (i.e. the ones with most zeros) should be
+            // the ones on the right hand side.  Therefore we presort the
+            // columns of the matrix:
+            typedef std::pair<unsigned,unsigned> uintpair;
+            std::vector<uintpair> c_zeros;  // number of zeros in column
+            for (unsigned c=0; c<col; ++c) {
+                unsigned acc = 0;
+                for (unsigned r=0; r<row; ++r)
+                    if (m[r*col+c].is_zero())
+                        ++acc;
+                c_zeros.push_back(uintpair(acc,c));
+            }
+            sort(c_zeros.begin(),c_zeros.end());
+            std::vector<unsigned> pre_sort;
+            for (std::vector<uintpair>::iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
+                pre_sort.push_back(i->second);
+            int sign = permutation_sign(pre_sort);
+            exvector result(row*col);  // represents sorted matrix
+            unsigned c = 0;
+            for (std::vector<unsigned>::iterator i=pre_sort.begin();
+                 i!=pre_sort.end();
+                 ++i,++c) {
+                for (unsigned r=0; r<row; ++r)
+                    result[r*col+c] = m[r*col+(*i)];
+            }
+            
+            if (normal_flag)
+                return (sign*matrix(row,col,result).determinant_minor()).normal();
+            else
+                return sign*matrix(row,col,result).determinant_minor();
+        }
     }
-    
-    if (normal_flag)
-        return sign*matrix(row,col,result).determinant_minor().normal();
-    return sign*matrix(row,col,result).determinant_minor();
 }
 
 
@@ -538,7 +572,6 @@ ex matrix::trace(void) const
 {
     if (row != col)
         throw (std::logic_error("matrix::trace(): matrix not square"));
-    GINAC_ASSERT(row*col==m.capacity());
     
     ex tr;
     for (unsigned r=0; r<col; ++r)
@@ -613,15 +646,19 @@ matrix matrix::inverse(void) const
     if (row != col)
         throw (std::logic_error("matrix::inverse(): matrix not square"));
     
+    // NOTE: the Gauss-Jordan elimination used here can in principle be
+    // replaced this by two clever calls to gauss_elimination() and some to
+    // transpose().  Wouldn't be more efficient (maybe less?), just more
+    // orthogonal.
     matrix tmp(row,col);
     // set tmp to the unit matrix
     for (unsigned i=0; i<col; ++i)
         tmp.m[i*col+i] = _ex1();
-
+    
     // create a copy of this matrix
     matrix cpy(*this);
     for (unsigned r1=0; r1<row; ++r1) {
-        int indx = cpy.pivot(r1);
+        int indx = cpy.pivot(r1, r1);
         if (indx == -1) {
             throw (std::runtime_error("matrix::inverse(): singular matrix"));
         }
@@ -639,212 +676,121 @@ matrix matrix::inverse(void) const
                 ex a2 = cpy.m[r2*col+r1];
                 for (unsigned c=0; c<col; ++c) {
                     cpy.m[r2*col+c] -= a2 * cpy.m[r1*col+c];
+                    if (!cpy.m[r2*col+c].info(info_flags::numeric))
+                        cpy.m[r2*col+c] = cpy.m[r2*col+c].normal();
                     tmp.m[r2*col+c] -= a2 * tmp.m[r1*col+c];
+                    if (!tmp.m[r2*col+c].info(info_flags::numeric))
+                        tmp.m[r2*col+c] = tmp.m[r2*col+c].normal();
                 }
             }
         }
     }
-    return tmp;
-}
-
-// superfluous helper function, to be removed:
-void matrix::swap(unsigned r1, unsigned c1, unsigned r2 ,unsigned c2)
-{
-    ensure_if_modifiable();
     
-    ex tmp = (*this)(r1,c1);
-    set(r1,c1,(*this)(r2,c2));
-    set(r2,c2,tmp);
+    return tmp;
 }
 
 
-/** Solve a set of equations for an m x n matrix by fraction-free Gaussian
- *  elimination.  Based on algorithm 9.1 from 'Algorithms for Computer Algebra'
- *  by Keith O. Geddes et al.
+/** Solve a linear system consisting of a m x n matrix and a m x p right hand
+ *  side by applying an elimination scheme to the augmented matrix.
  *
- *  @param vars n x p matrix
+ *  @param vars n x p matrix, all elements must be symbols 
  *  @param rhs m x p matrix
+ *  @return n x p solution matrix
  *  @exception logic_error (incompatible matrices)
- *  @exception runtime_error (singular matrix) */
-matrix matrix::fraction_free_elim(const matrix & vars,
-                                  const matrix & rhs) const
+ *  @exception invalid_argument (1st argument must be matrix of symbols)
+ *  @exception runtime_error (inconsistent linear system)
+ *  @see       solve_algo */
+matrix matrix::solve(const matrix & vars,
+                     const matrix & rhs,
+                     unsigned algo) const
 {
-    // FIXME: use implementation of matrix::fraction_free_elimination
-    if ((row != rhs.row) || (col != vars.row) || (rhs.col != vars.col))
-        throw (std::logic_error("matrix::fraction_free_elim(): incompatible matrices"));
-    
-    matrix a(*this);  // make a copy of the matrix
-    matrix b(rhs);    // make a copy of the rhs vector
+    const unsigned m = this->rows();
+    const unsigned n = this->cols();
+    const unsigned p = rhs.cols();
     
-    // given an m x n matrix a, reduce it to upper echelon form
-    unsigned m = a.row;
-    unsigned n = a.col;
-    int sign = 1;
-    ex divisor = 1;
-    unsigned r = 0;
+    // syntax checks    
+    if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
+        throw (std::logic_error("matrix::solve(): incompatible matrices"));
+    for (unsigned ro=0; ro<n; ++ro)
+        for (unsigned co=0; co<p; ++co)
+            if (!vars(ro,co).info(info_flags::symbol))
+                throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
     
-    // eliminate below row r, with pivot in column k
-    for (unsigned k=0; (k<n)&&(r<m); ++k) {
-        // find a nonzero pivot
-        unsigned p;
-        for (p=r; (p<m)&&(a(p,k).is_zero()); ++p) {}
-        // pivot is in row p
-        if (p<m) {
-            if (p!=r) {
-                // swap rows p and r
-                for (unsigned j=k; j<n; ++j)
-                    a.swap(p,j,r,j);
-                b.swap(p,0,r,0);
-                // keep track of sign changes due to row exchange
-                sign *= -1;
-            }
-            for (unsigned i=r+1; i<m; ++i) {
-                for (unsigned j=k+1; j<n; ++j) {
-                    a.set(i,j,(a(r,k)*a(i,j)
-                              -a(r,j)*a(i,k))/divisor);
-                    a.set(i,j,a(i,j).normal());
-                }
-                b.set(i,0,(a(r,k)*b(i,0)
-                          -b(r,0)*a(i,k))/divisor);
-                b.set(i,0,b(i,0).normal());
-                a.set(i,k,0);
-            }
-            divisor = a(r,k);
-            ++r;
-        }
+    // build the augmented matrix of *this with rhs attached to the right
+    matrix aug(m,n+p);
+    for (unsigned r=0; r<m; ++r) {
+        for (unsigned c=0; c<n; ++c)
+            aug.m[r*(n+p)+c] = this->m[r*n+c];
+        for (unsigned c=0; c<p; ++c)
+            aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
     }
     
-#ifdef DO_GINAC_ASSERT
-    // test if we really have an upper echelon matrix
-    int zero_in_last_row = -1;
-    for (unsigned r=0; r<m; ++r) {
-        int zero_in_this_row=0;
-        for (unsigned c=0; c<n; ++c) {
-            if (a(r,c).is_zero())
-               zero_in_this_row++;
-            else
-                break;
-        }
-        GINAC_ASSERT((zero_in_this_row>zero_in_last_row)||(zero_in_this_row=n));
-        zero_in_last_row = zero_in_this_row;
+    // Gather some statistical information about the augmented matrix:
+    bool numeric_flag = true;
+    for (exvector::const_iterator r=aug.m.begin(); r!=aug.m.end(); ++r) {
+        if (!(*r).info(info_flags::numeric))
+            numeric_flag = false;
     }
-#endif // def DO_GINAC_ASSERT
     
-    // assemble solution
-    matrix sol(n,1);
-    unsigned last_assigned_sol = n+1;
-    for (int r=m-1; r>=0; --r) {
-        unsigned first_non_zero = 1;
-        while ((first_non_zero<=n)&&(a(r,first_non_zero-1).is_zero()))
-            first_non_zero++;
-        if (first_non_zero>n) {
-            // row consists only of zeroes, corresponding rhs must be 0 as well
-            if (!b(r,0).is_zero()) {
-                throw (std::runtime_error("matrix::fraction_free_elim(): singular matrix"));
-            }
-        } else {
-            // assign solutions for vars between first_non_zero+1 and
-            // last_assigned_sol-1: free parameters
-            for (unsigned c=first_non_zero; c<last_assigned_sol-1; ++c)
-                sol.set(c,0,vars(c,0));
-            ex e = b(r,0);
-            for (unsigned c=first_non_zero; c<n; ++c)
-                e -= a(r,c)*sol(c,0);
-            sol.set(first_non_zero-1,0,
-                    (e/a(r,first_non_zero-1)).normal());
-            last_assigned_sol = first_non_zero;
-        }
+    // Here is the heuristics in case this routine has to decide:
+    if (algo == solve_algo::automatic) {
+        // Bareiss (fraction-free) elimination is generally a good guess:
+        algo = solve_algo::bareiss;
+        // For m<3, Bareiss elimination is equivalent to division free
+        // elimination but has more logistic overhead
+        if (m<3)
+            algo = solve_algo::divfree;
+        // This overrides any prior decisions.
+        if (numeric_flag)
+            algo = solve_algo::gauss;
     }
-    // assign solutions for vars between 1 and
-    // last_assigned_sol-1: free parameters
-    for (unsigned c=0; c<last_assigned_sol-1; ++c)
-        sol.set(c,0,vars(c,0));
     
-#ifdef DO_GINAC_ASSERT
-    // test solution with echelon matrix
-    for (unsigned r=0; r<m; ++r) {
-        ex e = 0;
-        for (unsigned c=0; c<n; ++c)
-            e += a(r,c)*sol(c,0);
-        if (!(e-b(r,0)).normal().is_zero()) {
-            cout << "e=" << e;
-            cout << "b(" << r <<",0)=" << b(r,0) << endl;
-            cout << "diff=" << (e-b(r,0)).normal() << endl;
-        }
-        GINAC_ASSERT((e-b(r,0)).normal().is_zero());
+    // Eliminate the augmented matrix:
+    switch(algo) {
+        case solve_algo::gauss:
+            aug.gauss_elimination();
+        case solve_algo::divfree:
+            aug.division_free_elimination();
+        case solve_algo::bareiss:
+        default:
+            aug.fraction_free_elimination();
     }
     
-    // test solution with original matrix
-    for (unsigned r=0; r<m; ++r) {
-        ex e = 0;
-        for (unsigned c=0; c<n; ++c)
-            e += this->m[r*cols()+c]*sol(c,0);
-        try {
-            if (!(e-rhs(r,0)).normal().is_zero()) {
-                cout << "e==" << e << endl;
-                e.printtree(cout);
-                ex en = e.normal();
-                cout << "e.normal()=" << en << endl;
-                en.printtree(cout);
-                cout << "rhs(" << r <<",0)=" << rhs(r,0) << endl;
-                cout << "diff=" << (e-rhs(r,0)).normal() << endl;
+    // assemble the solution matrix:
+    matrix sol(n,p);
+    for (unsigned co=0; co<p; ++co) {
+        unsigned last_assigned_sol = n+1;
+        for (int r=m-1; r>=0; --r) {
+            unsigned fnz = 1;    // first non-zero in row
+            while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
+                ++fnz;
+            if (fnz>n) {
+                // row consists only of zeros, corresponding rhs must be 0, too
+                if (!aug.m[r*(n+p)+n+co].is_zero()) {
+                    throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
+                }
+            } else {
+                // assign solutions for vars between fnz+1 and
+                // last_assigned_sol-1: free parameters
+                for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
+                    sol.set(c,co,vars.m[c*p+co]);
+                ex e = aug.m[r*(n+p)+n+co];
+                for (unsigned c=fnz; c<n; ++c)
+                    e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
+                sol.set(fnz-1,co,
+                        (e/(aug.m[r*(n+p)+(fnz-1)])).normal());
+                last_assigned_sol = fnz;
             }
-        } catch (...) {
-            ex xxx = e - rhs(r,0);
-            cerr << "xxx=" << xxx << endl << endl;
         }
-        GINAC_ASSERT((e-rhs(r,0)).normal().is_zero());
+        // assign solutions for vars between 1 and
+        // last_assigned_sol-1: free parameters
+        for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
+            sol.set(ro,co,vars(ro,co));
     }
-#endif // def DO_GINAC_ASSERT
     
     return sol;
 }
 
-/** Solve a set of equations for an m x n matrix.
- *
- *  @param vars n x p matrix
- *  @param rhs m x p matrix
- *  @exception logic_error (incompatible matrices)
- *  @exception runtime_error (singular matrix) */
-matrix matrix::solve(const matrix & vars,
-                     const matrix & rhs) const
-{
-    if ((row != rhs.row) || (col != vars.row) || (rhs.col != vars.col))
-        throw (std::logic_error("matrix::solve(): incompatible matrices"));
-    
-    throw (std::runtime_error("FIXME: need implementation."));
-}
-
-/** Old and obsolete interface: */
-matrix matrix::old_solve(const matrix & v) const
-{
-    if ((v.row != col) || (col != v.row))
-        throw (std::logic_error("matrix::solve(): incompatible matrices"));
-    
-    // build the augmented matrix of *this with v attached to the right
-    matrix tmp(row,col+v.col);
-    for (unsigned r=0; r<row; ++r) {
-        for (unsigned c=0; c<col; ++c)
-            tmp.m[r*tmp.col+c] = this->m[r*col+c];
-        for (unsigned c=0; c<v.col; ++c)
-            tmp.m[r*tmp.col+c+col] = v.m[r*v.col+c];
-    }
-    // cout << "augmented: " << tmp << endl;
-    tmp.gauss_elimination();
-    // cout << "degaussed: " << tmp << endl;
-    // assemble the solution matrix
-    exvector sol(v.row*v.col);
-    for (unsigned c=0; c<v.col; ++c) {
-        for (unsigned r=row; r>0; --r) {
-            for (unsigned i=r; i<col; ++i)
-                sol[(r-1)*v.col+c] -= tmp.m[(r-1)*tmp.col+i]*sol[i*v.col+c];
-            sol[(r-1)*v.col+c] += tmp.m[(r-1)*tmp.col+col+c];
-            sol[(r-1)*v.col+c] = (sol[(r-1)*v.col+c]/tmp.m[(r-1)*tmp.col+(r-1)]).normal();
-        }
-    }
-    return matrix(v.row, v.col, sol);
-}
-
 
 // protected
 
@@ -861,11 +807,12 @@ matrix matrix::old_solve(const matrix & v) const
 ex matrix::determinant_minor(void) const
 {
     // for small matrices the algorithm does not make any sense:
-    if (this->row==1)
-        return m[0];
-    if (this->row==2)
+    const unsigned n = this->cols();
+    if (n==1)
+        return m[0].expand();
+    if (n==2)
         return (m[0]*m[3]-m[2]*m[1]).expand();
-    if (this->row==3)
+    if (n==3)
         return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
                 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
                 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
@@ -873,8 +820,8 @@ ex matrix::determinant_minor(void) const
     // This algorithm can best be understood by looking at a naive
     // implementation of Laplace-expansion, like this one:
     // ex det;
-    // matrix minorM(this->row-1,this->col-1);
-    // for (unsigned r1=0; r1<this->row; ++r1) {
+    // matrix minorM(this->rows()-1,this->cols()-1);
+    // for (unsigned r1=0; r1<this->rows(); ++r1) {
     //     // shortcut if element(r1,0) vanishes
     //     if (m[r1*col].is_zero())
     //         continue;
@@ -904,10 +851,10 @@ ex matrix::determinant_minor(void) const
     
     // Unique flipper counter for partitioning into minors
     std::vector<unsigned> Pkey;
-    Pkey.reserve(this->col);
+    Pkey.reserve(n);
     // key for minor determinant (a subpartition of Pkey)
     std::vector<unsigned> Mkey;
-    Mkey.reserve(this->col-1);
+    Mkey.reserve(n-1);
     // we store our subminors in maps, keys being the rows they arise from
     typedef std::map<std::vector<unsigned>,class ex> Rmap;
     typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
@@ -915,34 +862,34 @@ ex matrix::determinant_minor(void) const
     Rmap B;
     ex det;
     // initialize A with last column:
-    for (unsigned r=0; r<this->col; ++r) {
+    for (unsigned r=0; r<n; ++r) {
         Pkey.erase(Pkey.begin(),Pkey.end());
         Pkey.push_back(r);
-        A.insert(Rmap_value(Pkey,m[this->col*r+this->col-1]));
+        A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
     }
     // proceed from right to left through matrix
-    for (int c=this->col-2; c>=0; --c) {
+    for (int c=n-2; c>=0; --c) {
         Pkey.erase(Pkey.begin(),Pkey.end());  // don't change capacity
         Mkey.erase(Mkey.begin(),Mkey.end());
-        for (unsigned i=0; i<this->col-c; ++i)
+        for (unsigned i=0; i<n-c; ++i)
             Pkey.push_back(i);
         unsigned fc = 0;  // controls logic for our strange flipper counter
         do {
             det = _ex0();
-            for (unsigned r=0; r<this->col-c; ++r) {
+            for (unsigned r=0; r<n-c; ++r) {
                 // maybe there is nothing to do?
-                if (m[Pkey[r]*this->col+c].is_zero())
+                if (m[Pkey[r]*n+c].is_zero())
                     continue;
                 // create the sorted key for all possible minors
                 Mkey.erase(Mkey.begin(),Mkey.end());
-                for (unsigned i=0; i<this->col-c; ++i)
+                for (unsigned i=0; i<n-c; ++i)
                     if (i!=r)
                         Mkey.push_back(Pkey[i]);
                 // Fetch the minors and compute the new determinant
                 if (r%2)
-                    det -= m[Pkey[r]*this->col+c]*A[Mkey];
+                    det -= m[Pkey[r]*n+c]*A[Mkey];
                 else
-                    det += m[Pkey[r]*this->col+c]*A[Mkey];
+                    det += m[Pkey[r]*n+c]*A[Mkey];
             }
             // prevent build-up of deep nesting of expressions saves time:
             det = det.expand();
@@ -950,13 +897,13 @@ ex matrix::determinant_minor(void) const
             if (!det.is_zero())
                 B.insert(Rmap_value(Pkey,det));
             // increment our strange flipper counter
-            for (fc=this->col-c; fc>0; --fc) {
+            for (fc=n-c; fc>0; --fc) {
                 ++Pkey[fc-1];
                 if (Pkey[fc-1]<fc+c)
                     break;
             }
-            if (fc<this->col-c)
-                for (unsigned j=fc; j<this->col-c; ++j)
+            if (fc<n-c)
+                for (unsigned j=fc; j<n-c; ++j)
                     Pkey[j] = Pkey[j-1]+1;
         } while(fc);
         // next column, so change the role of A and B:
@@ -968,28 +915,51 @@ ex matrix::determinant_minor(void) const
 }
 
 
-/** Perform the steps of an ordinary Gaussian elimination to bring the matrix
- *  into an upper echelon form.
+/** Perform the steps of an ordinary Gaussian elimination to bring the m x n
+ *  matrix into an upper echelon form.  The algorithm is ok for matrices
+ *  with numeric coefficients but quite unsuited for symbolic matrices.
  *
+ *  @param det may be set to true to save a lot of space if one is only
+ *  interested in the diagonal elements (i.e. for calculating determinants).
+ *  The others are set to zero in this case.
  *  @return sign is 1 if an even number of rows was swapped, -1 if an odd
  *  number of rows was swapped and 0 if the matrix is singular. */
-int matrix::gauss_elimination(void)
+int matrix::gauss_elimination(const bool det)
 {
     ensure_if_modifiable();
+    const unsigned m = this->rows();
+    const unsigned n = this->cols();
+    GINAC_ASSERT(!det || n==m);
     int sign = 1;
-    ex piv;
-    for (unsigned r1=0; r1<row-1; ++r1) {
-        int indx = pivot(r1);
-        if (indx == -1)
-            return 0;  // Note: leaves *this in a messy state.
-        if (indx > 0)
-            sign = -sign;
-        for (unsigned r2=r1+1; r2<row; ++r2) {
-            piv = this->m[r2*col+r1] / this->m[r1*col+r1];
-            for (unsigned c=r1+1; c<col; ++c)
-                this->m[r2*col+c] -= piv * this->m[r1*col+c];
-            for (unsigned c=0; c<=r1; ++c)
-                this->m[r2*col+c] = _ex0();
+    
+    unsigned r0 = 0;
+    for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
+        int indx = pivot(r0, r1, true);
+        if (indx == -1) {
+            sign = 0;
+            if (det)
+                return 0;  // leaves *this in a messy state
+        }
+        if (indx>=0) {
+            if (indx > 0)
+                sign = -sign;
+            for (unsigned r2=r0+1; r2<m; ++r2) {
+                ex piv = this->m[r2*n+r1] / this->m[r0*n+r1];
+                for (unsigned c=r1+1; c<n; ++c) {
+                    this->m[r2*n+c] -= piv * this->m[r0*n+c];
+                    if (!this->m[r2*n+c].info(info_flags::numeric))
+                        this->m[r2*n+c] = this->m[r2*n+c].normal();
+                }
+                // fill up left hand side with zeros
+                for (unsigned c=0; c<=r1; ++c)
+                    this->m[r2*n+c] = _ex0();
+            }
+            if (det) {
+                // save space by deleting no longer needed elements
+                for (unsigned c=r0+1; c<n; ++c)
+                    this->m[r0*n+c] = _ex0();
+            }
+            ++r0;
         }
     }
     
@@ -997,26 +967,46 @@ int matrix::gauss_elimination(void)
 }
 
 
-/** Perform the steps of division free elimination to bring the matrix
+/** Perform the steps of division free elimination to bring the m x n matrix
  *  into an upper echelon form.
  *
+ *  @param det may be set to true to save a lot of space if one is only
+ *  interested in the diagonal elements (i.e. for calculating determinants).
+ *  The others are set to zero in this case.
  *  @return sign is 1 if an even number of rows was swapped, -1 if an odd
  *  number of rows was swapped and 0 if the matrix is singular. */
-int matrix::division_free_elimination(void)
+int matrix::division_free_elimination(const bool det)
 {
-    int sign = 1;
     ensure_if_modifiable();
-    for (unsigned r1=0; r1<row-1; ++r1) {
-        int indx = pivot(r1);
-        if (indx==-1)
-            return 0;  // Note: leaves *this in a messy state.
-        if (indx>0)
-            sign = -sign;
-        for (unsigned r2=r1+1; r2<row; ++r2) {
-            for (unsigned c=r1+1; c<col; ++c)
-                this->m[r2*col+c] = this->m[r1*col+r1]*this->m[r2*col+c] - this->m[r2*col+r1]*this->m[r1*col+c];
-            for (unsigned c=0; c<=r1; ++c)
-                this->m[r2*col+c] = _ex0();
+    const unsigned m = this->rows();
+    const unsigned n = this->cols();
+    GINAC_ASSERT(!det || n==m);
+    int sign = 1;
+    
+    unsigned r0 = 0;
+    for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
+        int indx = pivot(r0, r1, true);
+        if (indx==-1) {
+            sign = 0;
+            if (det)
+                return 0;  // leaves *this in a messy state
+        }
+        if (indx>=0) {
+            if (indx>0)
+                sign = -sign;
+            for (unsigned r2=r0+1; r2<m; ++r2) {
+                for (unsigned c=r1+1; c<n; ++c)
+                    this->m[r2*n+c] = (this->m[r0*n+r1]*this->m[r2*n+c] - this->m[r2*n+r1]*this->m[r0*n+c]).expand();
+                // fill up left hand side with zeros
+                for (unsigned c=0; c<=r1; ++c)
+                    this->m[r2*n+c] = _ex0();
+            }
+            if (det) {
+                // save space by deleting no longer needed elements
+                for (unsigned c=r0+1; c<n; ++c)
+                    this->m[r0*n+c] = _ex0();
+            }
+            ++r0;
         }
     }
     
@@ -1030,11 +1020,11 @@ int matrix::division_free_elimination(void)
  *  is possible, since we know the divisor at each step.
  *  
  *  @param det may be set to true to save a lot of space if one is only
- *  interested in the last element (i.e. for calculating determinants), the
+ *  interested in the last element (i.e. for calculating determinants). The
  *  others are set to zero in this case.
  *  @return sign is 1 if an even number of rows was swapped, -1 if an odd
  *  number of rows was swapped and 0 if the matrix is singular. */
-int matrix::fraction_free_elimination(bool det)
+int matrix::fraction_free_elimination(const bool det)
 {
     // Method:
     // (single-step fraction free elimination scheme, already known to Jordan)
@@ -1060,12 +1050,13 @@ int matrix::fraction_free_elimination(bool det)
     // and D{m[k+1](r,c)} by
     //     D{m[k-1](k-1,k-1)}.
     
-    GINAC_ASSERT(det || row==col);
     ensure_if_modifiable();
-    if (rows()==1)
-        return 1;
-    
+    const unsigned m = this->rows();
+    const unsigned n = this->cols();
+    GINAC_ASSERT(!det || n==m);
     int sign = 1;
+    if (m==1)
+        return 1;
     ex divisor_n = 1;
     ex divisor_d = 1;
     ex dividend_n;
@@ -1079,62 +1070,70 @@ int matrix::fraction_free_elimination(bool det)
     // need GCDs) since the elements of *this might be unnormalized, which
     // makes things more complicated than they need to be.
     matrix tmp_n(*this);
-    matrix tmp_d(row,col);  // for denominators, if needed
+    matrix tmp_d(m,n);  // for denominators, if needed
     lst srl;  // symbol replacement list
-    exvector::iterator it = m.begin();
+    exvector::iterator it = this->m.begin();
     exvector::iterator tmp_n_it = tmp_n.m.begin();
     exvector::iterator tmp_d_it = tmp_d.m.begin();
-    for (; it!= m.end(); ++it, ++tmp_n_it, ++tmp_d_it) {
+    for (; it!= this->m.end(); ++it, ++tmp_n_it, ++tmp_d_it) {
         (*tmp_n_it) = (*it).normal().to_rational(srl);
         (*tmp_d_it) = (*tmp_n_it).denom();
         (*tmp_n_it) = (*tmp_n_it).numer();
     }
     
-    for (unsigned r1=0; r1<row-1; ++r1) {
-        int indx = tmp_n.pivot(r1);
-        if (det && indx==-1)
-            return 0;  // FIXME: what to do if det is false, some day?
-        if (indx>0) {
-            sign = -sign;
-            // rows r1 and indx were swapped, so pivot matrix tmp_d:
-            for (unsigned c=0; c<col; ++c)
-                tmp_d.m[row*indx+c].swap(tmp_d.m[row*r1+c]);
+    unsigned r0 = 0;
+    for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
+        int indx = tmp_n.pivot(r0, r1, true);
+        if (indx==-1) {
+            sign = 0;
+            if (det)
+                return 0;
         }
-        if (r1>0) {
-            divisor_n = tmp_n.m[(r1-1)*col+(r1-1)].expand();
-            divisor_d = tmp_d.m[(r1-1)*col+(r1-1)].expand();
-            // save space by deleting no longer needed elements:
-            if (det) {
-                for (unsigned c=0; c<col; ++c) {
-                    tmp_n.m[(r1-1)*col+c] = 0;
-                    tmp_d.m[(r1-1)*col+c] = 1;
+        if (indx>=0) {
+            if (indx>0) {
+                sign = -sign;
+                // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
+                for (unsigned c=r1; c<n; ++c)
+                    tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
+            }
+            for (unsigned r2=r0+1; r2<m; ++r2) {
+                for (unsigned c=r1+1; c<n; ++c) {
+                    dividend_n = (tmp_n.m[r0*n+r1]*tmp_n.m[r2*n+c]*
+                                  tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]
+                                 -tmp_n.m[r2*n+r1]*tmp_n.m[r0*n+c]*
+                                  tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
+                    dividend_d = (tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]*
+                                  tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
+                    bool check = divide(dividend_n, divisor_n,
+                                        tmp_n.m[r2*n+c], true);
+                    check &= divide(dividend_d, divisor_d,
+                                    tmp_d.m[r2*n+c], true);
+                    GINAC_ASSERT(check);
                 }
+                // fill up left hand side with zeros
+                for (unsigned c=0; c<=r1; ++c)
+                    tmp_n.m[r2*n+c] = _ex0();
             }
-        }
-        for (unsigned r2=r1+1; r2<row; ++r2) {
-            for (unsigned c=r1+1; c<col; ++c) {
-                dividend_n = (tmp_n.m[r1*col+r1]*tmp_n.m[r2*col+c]*
-                              tmp_d.m[r2*col+r1]*tmp_d.m[r1*col+c]
-                             -tmp_n.m[r2*col+r1]*tmp_n.m[r1*col+c]*
-                              tmp_d.m[r1*col+r1]*tmp_d.m[r2*col+c]).expand();
-                dividend_d = (tmp_d.m[r2*col+r1]*tmp_d.m[r1*col+c]*
-                              tmp_d.m[r1*col+r1]*tmp_d.m[r2*col+c]).expand();
-                bool check = divide(dividend_n, divisor_n,
-                                    tmp_n.m[r2*col+c],true);
-                check &= divide(dividend_d, divisor_d,
-                                tmp_d.m[r2*col+c],true);
-                GINAC_ASSERT(check);
+            if ((r1<n-1)&&(r0<m-1)) {
+                // compute next iteration's divisor
+                divisor_n = tmp_n.m[r0*n+r1].expand();
+                divisor_d = tmp_d.m[r0*n+r1].expand();
+                if (det) {
+                    // save space by deleting no longer needed elements
+                    for (unsigned c=0; c<n; ++c) {
+                        tmp_n.m[r0*n+c] = _ex0();
+                        tmp_d.m[r0*n+c] = _ex1();
+                    }
+                }
             }
-            // fill up left hand side.
-            for (unsigned c=0; c<=r1; ++c)
-                tmp_n.m[r2*col+c] = _ex0();
+            ++r0;
         }
     }
     // repopulate *this matrix:
-    it = m.begin();
+    it = this->m.begin();
     tmp_n_it = tmp_n.m.begin();
     tmp_d_it = tmp_d.m.begin();
-    for (; it!= m.end(); ++it, ++tmp_n_it, ++tmp_d_it)
+    for (; it!= this->m.end(); ++it, ++tmp_n_it, ++tmp_d_it)
         (*it) = ((*tmp_n_it)/(*tmp_d_it)).subs(srl);
     
     return sign;
@@ -1147,68 +1146,72 @@ int matrix::fraction_free_elimination(bool det)
  *  where the element was found.  With (symbolic==true) it does the same thing
  *  with the first non-zero element.
  *
- *  @param ro is the row to be inspected
+ *  @param ro is the row from where to begin
+ *  @param co is the column to be inspected
  *  @param symbolic signal if we want the first non-zero element to be pivoted
  *  (true) or the one with the largest absolute value (false).
  *  @return 0 if no interchange occured, -1 if all are zero (usually signaling
  *  a degeneracy) and positive integer k means that rows ro and k were swapped.
  */
-int matrix::pivot(unsigned ro, bool symbolic)
+int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
 {
     unsigned k = ro;
-    
-    if (symbolic) {  // search first non-zero
-        for (unsigned r=ro; r<row; ++r) {
-            if (!m[r*col+ro].is_zero()) {
-                k = r;
-                break;
-            }
-        }
-    } else {  // search largest
-        numeric tmp(0);
-        numeric maxn(-1);
-        for (unsigned r=ro; r<row; ++r) {
-            GINAC_ASSERT(is_ex_of_type(m[r*col+ro],numeric));
-            if ((tmp = abs(ex_to_numeric(m[r*col+ro]))) > maxn &&
-                !tmp.is_zero()) {
-                maxn = tmp;
-                k = r;
+    if (symbolic) {
+        // search first non-zero element in column co beginning at row ro
+        while ((k<row) && (this->m[k*col+co].expand().is_zero()))
+            ++k;
+    } else {
+        // search largest element in column co beginning at row ro
+        GINAC_ASSERT(is_ex_of_type(this->m[k*col+co],numeric));
+        unsigned kmax = k+1;
+        numeric mmax = abs(ex_to_numeric(m[kmax*col+co]));
+        while (kmax<row) {
+            GINAC_ASSERT(is_ex_of_type(this->m[kmax*col+co],numeric));
+            numeric tmp = ex_to_numeric(this->m[kmax*col+co]);
+            if (abs(tmp) > mmax) {
+                mmax = tmp;
+                k = kmax;
             }
+            ++kmax;
         }
+        if (!mmax.is_zero())
+            k = kmax;
     }
-    if (m[k*col+ro].is_zero())
+    if (k==row)
+        // all elements in column co below row ro vanish
         return -1;
-    if (k!=ro) {  // swap rows
-        ensure_if_modifiable();
-        for (unsigned c=0; c<col; ++c) {
-            m[k*col+c].swap(m[ro*col+c]);
-        }
-        return k;
-    }
-    return 0;
+    if (k==ro)
+        // matrix needs no pivoting
+        return 0;
+    // matrix needs pivoting, so swap rows k and ro
+    ensure_if_modifiable();
+    for (unsigned c=0; c<col; ++c)
+        m[k*col+c].swap(m[ro*col+c]);
+    
+    return k;
 }
 
 /** Convert list of lists to matrix. */
 ex lst_to_matrix(const ex &l)
 {
-       if (!is_ex_of_type(l, lst))
-               throw(std::invalid_argument("argument to lst_to_matrix() must be a lst"));
-
-       // Find number of rows and columns
-       unsigned rows = l.nops(), cols = 0, i, j;
-       for (i=0; i<rows; i++)
-               if (l.op(i).nops() > cols)
-                       cols = l.op(i).nops();
-
-       // Allocate and fill matrix
-       matrix &m = *new matrix(rows, cols);
-       for (i=0; i<rows; i++)
-               for (j=0; j<cols; j++)
-                       if (l.op(i).nops() > j)
-                               m.set(i, j, l.op(i).op(j));
-                       else
-                               m.set(i, j, ex(0));
-       return m;
+    if (!is_ex_of_type(l, lst))
+        throw(std::invalid_argument("argument to lst_to_matrix() must be a lst"));
+    
+    // Find number of rows and columns
+    unsigned rows = l.nops(), cols = 0, i, j;
+    for (i=0; i<rows; i++)
+        if (l.op(i).nops() > cols)
+            cols = l.op(i).nops();
+    
+    // Allocate and fill matrix
+    matrix &m = *new matrix(rows, cols);
+    for (i=0; i<rows; i++)
+        for (j=0; j<cols; j++)
+            if (l.op(i).nops() > j)
+                m.set(i, j, l.op(i).op(j));
+            else
+                m.set(i, j, ex(0));
+    return m;
 }
 
 //////////