]> www.ginac.de Git - ginac.git/blobdiff - ginac/matrix.cpp
unit_matrix(r,c) can return an evaluated matrix object.
[ginac.git] / ginac / matrix.cpp
index db5916178ff3617d0513af5e74e7426f99069d5b..82c2255fb76c30440a46d7de495b516d50940289 100644 (file)
@@ -3,7 +3,7 @@
  *  Implementation of symbolic matrices */
 
 /*
- *  GiNaC Copyright (C) 1999-2003 Johannes Gutenberg University Mainz, Germany
+ *  GiNaC Copyright (C) 1999-2015 Johannes Gutenberg University Mainz, Germany
  *
  *  This program is free software; you can redistribute it and/or modify
  *  it under the terms of the GNU General Public License as published by
  *
  *  You should have received a copy of the GNU General Public License
  *  along with this program; if not, write to the Free Software
- *  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
+ *  Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
  */
 
-#include <string>
-#include <iostream>
-#include <sstream>
-#include <algorithm>
-#include <map>
-#include <stdexcept>
-
 #include "matrix.h"
 #include "numeric.h"
 #include "lst.h"
 #include "archive.h"
 #include "utils.h"
 
+#include <algorithm>
+#include <iostream>
+#include <map>
+#include <sstream>
+#include <stdexcept>
+#include <string>
+
 namespace GiNaC {
 
 GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(matrix, basic,
   print_func<print_context>(&matrix::do_print).
   print_func<print_latex>(&matrix::do_print_latex).
-  print_func<print_tree>(&basic::do_print_tree).
+  print_func<print_tree>(&matrix::do_print_tree).
   print_func<print_python_repr>(&matrix::do_print_python_repr))
 
 //////////
@@ -53,9 +53,9 @@ GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(matrix, basic,
 //////////
 
 /** Default ctor.  Initializes to 1 x 1-dimensional zero-matrix. */
-matrix::matrix() : inherited(TINFO_matrix), row(1), col(1)
+matrix::matrix() : row(1), col(1), m(1, _ex0)
 {
-       m.push_back(_ex0);
+       setflag(status_flags::not_shareable);
 }
 
 //////////
@@ -68,26 +68,28 @@ matrix::matrix() : inherited(TINFO_matrix), row(1), col(1)
  *
  *  @param r number of rows
  *  @param c number of cols */
-matrix::matrix(unsigned r, unsigned c)
-  : inherited(TINFO_matrix), row(r), col(c)
+matrix::matrix(unsigned r, unsigned c) : row(r), col(c), m(r*c, _ex0)
 {
-       m.resize(r*c, _ex0);
+       setflag(status_flags::not_shareable);
 }
 
 // protected
 
 /** Ctor from representation, for internal use only. */
 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
-  : inherited(TINFO_matrix), row(r), col(c), m(m2) {}
+  : row(r), col(c), m(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
  *  thrown away. */
 matrix::matrix(unsigned r, unsigned c, const lst & l)
-  : inherited(TINFO_matrix), row(r), col(c)
+  : row(r), col(c), m(r*c, _ex0)
 {
-       m.resize(r*c, _ex0);
+       setflag(status_flags::not_shareable);
 
        size_t i = 0;
        for (lst::const_iterator it = l.begin(); it != l.end(); ++it, ++i) {
@@ -103,19 +105,25 @@ matrix::matrix(unsigned r, unsigned c, const lst & l)
 // archiving
 //////////
 
-matrix::matrix(const archive_node &n, lst &sym_lst) : inherited(n, sym_lst)
+void matrix::read_archive(const archive_node &n, lst &sym_lst)
 {
+       inherited::read_archive(n, sym_lst);
+
        if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
                throw (std::runtime_error("unknown matrix dimensions in archive"));
        m.reserve(row * col);
-       for (unsigned int i=0; true; i++) {
+       // XXX: default ctor inserts a zero element, we need to erase it here.
+       m.pop_back();
+       archive_node::archive_node_cit first = n.find_first("m");
+       archive_node::archive_node_cit last = n.find_last("m");
+       ++last;
+       for (archive_node::archive_node_cit i=first; i != last; ++i) {
                ex e;
-               if (n.find_ex("m", e, sym_lst, i))
-                       m.push_back(e);
-               else
-                       break;
+               n.find_ex_by_loc(i, e, sym_lst);
+               m.push_back(e);
        }
 }
+GINAC_BIND_UNARCHIVER(matrix);
 
 void matrix::archive(archive_node &n) const
 {
@@ -129,8 +137,6 @@ void matrix::archive(archive_node &n) const
        }
 }
 
-DEFAULT_UNARCHIVE(matrix)
-
 //////////
 // functions overriding virtual functions from base classes
 //////////
@@ -216,7 +222,7 @@ ex matrix::eval(int level) const
                        m2[r*col+c] = m[r*col+c].eval(level);
        
        return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
-                                                                                          status_flags::evaluated);
+                                                  status_flags::evaluated);
 }
 
 ex matrix::subs(const exmap & mp, unsigned options) const
@@ -229,6 +235,52 @@ ex matrix::subs(const exmap & mp, unsigned options) const
        return matrix(row, col, m2).subs_one_level(mp, options);
 }
 
+/** Complex conjugate every matrix entry. */
+ex matrix::conjugate() const
+{
+       exvector * ev = 0;
+       for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i) {
+               ex x = i->conjugate();
+               if (ev) {
+                       ev->push_back(x);
+                       continue;
+               }
+               if (are_ex_trivially_equal(x, *i)) {
+                       continue;
+               }
+               ev = new exvector;
+               ev->reserve(m.size());
+               for (exvector::const_iterator j=m.begin(); j!=i; ++j) {
+                       ev->push_back(*j);
+               }
+               ev->push_back(x);
+       }
+       if (ev) {
+               ex result = matrix(row, col, *ev);
+               delete ev;
+               return result;
+       }
+       return *this;
+}
+
+ex matrix::real_part() const
+{
+       exvector v;
+       v.reserve(m.size());
+       for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
+               v.push_back(i->real_part());
+       return matrix(row, col, v);
+}
+
+ex matrix::imag_part() const
+{
+       exvector v;
+       v.reserve(m.size());
+       for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
+               v.push_back(i->imag_part());
+       return matrix(row, col, v);
+}
+
 // protected
 
 int matrix::compare_same_type(const basic & other) const
@@ -547,10 +599,11 @@ matrix matrix::mul(const matrix & other) const
        
        for (unsigned r1=0; r1<this->rows(); ++r1) {
                for (unsigned c=0; c<this->cols(); ++c) {
+                       // Quick test: can we shortcut?
                        if (m[r1*col+c].is_zero())
                                continue;
                        for (unsigned r2=0; r2<other.cols(); ++r2)
-                               prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]).expand();
+                               prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]);
                }
        }
        return matrix(row, other.col, prod);
@@ -614,12 +667,12 @@ matrix matrix::pow(const ex & expn) const
                        // that this is not entirely optimal but close to optimal and
                        // "better" algorithms are much harder to implement.  (See Knuth,
                        // TAoCP2, section "Evaluation of Powers" for a good discussion.)
-                       while (b!=_num1) {
+                       while (b!=*_num1_p) {
                                if (b.is_odd()) {
                                        C = C.mul(A);
                                        --b;
                                }
-                               b /= _num2;  // still integer.
+                               b /= *_num2_p;  // still integer.
                                A = A.mul(A);
                        }
                        return A.mul(C);
@@ -697,12 +750,12 @@ ex matrix::determinant(unsigned algo) const
        unsigned sparse_count = 0;  // counts non-zero elements
        exvector::const_iterator r = m.begin(), rend = m.end();
        while (r != rend) {
-               lst srl;  // symbol replacement list
+               if (!r->info(info_flags::numeric))
+                       numeric_flag = false;
+               exmap srl;  // symbol replacement list
                ex rtest = r->to_rational(srl);
                if (!rtest.is_zero())
                        ++sparse_count;
-               if (!rtest.info(info_flags::numeric))
-                       numeric_flag = false;
                if (!rtest.info(info_flags::crational_polynomial) &&
                         rtest.info(info_flags::rational_function))
                        normal_flag = true;
@@ -731,7 +784,7 @@ ex matrix::determinant(unsigned algo) const
                else
                        return m[0].expand();
        }
-       
+
        // Compute the determinant
        switch(algo) {
                case determinant_algo::gauss: {
@@ -827,7 +880,7 @@ ex matrix::trace() const
                tr += m[r*col+r];
        
        if (tr.info(info_flags::rational_function) &&
-               !tr.info(info_flags::crational_polynomial))
+          !tr.info(info_flags::crational_polynomial))
                return tr.normal();
        else
                return tr.expand();
@@ -835,7 +888,7 @@ ex matrix::trace() const
 
 
 /** Characteristic Polynomial.  Following mathematica notation the
- *  characteristic polynomial of a matrix M is defined as the determiant of
+ *  characteristic polynomial of a matrix M is defined as the determinant of
  *  (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
  *  as M.  Note that some CASs define it with a sign inside the determinant
  *  which gives rise to an overall sign if the dimension is odd.  This method
@@ -845,7 +898,7 @@ ex matrix::trace() const
  *  @return    characteristic polynomial as new expression
  *  @exception logic_error (matrix not square)
  *  @see       matrix::determinant() */
-ex matrix::charpoly(const symbol & lambda) const
+ex matrix::charpoly(const ex & lambda) const
 {
        if (row != col)
                throw (std::logic_error("matrix::charpoly(): matrix not square"));
@@ -865,13 +918,13 @@ ex matrix::charpoly(const symbol & lambda) const
 
                matrix B(*this);
                ex c = B.trace();
-               ex poly = power(lambda,row)-c*power(lambda,row-1);
+               ex poly = power(lambda, row) - c*power(lambda, row-1);
                for (unsigned i=1; i<row; ++i) {
                        for (unsigned j=0; j<row; ++j)
                                B.m[j*col+j] -= c;
                        B = this->mul(B);
                        c = B.trace() / ex(i+1);
-                       poly -= c*power(lambda,row-i-1);
+                       poly -= c*power(lambda, row-i-1);
                }
                if (row%2)
                        return -poly;
@@ -933,14 +986,15 @@ matrix matrix::inverse() const
  *
  *  @param vars n x p matrix, all elements must be symbols 
  *  @param rhs m x p matrix
+ *  @param algo selects the solving algorithm
  *  @return n x p solution matrix
  *  @exception logic_error (incompatible matrices)
  *  @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
+                     const matrix & rhs,
+                     unsigned algo) const
 {
        const unsigned m = this->rows();
        const unsigned n = this->cols();
@@ -1033,6 +1087,29 @@ matrix matrix::solve(const matrix & vars,
 }
 
 
+/** Compute the rank of this matrix. */
+unsigned matrix::rank() const
+{
+       // Method:
+       // Transform this matrix into upper echelon form and then count the
+       // number of non-zero rows.
+
+       GINAC_ASSERT(row*col==m.capacity());
+
+       // Actually, any elimination scheme will do since we are only
+       // interested in the echelon matrix' zeros.
+       matrix to_eliminate = *this;
+       to_eliminate.fraction_free_elimination();
+
+       unsigned r = row*col;  // index of last non-zero element
+       while (r--) {
+               if (!to_eliminate.m[r].is_zero())
+                       return 1+r/col;
+       }
+       return 0;
+}
+
+
 // protected
 
 /** Recursive determinant for small matrices having at least one symbolic
@@ -1041,7 +1118,7 @@ matrix matrix::solve(const matrix & vars,
  *  more than once.  According to W.M.Gentleman and S.C.Johnson this algorithm
  *  is better than elimination schemes for matrices of sparse multivariate
  *  polynomials and also for matrices of dense univariate polynomials if the
- *  matrix' dimesion is larger than 7.
+ *  matrix' dimension is larger than 7.
  *
  *  @return the determinant as a new expression (in expanded form)
  *  @see matrix::determinant() */
@@ -1148,7 +1225,7 @@ ex matrix::determinant_minor() const
                                        Pkey[j] = Pkey[j-1]+1;
                } while(fc);
                // next column, so change the role of A and B:
-               A = B;
+               A.swap(B);
                B.clear();
        }
        
@@ -1174,8 +1251,8 @@ int matrix::gauss_elimination(const bool det)
        int sign = 1;
        
        unsigned r0 = 0;
-       for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
-               int indx = pivot(r0, r1, true);
+       for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
+               int indx = pivot(r0, c0, true);
                if (indx == -1) {
                        sign = 0;
                        if (det)
@@ -1185,17 +1262,17 @@ int matrix::gauss_elimination(const bool det)
                        if (indx > 0)
                                sign = -sign;
                        for (unsigned r2=r0+1; r2<m; ++r2) {
-                               if (!this->m[r2*n+r1].is_zero()) {
+                               if (!this->m[r2*n+c0].is_zero()) {
                                        // yes, there is something to do in this row
-                                       ex piv = this->m[r2*n+r1] / this->m[r0*n+r1];
-                                       for (unsigned c=r1+1; c<n; ++c) {
+                                       ex piv = this->m[r2*n+c0] / this->m[r0*n+c0];
+                                       for (unsigned c=c0+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)
+                               for (unsigned c=r0; c<=c0; ++c)
                                        this->m[r2*n+c] = _ex0;
                        }
                        if (det) {
@@ -1206,7 +1283,12 @@ int matrix::gauss_elimination(const bool det)
                        ++r0;
                }
        }
-       
+       // clear remaining rows
+       for (unsigned r=r0+1; r<m; ++r) {
+               for (unsigned c=0; c<n; ++c)
+                       this->m[r*n+c] = _ex0;
+       }
+
        return sign;
 }
 
@@ -1228,8 +1310,8 @@ int matrix::division_free_elimination(const bool det)
        int sign = 1;
        
        unsigned r0 = 0;
-       for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
-               int indx = pivot(r0, r1, true);
+       for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
+               int indx = pivot(r0, c0, true);
                if (indx==-1) {
                        sign = 0;
                        if (det)
@@ -1239,10 +1321,10 @@ int matrix::division_free_elimination(const bool det)
                        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();
+                               for (unsigned c=c0+1; c<n; ++c)
+                                       this->m[r2*n+c] = (this->m[r0*n+c0]*this->m[r2*n+c] - this->m[r2*n+c0]*this->m[r0*n+c]).expand();
                                // fill up left hand side with zeros
-                               for (unsigned c=0; c<=r1; ++c)
+                               for (unsigned c=r0; c<=c0; ++c)
                                        this->m[r2*n+c] = _ex0;
                        }
                        if (det) {
@@ -1253,7 +1335,12 @@ int matrix::division_free_elimination(const bool det)
                        ++r0;
                }
        }
-       
+       // clear remaining rows
+       for (unsigned r=r0+1; r<m; ++r) {
+               for (unsigned c=0; c<n; ++c)
+                       this->m[r*n+c] = _ex0;
+       }
+
        return sign;
 }
 
@@ -1315,7 +1402,7 @@ int matrix::fraction_free_elimination(const bool det)
        // makes things more complicated than they need to be.
        matrix tmp_n(*this);
        matrix tmp_d(m,n);  // for denominators, if needed
-       lst srl;  // symbol replacement list
+       exmap srl;  // symbol replacement list
        exvector::const_iterator cit = this->m.begin(), citend = this->m.end();
        exvector::iterator tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
        while (cit != citend) {
@@ -1326,28 +1413,37 @@ int matrix::fraction_free_elimination(const bool det)
        }
        
        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) {
+       for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
+               // When trying to find a pivot, we should try a bit harder than expand().
+               // Searching the first non-zero element in-place here instead of calling
+               // pivot() allows us to do no more substitutions and back-substitutions
+               // than are actually necessary.
+               unsigned indx = r0;
+               while ((indx<m) &&
+                      (tmp_n[indx*n+c0].subs(srl, subs_options::no_pattern).expand().is_zero()))
+                       ++indx;
+               if (indx==m) {
+                       // all elements in column c0 below row r0 vanish
                        sign = 0;
                        if (det)
                                return 0;
-               }
-               if (indx>=0) {
-                       if (indx>0) {
+               } else {
+                       if (indx>r0) {
+                               // Matrix needs pivoting, swap rows r0 and indx of tmp_n and tmp_d.
                                sign = -sign;
-                               // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
-                               for (unsigned c=r1; c<n; ++c)
+                               for (unsigned c=c0; c<n; ++c) {
+                                       tmp_n.m[n*indx+c].swap(tmp_n.m[n*r0+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();
+                               for (unsigned c=c0+1; c<n; ++c) {
+                                       dividend_n = (tmp_n.m[r0*n+c0]*tmp_n.m[r2*n+c]*
+                                                     tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]
+                                                    -tmp_n.m[r2*n+c0]*tmp_n.m[r0*n+c]*
+                                                     tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
+                                       dividend_d = (tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]*
+                                                     tmp_d.m[r0*n+c0]*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,
@@ -1355,13 +1451,13 @@ int matrix::fraction_free_elimination(const bool det)
                                        GINAC_ASSERT(check);
                                }
                                // fill up left hand side with zeros
-                               for (unsigned c=0; c<=r1; ++c)
+                               for (unsigned c=r0; c<=c0; ++c)
                                        tmp_n.m[r2*n+c] = _ex0;
                        }
-                       if ((r1<n-1)&&(r0<m-1)) {
+                       if (c0<n && 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();
+                               divisor_n = tmp_n.m[r0*n+c0].expand();
+                               divisor_d = tmp_d.m[r0*n+c0].expand();
                                if (det) {
                                        // save space by deleting no longer needed elements
                                        for (unsigned c=0; c<n; ++c) {
@@ -1373,6 +1469,12 @@ int matrix::fraction_free_elimination(const bool det)
                        ++r0;
                }
        }
+       // clear remaining rows
+       for (unsigned r=r0+1; r<m; ++r) {
+               for (unsigned c=0; c<n; ++c)
+                       tmp_n.m[r*n+c] = _ex0;
+       }
+
        // repopulate *this matrix:
        exvector::iterator it = this->m.begin(), itend = this->m.end();
        tmp_n_it = tmp_n.m.begin();
@@ -1394,7 +1496,7 @@ int matrix::fraction_free_elimination(const bool det)
  *  @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
+ *  @return 0 if no interchange occurred, -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, unsigned co, bool symbolic)
@@ -1435,6 +1537,16 @@ int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
        return k;
 }
 
+/** Function to check that all elements of the matrix are zero.
+ */
+bool matrix::is_zero_matrix() const
+{
+       for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i) 
+               if(!(i->is_zero()))
+                       return false;
+       return true;
+}
+
 ex lst_to_matrix(const lst & l)
 {
        lst::const_iterator itr, itc;
@@ -1481,7 +1593,7 @@ ex diag_matrix(const lst & l)
 ex unit_matrix(unsigned r, unsigned c)
 {
        matrix &Id = *new matrix(r, c);
-       Id.setflag(status_flags::dynallocated);
+       Id.setflag(status_flags::dynallocated | status_flags::evaluated);
        for (unsigned i=0; i<r && i<c; i++)
                Id(i,i) = _ex1;
 
@@ -1525,4 +1637,52 @@ ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const
        return M;
 }
 
+ex reduced_matrix(const matrix& m, unsigned r, unsigned c)
+{
+       if (r+1>m.rows() || c+1>m.cols() || m.cols()<2 || m.rows()<2)
+               throw std::runtime_error("minor_matrix(): index out of bounds");
+
+       const unsigned rows = m.rows()-1;
+       const unsigned cols = m.cols()-1;
+       matrix &M = *new matrix(rows, cols);
+       M.setflag(status_flags::dynallocated | status_flags::evaluated);
+
+       unsigned ro = 0;
+       unsigned ro2 = 0;
+       while (ro2<rows) {
+               if (ro==r)
+                       ++ro;
+               unsigned co = 0;
+               unsigned co2 = 0;
+               while (co2<cols) {
+                       if (co==c)
+                               ++co;
+                       M(ro2,co2) = m(ro, co);
+                       ++co;
+                       ++co2;
+               }
+               ++ro;
+               ++ro2;
+       }
+
+       return M;
+}
+
+ex sub_matrix(const matrix&m, unsigned r, unsigned nr, unsigned c, unsigned nc)
+{
+       if (r+nr>m.rows() || c+nc>m.cols())
+               throw std::runtime_error("sub_matrix(): index out of bounds");
+
+       matrix &M = *new matrix(nr, nc);
+       M.setflag(status_flags::dynallocated | status_flags::evaluated);
+
+       for (unsigned ro=0; ro<nr; ++ro) {
+               for (unsigned co=0; co<nc; ++co) {
+                       M(ro,co) = m(ro+r,co+c);
+               }
+       }
+
+       return M;
+}
+
 } // namespace GiNaC