]> www.ginac.de Git - ginac.git/blob - ginac/matrix.cpp
- created AUTHORS and README files, updated INSTALL
[ginac.git] / ginac / matrix.cpp
1 /** @file matrix.cpp
2  *
3  *  Implementation of symbolic matrices */
4
5 #include <algorithm>
6 #include <stdexcept>
7
8 #include "ginac.h"
9
10 //////////
11 // default constructor, destructor, copy constructor, assignment operator
12 // and helpers:
13 //////////
14
15 // public
16
17 /** Default ctor.  Initializes to 1 x 1-dimensional zero-matrix. */
18 matrix::matrix()
19     : basic(TINFO_MATRIX), row(1), col(1)
20 {
21     debugmsg("matrix default constructor",LOGLEVEL_CONSTRUCT);
22     m.push_back(exZERO());
23 }
24
25 matrix::~matrix()
26 {
27     debugmsg("matrix destructor",LOGLEVEL_DESTRUCT);
28 }
29
30 matrix::matrix(matrix const & other)
31 {
32     debugmsg("matrix copy constructor",LOGLEVEL_CONSTRUCT);
33     copy(other);
34 }
35
36 matrix const & matrix::operator=(matrix const & other)
37 {
38     debugmsg("matrix operator=",LOGLEVEL_ASSIGNMENT);
39     if (this != &other) {
40         destroy(1);
41         copy(other);
42     }
43     return *this;
44 }
45
46 // protected
47
48 void matrix::copy(matrix const & other)
49 {
50     basic::copy(other);
51     row=other.row;
52     col=other.col;
53     m=other.m;  // use STL's vector copying
54 }
55
56 void matrix::destroy(bool call_parent)
57 {
58     if (call_parent) basic::destroy(call_parent);
59 }
60
61 //////////
62 // other constructors
63 //////////
64
65 // public
66
67 /** Very common ctor.  Initializes to r x c-dimensional zero-matrix.
68  *
69  *  @param r number of rows
70  *  @param c number of cols */
71 matrix::matrix(int r, int c)
72     : basic(TINFO_MATRIX), row(r), col(c)
73 {
74     debugmsg("matrix constructor from int,int",LOGLEVEL_CONSTRUCT);
75     m.resize(r*c, exZERO());
76 }
77
78 // protected
79
80 /** Ctor from representation, for internal use only. */
81 matrix::matrix(int r, int c, vector<ex> const & m2)
82     : basic(TINFO_MATRIX), row(r), col(c), m(m2)
83 {
84     debugmsg("matrix constructor from int,int,vector<ex>",LOGLEVEL_CONSTRUCT);
85 }
86
87 //////////
88 // functions overriding virtual functions from bases classes
89 //////////
90
91 // public
92
93 basic * matrix::duplicate() const
94 {
95     debugmsg("matrix duplicate",LOGLEVEL_DUPLICATE);
96     return new matrix(*this);
97 }
98
99 /** nops is defined to be rows x columns. */
100 int matrix::nops() const
101 {
102     return row*col;
103 }
104
105 /** returns matrix entry at position (i/col, i%col). */
106 ex & matrix::let_op(int const i)
107 {
108     return m[i];
109 }
110
111 /** expands the elements of a matrix entry by entry. */
112 ex matrix::expand(unsigned options) const
113 {
114     vector<ex> tmp(row*col);
115     for (int i=0; i<row*col; ++i) {
116         tmp[i]=m[i].expand(options);
117     }
118     return matrix(row, col, tmp);
119 }
120
121 /** Search ocurrences.  A matrix 'has' an expression if it is the expression
122  *  itself or one of the elements 'has' it. */
123 bool matrix::has(ex const & other) const
124 {
125     ASSERT(other.bp!=0);
126     
127     // tautology: it is the expression itself
128     if (is_equal(*other.bp)) return true;
129     
130     // search all the elements
131     for (vector<ex>::const_iterator r=m.begin(); r!=m.end(); ++r) {
132         if ((*r).has(other)) return true;
133     }
134     return false;
135 }
136
137 /** evaluate matrix entry by entry. */
138 ex matrix::eval(int level) const
139 {
140     debugmsg("matrix eval",LOGLEVEL_MEMBER_FUNCTION);
141     
142     // check if we have to do anything at all
143     if ((level==1)&&(flags & status_flags::evaluated)) {
144         return *this;
145     }
146     
147     // emergency break
148     if (level == -max_recursion_level) {
149         throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
150     }
151     
152     // eval() entry by entry
153     vector<ex> m2(row*col);
154     --level;    
155     for (int r=0; r<row; ++r) {
156         for (int c=0; c<col; ++c) {
157             m2[r*col+c] = m[r*col+c].eval(level);
158         }
159     }
160     
161     return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
162                                                status_flags::evaluated );
163 }
164
165 /** evaluate matrix numerically entry by entry. */
166 ex matrix::evalf(int level) const
167 {
168     debugmsg("matrix evalf",LOGLEVEL_MEMBER_FUNCTION);
169         
170     // check if we have to do anything at all
171     if (level==1) {
172         return *this;
173     }
174     
175     // emergency break
176     if (level == -max_recursion_level) {
177         throw (std::runtime_error("matrix::evalf(): recursion limit exceeded"));
178     }
179     
180     // evalf() entry by entry
181     vector<ex> m2(row*col);
182     --level;
183     for (int r=0; r<row; ++r) {
184         for (int c=0; c<col; ++c) {
185             m2[r*col+c] = m[r*col+c].evalf(level);
186         }
187     }
188     return matrix(row, col, m2);
189 }
190
191 // protected
192
193 int matrix::compare_same_type(basic const & other) const
194 {
195     ASSERT(is_exactly_of_type(other, matrix));
196     matrix const & o=static_cast<matrix &>(const_cast<basic &>(other));
197     
198     // compare number of rows
199     if (row != o.rows()) {
200         return row < o.rows() ? -1 : 1;
201     }
202     
203     // compare number of columns
204     if (col != o.cols()) {
205         return col < o.cols() ? -1 : 1;
206     }
207     
208     // equal number of rows and columns, compare individual elements
209     int cmpval;
210     for (int r=0; r<row; ++r) {
211         for (int c=0; c<col; ++c) {
212             cmpval=((*this)(r,c)).compare(o(r,c));
213             if (cmpval!=0) return cmpval;
214         }
215     }
216     // all elements are equal => matrices are equal;
217     return 0;
218 }
219
220 //////////
221 // non-virtual functions in this class
222 //////////
223
224 // public
225
226 /** Sum of matrices.
227  *
228  *  @exception logic_error (incompatible matrices) */
229 matrix matrix::add(matrix const & other) const
230 {
231     if (col != other.col || row != other.row) {
232         throw (std::logic_error("matrix::add(): incompatible matrices"));
233     }
234     
235     vector<ex> sum(this->m);
236     vector<ex>::iterator i;
237     vector<ex>::const_iterator ci;
238     for (i=sum.begin(), ci=other.m.begin();
239          i!=sum.end();
240          ++i, ++ci) {
241         (*i) += (*ci);
242     }
243     return matrix(row,col,sum);
244 }
245
246 /** Difference of matrices.
247  *
248  *  @exception logic_error (incompatible matrices) */
249 matrix matrix::sub(matrix const & other) const
250 {
251     if (col != other.col || row != other.row) {
252         throw (std::logic_error("matrix::sub(): incompatible matrices"));
253     }
254     
255     vector<ex> dif(this->m);
256     vector<ex>::iterator i;
257     vector<ex>::const_iterator ci;
258     for (i=dif.begin(), ci=other.m.begin();
259          i!=dif.end();
260          ++i, ++ci) {
261         (*i) -= (*ci);
262     }
263     return matrix(row,col,dif);
264 }
265
266 /** Product of matrices.
267  *
268  *  @exception logic_error (incompatible matrices) */
269 matrix matrix::mul(matrix const & other) const
270 {
271     if (col != other.row) {
272         throw (std::logic_error("matrix::mul(): incompatible matrices"));
273     }
274     
275     vector<ex> prod(row*other.col);
276     for (int i=0; i<row; ++i) {
277         for (int j=0; j<other.col; ++j) {
278             for (int l=0; l<col; ++l) {
279                 prod[i*other.col+j] += m[i*col+l] * other.m[l*other.col+j];
280             }
281         }
282     }
283     return matrix(row, other.col, prod);
284 }
285
286 /** operator() to access elements.
287  *
288  *  @param ro row of element
289  *  @param co column of element 
290  *  @exception range_error (index out of range) */
291 ex const & matrix::operator() (int ro, int co) const
292 {
293     if (ro<0 || ro>=row || co<0 || co>=col) {
294         throw (std::range_error("matrix::operator(): index out of range"));
295     }
296     
297     return m[ro*col+co];
298 }
299
300 /** Set individual elements manually.
301  *
302  *  @exception range_error (index out of range) */
303 matrix & matrix::set(int ro, int co, ex value)
304 {
305     if (ro<0 || ro>=row || co<0 || co>=col) {
306         throw (std::range_error("matrix::set(): index out of range"));
307     }
308     
309     ensure_if_modifiable();
310     m[ro*col+co]=value;
311     return *this;
312 }
313
314 /** Transposed of an m x n matrix, producing a new n x m matrix object that
315  *  represents the transposed. */
316 matrix matrix::transpose(void) const
317 {
318     vector<ex> trans(col*row);
319     
320     for (int r=0; r<col; ++r) {
321         for (int c=0; c<row; ++c) {
322             trans[r*row+c] = m[c*col+r];
323         }
324     }
325     return matrix(col,row,trans);
326 }
327
328 /* Determiant of purely numeric matrix, using pivoting. This routine is only
329  * called internally by matrix::determinant(). */
330 ex determinant_numeric(const matrix & M)
331 {
332     ASSERT(M.rows()==M.cols());  // cannot happen, just in case...
333     matrix tmp(M);
334     ex det=exONE();
335     ex piv;
336     
337     for (int r1=0; r1<M.rows(); ++r1) {
338         int indx = tmp.pivot(r1);
339         if (indx == -1) {
340             return exZERO();
341         }
342         if (indx != 0) {
343             det *= exMINUSONE();
344         }
345         det = det * tmp.m[r1*M.cols()+r1];
346         for (int r2=r1+1; r2<M.rows(); ++r2) {
347             piv = tmp.m[r2*M.cols()+r1] / tmp.m[r1*M.cols()+r1];
348             for (int c=r1+1; c<M.cols(); c++) {
349                 tmp.m[r2*M.cols()+c] -= piv * tmp.m[r1*M.cols()+c];
350             }
351         }
352     }
353     return det;
354 }
355
356 // Compute the sign of a permutation of a vector of things, used internally
357 // by determinant_symbolic_perm() where it is instantiated for int.
358 template <class T>
359 int permutation_sign(vector<T> s)
360 {
361     if (s.size() < 2)
362         return 0;
363     int sigma=1;
364     for (typename vector<T>::iterator i=s.begin(); i!=s.end()-1; ++i) {
365         for (typename vector<T>::iterator j=i+1; j!=s.end(); ++j) {
366             if (*i == *j)
367                 return 0;
368             if (*i > *j) {
369                 iter_swap(i,j);
370                 sigma = -sigma;
371             }
372         }
373     }
374     return sigma;
375 }
376
377 /** Determinant built by application of the full permutation group. This
378  *  routine is only called internally by matrix::determinant(). */
379 ex determinant_symbolic_perm(const matrix & M)
380 {
381     ASSERT(M.rows()==M.cols());  // cannot happen, just in case...
382     
383     if (M.rows()==1) {  // speed things up
384         return M(0,0);
385     }
386     
387     ex det;
388     ex term;
389     vector<int> sigma(M.cols());
390     for (int i=0; i<M.cols(); ++i) sigma[i]=i;
391     
392     do {
393         term = M(sigma[0],0);
394         for (int i=1; i<M.cols(); ++i) term *= M(sigma[i],i);
395         det += permutation_sign(sigma)*term;
396     } while (next_permutation(sigma.begin(), sigma.end()));
397     
398     return det;
399 }
400
401 /** Recursive determiant for small matrices having at least one symbolic entry.
402  *  This algorithm is also known as Laplace-expansion. This routine is only
403  *  called internally by matrix::determinant(). */
404 ex determinant_symbolic_minor(const matrix & M)
405 {
406     ASSERT(M.rows()==M.cols());  // cannot happen, just in case...
407     
408     if (M.rows()==1) {  // end of recursion
409         return M(0,0);
410     }
411     if (M.rows()==2) {  // speed things up
412         return (M(0,0)*M(1,1)-
413                 M(1,0)*M(0,1));
414     }
415     if (M.rows()==3) {  // speed things up even a little more
416         return ((M(2,1)*M(0,2)-M(2,2)*M(0,1))*M(1,0)+
417                 (M(1,2)*M(0,1)-M(1,1)*M(0,2))*M(2,0)+
418                 (M(2,2)*M(1,1)-M(2,1)*M(1,2))*M(0,0));
419     }
420     
421     ex det;
422     matrix minorM(M.rows()-1,M.cols()-1);
423     for (int r1=0; r1<M.rows(); ++r1) {
424         // assemble the minor matrix
425         for (int r=0; r<minorM.rows(); ++r) {
426             for (int c=0; c<minorM.cols(); ++c) {
427                 if (r<r1) {
428                     minorM.set(r,c,M(r,c+1));
429                 } else {
430                     minorM.set(r,c,M(r+1,c+1));
431                 }
432             }
433         }
434         // recurse down
435         if (r1%2) {
436             det -= M(r1,0) * determinant_symbolic_minor(minorM);
437         } else {
438             det += M(r1,0) * determinant_symbolic_minor(minorM);
439         }
440     }
441     return det;
442 }
443
444 /*  Leverrier algorithm for large matrices having at least one symbolic entry.
445  *  This routine is only called internally by matrix::determinant(). The
446  *  algorithm is deemed bad for symbolic matrices since it returns expressions
447  *  that are very hard to canonicalize. */
448 /*ex determinant_symbolic_leverrier(const matrix & M)
449  *{
450  *    ASSERT(M.rows()==M.cols());  // cannot happen, just in case...
451  *    
452  *    matrix B(M);
453  *    matrix I(M.row, M.col);
454  *    ex c=B.trace();
455  *    for (int i=1; i<M.row; ++i) {
456  *        for (int j=0; j<M.row; ++j)
457  *            I.m[j*M.col+j] = c;
458  *        B = M.mul(B.sub(I));
459  *        c = B.trace()/ex(i+1);
460  *    }
461  *    if (M.row%2) {
462  *        return c;
463  *    } else {
464  *        return -c;
465  *    }
466  *}*/
467
468 /** Determinant of square matrix.  This routine doesn't actually calculate the
469  *  determinant, it only implements some heuristics about which algorithm to
470  *  call.  When the parameter for normalization is explicitly turned off this
471  *  method does not normalize its result at the end, which might imply that
472  *  the symbolic 2x2 matrix [[a/(a-b),1],[b/(a-b),1]] is not immediatly
473  *  recognized to be unity.  (This is Mathematica's default behaviour, it
474  *  should be used with care.)
475  *
476  *  @param     normalized may be set to false if no normalization of the
477  *             result is desired (i.e. to force Mathematica behavior, Maple
478  *             does normalize the result).
479  *  @return    the determinant as a new expression
480  *  @exception logic_error (matrix not square) */
481 ex matrix::determinant(bool normalized) const
482 {
483     if (row != col) {
484         throw (std::logic_error("matrix::determinant(): matrix not square"));
485     }
486
487     // check, if there are non-numeric entries in the matrix:
488     for (vector<ex>::const_iterator r=m.begin(); r!=m.end(); ++r) {
489         if (!(*r).info(info_flags::numeric)) {
490             if (normalized) {
491                 return determinant_symbolic_minor(*this).normal();
492             } else {
493                 return determinant_symbolic_perm(*this);
494             }
495         }
496     }
497     // if it turns out that all elements are numeric
498     return determinant_numeric(*this);
499 }
500
501 /** Trace of a matrix.
502  *
503  *  @return    the sum of diagonal elements
504  *  @exception logic_error (matrix not square) */
505 ex matrix::trace(void) const
506 {
507     if (row != col) {
508         throw (std::logic_error("matrix::trace(): matrix not square"));
509     }
510     
511     ex tr;
512     for (int r=0; r<col; ++r) {
513         tr += m[r*col+r];
514     }
515     return tr;
516 }
517
518 /** Characteristic Polynomial.  The characteristic polynomial of a matrix M is
519  *  defined as the determiant of (M - lambda * 1) where 1 stands for the unit
520  *  matrix of the same dimension as M.  This method returns the characteristic
521  *  polynomial as a new expression.
522  *
523  *  @return    characteristic polynomial as new expression
524  *  @exception logic_error (matrix not square)
525  *  @see       matrix::determinant() */
526 ex matrix::charpoly(ex const & lambda) const
527 {
528     if (row != col) {
529         throw (std::logic_error("matrix::charpoly(): matrix not square"));
530     }
531     
532     matrix M(*this);
533     for (int r=0; r<col; ++r) {
534         M.m[r*col+r] -= lambda;
535     }
536     return (M.determinant());
537 }
538
539 /** Inverse of this matrix.
540  *
541  *  @return    the inverted matrix
542  *  @exception logic_error (matrix not square)
543  *  @exception runtime_error (singular matrix) */
544 matrix matrix::inverse(void) const
545 {
546     if (row != col) {
547         throw (std::logic_error("matrix::inverse(): matrix not square"));
548     }
549     
550     matrix tmp(row,col);
551     // set tmp to the unit matrix
552     for (int i=0; i<col; ++i) {
553         tmp.m[i*col+i] = exONE();
554     }
555     // create a copy of this matrix
556     matrix cpy(*this);
557     for (int r1=0; r1<row; ++r1) {
558         int indx = cpy.pivot(r1);
559         if (indx == -1) {
560             throw (std::runtime_error("matrix::inverse(): singular matrix"));
561         }
562         if (indx != 0) {  // swap rows r and indx of matrix tmp
563             for (int i=0; i<col; ++i) {
564                 tmp.m[r1*col+i].swap(tmp.m[indx*col+i]);
565             }
566         }
567         ex a1 = cpy.m[r1*col+r1];
568         for (int c=0; c<col; ++c) {
569             cpy.m[r1*col+c] /= a1;
570             tmp.m[r1*col+c] /= a1;
571         }
572         for (int r2=0; r2<row; ++r2) {
573             if (r2 != r1) {
574                 ex a2 = cpy.m[r2*col+r1];
575                 for (int c=0; c<col; ++c) {
576                     cpy.m[r2*col+c] -= a2 * cpy.m[r1*col+c];
577                     tmp.m[r2*col+c] -= a2 * tmp.m[r1*col+c];
578                 }
579             }
580         }
581     }
582     return tmp;
583 }
584
585 void matrix::ffe_swap(int r1, int c1, int r2 ,int c2)
586 {
587     ensure_if_modifiable();
588     
589     ex tmp=ffe_get(r1,c1);
590     ffe_set(r1,c1,ffe_get(r2,c2));
591     ffe_set(r2,c2,tmp);
592 }
593
594 void matrix::ffe_set(int r, int c, ex e)
595 {
596     set(r-1,c-1,e);
597 }
598
599 ex matrix::ffe_get(int r, int c) const
600 {
601     return operator()(r-1,c-1);
602 }
603
604 /** Solve a set of equations for an m x n matrix by fraction-free Gaussian
605  *  elimination. Based on algorithm 9.1 from 'Algorithms for Computer Algebra'
606  *  by Keith O. Geddes et al.
607  *
608  *  @param vars n x p matrix
609  *  @param rhs m x p matrix
610  *  @exception logic_error (incompatible matrices)
611  *  @exception runtime_error (singular matrix) */
612 matrix matrix::fraction_free_elim(matrix const & vars,
613                                   matrix const & rhs) const
614 {
615     if ((row != rhs.row) || (col != vars.row) || (rhs.col != vars.col)) {
616         throw (std::logic_error("matrix::solve(): incompatible matrices"));
617     }
618     
619     matrix a(*this); // make a copy of the matrix
620     matrix b(rhs);     // make a copy of the rhs vector
621     
622     // given an m x n matrix a, reduce it to upper echelon form
623     int m=a.row;
624     int n=a.col;
625     int sign=1;
626     ex divisor=1;
627     int r=1;
628     
629     // eliminate below row r, with pivot in column k
630     for (int k=1; (k<=n)&&(r<=m); ++k) {
631         // find a nonzero pivot
632         int p;
633         for (p=r; (p<=m)&&(a.ffe_get(p,k).is_equal(exZERO())); ++p) {}
634         // pivot is in row p
635         if (p<=m) {
636             if (p!=r) {
637                 // switch rows p and r
638                 for (int j=k; j<=n; ++j) {
639                     a.ffe_swap(p,j,r,j);
640                 }
641                 b.ffe_swap(p,1,r,1);
642                 // keep track of sign changes due to row exchange
643                 sign=-sign;
644             }
645             for (int i=r+1; i<=m; ++i) {
646                 for (int j=k+1; j<=n; ++j) {
647                     a.ffe_set(i,j,(a.ffe_get(r,k)*a.ffe_get(i,j)
648                                   -a.ffe_get(r,j)*a.ffe_get(i,k))/divisor);
649                     a.ffe_set(i,j,a.ffe_get(i,j).normal() /*.normal() */ );
650                 }
651                 b.ffe_set(i,1,(a.ffe_get(r,k)*b.ffe_get(i,1)
652                               -b.ffe_get(r,1)*a.ffe_get(i,k))/divisor);
653                 b.ffe_set(i,1,b.ffe_get(i,1).normal() /*.normal() */ );
654                 a.ffe_set(i,k,0);
655             }
656             divisor=a.ffe_get(r,k);
657             r++;
658         }
659     }
660     // optionally compute the determinant for square or augmented matrices
661     // if (r==m+1) { det=sign*divisor; } else { det=0; }
662     
663     /*
664     for (int r=1; r<=m; ++r) {
665         for (int c=1; c<=n; ++c) {
666             cout << a.ffe_get(r,c) << "\t";
667         }
668         cout << " | " <<  b.ffe_get(r,1) << endl;
669     }
670     */
671     
672 #ifdef DOASSERT
673     // test if we really have an upper echelon matrix
674     int zero_in_last_row=-1;
675     for (int r=1; r<=m; ++r) {
676         int zero_in_this_row=0;
677         for (int c=1; c<=n; ++c) {
678             if (a.ffe_get(r,c).is_equal(exZERO())) {
679                zero_in_this_row++;
680             } else {
681                 break;
682             }
683         }
684         ASSERT((zero_in_this_row>zero_in_last_row)||(zero_in_this_row=n));
685         zero_in_last_row=zero_in_this_row;
686     }
687 #endif // def DOASSERT
688     
689     // assemble solution
690     matrix sol(n,1);
691     int last_assigned_sol=n+1;
692     for (int r=m; r>0; --r) {
693         int first_non_zero=1;
694         while ((first_non_zero<=n)&&(a.ffe_get(r,first_non_zero).is_zero())) {
695             first_non_zero++;
696         }
697         if (first_non_zero>n) {
698             // row consists only of zeroes, corresponding rhs must be 0 as well
699             if (!b.ffe_get(r,1).is_zero()) {
700                 throw (std::runtime_error("matrix::fraction_free_elim(): singular matrix"));
701             }
702         } else {
703             // assign solutions for vars between first_non_zero+1 and
704             // last_assigned_sol-1: free parameters
705             for (int c=first_non_zero+1; c<=last_assigned_sol-1; ++c) {
706                 sol.ffe_set(c,1,vars.ffe_get(c,1));
707             }
708             ex e=b.ffe_get(r,1);
709             for (int c=first_non_zero+1; c<=n; ++c) {
710                 e=e-a.ffe_get(r,c)*sol.ffe_get(c,1);
711             }
712             sol.ffe_set(first_non_zero,1,
713                         (e/a.ffe_get(r,first_non_zero)).normal());
714             last_assigned_sol=first_non_zero;
715         }
716     }
717     // assign solutions for vars between 1 and
718     // last_assigned_sol-1: free parameters
719     for (int c=1; c<=last_assigned_sol-1; ++c) {
720         sol.ffe_set(c,1,vars.ffe_get(c,1));
721     }
722
723     /*
724     for (int c=1; c<=n; ++c) {
725         cout << vars.ffe_get(c,1) << "->" << sol.ffe_get(c,1) << endl;
726     }
727     */
728     
729 #ifdef DOASSERT
730     // test solution with echelon matrix
731     for (int r=1; r<=m; ++r) {
732         ex e=0;
733         for (int c=1; c<=n; ++c) {
734             e=e+a.ffe_get(r,c)*sol.ffe_get(c,1);
735         }
736         if (!(e-b.ffe_get(r,1)).normal().is_zero()) {
737             cout << "e=" << e;
738             cout << "b.ffe_get(" << r<<",1)=" << b.ffe_get(r,1) << endl;
739             cout << "diff=" << (e-b.ffe_get(r,1)).normal() << endl;
740         }
741         ASSERT((e-b.ffe_get(r,1)).normal().is_zero());
742     }
743
744     // test solution with original matrix
745     for (int r=1; r<=m; ++r) {
746         ex e=0;
747         for (int c=1; c<=n; ++c) {
748             e=e+ffe_get(r,c)*sol.ffe_get(c,1);
749         }
750         try {
751         if (!(e-rhs.ffe_get(r,1)).normal().is_zero()) {
752             cout << "e=" << e << endl;
753             e.printtree(cout);
754             ex en=e.normal();
755             cout << "e.normal()=" << en << endl;
756             en.printtree(cout);
757             cout << "rhs.ffe_get(" << r<<",1)=" << rhs.ffe_get(r,1) << endl;
758             cout << "diff=" << (e-rhs.ffe_get(r,1)).normal() << endl;
759         }
760         } catch (...) {
761             ex xxx=e-rhs.ffe_get(r,1);
762             cerr << "xxx=" << xxx << endl << endl;
763         }
764         ASSERT((e-rhs.ffe_get(r,1)).normal().is_zero());
765     }
766 #endif // def DOASSERT
767     
768     return sol;
769 }   
770     
771 /** Solve simultaneous set of equations. */
772 matrix matrix::solve(matrix const & v) const
773 {
774     if (!(row == col && col == v.row)) {
775         throw (std::logic_error("matrix::solve(): incompatible matrices"));
776     }
777     
778     // build the extended matrix of *this with v attached to the right
779     matrix tmp(row,col+v.col);
780     for (int r=0; r<row; ++r) {
781         for (int c=0; c<col; ++c) {
782             tmp.m[r*tmp.col+c] = m[r*col+c];
783         }
784         for (int c=0; c<v.col; ++c) {
785             tmp.m[r*tmp.col+c+col] = v.m[r*v.col+c];
786         }
787     }
788     for (int r1=0; r1<row; ++r1) {
789         int indx = tmp.pivot(r1);
790         if (indx == -1) {
791             throw (std::runtime_error("matrix::solve(): singular matrix"));
792         }
793         for (int c=r1; c<tmp.col; ++c) {
794             tmp.m[r1*tmp.col+c] /= tmp.m[r1*tmp.col+r1];
795         }
796         for (int r2=r1+1; r2<row; ++r2) {
797             for (int c=r1; c<tmp.col; ++c) {
798                 tmp.m[r2*tmp.col+c]
799                     -= tmp.m[r2*tmp.col+r1] * tmp.m[r1*tmp.col+c];
800             }
801         }
802     }
803     
804     // assemble the solution matrix
805     vector<ex> sol(v.row*v.col);
806     for (int c=0; c<v.col; ++c) {
807         for (int r=col-1; r>=0; --r) {
808             sol[r*v.col+c] = tmp[r*tmp.col+c];
809             for (int i=r+1; i<col; ++i) {
810                 sol[r*v.col+c]
811                     -= tmp[r*tmp.col+i] * sol[i*v.col+c];
812             }
813         }
814     }
815     return matrix(v.row, v.col, sol);
816 }
817
818 // protected
819
820 /** Partial pivoting method.
821  *  Usual pivoting returns the index to the element with the largest absolute
822  *  value and swaps the current row with the one where the element was found.
823  *  Here it does the same with the first non-zero element. (This works fine,
824  *  but may be far from optimal for numerics.) */
825 int matrix::pivot(int ro)
826 {
827     int k=ro;
828     
829     for (int r=ro; r<row; ++r) {
830         if (!m[r*col+ro].is_zero()) {
831             k = r;
832             break;
833         }
834     }
835     if (m[k*col+ro].is_zero()) {
836         return -1;
837     }
838     if (k!=ro) {  // swap rows
839         for (int c=0; c<col; ++c) {
840             m[k*col+c].swap(m[ro*col+c]);
841         }
842         return k;
843     }
844     return 0;
845 }
846
847 //////////
848 // global constants
849 //////////
850
851 const matrix some_matrix;
852 type_info const & typeid_matrix=typeid(some_matrix);