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