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