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