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