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"
46 GINAC_IMPLEMENT_REGISTERED_CLASS(matrix, basic)
49 // default constructor
52 /** Default ctor. Initializes to 1 x 1-dimensional zero-matrix. */
53 matrix::matrix() : inherited(TINFO_matrix), row(1), col(1)
64 /** Very common ctor. Initializes to r x c-dimensional zero-matrix.
66 * @param r number of rows
67 * @param c number of cols */
68 matrix::matrix(unsigned r, unsigned c)
69 : inherited(TINFO_matrix), row(r), col(c)
76 /** Ctor from representation, for internal use only. */
77 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
78 : inherited(TINFO_matrix), row(r), col(c), m(m2) {}
80 /** Construct matrix from (flat) list of elements. If the list has fewer
81 * elements than the matrix, the remaining matrix elements are set to zero.
82 * If the list has more elements than the matrix, the excessive elements are
84 matrix::matrix(unsigned r, unsigned c, const lst & l)
85 : inherited(TINFO_matrix), row(r), col(c)
90 for (lst::const_iterator it = l.begin(); it != l.end(); ++it, ++i) {
94 break; // matrix smaller than list: throw away excessive elements
103 matrix::matrix(const archive_node &n, lst &sym_lst) : inherited(n, sym_lst)
105 if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
106 throw (std::runtime_error("unknown matrix dimensions in archive"));
107 m.reserve(row * col);
108 for (unsigned int i=0; true; i++) {
110 if (n.find_ex("m", e, sym_lst, i))
117 void matrix::archive(archive_node &n) const
119 inherited::archive(n);
120 n.add_unsigned("row", row);
121 n.add_unsigned("col", col);
122 exvector::const_iterator i = m.begin(), iend = m.end();
129 DEFAULT_UNARCHIVE(matrix)
132 // functions overriding virtual functions from base classes
137 void matrix::print(const print_context & c, unsigned level) const
139 if (is_a<print_tree>(c)) {
141 inherited::print(c, level);
145 if (is_a<print_python_repr>(c))
146 c.s << class_name() << '(';
148 if (is_a<print_latex>(c))
149 c.s << "\\left(\\begin{array}{" << std::string(col,'c') << "}";
153 for (unsigned ro=0; ro<row; ++ro) {
154 if (!is_a<print_latex>(c))
156 for (unsigned co=0; co<col; ++co) {
157 m[ro*col+co].print(c);
159 if (is_a<print_latex>(c))
164 if (!is_a<print_latex>(c))
169 if (is_a<print_latex>(c))
176 if (is_a<print_latex>(c))
177 c.s << "\\end{array}\\right)";
181 if (is_a<print_python_repr>(c))
187 /** nops is defined to be rows x columns. */
188 size_t matrix::nops() const
190 return static_cast<size_t>(row) * static_cast<size_t>(col);
193 /** returns matrix entry at position (i/col, i%col). */
194 ex matrix::op(size_t i) const
196 GINAC_ASSERT(i<nops());
201 /** returns writable matrix entry at position (i/col, i%col). */
202 ex & matrix::let_op(size_t i)
204 GINAC_ASSERT(i<nops());
206 ensure_if_modifiable();
210 /** Evaluate matrix entry by entry. */
211 ex matrix::eval(int level) const
213 // check if we have to do anything at all
214 if ((level==1)&&(flags & status_flags::evaluated))
218 if (level == -max_recursion_level)
219 throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
221 // eval() entry by entry
222 exvector m2(row*col);
224 for (unsigned r=0; r<row; ++r)
225 for (unsigned c=0; c<col; ++c)
226 m2[r*col+c] = m[r*col+c].eval(level);
228 return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
229 status_flags::evaluated);
232 ex matrix::subs(const exmap & mp, unsigned options) const
234 exvector m2(row * col);
235 for (unsigned r=0; r<row; ++r)
236 for (unsigned c=0; c<col; ++c)
237 m2[r*col+c] = m[r*col+c].subs(mp, options);
239 return matrix(row, col, m2).subs_one_level(mp, options);
244 int matrix::compare_same_type(const basic & other) const
246 GINAC_ASSERT(is_exactly_a<matrix>(other));
247 const matrix &o = static_cast<const matrix &>(other);
249 // compare number of rows
251 return row < o.rows() ? -1 : 1;
253 // compare number of columns
255 return col < o.cols() ? -1 : 1;
257 // equal number of rows and columns, compare individual elements
259 for (unsigned r=0; r<row; ++r) {
260 for (unsigned c=0; c<col; ++c) {
261 cmpval = ((*this)(r,c)).compare(o(r,c));
262 if (cmpval!=0) return cmpval;
265 // all elements are equal => matrices are equal;
269 bool matrix::match_same_type(const basic & other) const
271 GINAC_ASSERT(is_exactly_a<matrix>(other));
272 const matrix & o = static_cast<const matrix &>(other);
274 // The number of rows and columns must be the same. This is necessary to
275 // prevent a 2x3 matrix from matching a 3x2 one.
276 return row == o.rows() && col == o.cols();
279 /** Automatic symbolic evaluation of an indexed matrix. */
280 ex matrix::eval_indexed(const basic & i) const
282 GINAC_ASSERT(is_a<indexed>(i));
283 GINAC_ASSERT(is_a<matrix>(i.op(0)));
285 bool all_indices_unsigned = static_cast<const indexed &>(i).all_index_values_are(info_flags::nonnegint);
290 // One index, must be one-dimensional vector
291 if (row != 1 && col != 1)
292 throw (std::runtime_error("matrix::eval_indexed(): vector must have exactly 1 index"));
294 const idx & i1 = ex_to<idx>(i.op(1));
299 if (!i1.get_dim().is_equal(row))
300 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
302 // Index numeric -> return vector element
303 if (all_indices_unsigned) {
304 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
306 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
307 return (*this)(n1, 0);
313 if (!i1.get_dim().is_equal(col))
314 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
316 // Index numeric -> return vector element
317 if (all_indices_unsigned) {
318 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
320 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
321 return (*this)(0, n1);
325 } else if (i.nops() == 3) {
328 const idx & i1 = ex_to<idx>(i.op(1));
329 const idx & i2 = ex_to<idx>(i.op(2));
331 if (!i1.get_dim().is_equal(row))
332 throw (std::runtime_error("matrix::eval_indexed(): dimension of first index must match number of rows"));
333 if (!i2.get_dim().is_equal(col))
334 throw (std::runtime_error("matrix::eval_indexed(): dimension of second index must match number of columns"));
336 // Pair of dummy indices -> compute trace
337 if (is_dummy_pair(i1, i2))
340 // Both indices numeric -> return matrix element
341 if (all_indices_unsigned) {
342 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int(), n2 = ex_to<numeric>(i2.get_value()).to_int();
344 throw (std::runtime_error("matrix::eval_indexed(): value of first index exceeds number of rows"));
346 throw (std::runtime_error("matrix::eval_indexed(): value of second index exceeds number of columns"));
347 return (*this)(n1, n2);
351 throw (std::runtime_error("matrix::eval_indexed(): matrix must have exactly 2 indices"));
356 /** Sum of two indexed matrices. */
357 ex matrix::add_indexed(const ex & self, const ex & other) const
359 GINAC_ASSERT(is_a<indexed>(self));
360 GINAC_ASSERT(is_a<matrix>(self.op(0)));
361 GINAC_ASSERT(is_a<indexed>(other));
362 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
364 // Only add two matrices
365 if (is_a<matrix>(other.op(0))) {
366 GINAC_ASSERT(other.nops() == 2 || other.nops() == 3);
368 const matrix &self_matrix = ex_to<matrix>(self.op(0));
369 const matrix &other_matrix = ex_to<matrix>(other.op(0));
371 if (self.nops() == 2 && other.nops() == 2) { // vector + vector
373 if (self_matrix.row == other_matrix.row)
374 return indexed(self_matrix.add(other_matrix), self.op(1));
375 else if (self_matrix.row == other_matrix.col)
376 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1));
378 } else if (self.nops() == 3 && other.nops() == 3) { // matrix + matrix
380 if (self.op(1).is_equal(other.op(1)) && self.op(2).is_equal(other.op(2)))
381 return indexed(self_matrix.add(other_matrix), self.op(1), self.op(2));
382 else if (self.op(1).is_equal(other.op(2)) && self.op(2).is_equal(other.op(1)))
383 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1), self.op(2));
388 // Don't know what to do, return unevaluated sum
392 /** Product of an indexed matrix with a number. */
393 ex matrix::scalar_mul_indexed(const ex & self, const numeric & other) const
395 GINAC_ASSERT(is_a<indexed>(self));
396 GINAC_ASSERT(is_a<matrix>(self.op(0)));
397 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
399 const matrix &self_matrix = ex_to<matrix>(self.op(0));
401 if (self.nops() == 2)
402 return indexed(self_matrix.mul(other), self.op(1));
403 else // self.nops() == 3
404 return indexed(self_matrix.mul(other), self.op(1), self.op(2));
407 /** Contraction of an indexed matrix with something else. */
408 bool matrix::contract_with(exvector::iterator self, exvector::iterator other, exvector & v) const
410 GINAC_ASSERT(is_a<indexed>(*self));
411 GINAC_ASSERT(is_a<indexed>(*other));
412 GINAC_ASSERT(self->nops() == 2 || self->nops() == 3);
413 GINAC_ASSERT(is_a<matrix>(self->op(0)));
415 // Only contract with other matrices
416 if (!is_a<matrix>(other->op(0)))
419 GINAC_ASSERT(other->nops() == 2 || other->nops() == 3);
421 const matrix &self_matrix = ex_to<matrix>(self->op(0));
422 const matrix &other_matrix = ex_to<matrix>(other->op(0));
424 if (self->nops() == 2) {
426 if (other->nops() == 2) { // vector * vector (scalar product)
428 if (self_matrix.col == 1) {
429 if (other_matrix.col == 1) {
430 // Column vector * column vector, transpose first vector
431 *self = self_matrix.transpose().mul(other_matrix)(0, 0);
433 // Column vector * row vector, swap factors
434 *self = other_matrix.mul(self_matrix)(0, 0);
437 if (other_matrix.col == 1) {
438 // Row vector * column vector, perfect
439 *self = self_matrix.mul(other_matrix)(0, 0);
441 // Row vector * row vector, transpose second vector
442 *self = self_matrix.mul(other_matrix.transpose())(0, 0);
448 } else { // vector * matrix
450 // B_i * A_ij = (B*A)_j (B is row vector)
451 if (is_dummy_pair(self->op(1), other->op(1))) {
452 if (self_matrix.row == 1)
453 *self = indexed(self_matrix.mul(other_matrix), other->op(2));
455 *self = indexed(self_matrix.transpose().mul(other_matrix), other->op(2));
460 // B_j * A_ij = (A*B)_i (B is column vector)
461 if (is_dummy_pair(self->op(1), other->op(2))) {
462 if (self_matrix.col == 1)
463 *self = indexed(other_matrix.mul(self_matrix), other->op(1));
465 *self = indexed(other_matrix.mul(self_matrix.transpose()), other->op(1));
471 } else if (other->nops() == 3) { // matrix * matrix
473 // A_ij * B_jk = (A*B)_ik
474 if (is_dummy_pair(self->op(2), other->op(1))) {
475 *self = indexed(self_matrix.mul(other_matrix), self->op(1), other->op(2));
480 // A_ij * B_kj = (A*Btrans)_ik
481 if (is_dummy_pair(self->op(2), other->op(2))) {
482 *self = indexed(self_matrix.mul(other_matrix.transpose()), self->op(1), other->op(1));
487 // A_ji * B_jk = (Atrans*B)_ik
488 if (is_dummy_pair(self->op(1), other->op(1))) {
489 *self = indexed(self_matrix.transpose().mul(other_matrix), self->op(2), other->op(2));
494 // A_ji * B_kj = (B*A)_ki
495 if (is_dummy_pair(self->op(1), other->op(2))) {
496 *self = indexed(other_matrix.mul(self_matrix), other->op(1), self->op(2));
507 // non-virtual functions in this class
514 * @exception logic_error (incompatible matrices) */
515 matrix matrix::add(const matrix & other) const
517 if (col != other.col || row != other.row)
518 throw std::logic_error("matrix::add(): incompatible matrices");
520 exvector sum(this->m);
521 exvector::iterator i = sum.begin(), end = sum.end();
522 exvector::const_iterator ci = other.m.begin();
526 return matrix(row,col,sum);
530 /** Difference of matrices.
532 * @exception logic_error (incompatible matrices) */
533 matrix matrix::sub(const matrix & other) const
535 if (col != other.col || row != other.row)
536 throw std::logic_error("matrix::sub(): incompatible matrices");
538 exvector dif(this->m);
539 exvector::iterator i = dif.begin(), end = dif.end();
540 exvector::const_iterator ci = other.m.begin();
544 return matrix(row,col,dif);
548 /** Product of matrices.
550 * @exception logic_error (incompatible matrices) */
551 matrix matrix::mul(const matrix & other) const
553 if (this->cols() != other.rows())
554 throw std::logic_error("matrix::mul(): incompatible matrices");
556 exvector prod(this->rows()*other.cols());
558 for (unsigned r1=0; r1<this->rows(); ++r1) {
559 for (unsigned c=0; c<this->cols(); ++c) {
560 if (m[r1*col+c].is_zero())
562 for (unsigned r2=0; r2<other.cols(); ++r2)
563 prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]).expand();
566 return matrix(row, other.col, prod);
570 /** Product of matrix and scalar. */
571 matrix matrix::mul(const numeric & other) const
573 exvector prod(row * col);
575 for (unsigned r=0; r<row; ++r)
576 for (unsigned c=0; c<col; ++c)
577 prod[r*col+c] = m[r*col+c] * other;
579 return matrix(row, col, prod);
583 /** Product of matrix and scalar expression. */
584 matrix matrix::mul_scalar(const ex & other) const
586 if (other.return_type() != return_types::commutative)
587 throw std::runtime_error("matrix::mul_scalar(): non-commutative scalar");
589 exvector prod(row * col);
591 for (unsigned r=0; r<row; ++r)
592 for (unsigned c=0; c<col; ++c)
593 prod[r*col+c] = m[r*col+c] * other;
595 return matrix(row, col, prod);
599 /** Power of a matrix. Currently handles integer exponents only. */
600 matrix matrix::pow(const ex & expn) const
603 throw (std::logic_error("matrix::pow(): matrix not square"));
605 if (is_exactly_a<numeric>(expn)) {
606 // Integer cases are computed by successive multiplication, using the
607 // obvious shortcut of storing temporaries, like A^4 == (A*A)*(A*A).
608 if (expn.info(info_flags::integer)) {
609 numeric b = ex_to<numeric>(expn);
611 if (expn.info(info_flags::negative)) {
618 for (unsigned r=0; r<row; ++r)
622 // This loop computes the representation of b in base 2 from right
623 // to left and multiplies the factors whenever needed. Note
624 // that this is not entirely optimal but close to optimal and
625 // "better" algorithms are much harder to implement. (See Knuth,
626 // TAoCP2, section "Evaluation of Powers" for a good discussion.)
632 b /= _num2; // still integer.
638 throw (std::runtime_error("matrix::pow(): don't know how to handle exponent"));
642 /** operator() to access elements for reading.
644 * @param ro row of element
645 * @param co column of element
646 * @exception range_error (index out of range) */
647 const ex & matrix::operator() (unsigned ro, unsigned co) const
649 if (ro>=row || co>=col)
650 throw (std::range_error("matrix::operator(): index out of range"));
656 /** operator() to access elements for writing.
658 * @param ro row of element
659 * @param co column of element
660 * @exception range_error (index out of range) */
661 ex & matrix::operator() (unsigned ro, unsigned co)
663 if (ro>=row || co>=col)
664 throw (std::range_error("matrix::operator(): index out of range"));
666 ensure_if_modifiable();
671 /** Transposed of an m x n matrix, producing a new n x m matrix object that
672 * represents the transposed. */
673 matrix matrix::transpose() const
675 exvector trans(this->cols()*this->rows());
677 for (unsigned r=0; r<this->cols(); ++r)
678 for (unsigned c=0; c<this->rows(); ++c)
679 trans[r*this->rows()+c] = m[c*this->cols()+r];
681 return matrix(this->cols(),this->rows(),trans);
684 /** Determinant of square matrix. This routine doesn't actually calculate the
685 * determinant, it only implements some heuristics about which algorithm to
686 * run. If all the elements of the matrix are elements of an integral domain
687 * the determinant is also in that integral domain and the result is expanded
688 * only. If one or more elements are from a quotient field the determinant is
689 * usually also in that quotient field and the result is normalized before it
690 * is returned. This implies that the determinant of the symbolic 2x2 matrix
691 * [[a/(a-b),1],[b/(a-b),1]] is returned as unity. (In this respect, it
692 * behaves like MapleV and unlike Mathematica.)
694 * @param algo allows to chose an algorithm
695 * @return the determinant as a new expression
696 * @exception logic_error (matrix not square)
697 * @see determinant_algo */
698 ex matrix::determinant(unsigned algo) const
701 throw (std::logic_error("matrix::determinant(): matrix not square"));
702 GINAC_ASSERT(row*col==m.capacity());
704 // Gather some statistical information about this matrix:
705 bool numeric_flag = true;
706 bool normal_flag = false;
707 unsigned sparse_count = 0; // counts non-zero elements
708 exvector::const_iterator r = m.begin(), rend = m.end();
710 lst srl; // symbol replacement list
711 ex rtest = r->to_rational(srl);
712 if (!rtest.is_zero())
714 if (!rtest.info(info_flags::numeric))
715 numeric_flag = false;
716 if (!rtest.info(info_flags::crational_polynomial) &&
717 rtest.info(info_flags::rational_function))
722 // Here is the heuristics in case this routine has to decide:
723 if (algo == determinant_algo::automatic) {
724 // Minor expansion is generally a good guess:
725 algo = determinant_algo::laplace;
726 // Does anybody know when a matrix is really sparse?
727 // Maybe <~row/2.236 nonzero elements average in a row?
728 if (row>3 && 5*sparse_count<=row*col)
729 algo = determinant_algo::bareiss;
730 // Purely numeric matrix can be handled by Gauss elimination.
731 // This overrides any prior decisions.
733 algo = determinant_algo::gauss;
736 // Trap the trivial case here, since some algorithms don't like it
738 // for consistency with non-trivial determinants...
740 return m[0].normal();
742 return m[0].expand();
745 // Compute the determinant
747 case determinant_algo::gauss: {
750 int sign = tmp.gauss_elimination(true);
751 for (unsigned d=0; d<row; ++d)
752 det *= tmp.m[d*col+d];
754 return (sign*det).normal();
756 return (sign*det).normal().expand();
758 case determinant_algo::bareiss: {
761 sign = tmp.fraction_free_elimination(true);
763 return (sign*tmp.m[row*col-1]).normal();
765 return (sign*tmp.m[row*col-1]).expand();
767 case determinant_algo::divfree: {
770 sign = tmp.division_free_elimination(true);
773 ex det = tmp.m[row*col-1];
774 // factor out accumulated bogus slag
775 for (unsigned d=0; d<row-2; ++d)
776 for (unsigned j=0; j<row-d-2; ++j)
777 det = (det/tmp.m[d*col+d]).normal();
780 case determinant_algo::laplace:
782 // This is the minor expansion scheme. We always develop such
783 // that the smallest minors (i.e, the trivial 1x1 ones) are on the
784 // rightmost column. For this to be efficient, empirical tests
785 // have shown that the emptiest columns (i.e. the ones with most
786 // zeros) should be the ones on the right hand side -- although
787 // this might seem counter-intuitive (and in contradiction to some
788 // literature like the FORM manual). Please go ahead and test it
789 // if you don't believe me! Therefore we presort the columns of
791 typedef std::pair<unsigned,unsigned> uintpair;
792 std::vector<uintpair> c_zeros; // number of zeros in column
793 for (unsigned c=0; c<col; ++c) {
795 for (unsigned r=0; r<row; ++r)
796 if (m[r*col+c].is_zero())
798 c_zeros.push_back(uintpair(acc,c));
800 std::sort(c_zeros.begin(),c_zeros.end());
801 std::vector<unsigned> pre_sort;
802 for (std::vector<uintpair>::const_iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
803 pre_sort.push_back(i->second);
804 std::vector<unsigned> pre_sort_test(pre_sort); // permutation_sign() modifies the vector so we make a copy here
805 int sign = permutation_sign(pre_sort_test.begin(), pre_sort_test.end());
806 exvector result(row*col); // represents sorted matrix
808 for (std::vector<unsigned>::const_iterator i=pre_sort.begin();
811 for (unsigned r=0; r<row; ++r)
812 result[r*col+c] = m[r*col+(*i)];
816 return (sign*matrix(row,col,result).determinant_minor()).normal();
818 return sign*matrix(row,col,result).determinant_minor();
824 /** Trace of a matrix. The result is normalized if it is in some quotient
825 * field and expanded only otherwise. This implies that the trace of the
826 * symbolic 2x2 matrix [[a/(a-b),x],[y,b/(b-a)]] is recognized to be unity.
828 * @return the sum of diagonal elements
829 * @exception logic_error (matrix not square) */
830 ex matrix::trace() const
833 throw (std::logic_error("matrix::trace(): matrix not square"));
836 for (unsigned r=0; r<col; ++r)
839 if (tr.info(info_flags::rational_function) &&
840 !tr.info(info_flags::crational_polynomial))
847 /** Characteristic Polynomial. Following mathematica notation the
848 * characteristic polynomial of a matrix M is defined as the determiant of
849 * (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
850 * as M. Note that some CASs define it with a sign inside the determinant
851 * which gives rise to an overall sign if the dimension is odd. This method
852 * returns the characteristic polynomial collected in powers of lambda as a
855 * @return characteristic polynomial as new expression
856 * @exception logic_error (matrix not square)
857 * @see matrix::determinant() */
858 ex matrix::charpoly(const symbol & lambda) const
861 throw (std::logic_error("matrix::charpoly(): matrix not square"));
863 bool numeric_flag = true;
864 exvector::const_iterator r = m.begin(), rend = m.end();
865 while (r!=rend && numeric_flag==true) {
866 if (!r->info(info_flags::numeric))
867 numeric_flag = false;
871 // The pure numeric case is traditionally rather common. Hence, it is
872 // trapped and we use Leverrier's algorithm which goes as row^3 for
873 // every coefficient. The expensive part is the matrix multiplication.
878 ex poly = power(lambda,row)-c*power(lambda,row-1);
879 for (unsigned i=1; i<row; ++i) {
880 for (unsigned j=0; j<row; ++j)
883 c = B.trace() / ex(i+1);
884 poly -= c*power(lambda,row-i-1);
894 for (unsigned r=0; r<col; ++r)
895 M.m[r*col+r] -= lambda;
897 return M.determinant().collect(lambda);
902 /** Inverse of this matrix.
904 * @return the inverted matrix
905 * @exception logic_error (matrix not square)
906 * @exception runtime_error (singular matrix) */
907 matrix matrix::inverse() const
910 throw (std::logic_error("matrix::inverse(): matrix not square"));
912 // This routine actually doesn't do anything fancy at all. We compute the
913 // inverse of the matrix A by solving the system A * A^{-1} == Id.
915 // First populate the identity matrix supposed to become the right hand side.
916 matrix identity(row,col);
917 for (unsigned i=0; i<row; ++i)
918 identity(i,i) = _ex1;
920 // Populate a dummy matrix of variables, just because of compatibility with
921 // matrix::solve() which wants this (for compatibility with under-determined
922 // systems of equations).
923 matrix vars(row,col);
924 for (unsigned r=0; r<row; ++r)
925 for (unsigned c=0; c<col; ++c)
926 vars(r,c) = symbol();
930 sol = this->solve(vars,identity);
931 } catch (const std::runtime_error & e) {
932 if (e.what()==std::string("matrix::solve(): inconsistent linear system"))
933 throw (std::runtime_error("matrix::inverse(): singular matrix"));
941 /** Solve a linear system consisting of a m x n matrix and a m x p right hand
942 * side by applying an elimination scheme to the augmented matrix.
944 * @param vars n x p matrix, all elements must be symbols
945 * @param rhs m x p matrix
946 * @return n x p solution matrix
947 * @exception logic_error (incompatible matrices)
948 * @exception invalid_argument (1st argument must be matrix of symbols)
949 * @exception runtime_error (inconsistent linear system)
951 matrix matrix::solve(const matrix & vars,
955 const unsigned m = this->rows();
956 const unsigned n = this->cols();
957 const unsigned p = rhs.cols();
960 if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
961 throw (std::logic_error("matrix::solve(): incompatible matrices"));
962 for (unsigned ro=0; ro<n; ++ro)
963 for (unsigned co=0; co<p; ++co)
964 if (!vars(ro,co).info(info_flags::symbol))
965 throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
967 // build the augmented matrix of *this with rhs attached to the right
969 for (unsigned r=0; r<m; ++r) {
970 for (unsigned c=0; c<n; ++c)
971 aug.m[r*(n+p)+c] = this->m[r*n+c];
972 for (unsigned c=0; c<p; ++c)
973 aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
976 // Gather some statistical information about the augmented matrix:
977 bool numeric_flag = true;
978 exvector::const_iterator r = aug.m.begin(), rend = aug.m.end();
979 while (r!=rend && numeric_flag==true) {
980 if (!r->info(info_flags::numeric))
981 numeric_flag = false;
985 // Here is the heuristics in case this routine has to decide:
986 if (algo == solve_algo::automatic) {
987 // Bareiss (fraction-free) elimination is generally a good guess:
988 algo = solve_algo::bareiss;
989 // For m<3, Bareiss elimination is equivalent to division free
990 // elimination but has more logistic overhead
992 algo = solve_algo::divfree;
993 // This overrides any prior decisions.
995 algo = solve_algo::gauss;
998 // Eliminate the augmented matrix:
1000 case solve_algo::gauss:
1001 aug.gauss_elimination();
1003 case solve_algo::divfree:
1004 aug.division_free_elimination();
1006 case solve_algo::bareiss:
1008 aug.fraction_free_elimination();
1011 // assemble the solution matrix:
1013 for (unsigned co=0; co<p; ++co) {
1014 unsigned last_assigned_sol = n+1;
1015 for (int r=m-1; r>=0; --r) {
1016 unsigned fnz = 1; // first non-zero in row
1017 while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
1020 // row consists only of zeros, corresponding rhs must be 0, too
1021 if (!aug.m[r*(n+p)+n+co].is_zero()) {
1022 throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
1025 // assign solutions for vars between fnz+1 and
1026 // last_assigned_sol-1: free parameters
1027 for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
1028 sol(c,co) = vars.m[c*p+co];
1029 ex e = aug.m[r*(n+p)+n+co];
1030 for (unsigned c=fnz; c<n; ++c)
1031 e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
1032 sol(fnz-1,co) = (e/(aug.m[r*(n+p)+(fnz-1)])).normal();
1033 last_assigned_sol = fnz;
1036 // assign solutions for vars between 1 and
1037 // last_assigned_sol-1: free parameters
1038 for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
1039 sol(ro,co) = vars(ro,co);
1048 /** Recursive determinant for small matrices having at least one symbolic
1049 * entry. The basic algorithm, known as Laplace-expansion, is enhanced by
1050 * some bookkeeping to avoid calculation of the same submatrices ("minors")
1051 * more than once. According to W.M.Gentleman and S.C.Johnson this algorithm
1052 * is better than elimination schemes for matrices of sparse multivariate
1053 * polynomials and also for matrices of dense univariate polynomials if the
1054 * matrix' dimesion is larger than 7.
1056 * @return the determinant as a new expression (in expanded form)
1057 * @see matrix::determinant() */
1058 ex matrix::determinant_minor() const
1060 // for small matrices the algorithm does not make any sense:
1061 const unsigned n = this->cols();
1063 return m[0].expand();
1065 return (m[0]*m[3]-m[2]*m[1]).expand();
1067 return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
1068 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
1069 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
1071 // This algorithm can best be understood by looking at a naive
1072 // implementation of Laplace-expansion, like this one:
1074 // matrix minorM(this->rows()-1,this->cols()-1);
1075 // for (unsigned r1=0; r1<this->rows(); ++r1) {
1076 // // shortcut if element(r1,0) vanishes
1077 // if (m[r1*col].is_zero())
1079 // // assemble the minor matrix
1080 // for (unsigned r=0; r<minorM.rows(); ++r) {
1081 // for (unsigned c=0; c<minorM.cols(); ++c) {
1083 // minorM(r,c) = m[r*col+c+1];
1085 // minorM(r,c) = m[(r+1)*col+c+1];
1088 // // recurse down and care for sign:
1090 // det -= m[r1*col] * minorM.determinant_minor();
1092 // det += m[r1*col] * minorM.determinant_minor();
1094 // return det.expand();
1095 // What happens is that while proceeding down many of the minors are
1096 // computed more than once. In particular, there are binomial(n,k)
1097 // kxk minors and each one is computed factorial(n-k) times. Therefore
1098 // it is reasonable to store the results of the minors. We proceed from
1099 // right to left. At each column c we only need to retrieve the minors
1100 // calculated in step c-1. We therefore only have to store at most
1101 // 2*binomial(n,n/2) minors.
1103 // Unique flipper counter for partitioning into minors
1104 std::vector<unsigned> Pkey;
1106 // key for minor determinant (a subpartition of Pkey)
1107 std::vector<unsigned> Mkey;
1109 // we store our subminors in maps, keys being the rows they arise from
1110 typedef std::map<std::vector<unsigned>,class ex> Rmap;
1111 typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
1115 // initialize A with last column:
1116 for (unsigned r=0; r<n; ++r) {
1117 Pkey.erase(Pkey.begin(),Pkey.end());
1119 A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
1121 // proceed from right to left through matrix
1122 for (int c=n-2; c>=0; --c) {
1123 Pkey.erase(Pkey.begin(),Pkey.end()); // don't change capacity
1124 Mkey.erase(Mkey.begin(),Mkey.end());
1125 for (unsigned i=0; i<n-c; ++i)
1127 unsigned fc = 0; // controls logic for our strange flipper counter
1130 for (unsigned r=0; r<n-c; ++r) {
1131 // maybe there is nothing to do?
1132 if (m[Pkey[r]*n+c].is_zero())
1134 // create the sorted key for all possible minors
1135 Mkey.erase(Mkey.begin(),Mkey.end());
1136 for (unsigned i=0; i<n-c; ++i)
1138 Mkey.push_back(Pkey[i]);
1139 // Fetch the minors and compute the new determinant
1141 det -= m[Pkey[r]*n+c]*A[Mkey];
1143 det += m[Pkey[r]*n+c]*A[Mkey];
1145 // prevent build-up of deep nesting of expressions saves time:
1147 // store the new determinant at its place in B:
1149 B.insert(Rmap_value(Pkey,det));
1150 // increment our strange flipper counter
1151 for (fc=n-c; fc>0; --fc) {
1153 if (Pkey[fc-1]<fc+c)
1157 for (unsigned j=fc; j<n-c; ++j)
1158 Pkey[j] = Pkey[j-1]+1;
1160 // next column, so change the role of A and B:
1169 /** Perform the steps of an ordinary Gaussian elimination to bring the m x n
1170 * matrix into an upper echelon form. The algorithm is ok for matrices
1171 * with numeric coefficients but quite unsuited for symbolic matrices.
1173 * @param det may be set to true to save a lot of space if one is only
1174 * interested in the diagonal elements (i.e. for calculating determinants).
1175 * The others are set to zero in this case.
1176 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1177 * number of rows was swapped and 0 if the matrix is singular. */
1178 int matrix::gauss_elimination(const bool det)
1180 ensure_if_modifiable();
1181 const unsigned m = this->rows();
1182 const unsigned n = this->cols();
1183 GINAC_ASSERT(!det || n==m);
1187 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
1188 int indx = pivot(r0, r1, true);
1192 return 0; // leaves *this in a messy state
1197 for (unsigned r2=r0+1; r2<m; ++r2) {
1198 if (!this->m[r2*n+r1].is_zero()) {
1199 // yes, there is something to do in this row
1200 ex piv = this->m[r2*n+r1] / this->m[r0*n+r1];
1201 for (unsigned c=r1+1; c<n; ++c) {
1202 this->m[r2*n+c] -= piv * this->m[r0*n+c];
1203 if (!this->m[r2*n+c].info(info_flags::numeric))
1204 this->m[r2*n+c] = this->m[r2*n+c].normal();
1207 // fill up left hand side with zeros
1208 for (unsigned c=0; c<=r1; ++c)
1209 this->m[r2*n+c] = _ex0;
1212 // save space by deleting no longer needed elements
1213 for (unsigned c=r0+1; c<n; ++c)
1214 this->m[r0*n+c] = _ex0;
1224 /** Perform the steps of division free elimination to bring the m x n matrix
1225 * into an upper echelon form.
1227 * @param det may be set to true to save a lot of space if one is only
1228 * interested in the diagonal elements (i.e. for calculating determinants).
1229 * The others are set to zero in this case.
1230 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1231 * number of rows was swapped and 0 if the matrix is singular. */
1232 int matrix::division_free_elimination(const bool det)
1234 ensure_if_modifiable();
1235 const unsigned m = this->rows();
1236 const unsigned n = this->cols();
1237 GINAC_ASSERT(!det || n==m);
1241 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
1242 int indx = pivot(r0, r1, true);
1246 return 0; // leaves *this in a messy state
1251 for (unsigned r2=r0+1; r2<m; ++r2) {
1252 for (unsigned c=r1+1; c<n; ++c)
1253 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();
1254 // fill up left hand side with zeros
1255 for (unsigned c=0; c<=r1; ++c)
1256 this->m[r2*n+c] = _ex0;
1259 // save space by deleting no longer needed elements
1260 for (unsigned c=r0+1; c<n; ++c)
1261 this->m[r0*n+c] = _ex0;
1271 /** Perform the steps of Bareiss' one-step fraction free elimination to bring
1272 * the matrix into an upper echelon form. Fraction free elimination means
1273 * that divide is used straightforwardly, without computing GCDs first. This
1274 * is possible, since we know the divisor at each step.
1276 * @param det may be set to true to save a lot of space if one is only
1277 * interested in the last element (i.e. for calculating determinants). The
1278 * others are set to zero in this case.
1279 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1280 * number of rows was swapped and 0 if the matrix is singular. */
1281 int matrix::fraction_free_elimination(const bool det)
1284 // (single-step fraction free elimination scheme, already known to Jordan)
1286 // Usual division-free elimination sets m[0](r,c) = m(r,c) and then sets
1287 // m[k+1](r,c) = m[k](k,k) * m[k](r,c) - m[k](r,k) * m[k](k,c).
1289 // Bareiss (fraction-free) elimination in addition divides that element
1290 // by m[k-1](k-1,k-1) for k>1, where it can be shown by means of the
1291 // Sylvester identity that this really divides m[k+1](r,c).
1293 // We also allow rational functions where the original prove still holds.
1294 // However, we must care for numerator and denominator separately and
1295 // "manually" work in the integral domains because of subtle cancellations
1296 // (see below). This blows up the bookkeeping a bit and the formula has
1297 // to be modified to expand like this (N{x} stands for numerator of x,
1298 // D{x} for denominator of x):
1299 // 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)}
1300 // -N{m[k](r,k)}*N{m[k](k,c)}*D{m[k](k,k)}*D{m[k](r,c)}
1301 // 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)}
1302 // where for k>1 we now divide N{m[k+1](r,c)} by
1303 // N{m[k-1](k-1,k-1)}
1304 // and D{m[k+1](r,c)} by
1305 // D{m[k-1](k-1,k-1)}.
1307 ensure_if_modifiable();
1308 const unsigned m = this->rows();
1309 const unsigned n = this->cols();
1310 GINAC_ASSERT(!det || n==m);
1319 // We populate temporary matrices to subsequently operate on. There is
1320 // one holding numerators and another holding denominators of entries.
1321 // This is a must since the evaluator (or even earlier mul's constructor)
1322 // might cancel some trivial element which causes divide() to fail. The
1323 // elements are normalized first (yes, even though this algorithm doesn't
1324 // need GCDs) since the elements of *this might be unnormalized, which
1325 // makes things more complicated than they need to be.
1326 matrix tmp_n(*this);
1327 matrix tmp_d(m,n); // for denominators, if needed
1328 lst srl; // symbol replacement list
1329 exvector::const_iterator cit = this->m.begin(), citend = this->m.end();
1330 exvector::iterator tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
1331 while (cit != citend) {
1332 ex nd = cit->normal().to_rational(srl).numer_denom();
1334 *tmp_n_it++ = nd.op(0);
1335 *tmp_d_it++ = nd.op(1);
1339 for (unsigned r1=0; (r1<n-1)&&(r0<m-1); ++r1) {
1340 int indx = tmp_n.pivot(r0, r1, true);
1349 // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
1350 for (unsigned c=r1; c<n; ++c)
1351 tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
1353 for (unsigned r2=r0+1; r2<m; ++r2) {
1354 for (unsigned c=r1+1; c<n; ++c) {
1355 dividend_n = (tmp_n.m[r0*n+r1]*tmp_n.m[r2*n+c]*
1356 tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]
1357 -tmp_n.m[r2*n+r1]*tmp_n.m[r0*n+c]*
1358 tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
1359 dividend_d = (tmp_d.m[r2*n+r1]*tmp_d.m[r0*n+c]*
1360 tmp_d.m[r0*n+r1]*tmp_d.m[r2*n+c]).expand();
1361 bool check = divide(dividend_n, divisor_n,
1362 tmp_n.m[r2*n+c], true);
1363 check &= divide(dividend_d, divisor_d,
1364 tmp_d.m[r2*n+c], true);
1365 GINAC_ASSERT(check);
1367 // fill up left hand side with zeros
1368 for (unsigned c=0; c<=r1; ++c)
1369 tmp_n.m[r2*n+c] = _ex0;
1371 if ((r1<n-1)&&(r0<m-1)) {
1372 // compute next iteration's divisor
1373 divisor_n = tmp_n.m[r0*n+r1].expand();
1374 divisor_d = tmp_d.m[r0*n+r1].expand();
1376 // save space by deleting no longer needed elements
1377 for (unsigned c=0; c<n; ++c) {
1378 tmp_n.m[r0*n+c] = _ex0;
1379 tmp_d.m[r0*n+c] = _ex1;
1386 // repopulate *this matrix:
1387 exvector::iterator it = this->m.begin(), itend = this->m.end();
1388 tmp_n_it = tmp_n.m.begin();
1389 tmp_d_it = tmp_d.m.begin();
1391 *it++ = ((*tmp_n_it++)/(*tmp_d_it++)).subs(srl);
1397 /** Partial pivoting method for matrix elimination schemes.
1398 * Usual pivoting (symbolic==false) returns the index to the element with the
1399 * largest absolute value in column ro and swaps the current row with the one
1400 * where the element was found. With (symbolic==true) it does the same thing
1401 * with the first non-zero element.
1403 * @param ro is the row from where to begin
1404 * @param co is the column to be inspected
1405 * @param symbolic signal if we want the first non-zero element to be pivoted
1406 * (true) or the one with the largest absolute value (false).
1407 * @return 0 if no interchange occured, -1 if all are zero (usually signaling
1408 * a degeneracy) and positive integer k means that rows ro and k were swapped.
1410 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
1414 // search first non-zero element in column co beginning at row ro
1415 while ((k<row) && (this->m[k*col+co].expand().is_zero()))
1418 // search largest element in column co beginning at row ro
1419 GINAC_ASSERT(is_exactly_a<numeric>(this->m[k*col+co]));
1420 unsigned kmax = k+1;
1421 numeric mmax = abs(ex_to<numeric>(m[kmax*col+co]));
1423 GINAC_ASSERT(is_exactly_a<numeric>(this->m[kmax*col+co]));
1424 numeric tmp = ex_to<numeric>(this->m[kmax*col+co]);
1425 if (abs(tmp) > mmax) {
1431 if (!mmax.is_zero())
1435 // all elements in column co below row ro vanish
1438 // matrix needs no pivoting
1440 // matrix needs pivoting, so swap rows k and ro
1441 ensure_if_modifiable();
1442 for (unsigned c=0; c<col; ++c)
1443 this->m[k*col+c].swap(this->m[ro*col+c]);
1448 ex lst_to_matrix(const lst & l)
1450 lst::const_iterator itr, itc;
1452 // Find number of rows and columns
1453 size_t rows = l.nops(), cols = 0;
1454 for (itr = l.begin(); itr != l.end(); ++itr) {
1455 if (!is_a<lst>(*itr))
1456 throw (std::invalid_argument("lst_to_matrix: argument must be a list of lists"));
1457 if (itr->nops() > cols)
1461 // Allocate and fill matrix
1462 matrix &M = *new matrix(rows, cols);
1463 M.setflag(status_flags::dynallocated);
1466 for (itr = l.begin(), i = 0; itr != l.end(); ++itr, ++i) {
1468 for (itc = ex_to<lst>(*itr).begin(), j = 0; itc != ex_to<lst>(*itr).end(); ++itc, ++j)
1475 ex diag_matrix(const lst & l)
1477 lst::const_iterator it;
1478 size_t dim = l.nops();
1480 // Allocate and fill matrix
1481 matrix &M = *new matrix(dim, dim);
1482 M.setflag(status_flags::dynallocated);
1485 for (it = l.begin(), i = 0; it != l.end(); ++it, ++i)
1491 ex unit_matrix(unsigned r, unsigned c)
1493 matrix &Id = *new matrix(r, c);
1494 Id.setflag(status_flags::dynallocated);
1495 for (unsigned i=0; i<r && i<c; i++)
1501 ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const std::string & tex_base_name)
1503 matrix &M = *new matrix(r, c);
1504 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1506 bool long_format = (r > 10 || c > 10);
1507 bool single_row = (r == 1 || c == 1);
1509 for (unsigned i=0; i<r; i++) {
1510 for (unsigned j=0; j<c; j++) {
1511 std::ostringstream s1, s2;
1513 s2 << tex_base_name << "_{";
1524 s1 << '_' << i << '_' << j;
1525 s2 << i << ';' << j << "}";
1528 s2 << i << j << '}';
1531 M(i, j) = symbol(s1.str(), s2.str());
1538 } // namespace GiNaC