]> www.ginac.de Git - ginac.git/blobdiff - ginac/matrix.cpp
Add optional matrix::rank() algorighm selection.
[ginac.git] / ginac / matrix.cpp
index 871c3f15d1c5f433238b6f5a3b354d1f17ef1ea2..6464ed8d52cd90731511ddf04e976d915fb734b9 100644 (file)
@@ -3,7 +3,7 @@
  *  Implementation of symbolic matrices */
 
 /*
- *  GiNaC Copyright (C) 1999-2015 Johannes Gutenberg University Mainz, Germany
+ *  GiNaC Copyright (C) 1999-2018 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
@@ -73,20 +73,6 @@ matrix::matrix(unsigned r, unsigned c) : row(r), col(c), m(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)
-  : 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
@@ -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
 //////////
@@ -207,28 +227,6 @@ ex & matrix::let_op(size_t i)
        return m[i];
 }
 
-/** Evaluate matrix entry by entry. */
-ex matrix::eval(int level) const
-{
-       // check if we have to do anything at all
-       if ((level==1)&&(flags & status_flags::evaluated))
-               return *this;
-       
-       // emergency break
-       if (level == -max_recursion_level)
-               throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
-       
-       // eval() entry by entry
-       exvector m2(row*col);
-       --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, std::move(m2)))->setflag(status_flags::dynallocated |
-                                                             status_flags::evaluated);
-}
-
 ex matrix::subs(const exmap & mp, unsigned options) const
 {
        exvector m2(row * col);
@@ -939,12 +937,19 @@ ex matrix::charpoly(const ex & lambda) const
 }
 
 
+/** Inverse of this matrix, with automatic algorithm selection. */
+matrix matrix::inverse() const
+{
+       return inverse(solve_algo::automatic);
+}
+
 /** Inverse of this matrix.
  *
+ *  @param algo selects the algorithm (one of solve_algo)
  *  @return    the inverted matrix
  *  @exception logic_error (matrix not square)
  *  @exception runtime_error (singular matrix) */
-matrix matrix::inverse() const
+matrix matrix::inverse(unsigned algo) const
 {
        if (row != col)
                throw (std::logic_error("matrix::inverse(): matrix not square"));
@@ -967,7 +972,7 @@ matrix matrix::inverse() const
        
        matrix sol(row,col);
        try {
-               sol = this->solve(vars,identity);
+               sol = this->solve(vars, identity, algo);
        } catch (const std::runtime_error & e) {
            if (e.what()==std::string("matrix::solve(): inconsistent linear system"))
                        throw (std::runtime_error("matrix::inverse(): singular matrix"));
@@ -981,7 +986,7 @@ matrix matrix::inverse() const
 /** 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, all elements must be symbols 
+ *  @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
@@ -997,8 +1002,8 @@ matrix matrix::solve(const matrix & vars,
        const unsigned n = this->cols();
        const unsigned p = rhs.cols();
        
-       // syntax checks    
-       if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
+       // syntax checks
+       if ((rhs.rows() != m) || (vars.rows() != n) || (vars.cols() != p))
                throw (std::logic_error("matrix::solve(): incompatible matrices"));
        for (unsigned ro=0; ro<n; ++ro)
                for (unsigned co=0; co<p; ++co)
@@ -1013,41 +1018,9 @@ matrix matrix::solve(const matrix & vars,
                for (unsigned c=0; c<p; ++c)
                        aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
        }
-       
-       // Gather some statistical information about the augmented matrix:
-       bool numeric_flag = true;
-       for (auto & r : aug.m) {
-               if (!r.info(info_flags::numeric)) {
-                       numeric_flag = false;
-                       break;
-               }
-       }
-       
-       // 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;
-       }
-       
+
        // Eliminate the augmented matrix:
-       switch(algo) {
-               case solve_algo::gauss:
-                       aug.gauss_elimination();
-                       break;
-               case solve_algo::divfree:
-                       aug.division_free_elimination();
-                       break;
-               case solve_algo::bareiss:
-               default:
-                       aug.fraction_free_elimination();
-       }
+       auto colid = aug.echelon_form(algo, n);
        
        // assemble the solution matrix:
        matrix sol(n,p);
@@ -1055,48 +1028,51 @@ matrix matrix::solve(const matrix & vars,
                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()))
+                       while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].normal().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()) {
+                               if (!aug.m[r*(n+p)+n+co].normal().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(c,co) = vars.m[c*p+co];
+                                       sol(colid[c],co) = vars.m[colid[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(fnz-1,co) = (e/(aug.m[r*(n+p)+(fnz-1)])).normal();
+                                       e -= aug.m[r*(n+p)+c]*sol.m[colid[c]*p+co];
+                               sol(colid[fnz-1],co) = (e/(aug.m[r*(n+p)+fnz-1])).normal();
                                last_assigned_sol = fnz;
                        }
                }
                // assign solutions for vars between 1 and
                // last_assigned_sol-1: free parameters
                for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
-                       sol(ro,co) = vars(ro,co);
+                       sol(colid[ro],co) = vars(colid[ro],co);
        }
        
        return sol;
 }
 
-
 /** Compute the rank of this matrix. */
 unsigned matrix::rank() const
+{
+       return rank(solve_algo::automatic);
+}
+
+/** Compute the rank of this matrix using the given algorithm,
+ *  which should be a member of enum solve_algo. */
+unsigned matrix::rank(unsigned solve_algo) 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();
+       to_eliminate.echelon_form(solve_algo, col);
 
        unsigned r = row*col;  // index of last non-zero element
        while (r--) {
@@ -1221,14 +1197,59 @@ ex matrix::determinant_minor() const
                                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:
-               A.swap(B);
-               B.clear();
+               // next column, clear B and change the role of A and B:
+               A = std::move(B);
        }
        
        return det;
 }
 
+std::vector<unsigned>
+matrix::echelon_form(unsigned algo, int n)
+{
+       // Here is the heuristics in case this routine has to decide:
+       if (algo == solve_algo::automatic) {
+               // Gather some statistical information about the augmented matrix:
+               bool numeric_flag = true;
+               for (auto & r : m) {
+                       if (!r.info(info_flags::numeric)) {
+                               numeric_flag = false;
+                               break;
+                       }
+               }
+               // Bareiss (fraction-free) elimination is generally a good guess:
+               algo = solve_algo::bareiss;
+               // For row<3, Bareiss elimination is equivalent to division free
+               // elimination but has more logistic overhead
+               if (row<3)
+                       algo = solve_algo::divfree;
+               // This overrides any prior decisions.
+               if (numeric_flag)
+                       algo = solve_algo::gauss;
+       }
+       // Eliminate the augmented matrix:
+       std::vector<unsigned> colid(col);
+       for (unsigned c = 0; c < col; c++) {
+               colid[c] = c;
+       }
+       switch(algo) {
+               case solve_algo::gauss:
+                       gauss_elimination();
+                       break;
+               case solve_algo::divfree:
+                       division_free_elimination();
+                       break;
+               case solve_algo::bareiss:
+                       fraction_free_elimination();
+                       break;
+               case solve_algo::markowitz:
+                       colid = markowitz_elimination(n);
+                       break;
+               default:
+                       throw std::invalid_argument("matrix::echelon_form(): 'algo' is not one of the solve_algo enum");
+       }
+       return colid;
+}
 
 /** 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
@@ -1289,6 +1310,119 @@ int matrix::gauss_elimination(const bool det)
        return sign;
 }
 
+/* Perform Markowitz-ordered Gaussian elimination (with full
+ * pivoting) on a matrix, constraining the choice of pivots to
+ * the first n columns (this simplifies handling of augmented
+ * matrices). Return the column id vector v, such that v[column]
+ * is the original number of the column before shuffling (v[i]==i
+ * for i >= n). */
+std::vector<unsigned>
+matrix::markowitz_elimination(unsigned n)
+{
+       GINAC_ASSERT(n <= col);
+       std::vector<int> rowcnt(row, 0);
+       std::vector<int> colcnt(col, 0);
+       // Normalize everything before start. We'll keep all the
+       // cells normalized throughout the algorithm to properly
+       // handle unnormal zeros.
+       for (unsigned r = 0; r < row; r++) {
+               for (unsigned c = 0; c < col; c++) {
+                       if (!m[r*col + c].is_zero()) {
+                               m[r*col + c] = m[r*col + c].normal();
+                               rowcnt[r]++;
+                               colcnt[c]++;
+                       }
+               }
+       }
+       std::vector<unsigned> colid(col);
+       for (unsigned c = 0; c < col; c++) {
+               colid[c] = c;
+       }
+       exvector ab(row);
+       for (unsigned k = 0; (k < col) && (k < row - 1); k++) {
+               // Find the pivot that minimizes (rowcnt[r]-1)*(colcnt[c]-1).
+               unsigned pivot_r = row + 1;
+               unsigned pivot_c = col + 1;
+               int pivot_m = row*col;
+               for (unsigned r = k; r < row; r++) {
+                       for (unsigned c = k; c < n; c++) {
+                               const ex &mrc = m[r*col + c];
+                               if (mrc.is_zero())
+                                       continue;
+                               GINAC_ASSERT(rowcnt[r] > 0);
+                               GINAC_ASSERT(colcnt[c] > 0);
+                               int measure = (rowcnt[r] - 1)*(colcnt[c] - 1);
+                               if (measure < pivot_m) {
+                                       pivot_m = measure;
+                                       pivot_r = r;
+                                       pivot_c = c;
+                               }
+                       }
+               }
+               if (pivot_m == row*col) {
+                       // The rest of the matrix is zero.
+                       break;
+               }
+               GINAC_ASSERT(k <= pivot_r && pivot_r < row);
+               GINAC_ASSERT(k <= pivot_c && pivot_c < col);
+               // Swap the pivot into (k, k).
+               if (pivot_c != k) {
+                       for (unsigned r = 0; r < row; r++) {
+                               m[r*col + pivot_c].swap(m[r*col + k]);
+                       }
+                       std::swap(colid[pivot_c], colid[k]);
+                       std::swap(colcnt[pivot_c], colcnt[k]);
+               }
+               if (pivot_r != k) {
+                       for (unsigned c = k; c < col; c++) {
+                               m[pivot_r*col + c].swap(m[k*col + c]);
+                       }
+                       std::swap(rowcnt[pivot_r], rowcnt[k]);
+               }
+               // No normalization before is_zero() here, because
+               // we maintain the matrix normalized throughout the
+               // algorithm.
+               ex a = m[k*col + k];
+               GINAC_ASSERT(!a.is_zero());
+               // Subtract the pivot row KJI-style (so: loop by pivot, then
+               // column, then row) to maximally exploit pivot row zeros (at
+               // the expense of the pivot column zeros). The speedup compared
+               // to the usual KIJ order is not really significant though...
+               for (unsigned r = k + 1; r < row; r++) {
+                       const ex &b = m[r*col + k];
+                       if (!b.is_zero()) {
+                               ab[r] = b/a;
+                               rowcnt[r]--;
+                       }
+               }
+               colcnt[k] = rowcnt[k] = 0;
+               for (unsigned c = k + 1; c < col; c++) {
+                       const ex &mr0c = m[k*col + c];
+                       if (mr0c.is_zero())
+                               continue;
+                       colcnt[c]--;
+                       for (unsigned r = k + 1; r < row; r++) {
+                               if (ab[r].is_zero())
+                                       continue;
+                               bool waszero = m[r*col + c].is_zero();
+                               m[r*col + c] = (m[r*col + c] - ab[r]*mr0c).normal();
+                               bool iszero = m[r*col + c].is_zero();
+                               if (waszero && !iszero) {
+                                       rowcnt[r]++;
+                                       colcnt[c]++;
+                               }
+                               if (!waszero && iszero) {
+                                       rowcnt[r]--;
+                                       colcnt[c]--;
+                               }
+                       }
+               }
+               for (unsigned r = k + 1; r < row; r++) {
+                       ab[r] = m[r*col + k] = _ex0;
+               }
+       }
+       return colid;
+}
 
 /** Perform the steps of division free elimination to bring the m x n matrix
  *  into an upper echelon form.
@@ -1553,8 +1687,7 @@ ex lst_to_matrix(const lst & l)
        }
 
        // Allocate and fill matrix
-       matrix &M = *new matrix(rows, cols);
-       M.setflag(status_flags::dynallocated);
+       matrix & M = dynallocate<matrix>(rows, cols);
 
        unsigned i = 0;
        for (auto & itr : l) {
@@ -1574,8 +1707,23 @@ ex diag_matrix(const lst & l)
        size_t dim = l.nops();
 
        // Allocate and fill matrix
-       matrix &M = *new matrix(dim, dim);
-       M.setflag(status_flags::dynallocated);
+       matrix & M = dynallocate<matrix>(dim, dim);
+
+       unsigned i = 0;
+       for (auto & it : l) {
+               M(i, i) = it;
+               ++i;
+       }
+
+       return M;
+}
+
+ex diag_matrix(std::initializer_list<ex> l)
+{
+       size_t dim = l.size();
+
+       // Allocate and fill matrix
+       matrix & M = dynallocate<matrix>(dim, dim);
 
        unsigned i = 0;
        for (auto & it : l) {
@@ -1588,8 +1736,8 @@ ex diag_matrix(const lst & l)
 
 ex unit_matrix(unsigned r, unsigned c)
 {
-       matrix &Id = *new matrix(r, c);
-       Id.setflag(status_flags::dynallocated);
+       matrix & Id = dynallocate<matrix>(r, c);
+       Id.setflag(status_flags::evaluated);
        for (unsigned i=0; i<r && i<c; i++)
                Id(i,i) = _ex1;
 
@@ -1598,8 +1746,8 @@ ex unit_matrix(unsigned r, unsigned c)
 
 ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const std::string & tex_base_name)
 {
-       matrix &M = *new matrix(r, c);
-       M.setflag(status_flags::dynallocated | status_flags::evaluated);
+       matrix & M = dynallocate<matrix>(r, c);
+       M.setflag(status_flags::evaluated);
 
        bool long_format = (r > 10 || c > 10);
        bool single_row = (r == 1 || c == 1);
@@ -1640,8 +1788,8 @@ ex reduced_matrix(const matrix& m, unsigned r, unsigned c)
 
        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);
+       matrix & M = dynallocate<matrix>(rows, cols);
+       M.setflag(status_flags::evaluated);
 
        unsigned ro = 0;
        unsigned ro2 = 0;
@@ -1669,8 +1817,8 @@ 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);
+       matrix & M = dynallocate<matrix>(nr, nc);
+       M.setflag(status_flags::evaluated);
 
        for (unsigned ro=0; ro<nr; ++ro) {
                for (unsigned co=0; co<nc; ++co) {