3 * Implementation of symbolic matrices */
6 * GiNaC Copyright (C) 1999-2006 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., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 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>(&matrix::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(&matrix::tinfo_static), 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(&matrix::tinfo_static), 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(&matrix::tinfo_static), 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(&matrix::tinfo_static), 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);
237 /** Complex conjugate every matrix entry. */
238 ex matrix::conjugate() const
241 for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i) {
242 ex x = i->conjugate();
247 if (are_ex_trivially_equal(x, *i)) {
251 ev->reserve(m.size());
252 for (exvector::const_iterator j=m.begin(); j!=i; ++j) {
258 ex result = matrix(row, col, *ev);
265 ex matrix::real_part() const
269 for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
270 v.push_back(i->real_part());
271 return matrix(row, col, v);
274 ex matrix::imag_part() const
278 for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
279 v.push_back(i->imag_part());
280 return matrix(row, col, v);
285 int matrix::compare_same_type(const basic & other) const
287 GINAC_ASSERT(is_exactly_a<matrix>(other));
288 const matrix &o = static_cast<const matrix &>(other);
290 // compare number of rows
292 return row < o.rows() ? -1 : 1;
294 // compare number of columns
296 return col < o.cols() ? -1 : 1;
298 // equal number of rows and columns, compare individual elements
300 for (unsigned r=0; r<row; ++r) {
301 for (unsigned c=0; c<col; ++c) {
302 cmpval = ((*this)(r,c)).compare(o(r,c));
303 if (cmpval!=0) return cmpval;
306 // all elements are equal => matrices are equal;
310 bool matrix::match_same_type(const basic & other) const
312 GINAC_ASSERT(is_exactly_a<matrix>(other));
313 const matrix & o = static_cast<const matrix &>(other);
315 // The number of rows and columns must be the same. This is necessary to
316 // prevent a 2x3 matrix from matching a 3x2 one.
317 return row == o.rows() && col == o.cols();
320 /** Automatic symbolic evaluation of an indexed matrix. */
321 ex matrix::eval_indexed(const basic & i) const
323 GINAC_ASSERT(is_a<indexed>(i));
324 GINAC_ASSERT(is_a<matrix>(i.op(0)));
326 bool all_indices_unsigned = static_cast<const indexed &>(i).all_index_values_are(info_flags::nonnegint);
331 // One index, must be one-dimensional vector
332 if (row != 1 && col != 1)
333 throw (std::runtime_error("matrix::eval_indexed(): vector must have exactly 1 index"));
335 const idx & i1 = ex_to<idx>(i.op(1));
340 if (!i1.get_dim().is_equal(row))
341 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
343 // Index numeric -> return vector element
344 if (all_indices_unsigned) {
345 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
347 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
348 return (*this)(n1, 0);
354 if (!i1.get_dim().is_equal(col))
355 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
357 // Index numeric -> return vector element
358 if (all_indices_unsigned) {
359 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
361 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
362 return (*this)(0, n1);
366 } else if (i.nops() == 3) {
369 const idx & i1 = ex_to<idx>(i.op(1));
370 const idx & i2 = ex_to<idx>(i.op(2));
372 if (!i1.get_dim().is_equal(row))
373 throw (std::runtime_error("matrix::eval_indexed(): dimension of first index must match number of rows"));
374 if (!i2.get_dim().is_equal(col))
375 throw (std::runtime_error("matrix::eval_indexed(): dimension of second index must match number of columns"));
377 // Pair of dummy indices -> compute trace
378 if (is_dummy_pair(i1, i2))
381 // Both indices numeric -> return matrix element
382 if (all_indices_unsigned) {
383 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int(), n2 = ex_to<numeric>(i2.get_value()).to_int();
385 throw (std::runtime_error("matrix::eval_indexed(): value of first index exceeds number of rows"));
387 throw (std::runtime_error("matrix::eval_indexed(): value of second index exceeds number of columns"));
388 return (*this)(n1, n2);
392 throw (std::runtime_error("matrix::eval_indexed(): matrix must have exactly 2 indices"));
397 /** Sum of two indexed matrices. */
398 ex matrix::add_indexed(const ex & self, const ex & other) const
400 GINAC_ASSERT(is_a<indexed>(self));
401 GINAC_ASSERT(is_a<matrix>(self.op(0)));
402 GINAC_ASSERT(is_a<indexed>(other));
403 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
405 // Only add two matrices
406 if (is_a<matrix>(other.op(0))) {
407 GINAC_ASSERT(other.nops() == 2 || other.nops() == 3);
409 const matrix &self_matrix = ex_to<matrix>(self.op(0));
410 const matrix &other_matrix = ex_to<matrix>(other.op(0));
412 if (self.nops() == 2 && other.nops() == 2) { // vector + vector
414 if (self_matrix.row == other_matrix.row)
415 return indexed(self_matrix.add(other_matrix), self.op(1));
416 else if (self_matrix.row == other_matrix.col)
417 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1));
419 } else if (self.nops() == 3 && other.nops() == 3) { // matrix + matrix
421 if (self.op(1).is_equal(other.op(1)) && self.op(2).is_equal(other.op(2)))
422 return indexed(self_matrix.add(other_matrix), self.op(1), self.op(2));
423 else if (self.op(1).is_equal(other.op(2)) && self.op(2).is_equal(other.op(1)))
424 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1), self.op(2));
429 // Don't know what to do, return unevaluated sum
433 /** Product of an indexed matrix with a number. */
434 ex matrix::scalar_mul_indexed(const ex & self, const numeric & other) const
436 GINAC_ASSERT(is_a<indexed>(self));
437 GINAC_ASSERT(is_a<matrix>(self.op(0)));
438 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
440 const matrix &self_matrix = ex_to<matrix>(self.op(0));
442 if (self.nops() == 2)
443 return indexed(self_matrix.mul(other), self.op(1));
444 else // self.nops() == 3
445 return indexed(self_matrix.mul(other), self.op(1), self.op(2));
448 /** Contraction of an indexed matrix with something else. */
449 bool matrix::contract_with(exvector::iterator self, exvector::iterator other, exvector & v) const
451 GINAC_ASSERT(is_a<indexed>(*self));
452 GINAC_ASSERT(is_a<indexed>(*other));
453 GINAC_ASSERT(self->nops() == 2 || self->nops() == 3);
454 GINAC_ASSERT(is_a<matrix>(self->op(0)));
456 // Only contract with other matrices
457 if (!is_a<matrix>(other->op(0)))
460 GINAC_ASSERT(other->nops() == 2 || other->nops() == 3);
462 const matrix &self_matrix = ex_to<matrix>(self->op(0));
463 const matrix &other_matrix = ex_to<matrix>(other->op(0));
465 if (self->nops() == 2) {
467 if (other->nops() == 2) { // vector * vector (scalar product)
469 if (self_matrix.col == 1) {
470 if (other_matrix.col == 1) {
471 // Column vector * column vector, transpose first vector
472 *self = self_matrix.transpose().mul(other_matrix)(0, 0);
474 // Column vector * row vector, swap factors
475 *self = other_matrix.mul(self_matrix)(0, 0);
478 if (other_matrix.col == 1) {
479 // Row vector * column vector, perfect
480 *self = self_matrix.mul(other_matrix)(0, 0);
482 // Row vector * row vector, transpose second vector
483 *self = self_matrix.mul(other_matrix.transpose())(0, 0);
489 } else { // vector * matrix
491 // B_i * A_ij = (B*A)_j (B is row vector)
492 if (is_dummy_pair(self->op(1), other->op(1))) {
493 if (self_matrix.row == 1)
494 *self = indexed(self_matrix.mul(other_matrix), other->op(2));
496 *self = indexed(self_matrix.transpose().mul(other_matrix), other->op(2));
501 // B_j * A_ij = (A*B)_i (B is column vector)
502 if (is_dummy_pair(self->op(1), other->op(2))) {
503 if (self_matrix.col == 1)
504 *self = indexed(other_matrix.mul(self_matrix), other->op(1));
506 *self = indexed(other_matrix.mul(self_matrix.transpose()), other->op(1));
512 } else if (other->nops() == 3) { // matrix * matrix
514 // A_ij * B_jk = (A*B)_ik
515 if (is_dummy_pair(self->op(2), other->op(1))) {
516 *self = indexed(self_matrix.mul(other_matrix), self->op(1), other->op(2));
521 // A_ij * B_kj = (A*Btrans)_ik
522 if (is_dummy_pair(self->op(2), other->op(2))) {
523 *self = indexed(self_matrix.mul(other_matrix.transpose()), self->op(1), other->op(1));
528 // A_ji * B_jk = (Atrans*B)_ik
529 if (is_dummy_pair(self->op(1), other->op(1))) {
530 *self = indexed(self_matrix.transpose().mul(other_matrix), self->op(2), other->op(2));
535 // A_ji * B_kj = (B*A)_ki
536 if (is_dummy_pair(self->op(1), other->op(2))) {
537 *self = indexed(other_matrix.mul(self_matrix), other->op(1), self->op(2));
548 // non-virtual functions in this class
555 * @exception logic_error (incompatible matrices) */
556 matrix matrix::add(const matrix & other) const
558 if (col != other.col || row != other.row)
559 throw std::logic_error("matrix::add(): incompatible matrices");
561 exvector sum(this->m);
562 exvector::iterator i = sum.begin(), end = sum.end();
563 exvector::const_iterator ci = other.m.begin();
567 return matrix(row,col,sum);
571 /** Difference of matrices.
573 * @exception logic_error (incompatible matrices) */
574 matrix matrix::sub(const matrix & other) const
576 if (col != other.col || row != other.row)
577 throw std::logic_error("matrix::sub(): incompatible matrices");
579 exvector dif(this->m);
580 exvector::iterator i = dif.begin(), end = dif.end();
581 exvector::const_iterator ci = other.m.begin();
585 return matrix(row,col,dif);
589 /** Product of matrices.
591 * @exception logic_error (incompatible matrices) */
592 matrix matrix::mul(const matrix & other) const
594 if (this->cols() != other.rows())
595 throw std::logic_error("matrix::mul(): incompatible matrices");
597 exvector prod(this->rows()*other.cols());
599 for (unsigned r1=0; r1<this->rows(); ++r1) {
600 for (unsigned c=0; c<this->cols(); ++c) {
601 // Quick test: can we shortcut?
602 if (m[r1*col+c].is_zero())
604 for (unsigned r2=0; r2<other.cols(); ++r2)
605 prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]);
608 return matrix(row, other.col, prod);
612 /** Product of matrix and scalar. */
613 matrix matrix::mul(const numeric & other) const
615 exvector prod(row * col);
617 for (unsigned r=0; r<row; ++r)
618 for (unsigned c=0; c<col; ++c)
619 prod[r*col+c] = m[r*col+c] * other;
621 return matrix(row, col, prod);
625 /** Product of matrix and scalar expression. */
626 matrix matrix::mul_scalar(const ex & other) const
628 if (other.return_type() != return_types::commutative)
629 throw std::runtime_error("matrix::mul_scalar(): non-commutative scalar");
631 exvector prod(row * col);
633 for (unsigned r=0; r<row; ++r)
634 for (unsigned c=0; c<col; ++c)
635 prod[r*col+c] = m[r*col+c] * other;
637 return matrix(row, col, prod);
641 /** Power of a matrix. Currently handles integer exponents only. */
642 matrix matrix::pow(const ex & expn) const
645 throw (std::logic_error("matrix::pow(): matrix not square"));
647 if (is_exactly_a<numeric>(expn)) {
648 // Integer cases are computed by successive multiplication, using the
649 // obvious shortcut of storing temporaries, like A^4 == (A*A)*(A*A).
650 if (expn.info(info_flags::integer)) {
651 numeric b = ex_to<numeric>(expn);
653 if (expn.info(info_flags::negative)) {
660 for (unsigned r=0; r<row; ++r)
664 // This loop computes the representation of b in base 2 from right
665 // to left and multiplies the factors whenever needed. Note
666 // that this is not entirely optimal but close to optimal and
667 // "better" algorithms are much harder to implement. (See Knuth,
668 // TAoCP2, section "Evaluation of Powers" for a good discussion.)
669 while (b!=*_num1_p) {
674 b /= *_num2_p; // still integer.
680 throw (std::runtime_error("matrix::pow(): don't know how to handle exponent"));
684 /** operator() to access elements for reading.
686 * @param ro row of element
687 * @param co column of element
688 * @exception range_error (index out of range) */
689 const ex & matrix::operator() (unsigned ro, unsigned co) const
691 if (ro>=row || co>=col)
692 throw (std::range_error("matrix::operator(): index out of range"));
698 /** operator() to access elements for writing.
700 * @param ro row of element
701 * @param co column of element
702 * @exception range_error (index out of range) */
703 ex & matrix::operator() (unsigned ro, unsigned co)
705 if (ro>=row || co>=col)
706 throw (std::range_error("matrix::operator(): index out of range"));
708 ensure_if_modifiable();
713 /** Transposed of an m x n matrix, producing a new n x m matrix object that
714 * represents the transposed. */
715 matrix matrix::transpose() const
717 exvector trans(this->cols()*this->rows());
719 for (unsigned r=0; r<this->cols(); ++r)
720 for (unsigned c=0; c<this->rows(); ++c)
721 trans[r*this->rows()+c] = m[c*this->cols()+r];
723 return matrix(this->cols(),this->rows(),trans);
726 /** Determinant of square matrix. This routine doesn't actually calculate the
727 * determinant, it only implements some heuristics about which algorithm to
728 * run. If all the elements of the matrix are elements of an integral domain
729 * the determinant is also in that integral domain and the result is expanded
730 * only. If one or more elements are from a quotient field the determinant is
731 * usually also in that quotient field and the result is normalized before it
732 * is returned. This implies that the determinant of the symbolic 2x2 matrix
733 * [[a/(a-b),1],[b/(a-b),1]] is returned as unity. (In this respect, it
734 * behaves like MapleV and unlike Mathematica.)
736 * @param algo allows to chose an algorithm
737 * @return the determinant as a new expression
738 * @exception logic_error (matrix not square)
739 * @see determinant_algo */
740 ex matrix::determinant(unsigned algo) const
743 throw (std::logic_error("matrix::determinant(): matrix not square"));
744 GINAC_ASSERT(row*col==m.capacity());
746 // Gather some statistical information about this matrix:
747 bool numeric_flag = true;
748 bool normal_flag = false;
749 unsigned sparse_count = 0; // counts non-zero elements
750 exvector::const_iterator r = m.begin(), rend = m.end();
752 if (!r->info(info_flags::numeric))
753 numeric_flag = false;
754 exmap srl; // symbol replacement list
755 ex rtest = r->to_rational(srl);
756 if (!rtest.is_zero())
758 if (!rtest.info(info_flags::crational_polynomial) &&
759 rtest.info(info_flags::rational_function))
764 // Here is the heuristics in case this routine has to decide:
765 if (algo == determinant_algo::automatic) {
766 // Minor expansion is generally a good guess:
767 algo = determinant_algo::laplace;
768 // Does anybody know when a matrix is really sparse?
769 // Maybe <~row/2.236 nonzero elements average in a row?
770 if (row>3 && 5*sparse_count<=row*col)
771 algo = determinant_algo::bareiss;
772 // Purely numeric matrix can be handled by Gauss elimination.
773 // This overrides any prior decisions.
775 algo = determinant_algo::gauss;
778 // Trap the trivial case here, since some algorithms don't like it
780 // for consistency with non-trivial determinants...
782 return m[0].normal();
784 return m[0].expand();
787 // Compute the determinant
789 case determinant_algo::gauss: {
792 int sign = tmp.gauss_elimination(true);
793 for (unsigned d=0; d<row; ++d)
794 det *= tmp.m[d*col+d];
796 return (sign*det).normal();
798 return (sign*det).normal().expand();
800 case determinant_algo::bareiss: {
803 sign = tmp.fraction_free_elimination(true);
805 return (sign*tmp.m[row*col-1]).normal();
807 return (sign*tmp.m[row*col-1]).expand();
809 case determinant_algo::divfree: {
812 sign = tmp.division_free_elimination(true);
815 ex det = tmp.m[row*col-1];
816 // factor out accumulated bogus slag
817 for (unsigned d=0; d<row-2; ++d)
818 for (unsigned j=0; j<row-d-2; ++j)
819 det = (det/tmp.m[d*col+d]).normal();
822 case determinant_algo::laplace:
824 // This is the minor expansion scheme. We always develop such
825 // that the smallest minors (i.e, the trivial 1x1 ones) are on the
826 // rightmost column. For this to be efficient, empirical tests
827 // have shown that the emptiest columns (i.e. the ones with most
828 // zeros) should be the ones on the right hand side -- although
829 // this might seem counter-intuitive (and in contradiction to some
830 // literature like the FORM manual). Please go ahead and test it
831 // if you don't believe me! Therefore we presort the columns of
833 typedef std::pair<unsigned,unsigned> uintpair;
834 std::vector<uintpair> c_zeros; // number of zeros in column
835 for (unsigned c=0; c<col; ++c) {
837 for (unsigned r=0; r<row; ++r)
838 if (m[r*col+c].is_zero())
840 c_zeros.push_back(uintpair(acc,c));
842 std::sort(c_zeros.begin(),c_zeros.end());
843 std::vector<unsigned> pre_sort;
844 for (std::vector<uintpair>::const_iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
845 pre_sort.push_back(i->second);
846 std::vector<unsigned> pre_sort_test(pre_sort); // permutation_sign() modifies the vector so we make a copy here
847 int sign = permutation_sign(pre_sort_test.begin(), pre_sort_test.end());
848 exvector result(row*col); // represents sorted matrix
850 for (std::vector<unsigned>::const_iterator i=pre_sort.begin();
853 for (unsigned r=0; r<row; ++r)
854 result[r*col+c] = m[r*col+(*i)];
858 return (sign*matrix(row,col,result).determinant_minor()).normal();
860 return sign*matrix(row,col,result).determinant_minor();
866 /** Trace of a matrix. The result is normalized if it is in some quotient
867 * field and expanded only otherwise. This implies that the trace of the
868 * symbolic 2x2 matrix [[a/(a-b),x],[y,b/(b-a)]] is recognized to be unity.
870 * @return the sum of diagonal elements
871 * @exception logic_error (matrix not square) */
872 ex matrix::trace() const
875 throw (std::logic_error("matrix::trace(): matrix not square"));
878 for (unsigned r=0; r<col; ++r)
881 if (tr.info(info_flags::rational_function) &&
882 !tr.info(info_flags::crational_polynomial))
889 /** Characteristic Polynomial. Following mathematica notation the
890 * characteristic polynomial of a matrix M is defined as the determiant of
891 * (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
892 * as M. Note that some CASs define it with a sign inside the determinant
893 * which gives rise to an overall sign if the dimension is odd. This method
894 * returns the characteristic polynomial collected in powers of lambda as a
897 * @return characteristic polynomial as new expression
898 * @exception logic_error (matrix not square)
899 * @see matrix::determinant() */
900 ex matrix::charpoly(const ex & lambda) const
903 throw (std::logic_error("matrix::charpoly(): matrix not square"));
905 bool numeric_flag = true;
906 exvector::const_iterator r = m.begin(), rend = m.end();
907 while (r!=rend && numeric_flag==true) {
908 if (!r->info(info_flags::numeric))
909 numeric_flag = false;
913 // The pure numeric case is traditionally rather common. Hence, it is
914 // trapped and we use Leverrier's algorithm which goes as row^3 for
915 // every coefficient. The expensive part is the matrix multiplication.
920 ex poly = power(lambda, row) - c*power(lambda, row-1);
921 for (unsigned i=1; i<row; ++i) {
922 for (unsigned j=0; j<row; ++j)
925 c = B.trace() / ex(i+1);
926 poly -= c*power(lambda, row-i-1);
936 for (unsigned r=0; r<col; ++r)
937 M.m[r*col+r] -= lambda;
939 return M.determinant().collect(lambda);
944 /** Inverse of this matrix.
946 * @return the inverted matrix
947 * @exception logic_error (matrix not square)
948 * @exception runtime_error (singular matrix) */
949 matrix matrix::inverse() const
952 throw (std::logic_error("matrix::inverse(): matrix not square"));
954 // This routine actually doesn't do anything fancy at all. We compute the
955 // inverse of the matrix A by solving the system A * A^{-1} == Id.
957 // First populate the identity matrix supposed to become the right hand side.
958 matrix identity(row,col);
959 for (unsigned i=0; i<row; ++i)
960 identity(i,i) = _ex1;
962 // Populate a dummy matrix of variables, just because of compatibility with
963 // matrix::solve() which wants this (for compatibility with under-determined
964 // systems of equations).
965 matrix vars(row,col);
966 for (unsigned r=0; r<row; ++r)
967 for (unsigned c=0; c<col; ++c)
968 vars(r,c) = symbol();
972 sol = this->solve(vars,identity);
973 } catch (const std::runtime_error & e) {
974 if (e.what()==std::string("matrix::solve(): inconsistent linear system"))
975 throw (std::runtime_error("matrix::inverse(): singular matrix"));
983 /** Solve a linear system consisting of a m x n matrix and a m x p right hand
984 * side by applying an elimination scheme to the augmented matrix.
986 * @param vars n x p matrix, all elements must be symbols
987 * @param rhs m x p matrix
988 * @param algo selects the solving algorithm
989 * @return n x p solution matrix
990 * @exception logic_error (incompatible matrices)
991 * @exception invalid_argument (1st argument must be matrix of symbols)
992 * @exception runtime_error (inconsistent linear system)
994 matrix matrix::solve(const matrix & vars,
998 const unsigned m = this->rows();
999 const unsigned n = this->cols();
1000 const unsigned p = rhs.cols();
1003 if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
1004 throw (std::logic_error("matrix::solve(): incompatible matrices"));
1005 for (unsigned ro=0; ro<n; ++ro)
1006 for (unsigned co=0; co<p; ++co)
1007 if (!vars(ro,co).info(info_flags::symbol))
1008 throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
1010 // build the augmented matrix of *this with rhs attached to the right
1012 for (unsigned r=0; r<m; ++r) {
1013 for (unsigned c=0; c<n; ++c)
1014 aug.m[r*(n+p)+c] = this->m[r*n+c];
1015 for (unsigned c=0; c<p; ++c)
1016 aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
1019 // Gather some statistical information about the augmented matrix:
1020 bool numeric_flag = true;
1021 exvector::const_iterator r = aug.m.begin(), rend = aug.m.end();
1022 while (r!=rend && numeric_flag==true) {
1023 if (!r->info(info_flags::numeric))
1024 numeric_flag = false;
1028 // Here is the heuristics in case this routine has to decide:
1029 if (algo == solve_algo::automatic) {
1030 // Bareiss (fraction-free) elimination is generally a good guess:
1031 algo = solve_algo::bareiss;
1032 // For m<3, Bareiss elimination is equivalent to division free
1033 // elimination but has more logistic overhead
1035 algo = solve_algo::divfree;
1036 // This overrides any prior decisions.
1038 algo = solve_algo::gauss;
1041 // Eliminate the augmented matrix:
1043 case solve_algo::gauss:
1044 aug.gauss_elimination();
1046 case solve_algo::divfree:
1047 aug.division_free_elimination();
1049 case solve_algo::bareiss:
1051 aug.fraction_free_elimination();
1054 // assemble the solution matrix:
1056 for (unsigned co=0; co<p; ++co) {
1057 unsigned last_assigned_sol = n+1;
1058 for (int r=m-1; r>=0; --r) {
1059 unsigned fnz = 1; // first non-zero in row
1060 while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
1063 // row consists only of zeros, corresponding rhs must be 0, too
1064 if (!aug.m[r*(n+p)+n+co].is_zero()) {
1065 throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
1068 // assign solutions for vars between fnz+1 and
1069 // last_assigned_sol-1: free parameters
1070 for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
1071 sol(c,co) = vars.m[c*p+co];
1072 ex e = aug.m[r*(n+p)+n+co];
1073 for (unsigned c=fnz; c<n; ++c)
1074 e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
1075 sol(fnz-1,co) = (e/(aug.m[r*(n+p)+(fnz-1)])).normal();
1076 last_assigned_sol = fnz;
1079 // assign solutions for vars between 1 and
1080 // last_assigned_sol-1: free parameters
1081 for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
1082 sol(ro,co) = vars(ro,co);
1089 /** Compute the rank of this matrix. */
1090 unsigned matrix::rank() const
1093 // Transform this matrix into upper echelon form and then count the
1094 // number of non-zero rows.
1096 GINAC_ASSERT(row*col==m.capacity());
1098 // Actually, any elimination scheme will do since we are only
1099 // interested in the echelon matrix' zeros.
1100 matrix to_eliminate = *this;
1101 to_eliminate.fraction_free_elimination();
1103 unsigned r = row*col; // index of last non-zero element
1105 if (!to_eliminate.m[r].is_zero())
1114 /** Recursive determinant for small matrices having at least one symbolic
1115 * entry. The basic algorithm, known as Laplace-expansion, is enhanced by
1116 * some bookkeeping to avoid calculation of the same submatrices ("minors")
1117 * more than once. According to W.M.Gentleman and S.C.Johnson this algorithm
1118 * is better than elimination schemes for matrices of sparse multivariate
1119 * polynomials and also for matrices of dense univariate polynomials if the
1120 * matrix' dimesion is larger than 7.
1122 * @return the determinant as a new expression (in expanded form)
1123 * @see matrix::determinant() */
1124 ex matrix::determinant_minor() const
1126 // for small matrices the algorithm does not make any sense:
1127 const unsigned n = this->cols();
1129 return m[0].expand();
1131 return (m[0]*m[3]-m[2]*m[1]).expand();
1133 return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
1134 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
1135 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
1137 // This algorithm can best be understood by looking at a naive
1138 // implementation of Laplace-expansion, like this one:
1140 // matrix minorM(this->rows()-1,this->cols()-1);
1141 // for (unsigned r1=0; r1<this->rows(); ++r1) {
1142 // // shortcut if element(r1,0) vanishes
1143 // if (m[r1*col].is_zero())
1145 // // assemble the minor matrix
1146 // for (unsigned r=0; r<minorM.rows(); ++r) {
1147 // for (unsigned c=0; c<minorM.cols(); ++c) {
1149 // minorM(r,c) = m[r*col+c+1];
1151 // minorM(r,c) = m[(r+1)*col+c+1];
1154 // // recurse down and care for sign:
1156 // det -= m[r1*col] * minorM.determinant_minor();
1158 // det += m[r1*col] * minorM.determinant_minor();
1160 // return det.expand();
1161 // What happens is that while proceeding down many of the minors are
1162 // computed more than once. In particular, there are binomial(n,k)
1163 // kxk minors and each one is computed factorial(n-k) times. Therefore
1164 // it is reasonable to store the results of the minors. We proceed from
1165 // right to left. At each column c we only need to retrieve the minors
1166 // calculated in step c-1. We therefore only have to store at most
1167 // 2*binomial(n,n/2) minors.
1169 // Unique flipper counter for partitioning into minors
1170 std::vector<unsigned> Pkey;
1172 // key for minor determinant (a subpartition of Pkey)
1173 std::vector<unsigned> Mkey;
1175 // we store our subminors in maps, keys being the rows they arise from
1176 typedef std::map<std::vector<unsigned>,class ex> Rmap;
1177 typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
1181 // initialize A with last column:
1182 for (unsigned r=0; r<n; ++r) {
1183 Pkey.erase(Pkey.begin(),Pkey.end());
1185 A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
1187 // proceed from right to left through matrix
1188 for (int c=n-2; c>=0; --c) {
1189 Pkey.erase(Pkey.begin(),Pkey.end()); // don't change capacity
1190 Mkey.erase(Mkey.begin(),Mkey.end());
1191 for (unsigned i=0; i<n-c; ++i)
1193 unsigned fc = 0; // controls logic for our strange flipper counter
1196 for (unsigned r=0; r<n-c; ++r) {
1197 // maybe there is nothing to do?
1198 if (m[Pkey[r]*n+c].is_zero())
1200 // create the sorted key for all possible minors
1201 Mkey.erase(Mkey.begin(),Mkey.end());
1202 for (unsigned i=0; i<n-c; ++i)
1204 Mkey.push_back(Pkey[i]);
1205 // Fetch the minors and compute the new determinant
1207 det -= m[Pkey[r]*n+c]*A[Mkey];
1209 det += m[Pkey[r]*n+c]*A[Mkey];
1211 // prevent build-up of deep nesting of expressions saves time:
1213 // store the new determinant at its place in B:
1215 B.insert(Rmap_value(Pkey,det));
1216 // increment our strange flipper counter
1217 for (fc=n-c; fc>0; --fc) {
1219 if (Pkey[fc-1]<fc+c)
1223 for (unsigned j=fc; j<n-c; ++j)
1224 Pkey[j] = Pkey[j-1]+1;
1226 // next column, so change the role of A and B:
1235 /** Perform the steps of an ordinary Gaussian elimination to bring the m x n
1236 * matrix into an upper echelon form. The algorithm is ok for matrices
1237 * with numeric coefficients but quite unsuited for symbolic matrices.
1239 * @param det may be set to true to save a lot of space if one is only
1240 * interested in the diagonal elements (i.e. for calculating determinants).
1241 * The others are set to zero in this case.
1242 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1243 * number of rows was swapped and 0 if the matrix is singular. */
1244 int matrix::gauss_elimination(const bool det)
1246 ensure_if_modifiable();
1247 const unsigned m = this->rows();
1248 const unsigned n = this->cols();
1249 GINAC_ASSERT(!det || n==m);
1253 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1254 int indx = pivot(r0, c0, true);
1258 return 0; // leaves *this in a messy state
1263 for (unsigned r2=r0+1; r2<m; ++r2) {
1264 if (!this->m[r2*n+c0].is_zero()) {
1265 // yes, there is something to do in this row
1266 ex piv = this->m[r2*n+c0] / this->m[r0*n+c0];
1267 for (unsigned c=c0+1; c<n; ++c) {
1268 this->m[r2*n+c] -= piv * this->m[r0*n+c];
1269 if (!this->m[r2*n+c].info(info_flags::numeric))
1270 this->m[r2*n+c] = this->m[r2*n+c].normal();
1273 // fill up left hand side with zeros
1274 for (unsigned c=r0; c<=c0; ++c)
1275 this->m[r2*n+c] = _ex0;
1278 // save space by deleting no longer needed elements
1279 for (unsigned c=r0+1; c<n; ++c)
1280 this->m[r0*n+c] = _ex0;
1285 // clear remaining rows
1286 for (unsigned r=r0+1; r<m; ++r) {
1287 for (unsigned c=0; c<n; ++c)
1288 this->m[r*n+c] = _ex0;
1295 /** Perform the steps of division free elimination to bring the m x n matrix
1296 * into an upper echelon form.
1298 * @param det may be set to true to save a lot of space if one is only
1299 * interested in the diagonal elements (i.e. for calculating determinants).
1300 * The others are set to zero in this case.
1301 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1302 * number of rows was swapped and 0 if the matrix is singular. */
1303 int matrix::division_free_elimination(const bool det)
1305 ensure_if_modifiable();
1306 const unsigned m = this->rows();
1307 const unsigned n = this->cols();
1308 GINAC_ASSERT(!det || n==m);
1312 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1313 int indx = pivot(r0, c0, true);
1317 return 0; // leaves *this in a messy state
1322 for (unsigned r2=r0+1; r2<m; ++r2) {
1323 for (unsigned c=c0+1; c<n; ++c)
1324 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();
1325 // fill up left hand side with zeros
1326 for (unsigned c=r0; c<=c0; ++c)
1327 this->m[r2*n+c] = _ex0;
1330 // save space by deleting no longer needed elements
1331 for (unsigned c=r0+1; c<n; ++c)
1332 this->m[r0*n+c] = _ex0;
1337 // clear remaining rows
1338 for (unsigned r=r0+1; r<m; ++r) {
1339 for (unsigned c=0; c<n; ++c)
1340 this->m[r*n+c] = _ex0;
1347 /** Perform the steps of Bareiss' one-step fraction free elimination to bring
1348 * the matrix into an upper echelon form. Fraction free elimination means
1349 * that divide is used straightforwardly, without computing GCDs first. This
1350 * is possible, since we know the divisor at each step.
1352 * @param det may be set to true to save a lot of space if one is only
1353 * interested in the last element (i.e. for calculating determinants). The
1354 * others are set to zero in this case.
1355 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1356 * number of rows was swapped and 0 if the matrix is singular. */
1357 int matrix::fraction_free_elimination(const bool det)
1360 // (single-step fraction free elimination scheme, already known to Jordan)
1362 // Usual division-free elimination sets m[0](r,c) = m(r,c) and then sets
1363 // m[k+1](r,c) = m[k](k,k) * m[k](r,c) - m[k](r,k) * m[k](k,c).
1365 // Bareiss (fraction-free) elimination in addition divides that element
1366 // by m[k-1](k-1,k-1) for k>1, where it can be shown by means of the
1367 // Sylvester identity that this really divides m[k+1](r,c).
1369 // We also allow rational functions where the original prove still holds.
1370 // However, we must care for numerator and denominator separately and
1371 // "manually" work in the integral domains because of subtle cancellations
1372 // (see below). This blows up the bookkeeping a bit and the formula has
1373 // to be modified to expand like this (N{x} stands for numerator of x,
1374 // D{x} for denominator of x):
1375 // 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)}
1376 // -N{m[k](r,k)}*N{m[k](k,c)}*D{m[k](k,k)}*D{m[k](r,c)}
1377 // 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)}
1378 // where for k>1 we now divide N{m[k+1](r,c)} by
1379 // N{m[k-1](k-1,k-1)}
1380 // and D{m[k+1](r,c)} by
1381 // D{m[k-1](k-1,k-1)}.
1383 ensure_if_modifiable();
1384 const unsigned m = this->rows();
1385 const unsigned n = this->cols();
1386 GINAC_ASSERT(!det || n==m);
1395 // We populate temporary matrices to subsequently operate on. There is
1396 // one holding numerators and another holding denominators of entries.
1397 // This is a must since the evaluator (or even earlier mul's constructor)
1398 // might cancel some trivial element which causes divide() to fail. The
1399 // elements are normalized first (yes, even though this algorithm doesn't
1400 // need GCDs) since the elements of *this might be unnormalized, which
1401 // makes things more complicated than they need to be.
1402 matrix tmp_n(*this);
1403 matrix tmp_d(m,n); // for denominators, if needed
1404 exmap srl; // symbol replacement list
1405 exvector::const_iterator cit = this->m.begin(), citend = this->m.end();
1406 exvector::iterator tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
1407 while (cit != citend) {
1408 ex nd = cit->normal().to_rational(srl).numer_denom();
1410 *tmp_n_it++ = nd.op(0);
1411 *tmp_d_it++ = nd.op(1);
1415 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1416 int indx = tmp_n.pivot(r0, c0, true);
1425 // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
1426 for (unsigned c=c0; c<n; ++c)
1427 tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
1429 for (unsigned r2=r0+1; r2<m; ++r2) {
1430 for (unsigned c=c0+1; c<n; ++c) {
1431 dividend_n = (tmp_n.m[r0*n+c0]*tmp_n.m[r2*n+c]*
1432 tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]
1433 -tmp_n.m[r2*n+c0]*tmp_n.m[r0*n+c]*
1434 tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
1435 dividend_d = (tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]*
1436 tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
1437 bool check = divide(dividend_n, divisor_n,
1438 tmp_n.m[r2*n+c], true);
1439 check &= divide(dividend_d, divisor_d,
1440 tmp_d.m[r2*n+c], true);
1441 GINAC_ASSERT(check);
1443 // fill up left hand side with zeros
1444 for (unsigned c=r0; c<=c0; ++c)
1445 tmp_n.m[r2*n+c] = _ex0;
1447 if (c0<n && r0<m-1) {
1448 // compute next iteration's divisor
1449 divisor_n = tmp_n.m[r0*n+c0].expand();
1450 divisor_d = tmp_d.m[r0*n+c0].expand();
1452 // save space by deleting no longer needed elements
1453 for (unsigned c=0; c<n; ++c) {
1454 tmp_n.m[r0*n+c] = _ex0;
1455 tmp_d.m[r0*n+c] = _ex1;
1462 // clear remaining rows
1463 for (unsigned r=r0+1; r<m; ++r) {
1464 for (unsigned c=0; c<n; ++c)
1465 tmp_n.m[r*n+c] = _ex0;
1468 // repopulate *this matrix:
1469 exvector::iterator it = this->m.begin(), itend = this->m.end();
1470 tmp_n_it = tmp_n.m.begin();
1471 tmp_d_it = tmp_d.m.begin();
1473 *it++ = ((*tmp_n_it++)/(*tmp_d_it++)).subs(srl, subs_options::no_pattern);
1479 /** Partial pivoting method for matrix elimination schemes.
1480 * Usual pivoting (symbolic==false) returns the index to the element with the
1481 * largest absolute value in column ro and swaps the current row with the one
1482 * where the element was found. With (symbolic==true) it does the same thing
1483 * with the first non-zero element.
1485 * @param ro is the row from where to begin
1486 * @param co is the column to be inspected
1487 * @param symbolic signal if we want the first non-zero element to be pivoted
1488 * (true) or the one with the largest absolute value (false).
1489 * @return 0 if no interchange occured, -1 if all are zero (usually signaling
1490 * a degeneracy) and positive integer k means that rows ro and k were swapped.
1492 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
1496 // search first non-zero element in column co beginning at row ro
1497 while ((k<row) && (this->m[k*col+co].expand().is_zero()))
1500 // search largest element in column co beginning at row ro
1501 GINAC_ASSERT(is_exactly_a<numeric>(this->m[k*col+co]));
1502 unsigned kmax = k+1;
1503 numeric mmax = abs(ex_to<numeric>(m[kmax*col+co]));
1505 GINAC_ASSERT(is_exactly_a<numeric>(this->m[kmax*col+co]));
1506 numeric tmp = ex_to<numeric>(this->m[kmax*col+co]);
1507 if (abs(tmp) > mmax) {
1513 if (!mmax.is_zero())
1517 // all elements in column co below row ro vanish
1520 // matrix needs no pivoting
1522 // matrix needs pivoting, so swap rows k and ro
1523 ensure_if_modifiable();
1524 for (unsigned c=0; c<col; ++c)
1525 this->m[k*col+c].swap(this->m[ro*col+c]);
1530 /** Function to check that all elements of the matrix are zero.
1532 bool matrix::is_zero_matrix() const
1534 for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
1540 ex lst_to_matrix(const lst & l)
1542 lst::const_iterator itr, itc;
1544 // Find number of rows and columns
1545 size_t rows = l.nops(), cols = 0;
1546 for (itr = l.begin(); itr != l.end(); ++itr) {
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)
1553 // Allocate and fill matrix
1554 matrix &M = *new matrix(rows, cols);
1555 M.setflag(status_flags::dynallocated);
1558 for (itr = l.begin(), i = 0; itr != l.end(); ++itr, ++i) {
1560 for (itc = ex_to<lst>(*itr).begin(), j = 0; itc != ex_to<lst>(*itr).end(); ++itc, ++j)
1567 ex diag_matrix(const lst & l)
1569 lst::const_iterator it;
1570 size_t dim = l.nops();
1572 // Allocate and fill matrix
1573 matrix &M = *new matrix(dim, dim);
1574 M.setflag(status_flags::dynallocated);
1577 for (it = l.begin(), i = 0; it != l.end(); ++it, ++i)
1583 ex unit_matrix(unsigned r, unsigned c)
1585 matrix &Id = *new matrix(r, c);
1586 Id.setflag(status_flags::dynallocated);
1587 for (unsigned i=0; i<r && i<c; i++)
1593 ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const std::string & tex_base_name)
1595 matrix &M = *new matrix(r, c);
1596 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1598 bool long_format = (r > 10 || c > 10);
1599 bool single_row = (r == 1 || c == 1);
1601 for (unsigned i=0; i<r; i++) {
1602 for (unsigned j=0; j<c; j++) {
1603 std::ostringstream s1, s2;
1605 s2 << tex_base_name << "_{";
1616 s1 << '_' << i << '_' << j;
1617 s2 << i << ';' << j << "}";
1620 s2 << i << j << '}';
1623 M(i, j) = symbol(s1.str(), s2.str());
1630 ex reduced_matrix(const matrix& m, unsigned r, unsigned c)
1632 if (r+1>m.rows() || c+1>m.cols() || m.cols()<2 || m.rows()<2)
1633 throw std::runtime_error("minor_matrix(): index out of bounds");
1635 const unsigned rows = m.rows()-1;
1636 const unsigned cols = m.cols()-1;
1637 matrix &M = *new matrix(rows, cols);
1638 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1650 M(ro2,co2) = m(ro, co);
1661 ex sub_matrix(const matrix&m, unsigned r, unsigned nr, unsigned c, unsigned nc)
1663 if (r+nr>m.rows() || c+nc>m.cols())
1664 throw std::runtime_error("sub_matrix(): index out of bounds");
1666 matrix &M = *new matrix(nr, nc);
1667 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1669 for (unsigned ro=0; ro<nr; ++ro) {
1670 for (unsigned co=0; co<nc; ++co) {
1671 M(ro,co) = m(ro+r,co+c);
1678 } // namespace GiNaC