3 * Implementation of symbolic matrices */
6 * GiNaC Copyright (C) 1999-2001 Johannes Gutenberg University Mainz, Germany
8 * This program is free software; you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation; either version 2 of the License, or
11 * (at your option) any later version.
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
18 * You should have received a copy of the GNU General Public License
19 * along with this program; if not, write to the Free Software
20 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
39 GINAC_IMPLEMENT_REGISTERED_CLASS(matrix, basic)
42 // default ctor, dtor, copy ctor, assignment operator and helpers:
47 /** Default ctor. Initializes to 1 x 1-dimensional zero-matrix. */
48 matrix::matrix() : inherited(TINFO_matrix), row(1), col(1)
50 debugmsg("matrix default ctor",LOGLEVEL_CONSTRUCT);
56 /** For use by copy ctor and assignment operator. */
57 void matrix::copy(const matrix & other)
59 inherited::copy(other);
62 m = other.m; // STL's vector copying invoked here
65 void matrix::destroy(bool call_parent)
67 if (call_parent) inherited::destroy(call_parent);
76 /** Very common ctor. Initializes to r x c-dimensional zero-matrix.
78 * @param r number of rows
79 * @param c number of cols */
80 matrix::matrix(unsigned r, unsigned c)
81 : inherited(TINFO_matrix), row(r), col(c)
83 debugmsg("matrix ctor from unsigned,unsigned",LOGLEVEL_CONSTRUCT);
84 m.resize(r*c, _ex0());
89 /** Ctor from representation, for internal use only. */
90 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
91 : inherited(TINFO_matrix), row(r), col(c), m(m2)
93 debugmsg("matrix ctor from unsigned,unsigned,exvector",LOGLEVEL_CONSTRUCT);
100 /** Construct object from archive_node. */
101 matrix::matrix(const archive_node &n, const lst &sym_lst) : inherited(n, sym_lst)
103 debugmsg("matrix ctor from archive_node", LOGLEVEL_CONSTRUCT);
104 if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
105 throw (std::runtime_error("unknown matrix dimensions in archive"));
106 m.reserve(row * col);
107 for (unsigned int i=0; true; i++) {
109 if (n.find_ex("m", e, sym_lst, i))
116 /** Unarchive the object. */
117 ex matrix::unarchive(const archive_node &n, const lst &sym_lst)
119 return (new matrix(n, sym_lst))->setflag(status_flags::dynallocated);
122 /** Archive the object. */
123 void matrix::archive(archive_node &n) const
125 inherited::archive(n);
126 n.add_unsigned("row", row);
127 n.add_unsigned("col", col);
128 exvector::const_iterator i = m.begin(), iend = m.end();
136 // functions overriding virtual functions from bases classes
141 void matrix::print(std::ostream & os, unsigned upper_precedence) const
143 debugmsg("matrix print",LOGLEVEL_PRINT);
145 for (unsigned r=0; r<row-1; ++r) {
147 for (unsigned c=0; c<col-1; ++c)
148 os << m[r*col+c] << ",";
149 os << m[col*(r+1)-1] << "]], ";
152 for (unsigned c=0; c<col-1; ++c)
153 os << m[(row-1)*col+c] << ",";
154 os << m[row*col-1] << "]] ]]";
157 void matrix::printraw(std::ostream & os) const
159 debugmsg("matrix printraw",LOGLEVEL_PRINT);
160 os << "matrix(" << row << "," << col <<",";
161 for (unsigned r=0; r<row-1; ++r) {
163 for (unsigned c=0; c<col-1; ++c)
164 os << m[r*col+c] << ",";
165 os << m[col*(r-1)-1] << "),";
168 for (unsigned c=0; c<col-1; ++c)
169 os << m[(row-1)*col+c] << ",";
170 os << m[row*col-1] << "))";
173 /** nops is defined to be rows x columns. */
174 unsigned matrix::nops() const
179 /** returns matrix entry at position (i/col, i%col). */
180 ex matrix::op(int i) const
185 /** returns matrix entry at position (i/col, i%col). */
186 ex & matrix::let_op(int i)
189 GINAC_ASSERT(i<nops());
194 /** expands the elements of a matrix entry by entry. */
195 ex matrix::expand(unsigned options) const
197 exvector tmp(row*col);
198 for (unsigned i=0; i<row*col; ++i)
199 tmp[i] = m[i].expand(options);
201 return matrix(row, col, tmp);
204 /** Search ocurrences. A matrix 'has' an expression if it is the expression
205 * itself or one of the elements 'has' it. */
206 bool matrix::has(const ex & other) const
208 GINAC_ASSERT(other.bp!=0);
210 // tautology: it is the expression itself
211 if (is_equal(*other.bp)) return true;
213 // search all the elements
214 for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r)
215 if ((*r).has(other)) return true;
220 /** evaluate matrix entry by entry. */
221 ex matrix::eval(int level) const
223 debugmsg("matrix eval",LOGLEVEL_MEMBER_FUNCTION);
225 // check if we have to do anything at all
226 if ((level==1)&&(flags & status_flags::evaluated))
230 if (level == -max_recursion_level)
231 throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
233 // eval() entry by entry
234 exvector m2(row*col);
236 for (unsigned r=0; r<row; ++r)
237 for (unsigned c=0; c<col; ++c)
238 m2[r*col+c] = m[r*col+c].eval(level);
240 return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
241 status_flags::evaluated );
244 /** evaluate matrix numerically entry by entry. */
245 ex matrix::evalf(int level) const
247 debugmsg("matrix evalf",LOGLEVEL_MEMBER_FUNCTION);
249 // check if we have to do anything at all
254 if (level == -max_recursion_level) {
255 throw (std::runtime_error("matrix::evalf(): recursion limit exceeded"));
258 // evalf() entry by entry
259 exvector m2(row*col);
261 for (unsigned r=0; r<row; ++r)
262 for (unsigned c=0; c<col; ++c)
263 m2[r*col+c] = m[r*col+c].evalf(level);
265 return matrix(row, col, m2);
270 int matrix::compare_same_type(const basic & other) const
272 GINAC_ASSERT(is_exactly_of_type(other, matrix));
273 const matrix & o = static_cast<matrix &>(const_cast<basic &>(other));
275 // compare number of rows
277 return row < o.rows() ? -1 : 1;
279 // compare number of columns
281 return col < o.cols() ? -1 : 1;
283 // equal number of rows and columns, compare individual elements
285 for (unsigned r=0; r<row; ++r) {
286 for (unsigned c=0; c<col; ++c) {
287 cmpval = ((*this)(r,c)).compare(o(r,c));
288 if (cmpval!=0) return cmpval;
291 // all elements are equal => matrices are equal;
296 // non-virtual functions in this class
303 * @exception logic_error (incompatible matrices) */
304 matrix matrix::add(const matrix & other) const
306 if (col != other.col || row != other.row)
307 throw (std::logic_error("matrix::add(): incompatible matrices"));
309 exvector sum(this->m);
310 exvector::iterator i;
311 exvector::const_iterator ci;
312 for (i=sum.begin(), ci=other.m.begin(); i!=sum.end(); ++i, ++ci)
315 return matrix(row,col,sum);
319 /** Difference of matrices.
321 * @exception logic_error (incompatible matrices) */
322 matrix matrix::sub(const matrix & other) const
324 if (col != other.col || row != other.row)
325 throw (std::logic_error("matrix::sub(): incompatible matrices"));
327 exvector dif(this->m);
328 exvector::iterator i;
329 exvector::const_iterator ci;
330 for (i=dif.begin(), ci=other.m.begin(); i!=dif.end(); ++i, ++ci)
333 return matrix(row,col,dif);
337 /** Product of matrices.
339 * @exception logic_error (incompatible matrices) */
340 matrix matrix::mul(const matrix & other) const
342 if (this->cols() != other.rows())
343 throw (std::logic_error("matrix::mul(): incompatible matrices"));
345 exvector prod(this->rows()*other.cols());
347 for (unsigned r1=0; r1<this->rows(); ++r1) {
348 for (unsigned c=0; c<this->cols(); ++c) {
349 if (m[r1*col+c].is_zero())
351 for (unsigned r2=0; r2<other.cols(); ++r2)
352 prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]).expand();
355 return matrix(row, other.col, prod);
359 /** operator() to access elements.
361 * @param ro row of element
362 * @param co column of element
363 * @exception range_error (index out of range) */
364 const ex & matrix::operator() (unsigned ro, unsigned co) const
366 if (ro>=row || co>=col)
367 throw (std::range_error("matrix::operator(): index out of range"));
373 /** Set individual elements manually.
375 * @exception range_error (index out of range) */
376 matrix & matrix::set(unsigned ro, unsigned co, ex value)
378 if (ro>=row || co>=col)
379 throw (std::range_error("matrix::set(): index out of range"));
381 ensure_if_modifiable();
382 m[ro*col+co] = value;
387 /** Transposed of an m x n matrix, producing a new n x m matrix object that
388 * represents the transposed. */
389 matrix matrix::transpose(void) const
391 exvector trans(this->cols()*this->rows());
393 for (unsigned r=0; r<this->cols(); ++r)
394 for (unsigned c=0; c<this->rows(); ++c)
395 trans[r*this->rows()+c] = m[c*this->cols()+r];
397 return matrix(this->cols(),this->rows(),trans);
401 /** Determinant of square matrix. This routine doesn't actually calculate the
402 * determinant, it only implements some heuristics about which algorithm to
403 * run. If all the elements of the matrix are elements of an integral domain
404 * the determinant is also in that integral domain and the result is expanded
405 * only. If one or more elements are from a quotient field the determinant is
406 * usually also in that quotient field and the result is normalized before it
407 * is returned. This implies that the determinant of the symbolic 2x2 matrix
408 * [[a/(a-b),1],[b/(a-b),1]] is returned as unity. (In this respect, it
409 * behaves like MapleV and unlike Mathematica.)
411 * @param algo allows to chose an algorithm
412 * @return the determinant as a new expression
413 * @exception logic_error (matrix not square)
414 * @see determinant_algo */
415 ex matrix::determinant(unsigned algo) const
418 throw (std::logic_error("matrix::determinant(): matrix not square"));
419 GINAC_ASSERT(row*col==m.capacity());
421 // Gather some statistical information about this matrix:
422 bool numeric_flag = true;
423 bool normal_flag = false;
424 unsigned sparse_count = 0; // counts non-zero elements
425 for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r) {
426 lst srl; // symbol replacement list
427 ex rtest = (*r).to_rational(srl);
428 if (!rtest.is_zero())
430 if (!rtest.info(info_flags::numeric))
431 numeric_flag = false;
432 if (!rtest.info(info_flags::crational_polynomial) &&
433 rtest.info(info_flags::rational_function))
437 // Here is the heuristics in case this routine has to decide:
438 if (algo == determinant_algo::automatic) {
439 // Minor expansion is generally a good guess:
440 algo = determinant_algo::laplace;
441 // Does anybody know when a matrix is really sparse?
442 // Maybe <~row/2.236 nonzero elements average in a row?
443 if (row>3 && 5*sparse_count<=row*col)
444 algo = determinant_algo::bareiss;
445 // Purely numeric matrix can be handled by Gauss elimination.
446 // This overrides any prior decisions.
448 algo = determinant_algo::gauss;
451 // Trap the trivial case here, since some algorithms don't like it
453 // for consistency with non-trivial determinants...
455 return m[0].normal();
457 return m[0].expand();
460 // Compute the determinant
462 case determinant_algo::gauss: {
465 int sign = tmp.gauss_elimination(true);
466 for (unsigned d=0; d<row; ++d)
467 det *= tmp.m[d*col+d];
469 return (sign*det).normal();
471 return (sign*det).normal().expand();
473 case determinant_algo::bareiss: {
476 sign = tmp.fraction_free_elimination(true);
478 return (sign*tmp.m[row*col-1]).normal();
480 return (sign*tmp.m[row*col-1]).expand();
482 case determinant_algo::divfree: {
485 sign = tmp.division_free_elimination(true);
488 ex det = tmp.m[row*col-1];
489 // factor out accumulated bogus slag
490 for (unsigned d=0; d<row-2; ++d)
491 for (unsigned j=0; j<row-d-2; ++j)
492 det = (det/tmp.m[d*col+d]).normal();
495 case determinant_algo::laplace:
497 // This is the minor expansion scheme. We always develop such
498 // that the smallest minors (i.e, the trivial 1x1 ones) are on the
499 // rightmost column. For this to be efficient it turns out that
500 // the emptiest columns (i.e. the ones with most zeros) should be
501 // the ones on the right hand side. Therefore we presort the
502 // columns of the matrix:
503 typedef std::pair<unsigned,unsigned> uintpair;
504 std::vector<uintpair> c_zeros; // number of zeros in column
505 for (unsigned c=0; c<col; ++c) {
507 for (unsigned r=0; r<row; ++r)
508 if (m[r*col+c].is_zero())
510 c_zeros.push_back(uintpair(acc,c));
512 sort(c_zeros.begin(),c_zeros.end());
513 std::vector<unsigned> pre_sort;
514 for (std::vector<uintpair>::iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
515 pre_sort.push_back(i->second);
516 int sign = permutation_sign(pre_sort);
517 exvector result(row*col); // represents sorted matrix
519 for (std::vector<unsigned>::iterator i=pre_sort.begin();
522 for (unsigned r=0; r<row; ++r)
523 result[r*col+c] = m[r*col+(*i)];
527 return (sign*matrix(row,col,result).determinant_minor()).normal();
529 return sign*matrix(row,col,result).determinant_minor();
535 /** Trace of a matrix. The result is normalized if it is in some quotient
536 * field and expanded only otherwise. This implies that the trace of the
537 * symbolic 2x2 matrix [[a/(a-b),x],[y,b/(b-a)]] is recognized to be unity.
539 * @return the sum of diagonal elements
540 * @exception logic_error (matrix not square) */
541 ex matrix::trace(void) const
544 throw (std::logic_error("matrix::trace(): matrix not square"));
547 for (unsigned r=0; r<col; ++r)
550 if (tr.info(info_flags::rational_function) &&
551 !tr.info(info_flags::crational_polynomial))
558 /** Characteristic Polynomial. Following mathematica notation the
559 * characteristic polynomial of a matrix M is defined as the determiant of
560 * (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
561 * as M. Note that some CASs define it with a sign inside the determinant
562 * which gives rise to an overall sign if the dimension is odd. This method
563 * returns the characteristic polynomial collected in powers of lambda as a
566 * @return characteristic polynomial as new expression
567 * @exception logic_error (matrix not square)
568 * @see matrix::determinant() */
569 ex matrix::charpoly(const symbol & lambda) const
572 throw (std::logic_error("matrix::charpoly(): matrix not square"));
574 bool numeric_flag = true;
575 for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r) {
576 if (!(*r).info(info_flags::numeric)) {
577 numeric_flag = false;
581 // The pure numeric case is traditionally rather common. Hence, it is
582 // trapped and we use Leverrier's algorithm which goes as row^3 for
583 // every coefficient. The expensive part is the matrix multiplication.
587 ex poly = power(lambda,row)-c*power(lambda,row-1);
588 for (unsigned i=1; i<row; ++i) {
589 for (unsigned j=0; j<row; ++j)
592 c = B.trace()/ex(i+1);
593 poly -= c*power(lambda,row-i-1);
602 for (unsigned r=0; r<col; ++r)
603 M.m[r*col+r] -= lambda;
605 return M.determinant().collect(lambda);
609 /** Inverse of this matrix.
611 * @return the inverted matrix
612 * @exception logic_error (matrix not square)
613 * @exception runtime_error (singular matrix) */
614 matrix matrix::inverse(void) const
617 throw (std::logic_error("matrix::inverse(): matrix not square"));
619 // NOTE: the Gauss-Jordan elimination used here can in principle be
620 // replaced by two clever calls to gauss_elimination() and some to
621 // transpose(). Wouldn't be more efficient (maybe less?), just more
624 // set tmp to the unit matrix
625 for (unsigned i=0; i<col; ++i)
626 tmp.m[i*col+i] = _ex1();
628 // create a copy of this matrix
630 for (unsigned r1=0; r1<row; ++r1) {
631 int indx = cpy.pivot(r1, r1);
633 throw (std::runtime_error("matrix::inverse(): singular matrix"));
635 if (indx != 0) { // swap rows r and indx of matrix tmp
636 for (unsigned i=0; i<col; ++i)
637 tmp.m[r1*col+i].swap(tmp.m[indx*col+i]);
639 ex a1 = cpy.m[r1*col+r1];
640 for (unsigned c=0; c<col; ++c) {
641 cpy.m[r1*col+c] /= a1;
642 tmp.m[r1*col+c] /= a1;
644 for (unsigned r2=0; r2<row; ++r2) {
646 if (!cpy.m[r2*col+r1].is_zero()) {
647 ex a2 = cpy.m[r2*col+r1];
648 // yes, there is something to do in this column
649 for (unsigned c=0; c<col; ++c) {
650 cpy.m[r2*col+c] -= a2 * cpy.m[r1*col+c];
651 if (!cpy.m[r2*col+c].info(info_flags::numeric))
652 cpy.m[r2*col+c] = cpy.m[r2*col+c].normal();
653 tmp.m[r2*col+c] -= a2 * tmp.m[r1*col+c];
654 if (!tmp.m[r2*col+c].info(info_flags::numeric))
655 tmp.m[r2*col+c] = tmp.m[r2*col+c].normal();
666 /** Solve a linear system consisting of a m x n matrix and a m x p right hand
667 * side by applying an elimination scheme to the augmented matrix.
669 * @param vars n x p matrix, all elements must be symbols
670 * @param rhs m x p matrix
671 * @return n x p solution matrix
672 * @exception logic_error (incompatible matrices)
673 * @exception invalid_argument (1st argument must be matrix of symbols)
674 * @exception runtime_error (inconsistent linear system)
676 matrix matrix::solve(const matrix & vars,
680 const unsigned m = this->rows();
681 const unsigned n = this->cols();
682 const unsigned p = rhs.cols();
685 if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
686 throw (std::logic_error("matrix::solve(): incompatible matrices"));
687 for (unsigned ro=0; ro<n; ++ro)
688 for (unsigned co=0; co<p; ++co)
689 if (!vars(ro,co).info(info_flags::symbol))
690 throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
692 // build the augmented matrix of *this with rhs attached to the right
694 for (unsigned r=0; r<m; ++r) {
695 for (unsigned c=0; c<n; ++c)
696 aug.m[r*(n+p)+c] = this->m[r*n+c];
697 for (unsigned c=0; c<p; ++c)
698 aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
701 // Gather some statistical information about the augmented matrix:
702 bool numeric_flag = true;
703 for (exvector::const_iterator r=aug.m.begin(); r!=aug.m.end(); ++r) {
704 if (!(*r).info(info_flags::numeric))
705 numeric_flag = false;
708 // Here is the heuristics in case this routine has to decide:
709 if (algo == solve_algo::automatic) {
710 // Bareiss (fraction-free) elimination is generally a good guess:
711 algo = solve_algo::bareiss;
712 // For m<3, Bareiss elimination is equivalent to division free
713 // elimination but has more logistic overhead
715 algo = solve_algo::divfree;
716 // This overrides any prior decisions.
718 algo = solve_algo::gauss;
721 // Eliminate the augmented matrix:
723 case solve_algo::gauss:
724 aug.gauss_elimination();
725 case solve_algo::divfree:
726 aug.division_free_elimination();
727 case solve_algo::bareiss:
729 aug.fraction_free_elimination();
732 // assemble the solution matrix:
734 for (unsigned co=0; co<p; ++co) {
735 unsigned last_assigned_sol = n+1;
736 for (int r=m-1; r>=0; --r) {
737 unsigned fnz = 1; // first non-zero in row
738 while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
741 // row consists only of zeros, corresponding rhs must be 0, too
742 if (!aug.m[r*(n+p)+n+co].is_zero()) {
743 throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
746 // assign solutions for vars between fnz+1 and
747 // last_assigned_sol-1: free parameters
748 for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
749 sol.set(c,co,vars.m[c*p+co]);
750 ex e = aug.m[r*(n+p)+n+co];
751 for (unsigned c=fnz; c<n; ++c)
752 e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
754 (e/(aug.m[r*(n+p)+(fnz-1)])).normal());
755 last_assigned_sol = fnz;
758 // assign solutions for vars between 1 and
759 // last_assigned_sol-1: free parameters
760 for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
761 sol.set(ro,co,vars(ro,co));
770 /** Recursive determinant for small matrices having at least one symbolic
771 * entry. The basic algorithm, known as Laplace-expansion, is enhanced by
772 * some bookkeeping to avoid calculation of the same submatrices ("minors")
773 * more than once. According to W.M.Gentleman and S.C.Johnson this algorithm
774 * is better than elimination schemes for matrices of sparse multivariate
775 * polynomials and also for matrices of dense univariate polynomials if the
776 * matrix' dimesion is larger than 7.
778 * @return the determinant as a new expression (in expanded form)
779 * @see matrix::determinant() */
780 ex matrix::determinant_minor(void) const
782 // for small matrices the algorithm does not make any sense:
783 const unsigned n = this->cols();
785 return m[0].expand();
787 return (m[0]*m[3]-m[2]*m[1]).expand();
789 return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
790 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
791 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
793 // This algorithm can best be understood by looking at a naive
794 // implementation of Laplace-expansion, like this one:
796 // matrix minorM(this->rows()-1,this->cols()-1);
797 // for (unsigned r1=0; r1<this->rows(); ++r1) {
798 // // shortcut if element(r1,0) vanishes
799 // if (m[r1*col].is_zero())
801 // // assemble the minor matrix
802 // for (unsigned r=0; r<minorM.rows(); ++r) {
803 // for (unsigned c=0; c<minorM.cols(); ++c) {
805 // minorM.set(r,c,m[r*col+c+1]);
807 // minorM.set(r,c,m[(r+1)*col+c+1]);
810 // // recurse down and care for sign:
812 // det -= m[r1*col] * minorM.determinant_minor();
814 // det += m[r1*col] * minorM.determinant_minor();
816 // return det.expand();
817 // What happens is that while proceeding down many of the minors are
818 // computed more than once. In particular, there are binomial(n,k)
819 // kxk minors and each one is computed factorial(n-k) times. Therefore
820 // it is reasonable to store the results of the minors. We proceed from
821 // right to left. At each column c we only need to retrieve the minors
822 // calculated in step c-1. We therefore only have to store at most
823 // 2*binomial(n,n/2) minors.
825 // Unique flipper counter for partitioning into minors
826 std::vector<unsigned> Pkey;
828 // key for minor determinant (a subpartition of Pkey)
829 std::vector<unsigned> Mkey;
831 // we store our subminors in maps, keys being the rows they arise from
832 typedef std::map<std::vector<unsigned>,class ex> Rmap;
833 typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
837 // initialize A with last column:
838 for (unsigned r=0; r<n; ++r) {
839 Pkey.erase(Pkey.begin(),Pkey.end());
841 A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
843 // proceed from right to left through matrix
844 for (int c=n-2; c>=0; --c) {
845 Pkey.erase(Pkey.begin(),Pkey.end()); // don't change capacity
846 Mkey.erase(Mkey.begin(),Mkey.end());
847 for (unsigned i=0; i<n-c; ++i)
849 unsigned fc = 0; // controls logic for our strange flipper counter
852 for (unsigned r=0; r<n-c; ++r) {
853 // maybe there is nothing to do?
854 if (m[Pkey[r]*n+c].is_zero())
856 // create the sorted key for all possible minors
857 Mkey.erase(Mkey.begin(),Mkey.end());
858 for (unsigned i=0; i<n-c; ++i)
860 Mkey.push_back(Pkey[i]);
861 // Fetch the minors and compute the new determinant
863 det -= m[Pkey[r]*n+c]*A[Mkey];
865 det += m[Pkey[r]*n+c]*A[Mkey];
867 // prevent build-up of deep nesting of expressions saves time:
869 // store the new determinant at its place in B:
871 B.insert(Rmap_value(Pkey,det));
872 // increment our strange flipper counter
873 for (fc=n-c; fc>0; --fc) {
879 for (unsigned j=fc; j<n-c; ++j)
880 Pkey[j] = Pkey[j-1]+1;
882 // next column, so change the role of A and B:
891 /** Perform the steps of an ordinary Gaussian elimination to bring the m x n
892 * matrix into an upper echelon form. The algorithm is ok for matrices
893 * with numeric coefficients but quite unsuited for symbolic matrices.
895 * @param det may be set to true to save a lot of space if one is only
896 * interested in the diagonal elements (i.e. for calculating determinants).
897 * The others are set to zero in this case.
898 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
899 * number of rows was swapped and 0 if the matrix is singular. */
900 int matrix::gauss_elimination(const bool det)
902 ensure_if_modifiable();
903 const unsigned m = this->rows();
904 const unsigned n = this->cols();
905 GINAC_ASSERT(!det || n==m);
909 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
910 int indx = pivot(r0, r1, true);
914 return 0; // leaves *this in a messy state
919 for (unsigned r2=r0+1; r2<m; ++r2) {
920 if (!this->m[r2*n+r1].is_zero()) {
921 // yes, there is something to do in this row
922 ex piv = this->m[r2*n+r1] / this->m[r0*n+r1];
923 for (unsigned c=r1+1; c<n; ++c) {
924 this->m[r2*n+c] -= piv * this->m[r0*n+c];
925 if (!this->m[r2*n+c].info(info_flags::numeric))
926 this->m[r2*n+c] = this->m[r2*n+c].normal();
929 // fill up left hand side with zeros
930 for (unsigned c=0; c<=r1; ++c)
931 this->m[r2*n+c] = _ex0();
934 // save space by deleting no longer needed elements
935 for (unsigned c=r0+1; c<n; ++c)
936 this->m[r0*n+c] = _ex0();
946 /** Perform the steps of division free elimination to bring the m x n matrix
947 * into an upper echelon form.
949 * @param det may be set to true to save a lot of space if one is only
950 * interested in the diagonal elements (i.e. for calculating determinants).
951 * The others are set to zero in this case.
952 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
953 * number of rows was swapped and 0 if the matrix is singular. */
954 int matrix::division_free_elimination(const bool det)
956 ensure_if_modifiable();
957 const unsigned m = this->rows();
958 const unsigned n = this->cols();
959 GINAC_ASSERT(!det || n==m);
963 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
964 int indx = pivot(r0, r1, true);
968 return 0; // leaves *this in a messy state
973 for (unsigned r2=r0+1; r2<m; ++r2) {
974 for (unsigned c=r1+1; c<n; ++c)
975 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();
976 // fill up left hand side with zeros
977 for (unsigned c=0; c<=r1; ++c)
978 this->m[r2*n+c] = _ex0();
981 // save space by deleting no longer needed elements
982 for (unsigned c=r0+1; c<n; ++c)
983 this->m[r0*n+c] = _ex0();
993 /** Perform the steps of Bareiss' one-step fraction free elimination to bring
994 * the matrix into an upper echelon form. Fraction free elimination means
995 * that divide is used straightforwardly, without computing GCDs first. This
996 * is possible, since we know the divisor at each step.
998 * @param det may be set to true to save a lot of space if one is only
999 * interested in the last element (i.e. for calculating determinants). The
1000 * others are set to zero in this case.
1001 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1002 * number of rows was swapped and 0 if the matrix is singular. */
1003 int matrix::fraction_free_elimination(const bool det)
1006 // (single-step fraction free elimination scheme, already known to Jordan)
1008 // Usual division-free elimination sets m[0](r,c) = m(r,c) and then sets
1009 // m[k+1](r,c) = m[k](k,k) * m[k](r,c) - m[k](r,k) * m[k](k,c).
1011 // Bareiss (fraction-free) elimination in addition divides that element
1012 // by m[k-1](k-1,k-1) for k>1, where it can be shown by means of the
1013 // Sylvester determinant that this really divides m[k+1](r,c).
1015 // We also allow rational functions where the original prove still holds.
1016 // However, we must care for numerator and denominator separately and
1017 // "manually" work in the integral domains because of subtle cancellations
1018 // (see below). This blows up the bookkeeping a bit and the formula has
1019 // to be modified to expand like this (N{x} stands for numerator of x,
1020 // D{x} for denominator of x):
1021 // N{m[k+1](r,c)} = N{m[k](k,k)}*N{m[k](r,c)}*D{m[k](r,k)}*D{m[k](k,c)}
1022 // -N{m[k](r,k)}*N{m[k](k,c)}*D{m[k](k,k)}*D{m[k](r,c)}
1023 // D{m[k+1](r,c)} = D{m[k](k,k)}*D{m[k](r,c)}*D{m[k](r,k)}*D{m[k](k,c)}
1024 // where for k>1 we now divide N{m[k+1](r,c)} by
1025 // N{m[k-1](k-1,k-1)}
1026 // and D{m[k+1](r,c)} by
1027 // D{m[k-1](k-1,k-1)}.
1029 ensure_if_modifiable();
1030 const unsigned m = this->rows();
1031 const unsigned n = this->cols();
1032 GINAC_ASSERT(!det || n==m);
1041 // We populate temporary matrices to subsequently operate on. There is
1042 // one holding numerators and another holding denominators of entries.
1043 // This is a must since the evaluator (or even earlier mul's constructor)
1044 // might cancel some trivial element which causes divide() to fail. The
1045 // elements are normalized first (yes, even though this algorithm doesn't
1046 // need GCDs) since the elements of *this might be unnormalized, which
1047 // makes things more complicated than they need to be.
1048 matrix tmp_n(*this);
1049 matrix tmp_d(m,n); // for denominators, if needed
1050 lst srl; // symbol replacement list
1051 exvector::iterator it = this->m.begin();
1052 exvector::iterator tmp_n_it = tmp_n.m.begin();
1053 exvector::iterator tmp_d_it = tmp_d.m.begin();
1054 for (; it!= this->m.end(); ++it, ++tmp_n_it, ++tmp_d_it) {
1055 (*tmp_n_it) = (*it).normal().to_rational(srl);
1056 (*tmp_d_it) = (*tmp_n_it).denom();
1057 (*tmp_n_it) = (*tmp_n_it).numer();
1061 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
1062 int indx = tmp_n.pivot(r0, r1, true);
1071 // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
1072 for (unsigned c=r1; c<n; ++c)
1073 tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
1075 for (unsigned r2=r0+1; r2<m; ++r2) {
1076 for (unsigned c=r1+1; c<n; ++c) {
1077 dividend_n = (tmp_n.m[r0*n+r1]*tmp_n.m[r2*n+c]*
1078 tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]
1079 -tmp_n.m[r2*n+r1]*tmp_n.m[r0*n+c]*
1080 tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
1081 dividend_d = (tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]*
1082 tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
1083 bool check = divide(dividend_n, divisor_n,
1084 tmp_n.m[r2*n+c], true);
1085 check &= divide(dividend_d, divisor_d,
1086 tmp_d.m[r2*n+c], true);
1087 GINAC_ASSERT(check);
1089 // fill up left hand side with zeros
1090 for (unsigned c=0; c<=r1; ++c)
1091 tmp_n.m[r2*n+c] = _ex0();
1093 if ((r1<n-1)&&(r0<m-1)) {
1094 // compute next iteration's divisor
1095 divisor_n = tmp_n.m[r0*n+r1].expand();
1096 divisor_d = tmp_d.m[r0*n+r1].expand();
1098 // save space by deleting no longer needed elements
1099 for (unsigned c=0; c<n; ++c) {
1100 tmp_n.m[r0*n+c] = _ex0();
1101 tmp_d.m[r0*n+c] = _ex1();
1108 // repopulate *this matrix:
1109 it = this->m.begin();
1110 tmp_n_it = tmp_n.m.begin();
1111 tmp_d_it = tmp_d.m.begin();
1112 for (; it!= this->m.end(); ++it, ++tmp_n_it, ++tmp_d_it)
1113 (*it) = ((*tmp_n_it)/(*tmp_d_it)).subs(srl);
1119 /** Partial pivoting method for matrix elimination schemes.
1120 * Usual pivoting (symbolic==false) returns the index to the element with the
1121 * largest absolute value in column ro and swaps the current row with the one
1122 * where the element was found. With (symbolic==true) it does the same thing
1123 * with the first non-zero element.
1125 * @param ro is the row from where to begin
1126 * @param co is the column to be inspected
1127 * @param symbolic signal if we want the first non-zero element to be pivoted
1128 * (true) or the one with the largest absolute value (false).
1129 * @return 0 if no interchange occured, -1 if all are zero (usually signaling
1130 * a degeneracy) and positive integer k means that rows ro and k were swapped.
1132 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
1136 // search first non-zero element in column co beginning at row ro
1137 while ((k<row) && (this->m[k*col+co].expand().is_zero()))
1140 // search largest element in column co beginning at row ro
1141 GINAC_ASSERT(is_ex_of_type(this->m[k*col+co],numeric));
1142 unsigned kmax = k+1;
1143 numeric mmax = abs(ex_to_numeric(m[kmax*col+co]));
1145 GINAC_ASSERT(is_ex_of_type(this->m[kmax*col+co],numeric));
1146 numeric tmp = ex_to_numeric(this->m[kmax*col+co]);
1147 if (abs(tmp) > mmax) {
1153 if (!mmax.is_zero())
1157 // all elements in column co below row ro vanish
1160 // matrix needs no pivoting
1162 // matrix needs pivoting, so swap rows k and ro
1163 ensure_if_modifiable();
1164 for (unsigned c=0; c<col; ++c)
1165 this->m[k*col+c].swap(this->m[ro*col+c]);
1170 /** Convert list of lists to matrix. */
1171 ex lst_to_matrix(const ex &l)
1173 if (!is_ex_of_type(l, lst))
1174 throw(std::invalid_argument("argument to lst_to_matrix() must be a lst"));
1176 // Find number of rows and columns
1177 unsigned rows = l.nops(), cols = 0, i, j;
1178 for (i=0; i<rows; i++)
1179 if (l.op(i).nops() > cols)
1180 cols = l.op(i).nops();
1182 // Allocate and fill matrix
1183 matrix &m = *new matrix(rows, cols);
1184 for (i=0; i<rows; i++)
1185 for (j=0; j<cols; j++)
1186 if (l.op(i).nops() > j)
1187 m.set(i, j, l.op(i).op(j));
1193 } // namespace GiNaC