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