3 * Implementation of symbolic matrices */
6 * GiNaC Copyright (C) 1999-2003 Johannes Gutenberg University Mainz, Germany
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.
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.
18 * You should have received a copy of the GNU General Public License
19 * along with this program; if not, write to the Free Software
20 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
38 #include "operators.h"
45 GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(matrix, basic,
46 print_func<print_context>(&matrix::do_print).
47 print_func<print_latex>(&matrix::do_print_latex).
48 print_func<print_tree>(&basic::do_print_tree).
49 print_func<print_python_repr>(&matrix::do_print_python_repr))
52 // default constructor
55 /** Default ctor. Initializes to 1 x 1-dimensional zero-matrix. */
56 matrix::matrix() : inherited(TINFO_matrix), row(1), col(1), m(1, _ex0)
58 setflag(status_flags::not_shareable);
67 /** Very common ctor. Initializes to r x c-dimensional zero-matrix.
69 * @param r number of rows
70 * @param c number of cols */
71 matrix::matrix(unsigned r, unsigned c)
72 : inherited(TINFO_matrix), row(r), col(c), m(r*c, _ex0)
74 setflag(status_flags::not_shareable);
79 /** Ctor from representation, for internal use only. */
80 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
81 : inherited(TINFO_matrix), row(r), col(c), m(m2)
83 setflag(status_flags::not_shareable);
86 /** Construct matrix from (flat) list of elements. If the list has fewer
87 * elements than the matrix, the remaining matrix elements are set to zero.
88 * If the list has more elements than the matrix, the excessive elements are
90 matrix::matrix(unsigned r, unsigned c, const lst & l)
91 : inherited(TINFO_matrix), row(r), col(c), m(r*c, _ex0)
93 setflag(status_flags::not_shareable);
96 for (lst::const_iterator it = l.begin(); it != l.end(); ++it, ++i) {
100 break; // matrix smaller than list: throw away excessive elements
109 matrix::matrix(const archive_node &n, lst &sym_lst) : inherited(n, sym_lst)
111 setflag(status_flags::not_shareable);
113 if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
114 throw (std::runtime_error("unknown matrix dimensions in archive"));
115 m.reserve(row * col);
116 for (unsigned int i=0; true; i++) {
118 if (n.find_ex("m", e, sym_lst, i))
125 void matrix::archive(archive_node &n) const
127 inherited::archive(n);
128 n.add_unsigned("row", row);
129 n.add_unsigned("col", col);
130 exvector::const_iterator i = m.begin(), iend = m.end();
137 DEFAULT_UNARCHIVE(matrix)
140 // functions overriding virtual functions from base classes
145 void matrix::print_elements(const print_context & c, const char *row_start, const char *row_end, const char *row_sep, const char *col_sep) const
147 for (unsigned ro=0; ro<row; ++ro) {
149 for (unsigned co=0; co<col; ++co) {
150 m[ro*col+co].print(c);
161 void matrix::do_print(const print_context & c, unsigned level) const
164 print_elements(c, "[", "]", ",", ",");
168 void matrix::do_print_latex(const print_latex & c, unsigned level) const
170 c.s << "\\left(\\begin{array}{" << std::string(col,'c') << "}";
171 print_elements(c, "", "", "\\\\", "&");
172 c.s << "\\end{array}\\right)";
175 void matrix::do_print_python_repr(const print_python_repr & c, unsigned level) const
177 c.s << class_name() << '(';
178 print_elements(c, "[", "]", ",", ",");
182 /** nops is defined to be rows x columns. */
183 size_t matrix::nops() const
185 return static_cast<size_t>(row) * static_cast<size_t>(col);
188 /** returns matrix entry at position (i/col, i%col). */
189 ex matrix::op(size_t i) const
191 GINAC_ASSERT(i<nops());
196 /** returns writable matrix entry at position (i/col, i%col). */
197 ex & matrix::let_op(size_t i)
199 GINAC_ASSERT(i<nops());
201 ensure_if_modifiable();
205 /** Evaluate matrix entry by entry. */
206 ex matrix::eval(int level) const
208 // check if we have to do anything at all
209 if ((level==1)&&(flags & status_flags::evaluated))
213 if (level == -max_recursion_level)
214 throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
216 // eval() entry by entry
217 exvector m2(row*col);
219 for (unsigned r=0; r<row; ++r)
220 for (unsigned c=0; c<col; ++c)
221 m2[r*col+c] = m[r*col+c].eval(level);
223 return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
224 status_flags::evaluated);
227 ex matrix::subs(const exmap & mp, unsigned options) const
229 exvector m2(row * col);
230 for (unsigned r=0; r<row; ++r)
231 for (unsigned c=0; c<col; ++c)
232 m2[r*col+c] = m[r*col+c].subs(mp, options);
234 return matrix(row, col, m2).subs_one_level(mp, options);
239 int matrix::compare_same_type(const basic & other) const
241 GINAC_ASSERT(is_exactly_a<matrix>(other));
242 const matrix &o = static_cast<const matrix &>(other);
244 // compare number of rows
246 return row < o.rows() ? -1 : 1;
248 // compare number of columns
250 return col < o.cols() ? -1 : 1;
252 // equal number of rows and columns, compare individual elements
254 for (unsigned r=0; r<row; ++r) {
255 for (unsigned c=0; c<col; ++c) {
256 cmpval = ((*this)(r,c)).compare(o(r,c));
257 if (cmpval!=0) return cmpval;
260 // all elements are equal => matrices are equal;
264 bool matrix::match_same_type(const basic & other) const
266 GINAC_ASSERT(is_exactly_a<matrix>(other));
267 const matrix & o = static_cast<const matrix &>(other);
269 // The number of rows and columns must be the same. This is necessary to
270 // prevent a 2x3 matrix from matching a 3x2 one.
271 return row == o.rows() && col == o.cols();
274 /** Automatic symbolic evaluation of an indexed matrix. */
275 ex matrix::eval_indexed(const basic & i) const
277 GINAC_ASSERT(is_a<indexed>(i));
278 GINAC_ASSERT(is_a<matrix>(i.op(0)));
280 bool all_indices_unsigned = static_cast<const indexed &>(i).all_index_values_are(info_flags::nonnegint);
285 // One index, must be one-dimensional vector
286 if (row != 1 && col != 1)
287 throw (std::runtime_error("matrix::eval_indexed(): vector must have exactly 1 index"));
289 const idx & i1 = ex_to<idx>(i.op(1));
294 if (!i1.get_dim().is_equal(row))
295 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
297 // Index numeric -> return vector element
298 if (all_indices_unsigned) {
299 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
301 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
302 return (*this)(n1, 0);
308 if (!i1.get_dim().is_equal(col))
309 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
311 // Index numeric -> return vector element
312 if (all_indices_unsigned) {
313 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
315 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
316 return (*this)(0, n1);
320 } else if (i.nops() == 3) {
323 const idx & i1 = ex_to<idx>(i.op(1));
324 const idx & i2 = ex_to<idx>(i.op(2));
326 if (!i1.get_dim().is_equal(row))
327 throw (std::runtime_error("matrix::eval_indexed(): dimension of first index must match number of rows"));
328 if (!i2.get_dim().is_equal(col))
329 throw (std::runtime_error("matrix::eval_indexed(): dimension of second index must match number of columns"));
331 // Pair of dummy indices -> compute trace
332 if (is_dummy_pair(i1, i2))
335 // Both indices numeric -> return matrix element
336 if (all_indices_unsigned) {
337 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int(), n2 = ex_to<numeric>(i2.get_value()).to_int();
339 throw (std::runtime_error("matrix::eval_indexed(): value of first index exceeds number of rows"));
341 throw (std::runtime_error("matrix::eval_indexed(): value of second index exceeds number of columns"));
342 return (*this)(n1, n2);
346 throw (std::runtime_error("matrix::eval_indexed(): matrix must have exactly 2 indices"));
351 /** Sum of two indexed matrices. */
352 ex matrix::add_indexed(const ex & self, const ex & other) const
354 GINAC_ASSERT(is_a<indexed>(self));
355 GINAC_ASSERT(is_a<matrix>(self.op(0)));
356 GINAC_ASSERT(is_a<indexed>(other));
357 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
359 // Only add two matrices
360 if (is_a<matrix>(other.op(0))) {
361 GINAC_ASSERT(other.nops() == 2 || other.nops() == 3);
363 const matrix &self_matrix = ex_to<matrix>(self.op(0));
364 const matrix &other_matrix = ex_to<matrix>(other.op(0));
366 if (self.nops() == 2 && other.nops() == 2) { // vector + vector
368 if (self_matrix.row == other_matrix.row)
369 return indexed(self_matrix.add(other_matrix), self.op(1));
370 else if (self_matrix.row == other_matrix.col)
371 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1));
373 } else if (self.nops() == 3 && other.nops() == 3) { // matrix + matrix
375 if (self.op(1).is_equal(other.op(1)) && self.op(2).is_equal(other.op(2)))
376 return indexed(self_matrix.add(other_matrix), self.op(1), self.op(2));
377 else if (self.op(1).is_equal(other.op(2)) && self.op(2).is_equal(other.op(1)))
378 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1), self.op(2));
383 // Don't know what to do, return unevaluated sum
387 /** Product of an indexed matrix with a number. */
388 ex matrix::scalar_mul_indexed(const ex & self, const numeric & other) const
390 GINAC_ASSERT(is_a<indexed>(self));
391 GINAC_ASSERT(is_a<matrix>(self.op(0)));
392 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
394 const matrix &self_matrix = ex_to<matrix>(self.op(0));
396 if (self.nops() == 2)
397 return indexed(self_matrix.mul(other), self.op(1));
398 else // self.nops() == 3
399 return indexed(self_matrix.mul(other), self.op(1), self.op(2));
402 /** Contraction of an indexed matrix with something else. */
403 bool matrix::contract_with(exvector::iterator self, exvector::iterator other, exvector & v) const
405 GINAC_ASSERT(is_a<indexed>(*self));
406 GINAC_ASSERT(is_a<indexed>(*other));
407 GINAC_ASSERT(self->nops() == 2 || self->nops() == 3);
408 GINAC_ASSERT(is_a<matrix>(self->op(0)));
410 // Only contract with other matrices
411 if (!is_a<matrix>(other->op(0)))
414 GINAC_ASSERT(other->nops() == 2 || other->nops() == 3);
416 const matrix &self_matrix = ex_to<matrix>(self->op(0));
417 const matrix &other_matrix = ex_to<matrix>(other->op(0));
419 if (self->nops() == 2) {
421 if (other->nops() == 2) { // vector * vector (scalar product)
423 if (self_matrix.col == 1) {
424 if (other_matrix.col == 1) {
425 // Column vector * column vector, transpose first vector
426 *self = self_matrix.transpose().mul(other_matrix)(0, 0);
428 // Column vector * row vector, swap factors
429 *self = other_matrix.mul(self_matrix)(0, 0);
432 if (other_matrix.col == 1) {
433 // Row vector * column vector, perfect
434 *self = self_matrix.mul(other_matrix)(0, 0);
436 // Row vector * row vector, transpose second vector
437 *self = self_matrix.mul(other_matrix.transpose())(0, 0);
443 } else { // vector * matrix
445 // B_i * A_ij = (B*A)_j (B is row vector)
446 if (is_dummy_pair(self->op(1), other->op(1))) {
447 if (self_matrix.row == 1)
448 *self = indexed(self_matrix.mul(other_matrix), other->op(2));
450 *self = indexed(self_matrix.transpose().mul(other_matrix), other->op(2));
455 // B_j * A_ij = (A*B)_i (B is column vector)
456 if (is_dummy_pair(self->op(1), other->op(2))) {
457 if (self_matrix.col == 1)
458 *self = indexed(other_matrix.mul(self_matrix), other->op(1));
460 *self = indexed(other_matrix.mul(self_matrix.transpose()), other->op(1));
466 } else if (other->nops() == 3) { // matrix * matrix
468 // A_ij * B_jk = (A*B)_ik
469 if (is_dummy_pair(self->op(2), other->op(1))) {
470 *self = indexed(self_matrix.mul(other_matrix), self->op(1), other->op(2));
475 // A_ij * B_kj = (A*Btrans)_ik
476 if (is_dummy_pair(self->op(2), other->op(2))) {
477 *self = indexed(self_matrix.mul(other_matrix.transpose()), self->op(1), other->op(1));
482 // A_ji * B_jk = (Atrans*B)_ik
483 if (is_dummy_pair(self->op(1), other->op(1))) {
484 *self = indexed(self_matrix.transpose().mul(other_matrix), self->op(2), other->op(2));
489 // A_ji * B_kj = (B*A)_ki
490 if (is_dummy_pair(self->op(1), other->op(2))) {
491 *self = indexed(other_matrix.mul(self_matrix), other->op(1), self->op(2));
502 // non-virtual functions in this class
509 * @exception logic_error (incompatible matrices) */
510 matrix matrix::add(const matrix & other) const
512 if (col != other.col || row != other.row)
513 throw std::logic_error("matrix::add(): incompatible matrices");
515 exvector sum(this->m);
516 exvector::iterator i = sum.begin(), end = sum.end();
517 exvector::const_iterator ci = other.m.begin();
521 return matrix(row,col,sum);
525 /** Difference of matrices.
527 * @exception logic_error (incompatible matrices) */
528 matrix matrix::sub(const matrix & other) const
530 if (col != other.col || row != other.row)
531 throw std::logic_error("matrix::sub(): incompatible matrices");
533 exvector dif(this->m);
534 exvector::iterator i = dif.begin(), end = dif.end();
535 exvector::const_iterator ci = other.m.begin();
539 return matrix(row,col,dif);
543 /** Product of matrices.
545 * @exception logic_error (incompatible matrices) */
546 matrix matrix::mul(const matrix & other) const
548 if (this->cols() != other.rows())
549 throw std::logic_error("matrix::mul(): incompatible matrices");
551 exvector prod(this->rows()*other.cols());
553 for (unsigned r1=0; r1<this->rows(); ++r1) {
554 for (unsigned c=0; c<this->cols(); ++c) {
555 if (m[r1*col+c].is_zero())
557 for (unsigned r2=0; r2<other.cols(); ++r2)
558 prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]).expand();
561 return matrix(row, other.col, prod);
565 /** Product of matrix and scalar. */
566 matrix matrix::mul(const numeric & other) const
568 exvector prod(row * col);
570 for (unsigned r=0; r<row; ++r)
571 for (unsigned c=0; c<col; ++c)
572 prod[r*col+c] = m[r*col+c] * other;
574 return matrix(row, col, prod);
578 /** Product of matrix and scalar expression. */
579 matrix matrix::mul_scalar(const ex & other) const
581 if (other.return_type() != return_types::commutative)
582 throw std::runtime_error("matrix::mul_scalar(): non-commutative scalar");
584 exvector prod(row * col);
586 for (unsigned r=0; r<row; ++r)
587 for (unsigned c=0; c<col; ++c)
588 prod[r*col+c] = m[r*col+c] * other;
590 return matrix(row, col, prod);
594 /** Power of a matrix. Currently handles integer exponents only. */
595 matrix matrix::pow(const ex & expn) const
598 throw (std::logic_error("matrix::pow(): matrix not square"));
600 if (is_exactly_a<numeric>(expn)) {
601 // Integer cases are computed by successive multiplication, using the
602 // obvious shortcut of storing temporaries, like A^4 == (A*A)*(A*A).
603 if (expn.info(info_flags::integer)) {
604 numeric b = ex_to<numeric>(expn);
606 if (expn.info(info_flags::negative)) {
613 for (unsigned r=0; r<row; ++r)
617 // This loop computes the representation of b in base 2 from right
618 // to left and multiplies the factors whenever needed. Note
619 // that this is not entirely optimal but close to optimal and
620 // "better" algorithms are much harder to implement. (See Knuth,
621 // TAoCP2, section "Evaluation of Powers" for a good discussion.)
627 b /= _num2; // still integer.
633 throw (std::runtime_error("matrix::pow(): don't know how to handle exponent"));
637 /** operator() to access elements for reading.
639 * @param ro row of element
640 * @param co column of element
641 * @exception range_error (index out of range) */
642 const ex & matrix::operator() (unsigned ro, unsigned co) const
644 if (ro>=row || co>=col)
645 throw (std::range_error("matrix::operator(): index out of range"));
651 /** operator() to access elements for writing.
653 * @param ro row of element
654 * @param co column of element
655 * @exception range_error (index out of range) */
656 ex & matrix::operator() (unsigned ro, unsigned co)
658 if (ro>=row || co>=col)
659 throw (std::range_error("matrix::operator(): index out of range"));
661 ensure_if_modifiable();
666 /** Transposed of an m x n matrix, producing a new n x m matrix object that
667 * represents the transposed. */
668 matrix matrix::transpose() const
670 exvector trans(this->cols()*this->rows());
672 for (unsigned r=0; r<this->cols(); ++r)
673 for (unsigned c=0; c<this->rows(); ++c)
674 trans[r*this->rows()+c] = m[c*this->cols()+r];
676 return matrix(this->cols(),this->rows(),trans);
679 /** Determinant of square matrix. This routine doesn't actually calculate the
680 * determinant, it only implements some heuristics about which algorithm to
681 * run. If all the elements of the matrix are elements of an integral domain
682 * the determinant is also in that integral domain and the result is expanded
683 * only. If one or more elements are from a quotient field the determinant is
684 * usually also in that quotient field and the result is normalized before it
685 * is returned. This implies that the determinant of the symbolic 2x2 matrix
686 * [[a/(a-b),1],[b/(a-b),1]] is returned as unity. (In this respect, it
687 * behaves like MapleV and unlike Mathematica.)
689 * @param algo allows to chose an algorithm
690 * @return the determinant as a new expression
691 * @exception logic_error (matrix not square)
692 * @see determinant_algo */
693 ex matrix::determinant(unsigned algo) const
696 throw (std::logic_error("matrix::determinant(): matrix not square"));
697 GINAC_ASSERT(row*col==m.capacity());
699 // Gather some statistical information about this matrix:
700 bool numeric_flag = true;
701 bool normal_flag = false;
702 unsigned sparse_count = 0; // counts non-zero elements
703 exvector::const_iterator r = m.begin(), rend = m.end();
705 exmap srl; // symbol replacement list
706 ex rtest = r->to_rational(srl);
707 if (!rtest.is_zero())
709 if (!rtest.info(info_flags::numeric))
710 numeric_flag = false;
711 if (!rtest.info(info_flags::crational_polynomial) &&
712 rtest.info(info_flags::rational_function))
717 // Here is the heuristics in case this routine has to decide:
718 if (algo == determinant_algo::automatic) {
719 // Minor expansion is generally a good guess:
720 algo = determinant_algo::laplace;
721 // Does anybody know when a matrix is really sparse?
722 // Maybe <~row/2.236 nonzero elements average in a row?
723 if (row>3 && 5*sparse_count<=row*col)
724 algo = determinant_algo::bareiss;
725 // Purely numeric matrix can be handled by Gauss elimination.
726 // This overrides any prior decisions.
728 algo = determinant_algo::gauss;
731 // Trap the trivial case here, since some algorithms don't like it
733 // for consistency with non-trivial determinants...
735 return m[0].normal();
737 return m[0].expand();
740 // Compute the determinant
742 case determinant_algo::gauss: {
745 int sign = tmp.gauss_elimination(true);
746 for (unsigned d=0; d<row; ++d)
747 det *= tmp.m[d*col+d];
749 return (sign*det).normal();
751 return (sign*det).normal().expand();
753 case determinant_algo::bareiss: {
756 sign = tmp.fraction_free_elimination(true);
758 return (sign*tmp.m[row*col-1]).normal();
760 return (sign*tmp.m[row*col-1]).expand();
762 case determinant_algo::divfree: {
765 sign = tmp.division_free_elimination(true);
768 ex det = tmp.m[row*col-1];
769 // factor out accumulated bogus slag
770 for (unsigned d=0; d<row-2; ++d)
771 for (unsigned j=0; j<row-d-2; ++j)
772 det = (det/tmp.m[d*col+d]).normal();
775 case determinant_algo::laplace:
777 // This is the minor expansion scheme. We always develop such
778 // that the smallest minors (i.e, the trivial 1x1 ones) are on the
779 // rightmost column. For this to be efficient, empirical tests
780 // have shown that the emptiest columns (i.e. the ones with most
781 // zeros) should be the ones on the right hand side -- although
782 // this might seem counter-intuitive (and in contradiction to some
783 // literature like the FORM manual). Please go ahead and test it
784 // if you don't believe me! Therefore we presort the columns of
786 typedef std::pair<unsigned,unsigned> uintpair;
787 std::vector<uintpair> c_zeros; // number of zeros in column
788 for (unsigned c=0; c<col; ++c) {
790 for (unsigned r=0; r<row; ++r)
791 if (m[r*col+c].is_zero())
793 c_zeros.push_back(uintpair(acc,c));
795 std::sort(c_zeros.begin(),c_zeros.end());
796 std::vector<unsigned> pre_sort;
797 for (std::vector<uintpair>::const_iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
798 pre_sort.push_back(i->second);
799 std::vector<unsigned> pre_sort_test(pre_sort); // permutation_sign() modifies the vector so we make a copy here
800 int sign = permutation_sign(pre_sort_test.begin(), pre_sort_test.end());
801 exvector result(row*col); // represents sorted matrix
803 for (std::vector<unsigned>::const_iterator i=pre_sort.begin();
806 for (unsigned r=0; r<row; ++r)
807 result[r*col+c] = m[r*col+(*i)];
811 return (sign*matrix(row,col,result).determinant_minor()).normal();
813 return sign*matrix(row,col,result).determinant_minor();
819 /** Trace of a matrix. The result is normalized if it is in some quotient
820 * field and expanded only otherwise. This implies that the trace of the
821 * symbolic 2x2 matrix [[a/(a-b),x],[y,b/(b-a)]] is recognized to be unity.
823 * @return the sum of diagonal elements
824 * @exception logic_error (matrix not square) */
825 ex matrix::trace() const
828 throw (std::logic_error("matrix::trace(): matrix not square"));
831 for (unsigned r=0; r<col; ++r)
834 if (tr.info(info_flags::rational_function) &&
835 !tr.info(info_flags::crational_polynomial))
842 /** Characteristic Polynomial. Following mathematica notation the
843 * characteristic polynomial of a matrix M is defined as the determiant of
844 * (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
845 * as M. Note that some CASs define it with a sign inside the determinant
846 * which gives rise to an overall sign if the dimension is odd. This method
847 * returns the characteristic polynomial collected in powers of lambda as a
850 * @return characteristic polynomial as new expression
851 * @exception logic_error (matrix not square)
852 * @see matrix::determinant() */
853 ex matrix::charpoly(const ex & lambda) const
856 throw (std::logic_error("matrix::charpoly(): matrix not square"));
858 bool numeric_flag = true;
859 exvector::const_iterator r = m.begin(), rend = m.end();
860 while (r!=rend && numeric_flag==true) {
861 if (!r->info(info_flags::numeric))
862 numeric_flag = false;
866 // The pure numeric case is traditionally rather common. Hence, it is
867 // trapped and we use Leverrier's algorithm which goes as row^3 for
868 // every coefficient. The expensive part is the matrix multiplication.
873 ex poly = power(lambda, row) - c*power(lambda, row-1);
874 for (unsigned i=1; i<row; ++i) {
875 for (unsigned j=0; j<row; ++j)
878 c = B.trace() / ex(i+1);
879 poly -= c*power(lambda, row-i-1);
889 for (unsigned r=0; r<col; ++r)
890 M.m[r*col+r] -= lambda;
892 return M.determinant().collect(lambda);
897 /** Inverse of this matrix.
899 * @return the inverted matrix
900 * @exception logic_error (matrix not square)
901 * @exception runtime_error (singular matrix) */
902 matrix matrix::inverse() const
905 throw (std::logic_error("matrix::inverse(): matrix not square"));
907 // This routine actually doesn't do anything fancy at all. We compute the
908 // inverse of the matrix A by solving the system A * A^{-1} == Id.
910 // First populate the identity matrix supposed to become the right hand side.
911 matrix identity(row,col);
912 for (unsigned i=0; i<row; ++i)
913 identity(i,i) = _ex1;
915 // Populate a dummy matrix of variables, just because of compatibility with
916 // matrix::solve() which wants this (for compatibility with under-determined
917 // systems of equations).
918 matrix vars(row,col);
919 for (unsigned r=0; r<row; ++r)
920 for (unsigned c=0; c<col; ++c)
921 vars(r,c) = symbol();
925 sol = this->solve(vars,identity);
926 } catch (const std::runtime_error & e) {
927 if (e.what()==std::string("matrix::solve(): inconsistent linear system"))
928 throw (std::runtime_error("matrix::inverse(): singular matrix"));
936 /** Solve a linear system consisting of a m x n matrix and a m x p right hand
937 * side by applying an elimination scheme to the augmented matrix.
939 * @param vars n x p matrix, all elements must be symbols
940 * @param rhs m x p matrix
941 * @return n x p solution matrix
942 * @exception logic_error (incompatible matrices)
943 * @exception invalid_argument (1st argument must be matrix of symbols)
944 * @exception runtime_error (inconsistent linear system)
946 matrix matrix::solve(const matrix & vars,
950 const unsigned m = this->rows();
951 const unsigned n = this->cols();
952 const unsigned p = rhs.cols();
955 if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
956 throw (std::logic_error("matrix::solve(): incompatible matrices"));
957 for (unsigned ro=0; ro<n; ++ro)
958 for (unsigned co=0; co<p; ++co)
959 if (!vars(ro,co).info(info_flags::symbol))
960 throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
962 // build the augmented matrix of *this with rhs attached to the right
964 for (unsigned r=0; r<m; ++r) {
965 for (unsigned c=0; c<n; ++c)
966 aug.m[r*(n+p)+c] = this->m[r*n+c];
967 for (unsigned c=0; c<p; ++c)
968 aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
971 // Gather some statistical information about the augmented matrix:
972 bool numeric_flag = true;
973 exvector::const_iterator r = aug.m.begin(), rend = aug.m.end();
974 while (r!=rend && numeric_flag==true) {
975 if (!r->info(info_flags::numeric))
976 numeric_flag = false;
980 // Here is the heuristics in case this routine has to decide:
981 if (algo == solve_algo::automatic) {
982 // Bareiss (fraction-free) elimination is generally a good guess:
983 algo = solve_algo::bareiss;
984 // For m<3, Bareiss elimination is equivalent to division free
985 // elimination but has more logistic overhead
987 algo = solve_algo::divfree;
988 // This overrides any prior decisions.
990 algo = solve_algo::gauss;
993 // Eliminate the augmented matrix:
995 case solve_algo::gauss:
996 aug.gauss_elimination();
998 case solve_algo::divfree:
999 aug.division_free_elimination();
1001 case solve_algo::bareiss:
1003 aug.fraction_free_elimination();
1006 // assemble the solution matrix:
1008 for (unsigned co=0; co<p; ++co) {
1009 unsigned last_assigned_sol = n+1;
1010 for (int r=m-1; r>=0; --r) {
1011 unsigned fnz = 1; // first non-zero in row
1012 while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
1015 // row consists only of zeros, corresponding rhs must be 0, too
1016 if (!aug.m[r*(n+p)+n+co].is_zero()) {
1017 throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
1020 // assign solutions for vars between fnz+1 and
1021 // last_assigned_sol-1: free parameters
1022 for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
1023 sol(c,co) = vars.m[c*p+co];
1024 ex e = aug.m[r*(n+p)+n+co];
1025 for (unsigned c=fnz; c<n; ++c)
1026 e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
1027 sol(fnz-1,co) = (e/(aug.m[r*(n+p)+(fnz-1)])).normal();
1028 last_assigned_sol = fnz;
1031 // assign solutions for vars between 1 and
1032 // last_assigned_sol-1: free parameters
1033 for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
1034 sol(ro,co) = vars(ro,co);
1043 /** Recursive determinant for small matrices having at least one symbolic
1044 * entry. The basic algorithm, known as Laplace-expansion, is enhanced by
1045 * some bookkeeping to avoid calculation of the same submatrices ("minors")
1046 * more than once. According to W.M.Gentleman and S.C.Johnson this algorithm
1047 * is better than elimination schemes for matrices of sparse multivariate
1048 * polynomials and also for matrices of dense univariate polynomials if the
1049 * matrix' dimesion is larger than 7.
1051 * @return the determinant as a new expression (in expanded form)
1052 * @see matrix::determinant() */
1053 ex matrix::determinant_minor() const
1055 // for small matrices the algorithm does not make any sense:
1056 const unsigned n = this->cols();
1058 return m[0].expand();
1060 return (m[0]*m[3]-m[2]*m[1]).expand();
1062 return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
1063 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
1064 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
1066 // This algorithm can best be understood by looking at a naive
1067 // implementation of Laplace-expansion, like this one:
1069 // matrix minorM(this->rows()-1,this->cols()-1);
1070 // for (unsigned r1=0; r1<this->rows(); ++r1) {
1071 // // shortcut if element(r1,0) vanishes
1072 // if (m[r1*col].is_zero())
1074 // // assemble the minor matrix
1075 // for (unsigned r=0; r<minorM.rows(); ++r) {
1076 // for (unsigned c=0; c<minorM.cols(); ++c) {
1078 // minorM(r,c) = m[r*col+c+1];
1080 // minorM(r,c) = m[(r+1)*col+c+1];
1083 // // recurse down and care for sign:
1085 // det -= m[r1*col] * minorM.determinant_minor();
1087 // det += m[r1*col] * minorM.determinant_minor();
1089 // return det.expand();
1090 // What happens is that while proceeding down many of the minors are
1091 // computed more than once. In particular, there are binomial(n,k)
1092 // kxk minors and each one is computed factorial(n-k) times. Therefore
1093 // it is reasonable to store the results of the minors. We proceed from
1094 // right to left. At each column c we only need to retrieve the minors
1095 // calculated in step c-1. We therefore only have to store at most
1096 // 2*binomial(n,n/2) minors.
1098 // Unique flipper counter for partitioning into minors
1099 std::vector<unsigned> Pkey;
1101 // key for minor determinant (a subpartition of Pkey)
1102 std::vector<unsigned> Mkey;
1104 // we store our subminors in maps, keys being the rows they arise from
1105 typedef std::map<std::vector<unsigned>,class ex> Rmap;
1106 typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
1110 // initialize A with last column:
1111 for (unsigned r=0; r<n; ++r) {
1112 Pkey.erase(Pkey.begin(),Pkey.end());
1114 A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
1116 // proceed from right to left through matrix
1117 for (int c=n-2; c>=0; --c) {
1118 Pkey.erase(Pkey.begin(),Pkey.end()); // don't change capacity
1119 Mkey.erase(Mkey.begin(),Mkey.end());
1120 for (unsigned i=0; i<n-c; ++i)
1122 unsigned fc = 0; // controls logic for our strange flipper counter
1125 for (unsigned r=0; r<n-c; ++r) {
1126 // maybe there is nothing to do?
1127 if (m[Pkey[r]*n+c].is_zero())
1129 // create the sorted key for all possible minors
1130 Mkey.erase(Mkey.begin(),Mkey.end());
1131 for (unsigned i=0; i<n-c; ++i)
1133 Mkey.push_back(Pkey[i]);
1134 // Fetch the minors and compute the new determinant
1136 det -= m[Pkey[r]*n+c]*A[Mkey];
1138 det += m[Pkey[r]*n+c]*A[Mkey];
1140 // prevent build-up of deep nesting of expressions saves time:
1142 // store the new determinant at its place in B:
1144 B.insert(Rmap_value(Pkey,det));
1145 // increment our strange flipper counter
1146 for (fc=n-c; fc>0; --fc) {
1148 if (Pkey[fc-1]<fc+c)
1152 for (unsigned j=fc; j<n-c; ++j)
1153 Pkey[j] = Pkey[j-1]+1;
1155 // next column, so change the role of A and B:
1164 /** Perform the steps of an ordinary Gaussian elimination to bring the m x n
1165 * matrix into an upper echelon form. The algorithm is ok for matrices
1166 * with numeric coefficients but quite unsuited for symbolic matrices.
1168 * @param det may be set to true to save a lot of space if one is only
1169 * interested in the diagonal elements (i.e. for calculating determinants).
1170 * The others are set to zero in this case.
1171 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1172 * number of rows was swapped and 0 if the matrix is singular. */
1173 int matrix::gauss_elimination(const bool det)
1175 ensure_if_modifiable();
1176 const unsigned m = this->rows();
1177 const unsigned n = this->cols();
1178 GINAC_ASSERT(!det || n==m);
1182 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
1183 int indx = pivot(r0, r1, true);
1187 return 0; // leaves *this in a messy state
1192 for (unsigned r2=r0+1; r2<m; ++r2) {
1193 if (!this->m[r2*n+r1].is_zero()) {
1194 // yes, there is something to do in this row
1195 ex piv = this->m[r2*n+r1] / this->m[r0*n+r1];
1196 for (unsigned c=r1+1; c<n; ++c) {
1197 this->m[r2*n+c] -= piv * this->m[r0*n+c];
1198 if (!this->m[r2*n+c].info(info_flags::numeric))
1199 this->m[r2*n+c] = this->m[r2*n+c].normal();
1202 // fill up left hand side with zeros
1203 for (unsigned c=0; c<=r1; ++c)
1204 this->m[r2*n+c] = _ex0;
1207 // save space by deleting no longer needed elements
1208 for (unsigned c=r0+1; c<n; ++c)
1209 this->m[r0*n+c] = _ex0;
1219 /** Perform the steps of division free elimination to bring the m x n matrix
1220 * into an upper echelon form.
1222 * @param det may be set to true to save a lot of space if one is only
1223 * interested in the diagonal elements (i.e. for calculating determinants).
1224 * The others are set to zero in this case.
1225 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1226 * number of rows was swapped and 0 if the matrix is singular. */
1227 int matrix::division_free_elimination(const bool det)
1229 ensure_if_modifiable();
1230 const unsigned m = this->rows();
1231 const unsigned n = this->cols();
1232 GINAC_ASSERT(!det || n==m);
1236 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
1237 int indx = pivot(r0, r1, true);
1241 return 0; // leaves *this in a messy state
1246 for (unsigned r2=r0+1; r2<m; ++r2) {
1247 for (unsigned c=r1+1; c<n; ++c)
1248 this->m[r2*n+c] = (this->m[r0*n+r1]*this->m[r2*n+c] - this->m[r2*n+r1]*this->m[r0*n+c]).expand();
1249 // fill up left hand side with zeros
1250 for (unsigned c=0; c<=r1; ++c)
1251 this->m[r2*n+c] = _ex0;
1254 // save space by deleting no longer needed elements
1255 for (unsigned c=r0+1; c<n; ++c)
1256 this->m[r0*n+c] = _ex0;
1266 /** Perform the steps of Bareiss' one-step fraction free elimination to bring
1267 * the matrix into an upper echelon form. Fraction free elimination means
1268 * that divide is used straightforwardly, without computing GCDs first. This
1269 * is possible, since we know the divisor at each step.
1271 * @param det may be set to true to save a lot of space if one is only
1272 * interested in the last element (i.e. for calculating determinants). The
1273 * others are set to zero in this case.
1274 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1275 * number of rows was swapped and 0 if the matrix is singular. */
1276 int matrix::fraction_free_elimination(const bool det)
1279 // (single-step fraction free elimination scheme, already known to Jordan)
1281 // Usual division-free elimination sets m[0](r,c) = m(r,c) and then sets
1282 // m[k+1](r,c) = m[k](k,k) * m[k](r,c) - m[k](r,k) * m[k](k,c).
1284 // Bareiss (fraction-free) elimination in addition divides that element
1285 // by m[k-1](k-1,k-1) for k>1, where it can be shown by means of the
1286 // Sylvester identity that this really divides m[k+1](r,c).
1288 // We also allow rational functions where the original prove still holds.
1289 // However, we must care for numerator and denominator separately and
1290 // "manually" work in the integral domains because of subtle cancellations
1291 // (see below). This blows up the bookkeeping a bit and the formula has
1292 // to be modified to expand like this (N{x} stands for numerator of x,
1293 // D{x} for denominator of x):
1294 // N{m[k+1](r,c)} = N{m[k](k,k)}*N{m[k](r,c)}*D{m[k](r,k)}*D{m[k](k,c)}
1295 // -N{m[k](r,k)}*N{m[k](k,c)}*D{m[k](k,k)}*D{m[k](r,c)}
1296 // D{m[k+1](r,c)} = D{m[k](k,k)}*D{m[k](r,c)}*D{m[k](r,k)}*D{m[k](k,c)}
1297 // where for k>1 we now divide N{m[k+1](r,c)} by
1298 // N{m[k-1](k-1,k-1)}
1299 // and D{m[k+1](r,c)} by
1300 // D{m[k-1](k-1,k-1)}.
1302 ensure_if_modifiable();
1303 const unsigned m = this->rows();
1304 const unsigned n = this->cols();
1305 GINAC_ASSERT(!det || n==m);
1314 // We populate temporary matrices to subsequently operate on. There is
1315 // one holding numerators and another holding denominators of entries.
1316 // This is a must since the evaluator (or even earlier mul's constructor)
1317 // might cancel some trivial element which causes divide() to fail. The
1318 // elements are normalized first (yes, even though this algorithm doesn't
1319 // need GCDs) since the elements of *this might be unnormalized, which
1320 // makes things more complicated than they need to be.
1321 matrix tmp_n(*this);
1322 matrix tmp_d(m,n); // for denominators, if needed
1323 exmap srl; // symbol replacement list
1324 exvector::const_iterator cit = this->m.begin(), citend = this->m.end();
1325 exvector::iterator tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
1326 while (cit != citend) {
1327 ex nd = cit->normal().to_rational(srl).numer_denom();
1329 *tmp_n_it++ = nd.op(0);
1330 *tmp_d_it++ = nd.op(1);
1334 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
1335 int indx = tmp_n.pivot(r0, r1, true);
1344 // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
1345 for (unsigned c=r1; c<n; ++c)
1346 tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
1348 for (unsigned r2=r0+1; r2<m; ++r2) {
1349 for (unsigned c=r1+1; c<n; ++c) {
1350 dividend_n = (tmp_n.m[r0*n+r1]*tmp_n.m[r2*n+c]*
1351 tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]
1352 -tmp_n.m[r2*n+r1]*tmp_n.m[r0*n+c]*
1353 tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
1354 dividend_d = (tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]*
1355 tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
1356 bool check = divide(dividend_n, divisor_n,
1357 tmp_n.m[r2*n+c], true);
1358 check &= divide(dividend_d, divisor_d,
1359 tmp_d.m[r2*n+c], true);
1360 GINAC_ASSERT(check);
1362 // fill up left hand side with zeros
1363 for (unsigned c=0; c<=r1; ++c)
1364 tmp_n.m[r2*n+c] = _ex0;
1366 if ((r1<n-1)&&(r0<m-1)) {
1367 // compute next iteration's divisor
1368 divisor_n = tmp_n.m[r0*n+r1].expand();
1369 divisor_d = tmp_d.m[r0*n+r1].expand();
1371 // save space by deleting no longer needed elements
1372 for (unsigned c=0; c<n; ++c) {
1373 tmp_n.m[r0*n+c] = _ex0;
1374 tmp_d.m[r0*n+c] = _ex1;
1381 // repopulate *this matrix:
1382 exvector::iterator it = this->m.begin(), itend = this->m.end();
1383 tmp_n_it = tmp_n.m.begin();
1384 tmp_d_it = tmp_d.m.begin();
1386 *it++ = ((*tmp_n_it++)/(*tmp_d_it++)).subs(srl, subs_options::no_pattern);
1392 /** Partial pivoting method for matrix elimination schemes.
1393 * Usual pivoting (symbolic==false) returns the index to the element with the
1394 * largest absolute value in column ro and swaps the current row with the one
1395 * where the element was found. With (symbolic==true) it does the same thing
1396 * with the first non-zero element.
1398 * @param ro is the row from where to begin
1399 * @param co is the column to be inspected
1400 * @param symbolic signal if we want the first non-zero element to be pivoted
1401 * (true) or the one with the largest absolute value (false).
1402 * @return 0 if no interchange occured, -1 if all are zero (usually signaling
1403 * a degeneracy) and positive integer k means that rows ro and k were swapped.
1405 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
1409 // search first non-zero element in column co beginning at row ro
1410 while ((k<row) && (this->m[k*col+co].expand().is_zero()))
1413 // search largest element in column co beginning at row ro
1414 GINAC_ASSERT(is_exactly_a<numeric>(this->m[k*col+co]));
1415 unsigned kmax = k+1;
1416 numeric mmax = abs(ex_to<numeric>(m[kmax*col+co]));
1418 GINAC_ASSERT(is_exactly_a<numeric>(this->m[kmax*col+co]));
1419 numeric tmp = ex_to<numeric>(this->m[kmax*col+co]);
1420 if (abs(tmp) > mmax) {
1426 if (!mmax.is_zero())
1430 // all elements in column co below row ro vanish
1433 // matrix needs no pivoting
1435 // matrix needs pivoting, so swap rows k and ro
1436 ensure_if_modifiable();
1437 for (unsigned c=0; c<col; ++c)
1438 this->m[k*col+c].swap(this->m[ro*col+c]);
1443 ex lst_to_matrix(const lst & l)
1445 lst::const_iterator itr, itc;
1447 // Find number of rows and columns
1448 size_t rows = l.nops(), cols = 0;
1449 for (itr = l.begin(); itr != l.end(); ++itr) {
1450 if (!is_a<lst>(*itr))
1451 throw (std::invalid_argument("lst_to_matrix: argument must be a list of lists"));
1452 if (itr->nops() > cols)
1456 // Allocate and fill matrix
1457 matrix &M = *new matrix(rows, cols);
1458 M.setflag(status_flags::dynallocated);
1461 for (itr = l.begin(), i = 0; itr != l.end(); ++itr, ++i) {
1463 for (itc = ex_to<lst>(*itr).begin(), j = 0; itc != ex_to<lst>(*itr).end(); ++itc, ++j)
1470 ex diag_matrix(const lst & l)
1472 lst::const_iterator it;
1473 size_t dim = l.nops();
1475 // Allocate and fill matrix
1476 matrix &M = *new matrix(dim, dim);
1477 M.setflag(status_flags::dynallocated);
1480 for (it = l.begin(), i = 0; it != l.end(); ++it, ++i)
1486 ex unit_matrix(unsigned r, unsigned c)
1488 matrix &Id = *new matrix(r, c);
1489 Id.setflag(status_flags::dynallocated);
1490 for (unsigned i=0; i<r && i<c; i++)
1496 ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const std::string & tex_base_name)
1498 matrix &M = *new matrix(r, c);
1499 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1501 bool long_format = (r > 10 || c > 10);
1502 bool single_row = (r == 1 || c == 1);
1504 for (unsigned i=0; i<r; i++) {
1505 for (unsigned j=0; j<c; j++) {
1506 std::ostringstream s1, s2;
1508 s2 << tex_base_name << "_{";
1519 s1 << '_' << i << '_' << j;
1520 s2 << i << ';' << j << "}";
1523 s2 << i << j << '}';
1526 M(i, j) = symbol(s1.str(), s2.str());
1533 } // namespace GiNaC