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