]> www.ginac.de Git - ginac.git/blob - ginac/matrix.cpp
6240ce20e9e36a14557e38a9918841ecf63923c7
[ginac.git] / ginac / matrix.cpp
1 /** @file matrix.cpp
2  *
3  *  Implementation of symbolic matrices */
4
5 /*
6  *  GiNaC Copyright (C) 1999-2000 Johannes Gutenberg University Mainz, Germany
7  *
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.
12  *
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.
17  *
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
21  */
22
23 #include <algorithm>
24 #include <map>
25 #include <stdexcept>
26
27 #include "matrix.h"
28 #include "archive.h"
29 #include "numeric.h"
30 #include "lst.h"
31 #include "utils.h"
32 #include "debugmsg.h"
33 #include "power.h"
34 #include "symbol.h"
35 #include "normal.h"
36
37 #ifndef NO_NAMESPACE_GINAC
38 namespace GiNaC {
39 #endif // ndef NO_NAMESPACE_GINAC
40
41 GINAC_IMPLEMENT_REGISTERED_CLASS(matrix, basic)
42
43 //////////
44 // default constructor, destructor, copy constructor, assignment operator
45 // and helpers:
46 //////////
47
48 // public
49
50 /** Default ctor.  Initializes to 1 x 1-dimensional zero-matrix. */
51 matrix::matrix()
52     : inherited(TINFO_matrix), row(1), col(1)
53 {
54     debugmsg("matrix default constructor",LOGLEVEL_CONSTRUCT);
55     m.push_back(_ex0());
56 }
57
58 matrix::~matrix()
59 {
60     debugmsg("matrix destructor",LOGLEVEL_DESTRUCT);
61 }
62
63 matrix::matrix(const matrix & other)
64 {
65     debugmsg("matrix copy constructor",LOGLEVEL_CONSTRUCT);
66     copy(other);
67 }
68
69 const matrix & matrix::operator=(const matrix & other)
70 {
71     debugmsg("matrix operator=",LOGLEVEL_ASSIGNMENT);
72     if (this != &other) {
73         destroy(1);
74         copy(other);
75     }
76     return *this;
77 }
78
79 // protected
80
81 void matrix::copy(const matrix & other)
82 {
83     inherited::copy(other);
84     row = other.row;
85     col = other.col;
86     m = other.m;  // STL's vector copying invoked here
87 }
88
89 void matrix::destroy(bool call_parent)
90 {
91     if (call_parent) inherited::destroy(call_parent);
92 }
93
94 //////////
95 // other constructors
96 //////////
97
98 // public
99
100 /** Very common ctor.  Initializes to r x c-dimensional zero-matrix.
101  *
102  *  @param r number of rows
103  *  @param c number of cols */
104 matrix::matrix(unsigned r, unsigned c)
105     : inherited(TINFO_matrix), row(r), col(c)
106 {
107     debugmsg("matrix constructor from unsigned,unsigned",LOGLEVEL_CONSTRUCT);
108     m.resize(r*c, _ex0());
109 }
110
111 // protected
112
113 /** Ctor from representation, for internal use only. */
114 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
115     : inherited(TINFO_matrix), row(r), col(c), m(m2)
116 {
117     debugmsg("matrix constructor from unsigned,unsigned,exvector",LOGLEVEL_CONSTRUCT);
118 }
119
120 //////////
121 // archiving
122 //////////
123
124 /** Construct object from archive_node. */
125 matrix::matrix(const archive_node &n, const lst &sym_lst) : inherited(n, sym_lst)
126 {
127     debugmsg("matrix constructor from archive_node", LOGLEVEL_CONSTRUCT);
128     if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
129         throw (std::runtime_error("unknown matrix dimensions in archive"));
130     m.reserve(row * col);
131     for (unsigned int i=0; true; i++) {
132         ex e;
133         if (n.find_ex("m", e, sym_lst, i))
134             m.push_back(e);
135         else
136             break;
137     }
138 }
139
140 /** Unarchive the object. */
141 ex matrix::unarchive(const archive_node &n, const lst &sym_lst)
142 {
143     return (new matrix(n, sym_lst))->setflag(status_flags::dynallocated);
144 }
145
146 /** Archive the object. */
147 void matrix::archive(archive_node &n) const
148 {
149     inherited::archive(n);
150     n.add_unsigned("row", row);
151     n.add_unsigned("col", col);
152     exvector::const_iterator i = m.begin(), iend = m.end();
153     while (i != iend) {
154         n.add_ex("m", *i);
155         ++i;
156     }
157 }
158
159 //////////
160 // functions overriding virtual functions from bases classes
161 //////////
162
163 // public
164
165 basic * matrix::duplicate() const
166 {
167     debugmsg("matrix duplicate",LOGLEVEL_DUPLICATE);
168     return new matrix(*this);
169 }
170
171 void matrix::print(std::ostream & os, unsigned upper_precedence) const
172 {
173     debugmsg("matrix print",LOGLEVEL_PRINT);
174     os << "[[ ";
175     for (unsigned r=0; r<row-1; ++r) {
176         os << "[[";
177         for (unsigned c=0; c<col-1; ++c)
178             os << m[r*col+c] << ",";
179         os << m[col*(r+1)-1] << "]], ";
180     }
181     os << "[[";
182     for (unsigned c=0; c<col-1; ++c)
183         os << m[(row-1)*col+c] << ",";
184     os << m[row*col-1] << "]] ]]";
185 }
186
187 void matrix::printraw(std::ostream & os) const
188 {
189     debugmsg("matrix printraw",LOGLEVEL_PRINT);
190     os << "matrix(" << row << "," << col <<",";
191     for (unsigned r=0; r<row-1; ++r) {
192         os << "(";
193         for (unsigned c=0; c<col-1; ++c)
194             os << m[r*col+c] << ",";
195         os << m[col*(r-1)-1] << "),";
196     }
197     os << "(";
198     for (unsigned c=0; c<col-1; ++c)
199         os << m[(row-1)*col+c] << ",";
200     os << m[row*col-1] << "))";
201 }
202
203 /** nops is defined to be rows x columns. */
204 unsigned matrix::nops() const
205 {
206     return row*col;
207 }
208
209 /** returns matrix entry at position (i/col, i%col). */
210 ex matrix::op(int i) const
211 {
212     return m[i];
213 }
214
215 /** returns matrix entry at position (i/col, i%col). */
216 ex & matrix::let_op(int i)
217 {
218     GINAC_ASSERT(i>=0);
219     GINAC_ASSERT(i<nops());
220     
221     return m[i];
222 }
223
224 /** expands the elements of a matrix entry by entry. */
225 ex matrix::expand(unsigned options) const
226 {
227     exvector tmp(row*col);
228     for (unsigned i=0; i<row*col; ++i)
229         tmp[i] = m[i].expand(options);
230     
231     return matrix(row, col, tmp);
232 }
233
234 /** Search ocurrences.  A matrix 'has' an expression if it is the expression
235  *  itself or one of the elements 'has' it. */
236 bool matrix::has(const ex & other) const
237 {
238     GINAC_ASSERT(other.bp!=0);
239     
240     // tautology: it is the expression itself
241     if (is_equal(*other.bp)) return true;
242     
243     // search all the elements
244     for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r)
245         if ((*r).has(other)) return true;
246     
247     return false;
248 }
249
250 /** evaluate matrix entry by entry. */
251 ex matrix::eval(int level) const
252 {
253     debugmsg("matrix eval",LOGLEVEL_MEMBER_FUNCTION);
254     
255     // check if we have to do anything at all
256     if ((level==1)&&(flags & status_flags::evaluated))
257         return *this;
258     
259     // emergency break
260     if (level == -max_recursion_level)
261         throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
262     
263     // eval() entry by entry
264     exvector m2(row*col);
265     --level;
266     for (unsigned r=0; r<row; ++r)
267         for (unsigned c=0; c<col; ++c)
268             m2[r*col+c] = m[r*col+c].eval(level);
269     
270     return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
271                                                status_flags::evaluated );
272 }
273
274 /** evaluate matrix numerically entry by entry. */
275 ex matrix::evalf(int level) const
276 {
277     debugmsg("matrix evalf",LOGLEVEL_MEMBER_FUNCTION);
278         
279     // check if we have to do anything at all
280     if (level==1)
281         return *this;
282     
283     // emergency break
284     if (level == -max_recursion_level) {
285         throw (std::runtime_error("matrix::evalf(): recursion limit exceeded"));
286     }
287     
288     // evalf() entry by entry
289     exvector m2(row*col);
290     --level;
291     for (unsigned r=0; r<row; ++r)
292         for (unsigned c=0; c<col; ++c)
293             m2[r*col+c] = m[r*col+c].evalf(level);
294     
295     return matrix(row, col, m2);
296 }
297
298 // protected
299
300 int matrix::compare_same_type(const basic & other) const
301 {
302     GINAC_ASSERT(is_exactly_of_type(other, matrix));
303     const matrix & o = static_cast<matrix &>(const_cast<basic &>(other));
304     
305     // compare number of rows
306     if (row != o.rows())
307         return row < o.rows() ? -1 : 1;
308     
309     // compare number of columns
310     if (col != o.cols())
311         return col < o.cols() ? -1 : 1;
312     
313     // equal number of rows and columns, compare individual elements
314     int cmpval;
315     for (unsigned r=0; r<row; ++r) {
316         for (unsigned c=0; c<col; ++c) {
317             cmpval = ((*this)(r,c)).compare(o(r,c));
318             if (cmpval!=0) return cmpval;
319         }
320     }
321     // all elements are equal => matrices are equal;
322     return 0;
323 }
324
325 //////////
326 // non-virtual functions in this class
327 //////////
328
329 // public
330
331 /** Sum of matrices.
332  *
333  *  @exception logic_error (incompatible matrices) */
334 matrix matrix::add(const matrix & other) const
335 {
336     if (col != other.col || row != other.row)
337         throw (std::logic_error("matrix::add(): incompatible matrices"));
338     
339     exvector sum(this->m);
340     exvector::iterator i;
341     exvector::const_iterator ci;
342     for (i=sum.begin(), ci=other.m.begin(); i!=sum.end(); ++i, ++ci)
343         (*i) += (*ci);
344     
345     return matrix(row,col,sum);
346 }
347
348
349 /** Difference of matrices.
350  *
351  *  @exception logic_error (incompatible matrices) */
352 matrix matrix::sub(const matrix & other) const
353 {
354     if (col != other.col || row != other.row)
355         throw (std::logic_error("matrix::sub(): incompatible matrices"));
356     
357     exvector dif(this->m);
358     exvector::iterator i;
359     exvector::const_iterator ci;
360     for (i=dif.begin(), ci=other.m.begin(); i!=dif.end(); ++i, ++ci)
361         (*i) -= (*ci);
362     
363     return matrix(row,col,dif);
364 }
365
366
367 /** Product of matrices.
368  *
369  *  @exception logic_error (incompatible matrices) */
370 matrix matrix::mul(const matrix & other) const
371 {
372     if (this->cols() != other.rows())
373         throw (std::logic_error("matrix::mul(): incompatible matrices"));
374     
375     exvector prod(this->rows()*other.cols());
376     
377     for (unsigned r1=0; r1<this->rows(); ++r1) {
378         for (unsigned c=0; c<this->cols(); ++c) {
379             if (m[r1*col+c].is_zero())
380                 continue;
381             for (unsigned r2=0; r2<other.cols(); ++r2)
382                 prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]).expand();
383         }
384     }
385     return matrix(row, other.col, prod);
386 }
387
388
389 /** operator() to access elements.
390  *
391  *  @param ro row of element
392  *  @param co column of element
393  *  @exception range_error (index out of range) */
394 const ex & matrix::operator() (unsigned ro, unsigned co) const
395 {
396     if (ro<0 || ro>=row || co<0 || co>=col)
397         throw (std::range_error("matrix::operator(): index out of range"));
398
399     return m[ro*col+co];
400 }
401
402
403 /** Set individual elements manually.
404  *
405  *  @exception range_error (index out of range) */
406 matrix & matrix::set(unsigned ro, unsigned co, ex value)
407 {
408     if (ro<0 || ro>=row || co<0 || co>=col)
409         throw (std::range_error("matrix::set(): index out of range"));
410     
411     ensure_if_modifiable();
412     m[ro*col+co] = value;
413     return *this;
414 }
415
416
417 /** Transposed of an m x n matrix, producing a new n x m matrix object that
418  *  represents the transposed. */
419 matrix matrix::transpose(void) const
420 {
421     exvector trans(this->cols()*this->rows());
422     
423     for (unsigned r=0; r<this->cols(); ++r)
424         for (unsigned c=0; c<this->rows(); ++c)
425             trans[r*this->rows()+c] = m[c*this->cols()+r];
426     
427     return matrix(this->cols(),this->rows(),trans);
428 }
429
430
431 /** Determinant of square matrix.  This routine doesn't actually calculate the
432  *  determinant, it only implements some heuristics about which algorithm to
433  *  run.  If all the elements of the matrix are elements of an integral domain
434  *  the determinant is also in that integral domain and the result is expanded
435  *  only.  If one or more elements are from a quotient field the determinant is
436  *  usually also in that quotient field and the result is normalized before it
437  *  is returned.  This implies that the determinant of the symbolic 2x2 matrix
438  *  [[a/(a-b),1],[b/(a-b),1]] is returned as unity.  (In this respect, it
439  *  behaves like MapleV and unlike Mathematica.)
440  *
441  *  @param     algo allows to chose an algorithm
442  *  @return    the determinant as a new expression
443  *  @exception logic_error (matrix not square)
444  *  @see       determinant_algo */
445 ex matrix::determinant(unsigned algo) const
446 {
447     if (row!=col)
448         throw (std::logic_error("matrix::determinant(): matrix not square"));
449     GINAC_ASSERT(row*col==m.capacity());
450     
451     // Gather some statistical information about this matrix:
452     bool numeric_flag = true;
453     bool normal_flag = false;
454     unsigned sparse_count = 0;  // counts non-zero elements
455     for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r) {
456         lst srl;  // symbol replacement list
457         ex rtest = (*r).to_rational(srl);
458         if (!rtest.is_zero())
459             ++sparse_count;
460         if (!rtest.info(info_flags::numeric))
461             numeric_flag = false;
462         if (!rtest.info(info_flags::crational_polynomial) &&
463              rtest.info(info_flags::rational_function))
464             normal_flag = true;
465     }
466     
467     // Here is the heuristics in case this routine has to decide:
468     if (algo == determinant_algo::automatic) {
469         // Minor expansion is generally a good guess:
470         algo = determinant_algo::laplace;
471         // Does anybody know when a matrix is really sparse?
472         // Maybe <~row/2.236 nonzero elements average in a row?
473         if (row>3 && 5*sparse_count<=row*col)
474             algo = determinant_algo::bareiss;
475         // Purely numeric matrix can be handled by Gauss elimination.
476         // This overrides any prior decisions.
477         if (numeric_flag)
478             algo = determinant_algo::gauss;
479     }
480     
481     // Trap the trivial case here, since some algorithms don't like it
482     if (this->row==1) {
483         // for consistency with non-trivial determinants...
484         if (normal_flag)
485             return m[0].normal();
486         else
487             return m[0].expand();
488     }
489     
490     // Compute the determinant
491     switch(algo) {
492         case determinant_algo::gauss: {
493             ex det = 1;
494             matrix tmp(*this);
495             int sign = tmp.gauss_elimination(true);
496             for (unsigned d=0; d<row; ++d)
497                 det *= tmp.m[d*col+d];
498             if (normal_flag)
499                 return (sign*det).normal();
500             else
501                 return (sign*det).normal().expand();
502         }
503         case determinant_algo::bareiss: {
504             matrix tmp(*this);
505             int sign;
506             sign = tmp.fraction_free_elimination(true);
507             if (normal_flag)
508                 return (sign*tmp.m[row*col-1]).normal();
509             else
510                 return (sign*tmp.m[row*col-1]).expand();
511         }
512         case determinant_algo::divfree: {
513             matrix tmp(*this);
514             int sign;
515             sign = tmp.division_free_elimination(true);
516             if (sign==0)
517                 return _ex0();
518             ex det = tmp.m[row*col-1];
519             // factor out accumulated bogus slag
520             for (unsigned d=0; d<row-2; ++d)
521                 for (unsigned j=0; j<row-d-2; ++j)
522                     det = (det/tmp.m[d*col+d]).normal();
523             return (sign*det);
524         }
525         case determinant_algo::laplace:
526         default: {
527             // This is the minor expansion scheme.  We always develop such
528             // that the smallest minors (i.e, the trivial 1x1 ones) are on the
529             // rightmost column.  For this to be efficient it turns out that
530             // the emptiest columns (i.e. the ones with most zeros) should be
531             // the ones on the right hand side.  Therefore we presort the
532             // columns of the matrix:
533             typedef std::pair<unsigned,unsigned> uintpair;
534             std::vector<uintpair> c_zeros;  // number of zeros in column
535             for (unsigned c=0; c<col; ++c) {
536                 unsigned acc = 0;
537                 for (unsigned r=0; r<row; ++r)
538                     if (m[r*col+c].is_zero())
539                         ++acc;
540                 c_zeros.push_back(uintpair(acc,c));
541             }
542             sort(c_zeros.begin(),c_zeros.end());
543             std::vector<unsigned> pre_sort;
544             for (std::vector<uintpair>::iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
545                 pre_sort.push_back(i->second);
546             int sign = permutation_sign(pre_sort);
547             exvector result(row*col);  // represents sorted matrix
548             unsigned c = 0;
549             for (std::vector<unsigned>::iterator i=pre_sort.begin();
550                  i!=pre_sort.end();
551                  ++i,++c) {
552                 for (unsigned r=0; r<row; ++r)
553                     result[r*col+c] = m[r*col+(*i)];
554             }
555             
556             if (normal_flag)
557                 return (sign*matrix(row,col,result).determinant_minor()).normal();
558             else
559                 return sign*matrix(row,col,result).determinant_minor();
560         }
561     }
562 }
563
564
565 /** Trace of a matrix.  The result is normalized if it is in some quotient
566  *  field and expanded only otherwise.  This implies that the trace of the
567  *  symbolic 2x2 matrix [[a/(a-b),x],[y,b/(b-a)]] is recognized to be unity.
568  *
569  *  @return    the sum of diagonal elements
570  *  @exception logic_error (matrix not square) */
571 ex matrix::trace(void) const
572 {
573     if (row != col)
574         throw (std::logic_error("matrix::trace(): matrix not square"));
575     
576     ex tr;
577     for (unsigned r=0; r<col; ++r)
578         tr += m[r*col+r];
579     
580     if (tr.info(info_flags::rational_function) &&
581         !tr.info(info_flags::crational_polynomial))
582         return tr.normal();
583     else
584         return tr.expand();
585 }
586
587
588 /** Characteristic Polynomial.  Following mathematica notation the
589  *  characteristic polynomial of a matrix M is defined as the determiant of
590  *  (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
591  *  as M.  Note that some CASs define it with a sign inside the determinant
592  *  which gives rise to an overall sign if the dimension is odd.  This method
593  *  returns the characteristic polynomial collected in powers of lambda as a
594  *  new expression.
595  *
596  *  @return    characteristic polynomial as new expression
597  *  @exception logic_error (matrix not square)
598  *  @see       matrix::determinant() */
599 ex matrix::charpoly(const symbol & lambda) const
600 {
601     if (row != col)
602         throw (std::logic_error("matrix::charpoly(): matrix not square"));
603     
604     bool numeric_flag = true;
605     for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r) {
606         if (!(*r).info(info_flags::numeric)) {
607             numeric_flag = false;
608         }
609     }
610     
611     // The pure numeric case is traditionally rather common.  Hence, it is
612     // trapped and we use Leverrier's algorithm which goes as row^3 for
613     // every coefficient.  The expensive part is the matrix multiplication.
614     if (numeric_flag) {
615         matrix B(*this);
616         ex c = B.trace();
617         ex poly = power(lambda,row)-c*power(lambda,row-1);
618         for (unsigned i=1; i<row; ++i) {
619             for (unsigned j=0; j<row; ++j)
620                 B.m[j*col+j] -= c;
621             B = this->mul(B);
622             c = B.trace()/ex(i+1);
623             poly -= c*power(lambda,row-i-1);
624         }
625         if (row%2)
626             return -poly;
627         else
628             return poly;
629     }
630     
631     matrix M(*this);
632     for (unsigned r=0; r<col; ++r)
633         M.m[r*col+r] -= lambda;
634     
635     return M.determinant().collect(lambda);
636 }
637
638
639 /** Inverse of this matrix.
640  *
641  *  @return    the inverted matrix
642  *  @exception logic_error (matrix not square)
643  *  @exception runtime_error (singular matrix) */
644 matrix matrix::inverse(void) const
645 {
646     if (row != col)
647         throw (std::logic_error("matrix::inverse(): matrix not square"));
648     
649     // NOTE: the Gauss-Jordan elimination used here can in principle be
650     // replaced this by two clever calls to gauss_elimination() and some to
651     // transpose().  Wouldn't be more efficient (maybe less?), just more
652     // orthogonal.
653     matrix tmp(row,col);
654     // set tmp to the unit matrix
655     for (unsigned i=0; i<col; ++i)
656         tmp.m[i*col+i] = _ex1();
657     
658     // create a copy of this matrix
659     matrix cpy(*this);
660     for (unsigned r1=0; r1<row; ++r1) {
661         int indx = cpy.pivot(r1, r1);
662         if (indx == -1) {
663             throw (std::runtime_error("matrix::inverse(): singular matrix"));
664         }
665         if (indx != 0) {  // swap rows r and indx of matrix tmp
666             for (unsigned i=0; i<col; ++i)
667                 tmp.m[r1*col+i].swap(tmp.m[indx*col+i]);
668         }
669         ex a1 = cpy.m[r1*col+r1];
670         for (unsigned c=0; c<col; ++c) {
671             cpy.m[r1*col+c] /= a1;
672             tmp.m[r1*col+c] /= a1;
673         }
674         for (unsigned r2=0; r2<row; ++r2) {
675             if (r2 != r1) {
676                 ex a2 = cpy.m[r2*col+r1];
677                 for (unsigned c=0; c<col; ++c) {
678                     cpy.m[r2*col+c] -= a2 * cpy.m[r1*col+c];
679                     if (!cpy.m[r2*col+c].info(info_flags::numeric))
680                         cpy.m[r2*col+c] = cpy.m[r2*col+c].normal();
681                     tmp.m[r2*col+c] -= a2 * tmp.m[r1*col+c];
682                     if (!tmp.m[r2*col+c].info(info_flags::numeric))
683                         tmp.m[r2*col+c] = tmp.m[r2*col+c].normal();
684                 }
685             }
686         }
687     }
688     
689     return tmp;
690 }
691
692
693 /** Solve a linear system consisting of a m x n matrix and a m x p right hand
694  *  side by applying an elimination scheme to the augmented matrix.
695  *
696  *  @param vars n x p matrix, all elements must be symbols 
697  *  @param rhs m x p matrix
698  *  @return n x p solution matrix
699  *  @exception logic_error (incompatible matrices)
700  *  @exception invalid_argument (1st argument must be matrix of symbols)
701  *  @exception runtime_error (inconsistent linear system)
702  *  @see       solve_algo */
703 matrix matrix::solve(const matrix & vars,
704                      const matrix & rhs,
705                      unsigned algo) const
706 {
707     const unsigned m = this->rows();
708     const unsigned n = this->cols();
709     const unsigned p = rhs.cols();
710     
711     // syntax checks    
712     if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
713         throw (std::logic_error("matrix::solve(): incompatible matrices"));
714     for (unsigned ro=0; ro<n; ++ro)
715         for (unsigned co=0; co<p; ++co)
716             if (!vars(ro,co).info(info_flags::symbol))
717                 throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
718     
719     // build the augmented matrix of *this with rhs attached to the right
720     matrix aug(m,n+p);
721     for (unsigned r=0; r<m; ++r) {
722         for (unsigned c=0; c<n; ++c)
723             aug.m[r*(n+p)+c] = this->m[r*n+c];
724         for (unsigned c=0; c<p; ++c)
725             aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
726     }
727     
728     // Gather some statistical information about the augmented matrix:
729     bool numeric_flag = true;
730     for (exvector::const_iterator r=aug.m.begin(); r!=aug.m.end(); ++r) {
731         if (!(*r).info(info_flags::numeric))
732             numeric_flag = false;
733     }
734     
735     // Here is the heuristics in case this routine has to decide:
736     if (algo == solve_algo::automatic) {
737         // Bareiss (fraction-free) elimination is generally a good guess:
738         algo = solve_algo::bareiss;
739         // For m<3, Bareiss elimination is equivalent to division free
740         // elimination but has more logistic overhead
741         if (m<3)
742             algo = solve_algo::divfree;
743         // This overrides any prior decisions.
744         if (numeric_flag)
745             algo = solve_algo::gauss;
746     }
747     
748     // Eliminate the augmented matrix:
749     switch(algo) {
750         case solve_algo::gauss:
751             aug.gauss_elimination();
752         case solve_algo::divfree:
753             aug.division_free_elimination();
754         case solve_algo::bareiss:
755         default:
756             aug.fraction_free_elimination();
757     }
758     
759     // assemble the solution matrix:
760     matrix sol(n,p);
761     for (unsigned co=0; co<p; ++co) {
762         unsigned last_assigned_sol = n+1;
763         for (int r=m-1; r>=0; --r) {
764             unsigned fnz = 1;    // first non-zero in row
765             while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
766                 ++fnz;
767             if (fnz>n) {
768                 // row consists only of zeros, corresponding rhs must be 0, too
769                 if (!aug.m[r*(n+p)+n+co].is_zero()) {
770                     throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
771                 }
772             } else {
773                 // assign solutions for vars between fnz+1 and
774                 // last_assigned_sol-1: free parameters
775                 for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
776                     sol.set(c,co,vars.m[c*p+co]);
777                 ex e = aug.m[r*(n+p)+n+co];
778                 for (unsigned c=fnz; c<n; ++c)
779                     e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
780                 sol.set(fnz-1,co,
781                         (e/(aug.m[r*(n+p)+(fnz-1)])).normal());
782                 last_assigned_sol = fnz;
783             }
784         }
785         // assign solutions for vars between 1 and
786         // last_assigned_sol-1: free parameters
787         for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
788             sol.set(ro,co,vars(ro,co));
789     }
790     
791     return sol;
792 }
793
794
795 // protected
796
797 /** Recursive determinant for small matrices having at least one symbolic
798  *  entry.  The basic algorithm, known as Laplace-expansion, is enhanced by
799  *  some bookkeeping to avoid calculation of the same submatrices ("minors")
800  *  more than once.  According to W.M.Gentleman and S.C.Johnson this algorithm
801  *  is better than elimination schemes for matrices of sparse multivariate
802  *  polynomials and also for matrices of dense univariate polynomials if the
803  *  matrix' dimesion is larger than 7.
804  *
805  *  @return the determinant as a new expression (in expanded form)
806  *  @see matrix::determinant() */
807 ex matrix::determinant_minor(void) const
808 {
809     // for small matrices the algorithm does not make any sense:
810     const unsigned n = this->cols();
811     if (n==1)
812         return m[0].expand();
813     if (n==2)
814         return (m[0]*m[3]-m[2]*m[1]).expand();
815     if (n==3)
816         return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
817                 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
818                 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
819     
820     // This algorithm can best be understood by looking at a naive
821     // implementation of Laplace-expansion, like this one:
822     // ex det;
823     // matrix minorM(this->rows()-1,this->cols()-1);
824     // for (unsigned r1=0; r1<this->rows(); ++r1) {
825     //     // shortcut if element(r1,0) vanishes
826     //     if (m[r1*col].is_zero())
827     //         continue;
828     //     // assemble the minor matrix
829     //     for (unsigned r=0; r<minorM.rows(); ++r) {
830     //         for (unsigned c=0; c<minorM.cols(); ++c) {
831     //             if (r<r1)
832     //                 minorM.set(r,c,m[r*col+c+1]);
833     //             else
834     //                 minorM.set(r,c,m[(r+1)*col+c+1]);
835     //         }
836     //     }
837     //     // recurse down and care for sign:
838     //     if (r1%2)
839     //         det -= m[r1*col] * minorM.determinant_minor();
840     //     else
841     //         det += m[r1*col] * minorM.determinant_minor();
842     // }
843     // return det.expand();
844     // What happens is that while proceeding down many of the minors are
845     // computed more than once.  In particular, there are binomial(n,k)
846     // kxk minors and each one is computed factorial(n-k) times.  Therefore
847     // it is reasonable to store the results of the minors.  We proceed from
848     // right to left.  At each column c we only need to retrieve the minors
849     // calculated in step c-1.  We therefore only have to store at most 
850     // 2*binomial(n,n/2) minors.
851     
852     // Unique flipper counter for partitioning into minors
853     std::vector<unsigned> Pkey;
854     Pkey.reserve(n);
855     // key for minor determinant (a subpartition of Pkey)
856     std::vector<unsigned> Mkey;
857     Mkey.reserve(n-1);
858     // we store our subminors in maps, keys being the rows they arise from
859     typedef std::map<std::vector<unsigned>,class ex> Rmap;
860     typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
861     Rmap A;
862     Rmap B;
863     ex det;
864     // initialize A with last column:
865     for (unsigned r=0; r<n; ++r) {
866         Pkey.erase(Pkey.begin(),Pkey.end());
867         Pkey.push_back(r);
868         A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
869     }
870     // proceed from right to left through matrix
871     for (int c=n-2; c>=0; --c) {
872         Pkey.erase(Pkey.begin(),Pkey.end());  // don't change capacity
873         Mkey.erase(Mkey.begin(),Mkey.end());
874         for (unsigned i=0; i<n-c; ++i)
875             Pkey.push_back(i);
876         unsigned fc = 0;  // controls logic for our strange flipper counter
877         do {
878             det = _ex0();
879             for (unsigned r=0; r<n-c; ++r) {
880                 // maybe there is nothing to do?
881                 if (m[Pkey[r]*n+c].is_zero())
882                     continue;
883                 // create the sorted key for all possible minors
884                 Mkey.erase(Mkey.begin(),Mkey.end());
885                 for (unsigned i=0; i<n-c; ++i)
886                     if (i!=r)
887                         Mkey.push_back(Pkey[i]);
888                 // Fetch the minors and compute the new determinant
889                 if (r%2)
890                     det -= m[Pkey[r]*n+c]*A[Mkey];
891                 else
892                     det += m[Pkey[r]*n+c]*A[Mkey];
893             }
894             // prevent build-up of deep nesting of expressions saves time:
895             det = det.expand();
896             // store the new determinant at its place in B:
897             if (!det.is_zero())
898                 B.insert(Rmap_value(Pkey,det));
899             // increment our strange flipper counter
900             for (fc=n-c; fc>0; --fc) {
901                 ++Pkey[fc-1];
902                 if (Pkey[fc-1]<fc+c)
903                     break;
904             }
905             if (fc<n-c)
906                 for (unsigned j=fc; j<n-c; ++j)
907                     Pkey[j] = Pkey[j-1]+1;
908         } while(fc);
909         // next column, so change the role of A and B:
910         A = B;
911         B.clear();
912     }
913     
914     return det;
915 }
916
917
918 /** Perform the steps of an ordinary Gaussian elimination to bring the m x n
919  *  matrix into an upper echelon form.  The algorithm is ok for matrices
920  *  with numeric coefficients but quite unsuited for symbolic matrices.
921  *
922  *  @param det may be set to true to save a lot of space if one is only
923  *  interested in the diagonal elements (i.e. for calculating determinants).
924  *  The others are set to zero in this case.
925  *  @return sign is 1 if an even number of rows was swapped, -1 if an odd
926  *  number of rows was swapped and 0 if the matrix is singular. */
927 int matrix::gauss_elimination(const bool det)
928 {
929     ensure_if_modifiable();
930     const unsigned m = this->rows();
931     const unsigned n = this->cols();
932     GINAC_ASSERT(!det || n==m);
933     int sign = 1;
934     
935     unsigned r0 = 0;
936     for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
937         int indx = pivot(r0, r1, true);
938         if (indx == -1) {
939             sign = 0;
940             if (det)
941                 return 0;  // leaves *this in a messy state
942         }
943         if (indx>=0) {
944             if (indx > 0)
945                 sign = -sign;
946             for (unsigned r2=r0+1; r2<m; ++r2) {
947                 ex piv = this->m[r2*n+r1] / this->m[r0*n+r1];
948                 for (unsigned c=r1+1; c<n; ++c) {
949                     this->m[r2*n+c] -= piv * this->m[r0*n+c];
950                     if (!this->m[r2*n+c].info(info_flags::numeric))
951                         this->m[r2*n+c] = this->m[r2*n+c].normal();
952                 }
953                 // fill up left hand side with zeros
954                 for (unsigned c=0; c<=r1; ++c)
955                     this->m[r2*n+c] = _ex0();
956             }
957             if (det) {
958                 // save space by deleting no longer needed elements
959                 for (unsigned c=r0+1; c<n; ++c)
960                     this->m[r0*n+c] = _ex0();
961             }
962             ++r0;
963         }
964     }
965     
966     return sign;
967 }
968
969
970 /** Perform the steps of division free elimination to bring the m x n matrix
971  *  into an upper echelon form.
972  *
973  *  @param det may be set to true to save a lot of space if one is only
974  *  interested in the diagonal elements (i.e. for calculating determinants).
975  *  The others are set to zero in this case.
976  *  @return sign is 1 if an even number of rows was swapped, -1 if an odd
977  *  number of rows was swapped and 0 if the matrix is singular. */
978 int matrix::division_free_elimination(const bool det)
979 {
980     ensure_if_modifiable();
981     const unsigned m = this->rows();
982     const unsigned n = this->cols();
983     GINAC_ASSERT(!det || n==m);
984     int sign = 1;
985     
986     unsigned r0 = 0;
987     for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
988         int indx = pivot(r0, r1, true);
989         if (indx==-1) {
990             sign = 0;
991             if (det)
992                 return 0;  // leaves *this in a messy state
993         }
994         if (indx>=0) {
995             if (indx>0)
996                 sign = -sign;
997             for (unsigned r2=r0+1; r2<m; ++r2) {
998                 for (unsigned c=r1+1; c<n; ++c)
999                     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();
1000                 // fill up left hand side with zeros
1001                 for (unsigned c=0; c<=r1; ++c)
1002                     this->m[r2*n+c] = _ex0();
1003             }
1004             if (det) {
1005                 // save space by deleting no longer needed elements
1006                 for (unsigned c=r0+1; c<n; ++c)
1007                     this->m[r0*n+c] = _ex0();
1008             }
1009             ++r0;
1010         }
1011     }
1012     
1013     return sign;
1014 }
1015
1016
1017 /** Perform the steps of Bareiss' one-step fraction free elimination to bring
1018  *  the matrix into an upper echelon form.  Fraction free elimination means
1019  *  that divide is used straightforwardly, without computing GCDs first.  This
1020  *  is possible, since we know the divisor at each step.
1021  *  
1022  *  @param det may be set to true to save a lot of space if one is only
1023  *  interested in the last element (i.e. for calculating determinants). The
1024  *  others are set to zero in this case.
1025  *  @return sign is 1 if an even number of rows was swapped, -1 if an odd
1026  *  number of rows was swapped and 0 if the matrix is singular. */
1027 int matrix::fraction_free_elimination(const bool det)
1028 {
1029     // Method:
1030     // (single-step fraction free elimination scheme, already known to Jordan)
1031     //
1032     // Usual division-free elimination sets m[0](r,c) = m(r,c) and then sets
1033     //     m[k+1](r,c) = m[k](k,k) * m[k](r,c) - m[k](r,k) * m[k](k,c).
1034     //
1035     // Bareiss (fraction-free) elimination in addition divides that element
1036     // by m[k-1](k-1,k-1) for k>1, where it can be shown by means of the
1037     // Sylvester determinant that this really divides m[k+1](r,c).
1038     //
1039     // We also allow rational functions where the original prove still holds.
1040     // However, we must care for numerator and denominator separately and
1041     // "manually" work in the integral domains because of subtle cancellations
1042     // (see below).  This blows up the bookkeeping a bit and the formula has
1043     // to be modified to expand like this (N{x} stands for numerator of x,
1044     // D{x} for denominator of x):
1045     //     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)}
1046     //                     -N{m[k](r,k)}*N{m[k](k,c)}*D{m[k](k,k)}*D{m[k](r,c)}
1047     //     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)}
1048     // where for k>1 we now divide N{m[k+1](r,c)} by
1049     //     N{m[k-1](k-1,k-1)}
1050     // and D{m[k+1](r,c)} by
1051     //     D{m[k-1](k-1,k-1)}.
1052     
1053     ensure_if_modifiable();
1054     const unsigned m = this->rows();
1055     const unsigned n = this->cols();
1056     GINAC_ASSERT(!det || n==m);
1057     int sign = 1;
1058     if (m==1)
1059         return 1;
1060     ex divisor_n = 1;
1061     ex divisor_d = 1;
1062     ex dividend_n;
1063     ex dividend_d;
1064     
1065     // We populate temporary matrices to subsequently operate on.  There is
1066     // one holding numerators and another holding denominators of entries.
1067     // This is a must since the evaluator (or even earlier mul's constructor)
1068     // might cancel some trivial element which causes divide() to fail.  The
1069     // elements are normalized first (yes, even though this algorithm doesn't
1070     // need GCDs) since the elements of *this might be unnormalized, which
1071     // makes things more complicated than they need to be.
1072     matrix tmp_n(*this);
1073     matrix tmp_d(m,n);  // for denominators, if needed
1074     lst srl;  // symbol replacement list
1075     exvector::iterator it = this->m.begin();
1076     exvector::iterator tmp_n_it = tmp_n.m.begin();
1077     exvector::iterator tmp_d_it = tmp_d.m.begin();
1078     for (; it!= this->m.end(); ++it, ++tmp_n_it, ++tmp_d_it) {
1079         (*tmp_n_it) = (*it).normal().to_rational(srl);
1080         (*tmp_d_it) = (*tmp_n_it).denom();
1081         (*tmp_n_it) = (*tmp_n_it).numer();
1082     }
1083     
1084     unsigned r0 = 0;
1085     for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
1086         int indx = tmp_n.pivot(r0, r1, true);
1087         if (indx==-1) {
1088             sign = 0;
1089             if (det)
1090                 return 0;
1091         }
1092         if (indx>=0) {
1093             if (indx>0) {
1094                 sign = -sign;
1095                 // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
1096                 for (unsigned c=r1; c<n; ++c)
1097                     tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
1098             }
1099             for (unsigned r2=r0+1; r2<m; ++r2) {
1100                 for (unsigned c=r1+1; c<n; ++c) {
1101                     dividend_n = (tmp_n.m[r0*n+r1]*tmp_n.m[r2*n+c]*
1102                                   tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]
1103                                  -tmp_n.m[r2*n+r1]*tmp_n.m[r0*n+c]*
1104                                   tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
1105                     dividend_d = (tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]*
1106                                   tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
1107                     bool check = divide(dividend_n, divisor_n,
1108                                         tmp_n.m[r2*n+c], true);
1109                     check &= divide(dividend_d, divisor_d,
1110                                     tmp_d.m[r2*n+c], true);
1111                     GINAC_ASSERT(check);
1112                 }
1113                 // fill up left hand side with zeros
1114                 for (unsigned c=0; c<=r1; ++c)
1115                     tmp_n.m[r2*n+c] = _ex0();
1116             }
1117             if ((r1<n-1)&&(r0<m-1)) {
1118                 // compute next iteration's divisor
1119                 divisor_n = tmp_n.m[r0*n+r1].expand();
1120                 divisor_d = tmp_d.m[r0*n+r1].expand();
1121                 if (det) {
1122                     // save space by deleting no longer needed elements
1123                     for (unsigned c=0; c<n; ++c) {
1124                         tmp_n.m[r0*n+c] = _ex0();
1125                         tmp_d.m[r0*n+c] = _ex1();
1126                     }
1127                 }
1128             }
1129             ++r0;
1130         }
1131     }
1132     // repopulate *this matrix:
1133     it = this->m.begin();
1134     tmp_n_it = tmp_n.m.begin();
1135     tmp_d_it = tmp_d.m.begin();
1136     for (; it!= this->m.end(); ++it, ++tmp_n_it, ++tmp_d_it)
1137         (*it) = ((*tmp_n_it)/(*tmp_d_it)).subs(srl);
1138     
1139     return sign;
1140 }
1141
1142
1143 /** Partial pivoting method for matrix elimination schemes.
1144  *  Usual pivoting (symbolic==false) returns the index to the element with the
1145  *  largest absolute value in column ro and swaps the current row with the one
1146  *  where the element was found.  With (symbolic==true) it does the same thing
1147  *  with the first non-zero element.
1148  *
1149  *  @param ro is the row from where to begin
1150  *  @param co is the column to be inspected
1151  *  @param symbolic signal if we want the first non-zero element to be pivoted
1152  *  (true) or the one with the largest absolute value (false).
1153  *  @return 0 if no interchange occured, -1 if all are zero (usually signaling
1154  *  a degeneracy) and positive integer k means that rows ro and k were swapped.
1155  */
1156 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
1157 {
1158     unsigned k = ro;
1159     if (symbolic) {
1160         // search first non-zero element in column co beginning at row ro
1161         while ((k<row) && (this->m[k*col+co].expand().is_zero()))
1162             ++k;
1163     } else {
1164         // search largest element in column co beginning at row ro
1165         GINAC_ASSERT(is_ex_of_type(this->m[k*col+co],numeric));
1166         unsigned kmax = k+1;
1167         numeric mmax = abs(ex_to_numeric(m[kmax*col+co]));
1168         while (kmax<row) {
1169             GINAC_ASSERT(is_ex_of_type(this->m[kmax*col+co],numeric));
1170             numeric tmp = ex_to_numeric(this->m[kmax*col+co]);
1171             if (abs(tmp) > mmax) {
1172                 mmax = tmp;
1173                 k = kmax;
1174             }
1175             ++kmax;
1176         }
1177         if (!mmax.is_zero())
1178             k = kmax;
1179     }
1180     if (k==row)
1181         // all elements in column co below row ro vanish
1182         return -1;
1183     if (k==ro)
1184         // matrix needs no pivoting
1185         return 0;
1186     // matrix needs pivoting, so swap rows k and ro
1187     ensure_if_modifiable();
1188     for (unsigned c=0; c<col; ++c)
1189         m[k*col+c].swap(m[ro*col+c]);
1190     
1191     return k;
1192 }
1193
1194 /** Convert list of lists to matrix. */
1195 ex lst_to_matrix(const ex &l)
1196 {
1197     if (!is_ex_of_type(l, lst))
1198         throw(std::invalid_argument("argument to lst_to_matrix() must be a lst"));
1199     
1200     // Find number of rows and columns
1201     unsigned rows = l.nops(), cols = 0, i, j;
1202     for (i=0; i<rows; i++)
1203         if (l.op(i).nops() > cols)
1204             cols = l.op(i).nops();
1205     
1206     // Allocate and fill matrix
1207     matrix &m = *new matrix(rows, cols);
1208     for (i=0; i<rows; i++)
1209         for (j=0; j<cols; j++)
1210             if (l.op(i).nops() > j)
1211                 m.set(i, j, l.op(i).op(j));
1212             else
1213                 m.set(i, j, ex(0));
1214     return m;
1215 }
1216
1217 //////////
1218 // global constants
1219 //////////
1220
1221 const matrix some_matrix;
1222 const type_info & typeid_matrix=typeid(some_matrix);
1223
1224 #ifndef NO_NAMESPACE_GINAC
1225 } // namespace GiNaC
1226 #endif // ndef NO_NAMESPACE_GINAC