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