]> www.ginac.de Git - ginac.git/blob - ginac/matrix.cpp
- Fixed a logic error in numeric::info().
[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 "utils.h"
30 #include "debugmsg.h"
31 #include "numeric.h"
32
33 #ifndef NO_NAMESPACE_GINAC
34 namespace GiNaC {
35 #endif // ndef NO_NAMESPACE_GINAC
36
37 GINAC_IMPLEMENT_REGISTERED_CLASS(matrix, basic)
38
39 //////////
40 // default constructor, destructor, copy constructor, assignment operator
41 // and helpers:
42 //////////
43
44 // public
45
46 /** Default ctor.  Initializes to 1 x 1-dimensional zero-matrix. */
47 matrix::matrix()
48     : inherited(TINFO_matrix), row(1), col(1)
49 {
50     debugmsg("matrix default constructor",LOGLEVEL_CONSTRUCT);
51     m.push_back(_ex0());
52 }
53
54 matrix::~matrix()
55 {
56     debugmsg("matrix destructor",LOGLEVEL_DESTRUCT);
57 }
58
59 matrix::matrix(const matrix & other)
60 {
61     debugmsg("matrix copy constructor",LOGLEVEL_CONSTRUCT);
62     copy(other);
63 }
64
65 const matrix & matrix::operator=(const matrix & other)
66 {
67     debugmsg("matrix operator=",LOGLEVEL_ASSIGNMENT);
68     if (this != &other) {
69         destroy(1);
70         copy(other);
71     }
72     return *this;
73 }
74
75 // protected
76
77 void matrix::copy(const matrix & other)
78 {
79     inherited::copy(other);
80     row=other.row;
81     col=other.col;
82     m=other.m;  // use STL's vector copying
83 }
84
85 void matrix::destroy(bool call_parent)
86 {
87     if (call_parent) inherited::destroy(call_parent);
88 }
89
90 //////////
91 // other constructors
92 //////////
93
94 // public
95
96 /** Very common ctor.  Initializes to r x c-dimensional zero-matrix.
97  *
98  *  @param r number of rows
99  *  @param c number of cols */
100 matrix::matrix(unsigned r, unsigned c)
101     : inherited(TINFO_matrix), row(r), col(c)
102 {
103     debugmsg("matrix constructor from unsigned,unsigned",LOGLEVEL_CONSTRUCT);
104     m.resize(r*c, _ex0());
105 }
106
107 // protected
108
109 /** Ctor from representation, for internal use only. */
110 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
111     : inherited(TINFO_matrix), row(r), col(c), m(m2)
112 {
113     debugmsg("matrix constructor from unsigned,unsigned,exvector",LOGLEVEL_CONSTRUCT);
114 }
115
116 //////////
117 // archiving
118 //////////
119
120 /** Construct object from archive_node. */
121 matrix::matrix(const archive_node &n, const lst &sym_lst) : inherited(n, sym_lst)
122 {
123     debugmsg("matrix constructor from archive_node", LOGLEVEL_CONSTRUCT);
124     if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
125         throw (std::runtime_error("unknown matrix dimensions in archive"));
126     m.reserve(row * col);
127     for (unsigned int i=0; true; i++) {
128         ex e;
129         if (n.find_ex("m", e, sym_lst, i))
130             m.push_back(e);
131         else
132             break;
133     }
134 }
135
136 /** Unarchive the object. */
137 ex matrix::unarchive(const archive_node &n, const lst &sym_lst)
138 {
139     return (new matrix(n, sym_lst))->setflag(status_flags::dynallocated);
140 }
141
142 /** Archive the object. */
143 void matrix::archive(archive_node &n) const
144 {
145     inherited::archive(n);
146     n.add_unsigned("row", row);
147     n.add_unsigned("col", col);
148     exvector::const_iterator i = m.begin(), iend = m.end();
149     while (i != iend) {
150         n.add_ex("m", *i);
151         i++;
152     }
153 }
154
155 //////////
156 // functions overriding virtual functions from bases classes
157 //////////
158
159 // public
160
161 basic * matrix::duplicate() const
162 {
163     debugmsg("matrix duplicate",LOGLEVEL_DUPLICATE);
164     return new matrix(*this);
165 }
166
167 void matrix::print(ostream & os, unsigned upper_precedence) const
168 {
169     debugmsg("matrix print",LOGLEVEL_PRINT);
170     os << "[[ ";
171     for (unsigned r=0; r<row-1; ++r) {
172         os << "[[";
173         for (unsigned c=0; c<col-1; ++c) {
174             os << m[r*col+c] << ",";
175         }
176         os << m[col*(r+1)-1] << "]], ";
177     }
178     os << "[[";
179     for (unsigned c=0; c<col-1; ++c) {
180         os << m[(row-1)*col+c] << ",";
181     }
182     os << m[row*col-1] << "]] ]]";
183 }
184
185 void matrix::printraw(ostream & os) const
186 {
187     debugmsg("matrix printraw",LOGLEVEL_PRINT);
188     os << "matrix(" << row << "," << col <<",";
189     for (unsigned r=0; r<row-1; ++r) {
190         os << "(";
191         for (unsigned c=0; c<col-1; ++c) {
192             os << m[r*col+c] << ",";
193         }
194         os << m[col*(r-1)-1] << "),";
195     }
196     os << "(";
197     for (unsigned c=0; c<col-1; ++c) {
198         os << m[(row-1)*col+c] << ",";
199     }
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     return m[i];
219 }
220
221 /** expands the elements of a matrix entry by entry. */
222 ex matrix::expand(unsigned options) const
223 {
224     exvector tmp(row*col);
225     for (unsigned i=0; i<row*col; ++i) {
226         tmp[i]=m[i].expand(options);
227     }
228     return matrix(row, col, tmp);
229 }
230
231 /** Search ocurrences.  A matrix 'has' an expression if it is the expression
232  *  itself or one of the elements 'has' it. */
233 bool matrix::has(const ex & other) const
234 {
235     GINAC_ASSERT(other.bp!=0);
236     
237     // tautology: it is the expression itself
238     if (is_equal(*other.bp)) return true;
239     
240     // search all the elements
241     for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r) {
242         if ((*r).has(other)) return true;
243     }
244     return false;
245 }
246
247 /** evaluate matrix entry by entry. */
248 ex matrix::eval(int level) const
249 {
250     debugmsg("matrix eval",LOGLEVEL_MEMBER_FUNCTION);
251     
252     // check if we have to do anything at all
253     if ((level==1)&&(flags & status_flags::evaluated)) {
254         return *this;
255     }
256     
257     // emergency break
258     if (level == -max_recursion_level) {
259         throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
260     }
261     
262     // eval() entry by entry
263     exvector m2(row*col);
264     --level;    
265     for (unsigned r=0; r<row; ++r) {
266         for (unsigned c=0; c<col; ++c) {
267             m2[r*col+c] = m[r*col+c].eval(level);
268         }
269     }
270     
271     return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
272                                                status_flags::evaluated );
273 }
274
275 /** evaluate matrix numerically entry by entry. */
276 ex matrix::evalf(int level) const
277 {
278     debugmsg("matrix evalf",LOGLEVEL_MEMBER_FUNCTION);
279         
280     // check if we have to do anything at all
281     if (level==1) {
282         return *this;
283     }
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     
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 /** 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     
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 /** 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     
383     exvector prod(row*other.col);
384     for (unsigned i=0; i<row; ++i) {
385         for (unsigned j=0; j<other.col; ++j) {
386             for (unsigned l=0; l<col; ++l) {
387                 prod[i*other.col+j] += m[i*col+l] * other.m[l*other.col+j];
388             }
389         }
390     }
391     return matrix(row, other.col, prod);
392 }
393
394 /** operator() to access elements.
395  *
396  *  @param ro row of element
397  *  @param co column of element 
398  *  @exception range_error (index out of range) */
399 const ex & matrix::operator() (unsigned ro, unsigned co) const
400 {
401     if (ro<0 || ro>=row || co<0 || co>=col) {
402         throw (std::range_error("matrix::operator(): index out of range"));
403     }
404     
405     return m[ro*col+co];
406 }
407
408 /** Set individual elements manually.
409  *
410  *  @exception range_error (index out of range) */
411 matrix & matrix::set(unsigned ro, unsigned co, ex value)
412 {
413     if (ro<0 || ro>=row || co<0 || co>=col) {
414         throw (std::range_error("matrix::set(): index out of range"));
415     }
416     
417     ensure_if_modifiable();
418     m[ro*col+co] = value;
419     return *this;
420 }
421
422 /** Transposed of an m x n matrix, producing a new n x m matrix object that
423  *  represents the transposed. */
424 matrix matrix::transpose(void) const
425 {
426     exvector trans(col*row);
427     
428     for (unsigned r=0; r<col; ++r)
429         for (unsigned c=0; c<row; ++c)
430             trans[r*row+c] = m[c*col+r];
431     
432     return matrix(col,row,trans);
433 }
434
435 /*  Leverrier algorithm for large matrices having at least one symbolic entry.
436  *  This routine is only called internally by matrix::determinant(). The
437  *  algorithm is very bad for symbolic matrices since it returns expressions
438  *  that are quite hard to expand. */
439 /*ex matrix::determinant_symbolic_leverrier(const matrix & M)
440  *{
441  *    GINAC_ASSERT(M.rows()==M.cols());  // cannot happen, just in case...
442  *    
443  *    matrix B(M);
444  *    matrix I(M.row, M.col);
445  *    ex c=B.trace();
446  *    for (unsigned i=1; i<M.row; ++i) {
447  *        for (unsigned j=0; j<M.row; ++j)
448  *            I.m[j*M.col+j] = c;
449  *        B = M.mul(B.sub(I));
450  *        c = B.trace()/ex(i+1);
451  *    }
452  *    if (M.row%2) {
453  *        return c;
454  *    } else {
455  *        return -c;
456  *    }
457  *}*/
458
459 /** Determinant of square matrix.  This routine doesn't actually calculate the
460  *  determinant, it only implements some heuristics about which algorithm to
461  *  call.  When the parameter for normalization is explicitly turned off this
462  *  method does not normalize its result at the end, which might imply that
463  *  the symbolic 2x2 matrix [[a/(a-b),1],[b/(a-b),1]] is not immediatly
464  *  recognized to be unity.  (This is Mathematica's default behaviour, it
465  *  should be used with care.)
466  *
467  *  @param     normalized may be set to false if no normalization of the
468  *             result is desired (i.e. to force Mathematica behavior, Maple
469  *             does normalize the result).
470  *  @return    the determinant as a new expression
471  *  @exception logic_error (matrix not square) */
472 ex matrix::determinant(bool normalized) const
473 {
474     if (row != col) {
475         throw (std::logic_error("matrix::determinant(): matrix not square"));
476     }
477
478     // check, if there are non-numeric entries in the matrix:
479     for (exvector::const_iterator r=m.begin(); r!=m.end(); ++r) {
480         if (!(*r).info(info_flags::numeric)) {
481             if (normalized)
482                 // return determinant_symbolic_minor().normal();
483                 return determinant_symbolic_minor().normal();
484             else
485                 return determinant_symbolic_perm();
486         }
487     }
488     // if it turns out that all elements are numeric
489     return determinant_numeric();
490 }
491
492 /** Trace of a matrix.
493  *
494  *  @return    the sum of diagonal elements
495  *  @exception logic_error (matrix not square) */
496 ex matrix::trace(void) const
497 {
498     if (row != col) {
499         throw (std::logic_error("matrix::trace(): matrix not square"));
500     }
501     
502     ex tr;
503     for (unsigned r=0; r<col; ++r)
504         tr += m[r*col+r];
505
506     return tr;
507 }
508
509 /** Characteristic Polynomial.  The characteristic polynomial of a matrix M is
510  *  defined as the determiant of (M - lambda * 1) where 1 stands for the unit
511  *  matrix of the same dimension as M.  This method returns the characteristic
512  *  polynomial as a new expression.
513  *
514  *  @return    characteristic polynomial as new expression
515  *  @exception logic_error (matrix not square)
516  *  @see       matrix::determinant() */
517 ex matrix::charpoly(const ex & lambda) const
518 {
519     if (row != col) {
520         throw (std::logic_error("matrix::charpoly(): matrix not square"));
521     }
522     
523     matrix M(*this);
524     for (unsigned r=0; r<col; ++r)
525         M.m[r*col+r] -= lambda;
526     
527     return (M.determinant());
528 }
529
530 /** Inverse of this matrix.
531  *
532  *  @return    the inverted matrix
533  *  @exception logic_error (matrix not square)
534  *  @exception runtime_error (singular matrix) */
535 matrix matrix::inverse(void) const
536 {
537     if (row != col) {
538         throw (std::logic_error("matrix::inverse(): matrix not square"));
539     }
540     
541     matrix tmp(row,col);
542     // set tmp to the unit matrix
543     for (unsigned i=0; i<col; ++i)
544         tmp.m[i*col+i] = _ex1();
545
546     // create a copy of this matrix
547     matrix cpy(*this);
548     for (unsigned r1=0; r1<row; ++r1) {
549         int indx = cpy.pivot(r1);
550         if (indx == -1) {
551             throw (std::runtime_error("matrix::inverse(): singular matrix"));
552         }
553         if (indx != 0) {  // swap rows r and indx of matrix tmp
554             for (unsigned i=0; i<col; ++i) {
555                 tmp.m[r1*col+i].swap(tmp.m[indx*col+i]);
556             }
557         }
558         ex a1 = cpy.m[r1*col+r1];
559         for (unsigned c=0; c<col; ++c) {
560             cpy.m[r1*col+c] /= a1;
561             tmp.m[r1*col+c] /= a1;
562         }
563         for (unsigned r2=0; r2<row; ++r2) {
564             if (r2 != r1) {
565                 ex a2 = cpy.m[r2*col+r1];
566                 for (unsigned c=0; c<col; ++c) {
567                     cpy.m[r2*col+c] -= a2 * cpy.m[r1*col+c];
568                     tmp.m[r2*col+c] -= a2 * tmp.m[r1*col+c];
569                 }
570             }
571         }
572     }
573     return tmp;
574 }
575
576 // superfluous helper function
577 void matrix::ffe_swap(unsigned r1, unsigned c1, unsigned r2 ,unsigned c2)
578 {
579     ensure_if_modifiable();
580     
581     ex tmp = ffe_get(r1,c1);
582     ffe_set(r1,c1,ffe_get(r2,c2));
583     ffe_set(r2,c2,tmp);
584 }
585
586 // superfluous helper function
587 void matrix::ffe_set(unsigned r, unsigned c, ex e)
588 {
589     set(r-1,c-1,e);
590 }
591
592 // superfluous helper function
593 ex matrix::ffe_get(unsigned r, unsigned c) const
594 {
595     return operator()(r-1,c-1);
596 }
597
598 /** Solve a set of equations for an m x n matrix by fraction-free Gaussian
599  *  elimination.  Based on algorithm 9.1 from 'Algorithms for Computer Algebra'
600  *  by Keith O. Geddes et al.
601  *
602  *  @param vars n x p matrix
603  *  @param rhs m x p matrix
604  *  @exception logic_error (incompatible matrices)
605  *  @exception runtime_error (singular matrix) */
606 matrix matrix::fraction_free_elim(const matrix & vars,
607                                   const matrix & rhs) const
608 {
609     // FIXME: implement a Sasaki-Murao scheme which avoids division at all!
610     if ((row != rhs.row) || (col != vars.row) || (rhs.col != vars.col))
611         throw (std::logic_error("matrix::fraction_free_elim(): incompatible matrices"));
612     
613     matrix a(*this);  // make a copy of the matrix
614     matrix b(rhs);    // make a copy of the rhs vector
615     
616     // given an m x n matrix a, reduce it to upper echelon form
617     unsigned m = a.row;
618     unsigned n = a.col;
619     int sign = 1;
620     ex divisor = 1;
621     unsigned r = 1;
622     
623     // eliminate below row r, with pivot in column k
624     for (unsigned k=1; (k<=n)&&(r<=m); ++k) {
625         // find a nonzero pivot
626         unsigned p;
627         for (p=r; (p<=m)&&(a.ffe_get(p,k).is_equal(_ex0())); ++p) {}
628         // pivot is in row p
629         if (p<=m) {
630             if (p!=r) {
631                 // switch rows p and r
632                 for (unsigned j=k; j<=n; ++j)
633                     a.ffe_swap(p,j,r,j);
634                 b.ffe_swap(p,1,r,1);
635                 // keep track of sign changes due to row exchange
636                 sign = -sign;
637             }
638             for (unsigned i=r+1; i<=m; ++i) {
639                 for (unsigned j=k+1; j<=n; ++j) {
640                     a.ffe_set(i,j,(a.ffe_get(r,k)*a.ffe_get(i,j)
641                                   -a.ffe_get(r,j)*a.ffe_get(i,k))/divisor);
642                     a.ffe_set(i,j,a.ffe_get(i,j).normal() /*.normal() */ );
643                 }
644                 b.ffe_set(i,1,(a.ffe_get(r,k)*b.ffe_get(i,1)
645                               -b.ffe_get(r,1)*a.ffe_get(i,k))/divisor);
646                 b.ffe_set(i,1,b.ffe_get(i,1).normal() /*.normal() */ );
647                 a.ffe_set(i,k,0);
648             }
649             divisor = a.ffe_get(r,k);
650             r++;
651         }
652     }
653     // optionally compute the determinant for square or augmented matrices
654     // if (r==m+1) { det = sign*divisor; } else { det = 0; }
655     
656     /*
657     for (unsigned r=1; r<=m; ++r) {
658         for (unsigned c=1; c<=n; ++c) {
659             cout << a.ffe_get(r,c) << "\t";
660         }
661         cout << " | " <<  b.ffe_get(r,1) << endl;
662     }
663     */
664     
665 #ifdef DO_GINAC_ASSERT
666     // test if we really have an upper echelon matrix
667     int zero_in_last_row = -1;
668     for (unsigned r=1; r<=m; ++r) {
669         int zero_in_this_row=0;
670         for (unsigned c=1; c<=n; ++c) {
671             if (a.ffe_get(r,c).is_equal(_ex0()))
672                zero_in_this_row++;
673             else
674                 break;
675         }
676         GINAC_ASSERT((zero_in_this_row>zero_in_last_row)||(zero_in_this_row=n));
677         zero_in_last_row = zero_in_this_row;
678     }
679 #endif // def DO_GINAC_ASSERT
680     
681     /*
682     cout << "after" << endl;
683     cout << "a=" << a << endl;
684     cout << "b=" << b << endl;
685     */
686     
687     // assemble solution
688     matrix sol(n,1);
689     unsigned last_assigned_sol = n+1;
690     for (unsigned r=m; r>0; --r) {
691         unsigned first_non_zero = 1;
692         while ((first_non_zero<=n)&&(a.ffe_get(r,first_non_zero).is_zero()))
693             first_non_zero++;
694         if (first_non_zero>n) {
695             // row consists only of zeroes, corresponding rhs must be 0 as well
696             if (!b.ffe_get(r,1).is_zero()) {
697                 throw (std::runtime_error("matrix::fraction_free_elim(): singular matrix"));
698             }
699         } else {
700             // assign solutions for vars between first_non_zero+1 and
701             // last_assigned_sol-1: free parameters
702             for (unsigned c=first_non_zero+1; c<=last_assigned_sol-1; ++c) {
703                 sol.ffe_set(c,1,vars.ffe_get(c,1));
704             }
705             ex e = b.ffe_get(r,1);
706             for (unsigned c=first_non_zero+1; c<=n; ++c) {
707                 e=e-a.ffe_get(r,c)*sol.ffe_get(c,1);
708             }
709             sol.ffe_set(first_non_zero,1,
710                         (e/a.ffe_get(r,first_non_zero)).normal());
711             last_assigned_sol = first_non_zero;
712         }
713     }
714     // assign solutions for vars between 1 and
715     // last_assigned_sol-1: free parameters
716     for (unsigned c=1; c<=last_assigned_sol-1; ++c)
717         sol.ffe_set(c,1,vars.ffe_get(c,1));
718     
719 #ifdef DO_GINAC_ASSERT
720     // test solution with echelon matrix
721     for (unsigned r=1; r<=m; ++r) {
722         ex e = 0;
723         for (unsigned c=1; c<=n; ++c)
724             e = e+a.ffe_get(r,c)*sol.ffe_get(c,1);
725         if (!(e-b.ffe_get(r,1)).normal().is_zero()) {
726             cout << "e=" << e;
727             cout << "b.ffe_get(" << r<<",1)=" << b.ffe_get(r,1) << endl;
728             cout << "diff=" << (e-b.ffe_get(r,1)).normal() << endl;
729         }
730         GINAC_ASSERT((e-b.ffe_get(r,1)).normal().is_zero());
731     }
732     
733     // test solution with original matrix
734     for (unsigned r=1; r<=m; ++r) {
735         ex e = 0;
736         for (unsigned c=1; c<=n; ++c)
737             e = e+ffe_get(r,c)*sol.ffe_get(c,1);
738         try {
739             if (!(e-rhs.ffe_get(r,1)).normal().is_zero()) {
740                 cout << "e=" << e << endl;
741                 e.printtree(cout);
742                 ex en = e.normal();
743                 cout << "e.normal()=" << en << endl;
744                 en.printtree(cout);
745                 cout << "rhs.ffe_get(" << r<<",1)=" << rhs.ffe_get(r,1) << endl;
746                 cout << "diff=" << (e-rhs.ffe_get(r,1)).normal() << endl;
747             }
748         } catch (...) {
749             ex xxx = e - rhs.ffe_get(r,1);
750             cerr << "xxx=" << xxx << endl << endl;
751         }
752         GINAC_ASSERT((e-rhs.ffe_get(r,1)).normal().is_zero());
753     }
754 #endif // def DO_GINAC_ASSERT
755     
756     return sol;
757 }
758
759 /** Solve a set of equations for an m x n matrix.
760  *
761  *  @param vars n x p matrix
762  *  @param rhs m x p matrix
763  *  @exception logic_error (incompatible matrices)
764  *  @exception runtime_error (singular matrix) */
765 matrix matrix::solve(const matrix & vars,
766                      const matrix & rhs) const
767 {
768     if ((row != rhs.row) || (col != vars.row) || (rhs.col != vars.col))
769         throw (std::logic_error("matrix::solve(): incompatible matrices"));
770     
771     throw (std::runtime_error("FIXME: need implementation."));
772 }
773
774 /** Old and obsolete interface: */
775 matrix matrix::old_solve(const matrix & v) const
776 {
777     if ((v.row != col) || (col != v.row))
778         throw (std::logic_error("matrix::solve(): incompatible matrices"));
779     
780     // build the augmented matrix of *this with v attached to the right
781     matrix tmp(row,col+v.col);
782     for (unsigned r=0; r<row; ++r) {
783         for (unsigned c=0; c<col; ++c)
784             tmp.m[r*tmp.col+c] = this->m[r*col+c];
785         for (unsigned c=0; c<v.col; ++c)
786             tmp.m[r*tmp.col+c+col] = v.m[r*v.col+c];
787     }
788     // cout << "augmented: " << tmp << endl;
789     tmp.gauss_elimination();
790     // cout << "degaussed: " << tmp << endl;
791     // assemble the solution matrix
792     exvector sol(v.row*v.col);
793     for (unsigned c=0; c<v.col; ++c) {
794         for (unsigned r=row; r>0; --r) {
795             for (unsigned i=r; i<col; ++i)
796                 sol[(r-1)*v.col+c] -= tmp.m[(r-1)*tmp.col+i]*sol[i*v.col+c];
797             sol[(r-1)*v.col+c] += tmp.m[(r-1)*tmp.col+col+c];
798             sol[(r-1)*v.col+c] = (sol[(r-1)*v.col+c]/tmp.m[(r-1)*tmp.col+(r-1)]).normal();
799         }
800     }
801     return matrix(v.row, v.col, sol);
802 }
803
804 // protected
805
806 /** Determinant of purely numeric matrix, using pivoting.
807  *
808  *  @see       matrix::determinant() */
809 ex matrix::determinant_numeric(void) const
810 {
811     matrix tmp(*this);
812     ex det = _ex1();
813     ex piv;
814     
815     for (unsigned r1=0; r1<row; ++r1) {
816         int indx = tmp.pivot(r1);
817         if (indx == -1)
818             return _ex0();
819         if (indx != 0)
820             det *= _ex_1();
821         det = det * tmp.m[r1*col+r1];
822         for (unsigned r2=r1+1; r2<row; ++r2) {
823             piv = tmp.m[r2*col+r1] / tmp.m[r1*col+r1];
824             for (unsigned c=r1+1; c<col; c++) {
825                 tmp.m[r2*col+c] -= piv * tmp.m[r1*col+c];
826             }
827         }
828     }
829     return det;
830 }
831
832 /** Recursive determinant for small matrices having at least one symbolic
833  *  entry.  The basic algorithm, known as Laplace-expansion, is enhanced by
834  *  some bookkeeping to avoid calculation of the same submatrices ("minors")
835  *  more than once.  According to W.M.Gentleman and S.C.Johnson this algorithm
836  *  is better than elimination schemes for sparse multivariate polynomials and
837  *  also for dense univariate polynomials once the dimesion becomes larger
838  *  than 7.
839  *
840  *  @see matrix::determinant() */
841 ex matrix::determinant_symbolic_minor(void) const
842 {
843     // for small matrices the algorithm does not make sense:
844     if (this->row==1)
845         return m[0];
846     if (this->row==2)
847         return (m[0]*m[3]-m[2]*m[1]);
848     if (this->row==3)
849         return ((m[4]*m[8]-m[5]*m[7])*m[0]-
850                 (m[1]*m[8]-m[2]*m[7])*m[3]+
851                 (m[1]*m[5]-m[4]*m[2])*m[6]);
852     
853     // This algorithm can best be understood by looking at a naive
854     // implementation of Laplace-expansion, like this one:
855     // ex det;
856     // matrix minorM(this->row-1,this->col-1);
857     // for (unsigned r1=0; r1<this->row; ++r1) {
858     //     // shortcut if element(r1,0) vanishes
859     //     if (m[r1*col].is_zero())
860     //         continue;
861     //     // assemble the minor matrix
862     //     for (unsigned r=0; r<minorM.rows(); ++r) {
863     //         for (unsigned c=0; c<minorM.cols(); ++c) {
864     //             if (r<r1)
865     //                 minorM.set(r,c,m[r*col+c+1]);
866     //             else
867     //                 minorM.set(r,c,m[(r+1)*col+c+1]);
868     //         }
869     //     }
870     //     // recurse down and care for sign:
871     //     if (r1%2)
872     //         det -= m[r1*col] * minorM.determinant_symbolic_minor();
873     //     else
874     //         det += m[r1*col] * minorM.determinant_symbolic_minor();
875     // }
876     // return det;
877     // What happens is that while proceeding down many of the minors are
878     // computed more than once.  In particular, there are binomial(n,k)
879     // kxk minors and each one is computed factorial(n-k) times.  Therefore
880     // it is reasonable to store the results of the minors.  We proceed from
881     // right to left.  At each column c we only need to retrieve the minors
882     // calculated in step c-1.  We therefore only have to store at most 
883     // 2*binomial(n,n/2) minors.
884     
885     // we store our subminors in maps, keys being the rows they arise from
886     typedef map<vector<unsigned>,class ex> Rmap;
887     typedef map<vector<unsigned>,class ex>::value_type Rmap_value;
888     Rmap A, B;
889     ex det;
890     vector<unsigned> Pkey;    // Unique flipper counter for partitioning into minors
891     Pkey.reserve(this->col);
892     vector<unsigned> Mkey;    // key for minor determinant (a subpartition of Pkey)
893     Mkey.reserve(this->col-1);
894     // initialize A with last column:
895     for (unsigned r=0; r<this->col; ++r) {
896         Pkey.erase(Pkey.begin(),Pkey.end());
897         Pkey.push_back(r);
898         A.insert(Rmap_value(Pkey,m[this->col*r+this->col-1]));
899     }
900     // proceed from right to left through matrix
901     for (int c=this->col-2; c>=0; --c) {
902         Pkey.erase(Pkey.begin(),Pkey.end());  // don't change capacity
903         Mkey.erase(Mkey.begin(),Mkey.end());
904         for (unsigned i=0; i<this->col-c; ++i)
905             Pkey.push_back(i);
906         unsigned fc = 0;  // controls logic for our strange flipper counter
907         do {
908             A.insert(Rmap_value(Pkey,_ex0()));
909             det = _ex0();
910             for (unsigned r=0; r<this->col-c; ++r) {
911                 // maybe there is nothing to do?
912                 if (m[Pkey[r]*this->col+c].is_zero())
913                     continue;
914                 // create the sorted key for all possible minors
915                 Mkey.erase(Mkey.begin(),Mkey.end());
916                 for (unsigned i=0; i<this->col-c; ++i)
917                     if (i!=r)
918                         Mkey.push_back(Pkey[i]);
919                 // Fetch the minors and compute the new determinant
920                 if (r%2)
921                     det -= m[Pkey[r]*this->col+c]*A[Mkey];
922                 else
923                     det += m[Pkey[r]*this->col+c]*A[Mkey];
924             }
925             // Store the new determinant at its place in B:
926             B.insert(Rmap_value(Pkey,det));
927             // increment our strange flipper counter
928             for (fc=this->col-c; fc>0; --fc) {
929                 ++Pkey[fc-1];
930                 if (Pkey[fc-1]<fc+c)
931                     break;
932             }
933             if (fc<this->col-c)
934                 for (unsigned j=fc; j<this->col-c; ++j)
935                     Pkey[j] = Pkey[j-1]+1;
936         } while(fc);
937         // change the role of A and B:
938         A = B;
939         B.clear();
940     }
941     
942     return det;
943 }
944
945 /** Determinant built by application of the full permutation group.  This
946  *  routine is only called internally by matrix::determinant(). */
947 ex matrix::determinant_symbolic_perm(void) const
948 {
949     if (rows()==1)  // speed things up
950         return m[0];
951     
952     ex det;
953     ex term;
954     vector<unsigned> sigma(col);
955     for (unsigned i=0; i<col; ++i)
956         sigma[i]=i;
957     
958     do {
959         term = (*this)(sigma[0],0);
960         for (unsigned i=1; i<col; ++i)
961             term *= (*this)(sigma[i],i);
962         det += permutation_sign(sigma)*term;
963     } while (next_permutation(sigma.begin(), sigma.end()));
964     
965     return det;
966 }
967
968 /** Perform the steps of an ordinary Gaussian elimination to bring the matrix
969  *  into an upper echelon form.
970  *
971  *  @return sign is 1 if an even number of rows was swapped, -1 if an odd
972  *  number of rows was swapped and 0 if the matrix is singular. */
973 int matrix::gauss_elimination(void)
974 {
975     int sign = 1;
976     ensure_if_modifiable();
977     for (unsigned r1=0; r1<row-1; ++r1) {
978         int indx = pivot(r1);
979         if (indx == -1)
980             return 0;  // Note: leaves *this in a messy state.
981         if (indx > 0)
982             sign = -sign;
983         for (unsigned r2=r1+1; r2<row; ++r2) {
984             for (unsigned c=r1+1; c<col; ++c)
985                 this->m[r2*col+c] -= this->m[r2*col+r1]*this->m[r1*col+c]/this->m[r1*col+r1];
986             for (unsigned c=0; c<=r1; ++c)
987                 this->m[r2*col+c] = _ex0();
988         }
989     }
990     return sign;
991 }
992
993 /** Partial pivoting method.
994  *  Usual pivoting (symbolic==false) returns the index to the element with the
995  *  largest absolute value in column ro and swaps the current row with the one
996  *  where the element was found.  With (symbolic==true) it does the same thing
997  *  with the first non-zero element.
998  *
999  *  @param ro is the row to be inspected
1000  *  @param symbolic signal if we want the first non-zero element to be pivoted
1001  *  (true) or the one with the largest absolute value (false).
1002  *  @return 0 if no interchange occured, -1 if all are zero (usually signaling
1003  *  a degeneracy) and positive integer k means that rows ro and k were swapped.
1004  */
1005 int matrix::pivot(unsigned ro, bool symbolic)
1006 {
1007     unsigned k = ro;
1008     
1009     if (symbolic) {  // search first non-zero
1010         for (unsigned r=ro; r<row; ++r) {
1011             if (!m[r*col+ro].is_zero()) {
1012                 k = r;
1013                 break;
1014             }
1015         }
1016     } else {  // search largest
1017         numeric tmp(0);
1018         numeric maxn(-1);
1019         for (unsigned r=ro; r<row; ++r) {
1020             GINAC_ASSERT(is_ex_of_type(m[r*col+ro],numeric));
1021             if ((tmp = abs(ex_to_numeric(m[r*col+ro]))) > maxn &&
1022                 !tmp.is_zero()) {
1023                 maxn = tmp;
1024                 k = r;
1025             }
1026         }
1027     }
1028     if (m[k*col+ro].is_zero())
1029         return -1;
1030     if (k!=ro) {  // swap rows
1031         ensure_if_modifiable();
1032         for (unsigned c=0; c<col; ++c) {
1033             m[k*col+c].swap(m[ro*col+c]);
1034         }
1035         return k;
1036     }
1037     return 0;
1038 }
1039
1040 //////////
1041 // global constants
1042 //////////
1043
1044 const matrix some_matrix;
1045 const type_info & typeid_matrix=typeid(some_matrix);
1046
1047 #ifndef NO_NAMESPACE_GINAC
1048 } // namespace GiNaC
1049 #endif // ndef NO_NAMESPACE_GINAC