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