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
32 #include "relational.h"
41 #include "operators.h"
48 GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(matrix, basic,
49 print_func<print_context>(&matrix::do_print).
50 print_func<print_latex>(&matrix::do_print_latex).
51 print_func<print_tree>(&matrix::do_print_tree).
52 print_func<print_python_repr>(&matrix::do_print_python_repr))
55 // default constructor
58 /** Default ctor. Initializes to 1 x 1-dimensional zero-matrix. */
59 matrix::matrix() : inherited(&matrix::tinfo_static), row(1), col(1), m(1, _ex0)
61 setflag(status_flags::not_shareable);
70 /** Very common ctor. Initializes to r x c-dimensional zero-matrix.
72 * @param r number of rows
73 * @param c number of cols */
74 matrix::matrix(unsigned r, unsigned c)
75 : inherited(&matrix::tinfo_static), row(r), col(c), m(r*c, _ex0)
77 setflag(status_flags::not_shareable);
82 /** Ctor from representation, for internal use only. */
83 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
84 : inherited(&matrix::tinfo_static), row(r), col(c), m(m2)
86 setflag(status_flags::not_shareable);
89 /** Construct matrix from (flat) list of elements. If the list has fewer
90 * elements than the matrix, the remaining matrix elements are set to zero.
91 * If the list has more elements than the matrix, the excessive elements are
93 matrix::matrix(unsigned r, unsigned c, const lst & l)
94 : inherited(&matrix::tinfo_static), row(r), col(c), m(r*c, _ex0)
96 setflag(status_flags::not_shareable);
99 for (lst::const_iterator it = l.begin(); it != l.end(); ++it, ++i) {
103 break; // matrix smaller than list: throw away excessive elements
112 matrix::matrix(const archive_node &n, lst &sym_lst) : inherited(n, sym_lst)
114 setflag(status_flags::not_shareable);
116 if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
117 throw (std::runtime_error("unknown matrix dimensions in archive"));
118 m.reserve(row * col);
119 for (unsigned int i=0; true; i++) {
121 if (n.find_ex("m", e, sym_lst, i))
128 void matrix::archive(archive_node &n) const
130 inherited::archive(n);
131 n.add_unsigned("row", row);
132 n.add_unsigned("col", col);
133 exvector::const_iterator i = m.begin(), iend = m.end();
140 DEFAULT_UNARCHIVE(matrix)
143 // functions overriding virtual functions from base classes
148 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
150 for (unsigned ro=0; ro<row; ++ro) {
152 for (unsigned co=0; co<col; ++co) {
153 m[ro*col+co].print(c);
164 void matrix::do_print(const print_context & c, unsigned level) const
167 print_elements(c, "[", "]", ",", ",");
171 void matrix::do_print_latex(const print_latex & c, unsigned level) const
173 c.s << "\\left(\\begin{array}{" << std::string(col,'c') << "}";
174 print_elements(c, "", "", "\\\\", "&");
175 c.s << "\\end{array}\\right)";
178 void matrix::do_print_python_repr(const print_python_repr & c, unsigned level) const
180 c.s << class_name() << '(';
181 print_elements(c, "[", "]", ",", ",");
185 /** nops is defined to be rows x columns. */
186 size_t matrix::nops() const
188 return static_cast<size_t>(row) * static_cast<size_t>(col);
191 /** returns matrix entry at position (i/col, i%col). */
192 ex matrix::op(size_t i) const
194 GINAC_ASSERT(i<nops());
199 /** returns writable matrix entry at position (i/col, i%col). */
200 ex & matrix::let_op(size_t i)
202 GINAC_ASSERT(i<nops());
204 ensure_if_modifiable();
208 /** Evaluate matrix entry by entry. */
209 ex matrix::eval(int level) const
211 // check if we have to do anything at all
212 if ((level==1)&&(flags & status_flags::evaluated))
216 if (level == -max_recursion_level)
217 throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
219 // eval() entry by entry
220 exvector m2(row*col);
222 for (unsigned r=0; r<row; ++r)
223 for (unsigned c=0; c<col; ++c)
224 m2[r*col+c] = m[r*col+c].eval(level);
226 return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
227 status_flags::evaluated);
230 ex matrix::subs(const exmap & mp, unsigned options) const
232 exvector m2(row * col);
233 for (unsigned r=0; r<row; ++r)
234 for (unsigned c=0; c<col; ++c)
235 m2[r*col+c] = m[r*col+c].subs(mp, options);
237 return matrix(row, col, m2).subs_one_level(mp, options);
240 /** Complex conjugate every matrix entry. */
241 ex matrix::conjugate() const
244 for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i) {
245 ex x = i->conjugate();
250 if (are_ex_trivially_equal(x, *i)) {
254 ev->reserve(m.size());
255 for (exvector::const_iterator j=m.begin(); j!=i; ++j) {
261 ex result = matrix(row, col, *ev);
270 int matrix::compare_same_type(const basic & other) const
272 GINAC_ASSERT(is_exactly_a<matrix>(other));
273 const matrix &o = static_cast<const matrix &>(other);
275 // compare number of rows
277 return row < o.rows() ? -1 : 1;
279 // compare number of columns
281 return col < o.cols() ? -1 : 1;
283 // equal number of rows and columns, compare individual elements
285 for (unsigned r=0; r<row; ++r) {
286 for (unsigned c=0; c<col; ++c) {
287 cmpval = ((*this)(r,c)).compare(o(r,c));
288 if (cmpval!=0) return cmpval;
291 // all elements are equal => matrices are equal;
295 bool matrix::match_same_type(const basic & other) const
297 GINAC_ASSERT(is_exactly_a<matrix>(other));
298 const matrix & o = static_cast<const matrix &>(other);
300 // The number of rows and columns must be the same. This is necessary to
301 // prevent a 2x3 matrix from matching a 3x2 one.
302 return row == o.rows() && col == o.cols();
305 /** Automatic symbolic evaluation of an indexed matrix. */
306 ex matrix::eval_indexed(const basic & i) const
308 GINAC_ASSERT(is_a<indexed>(i));
309 GINAC_ASSERT(is_a<matrix>(i.op(0)));
311 bool all_indices_unsigned = static_cast<const indexed &>(i).all_index_values_are(info_flags::nonnegint);
316 // One index, must be one-dimensional vector
317 if (row != 1 && col != 1)
318 throw (std::runtime_error("matrix::eval_indexed(): vector must have exactly 1 index"));
320 const idx & i1 = ex_to<idx>(i.op(1));
325 if (!i1.get_dim().is_equal(row))
326 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
328 // Index numeric -> return vector element
329 if (all_indices_unsigned) {
330 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
332 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
333 return (*this)(n1, 0);
339 if (!i1.get_dim().is_equal(col))
340 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
342 // Index numeric -> return vector element
343 if (all_indices_unsigned) {
344 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
346 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
347 return (*this)(0, n1);
351 } else if (i.nops() == 3) {
354 const idx & i1 = ex_to<idx>(i.op(1));
355 const idx & i2 = ex_to<idx>(i.op(2));
357 if (!i1.get_dim().is_equal(row))
358 throw (std::runtime_error("matrix::eval_indexed(): dimension of first index must match number of rows"));
359 if (!i2.get_dim().is_equal(col))
360 throw (std::runtime_error("matrix::eval_indexed(): dimension of second index must match number of columns"));
362 // Pair of dummy indices -> compute trace
363 if (is_dummy_pair(i1, i2))
366 // Both indices numeric -> return matrix element
367 if (all_indices_unsigned) {
368 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int(), n2 = ex_to<numeric>(i2.get_value()).to_int();
370 throw (std::runtime_error("matrix::eval_indexed(): value of first index exceeds number of rows"));
372 throw (std::runtime_error("matrix::eval_indexed(): value of second index exceeds number of columns"));
373 return (*this)(n1, n2);
377 throw (std::runtime_error("matrix::eval_indexed(): matrix must have exactly 2 indices"));
382 /** Sum of two indexed matrices. */
383 ex matrix::add_indexed(const ex & self, const ex & other) const
385 GINAC_ASSERT(is_a<indexed>(self));
386 GINAC_ASSERT(is_a<matrix>(self.op(0)));
387 GINAC_ASSERT(is_a<indexed>(other));
388 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
390 // Only add two matrices
391 if (is_a<matrix>(other.op(0))) {
392 GINAC_ASSERT(other.nops() == 2 || other.nops() == 3);
394 const matrix &self_matrix = ex_to<matrix>(self.op(0));
395 const matrix &other_matrix = ex_to<matrix>(other.op(0));
397 if (self.nops() == 2 && other.nops() == 2) { // vector + vector
399 if (self_matrix.row == other_matrix.row)
400 return indexed(self_matrix.add(other_matrix), self.op(1));
401 else if (self_matrix.row == other_matrix.col)
402 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1));
404 } else if (self.nops() == 3 && other.nops() == 3) { // matrix + matrix
406 if (self.op(1).is_equal(other.op(1)) && self.op(2).is_equal(other.op(2)))
407 return indexed(self_matrix.add(other_matrix), self.op(1), self.op(2));
408 else if (self.op(1).is_equal(other.op(2)) && self.op(2).is_equal(other.op(1)))
409 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1), self.op(2));
414 // Don't know what to do, return unevaluated sum
418 /** Product of an indexed matrix with a number. */
419 ex matrix::scalar_mul_indexed(const ex & self, const numeric & other) const
421 GINAC_ASSERT(is_a<indexed>(self));
422 GINAC_ASSERT(is_a<matrix>(self.op(0)));
423 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
425 const matrix &self_matrix = ex_to<matrix>(self.op(0));
427 if (self.nops() == 2)
428 return indexed(self_matrix.mul(other), self.op(1));
429 else // self.nops() == 3
430 return indexed(self_matrix.mul(other), self.op(1), self.op(2));
433 /** Contraction of an indexed matrix with something else. */
434 bool matrix::contract_with(exvector::iterator self, exvector::iterator other, exvector & v) const
436 GINAC_ASSERT(is_a<indexed>(*self));
437 GINAC_ASSERT(is_a<indexed>(*other));
438 GINAC_ASSERT(self->nops() == 2 || self->nops() == 3);
439 GINAC_ASSERT(is_a<matrix>(self->op(0)));
441 // Only contract with other matrices
442 if (!is_a<matrix>(other->op(0)))
445 GINAC_ASSERT(other->nops() == 2 || other->nops() == 3);
447 const matrix &self_matrix = ex_to<matrix>(self->op(0));
448 const matrix &other_matrix = ex_to<matrix>(other->op(0));
450 if (self->nops() == 2) {
452 if (other->nops() == 2) { // vector * vector (scalar product)
454 if (self_matrix.col == 1) {
455 if (other_matrix.col == 1) {
456 // Column vector * column vector, transpose first vector
457 *self = self_matrix.transpose().mul(other_matrix)(0, 0);
459 // Column vector * row vector, swap factors
460 *self = other_matrix.mul(self_matrix)(0, 0);
463 if (other_matrix.col == 1) {
464 // Row vector * column vector, perfect
465 *self = self_matrix.mul(other_matrix)(0, 0);
467 // Row vector * row vector, transpose second vector
468 *self = self_matrix.mul(other_matrix.transpose())(0, 0);
474 } else { // vector * matrix
476 // B_i * A_ij = (B*A)_j (B is row vector)
477 if (is_dummy_pair(self->op(1), other->op(1))) {
478 if (self_matrix.row == 1)
479 *self = indexed(self_matrix.mul(other_matrix), other->op(2));
481 *self = indexed(self_matrix.transpose().mul(other_matrix), other->op(2));
486 // B_j * A_ij = (A*B)_i (B is column vector)
487 if (is_dummy_pair(self->op(1), other->op(2))) {
488 if (self_matrix.col == 1)
489 *self = indexed(other_matrix.mul(self_matrix), other->op(1));
491 *self = indexed(other_matrix.mul(self_matrix.transpose()), other->op(1));
497 } else if (other->nops() == 3) { // matrix * matrix
499 // A_ij * B_jk = (A*B)_ik
500 if (is_dummy_pair(self->op(2), other->op(1))) {
501 *self = indexed(self_matrix.mul(other_matrix), self->op(1), other->op(2));
506 // A_ij * B_kj = (A*Btrans)_ik
507 if (is_dummy_pair(self->op(2), other->op(2))) {
508 *self = indexed(self_matrix.mul(other_matrix.transpose()), self->op(1), other->op(1));
513 // A_ji * B_jk = (Atrans*B)_ik
514 if (is_dummy_pair(self->op(1), other->op(1))) {
515 *self = indexed(self_matrix.transpose().mul(other_matrix), self->op(2), other->op(2));
520 // A_ji * B_kj = (B*A)_ki
521 if (is_dummy_pair(self->op(1), other->op(2))) {
522 *self = indexed(other_matrix.mul(self_matrix), other->op(1), self->op(2));
533 // non-virtual functions in this class
540 * @exception logic_error (incompatible matrices) */
541 matrix matrix::add(const matrix & other) const
543 if (col != other.col || row != other.row)
544 throw std::logic_error("matrix::add(): incompatible matrices");
546 exvector sum(this->m);
547 exvector::iterator i = sum.begin(), end = sum.end();
548 exvector::const_iterator ci = other.m.begin();
552 return matrix(row,col,sum);
556 /** Difference of matrices.
558 * @exception logic_error (incompatible matrices) */
559 matrix matrix::sub(const matrix & other) const
561 if (col != other.col || row != other.row)
562 throw std::logic_error("matrix::sub(): incompatible matrices");
564 exvector dif(this->m);
565 exvector::iterator i = dif.begin(), end = dif.end();
566 exvector::const_iterator ci = other.m.begin();
570 return matrix(row,col,dif);
574 /** Product of matrices.
576 * @exception logic_error (incompatible matrices) */
577 matrix matrix::mul(const matrix & other) const
579 if (this->cols() != other.rows())
580 throw std::logic_error("matrix::mul(): incompatible matrices");
582 exvector prod(this->rows()*other.cols());
584 for (unsigned r1=0; r1<this->rows(); ++r1) {
585 for (unsigned c=0; c<this->cols(); ++c) {
586 // Quick test: can we shortcut?
587 if (m[r1*col+c].is_zero())
589 for (unsigned r2=0; r2<other.cols(); ++r2)
590 prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]);
593 return matrix(row, other.col, prod);
597 /** Product of matrix and scalar. */
598 matrix matrix::mul(const numeric & other) const
600 exvector prod(row * col);
602 for (unsigned r=0; r<row; ++r)
603 for (unsigned c=0; c<col; ++c)
604 prod[r*col+c] = m[r*col+c] * other;
606 return matrix(row, col, prod);
610 /** Product of matrix and scalar expression. */
611 matrix matrix::mul_scalar(const ex & other) const
613 if (other.return_type() != return_types::commutative)
614 throw std::runtime_error("matrix::mul_scalar(): non-commutative scalar");
616 exvector prod(row * col);
618 for (unsigned r=0; r<row; ++r)
619 for (unsigned c=0; c<col; ++c)
620 prod[r*col+c] = m[r*col+c] * other;
622 return matrix(row, col, prod);
626 /** Power of a matrix. Currently handles integer exponents only. */
627 matrix matrix::pow(const ex & expn) const
630 throw (std::logic_error("matrix::pow(): matrix not square"));
632 if (is_exactly_a<numeric>(expn)) {
633 // Integer cases are computed by successive multiplication, using the
634 // obvious shortcut of storing temporaries, like A^4 == (A*A)*(A*A).
635 if (expn.info(info_flags::integer)) {
636 numeric b = ex_to<numeric>(expn);
638 if (expn.info(info_flags::negative)) {
645 for (unsigned r=0; r<row; ++r)
649 // This loop computes the representation of b in base 2 from right
650 // to left and multiplies the factors whenever needed. Note
651 // that this is not entirely optimal but close to optimal and
652 // "better" algorithms are much harder to implement. (See Knuth,
653 // TAoCP2, section "Evaluation of Powers" for a good discussion.)
654 while (b!=*_num1_p) {
659 b /= *_num2_p; // still integer.
665 throw (std::runtime_error("matrix::pow(): don't know how to handle exponent"));
669 /** operator() to access elements for reading.
671 * @param ro row of element
672 * @param co column of element
673 * @exception range_error (index out of range) */
674 const ex & matrix::operator() (unsigned ro, unsigned co) const
676 if (ro>=row || co>=col)
677 throw (std::range_error("matrix::operator(): index out of range"));
683 /** operator() to access elements for writing.
685 * @param ro row of element
686 * @param co column of element
687 * @exception range_error (index out of range) */
688 ex & matrix::operator() (unsigned ro, unsigned co)
690 if (ro>=row || co>=col)
691 throw (std::range_error("matrix::operator(): index out of range"));
693 ensure_if_modifiable();
698 /** Transposed of an m x n matrix, producing a new n x m matrix object that
699 * represents the transposed. */
700 matrix matrix::transpose() const
702 exvector trans(this->cols()*this->rows());
704 for (unsigned r=0; r<this->cols(); ++r)
705 for (unsigned c=0; c<this->rows(); ++c)
706 trans[r*this->rows()+c] = m[c*this->cols()+r];
708 return matrix(this->cols(),this->rows(),trans);
711 /** Determinant of square matrix. This routine doesn't actually calculate the
712 * determinant, it only implements some heuristics about which algorithm to
713 * run. If all the elements of the matrix are elements of an integral domain
714 * the determinant is also in that integral domain and the result is expanded
715 * only. If one or more elements are from a quotient field the determinant is
716 * usually also in that quotient field and the result is normalized before it
717 * is returned. This implies that the determinant of the symbolic 2x2 matrix
718 * [[a/(a-b),1],[b/(a-b),1]] is returned as unity. (In this respect, it
719 * behaves like MapleV and unlike Mathematica.)
721 * @param algo allows to chose an algorithm
722 * @return the determinant as a new expression
723 * @exception logic_error (matrix not square)
724 * @see determinant_algo */
725 ex matrix::determinant(unsigned algo) const
728 throw (std::logic_error("matrix::determinant(): matrix not square"));
729 GINAC_ASSERT(row*col==m.capacity());
731 // Gather some statistical information about this matrix:
732 bool numeric_flag = true;
733 bool normal_flag = false;
734 unsigned sparse_count = 0; // counts non-zero elements
735 exvector::const_iterator r = m.begin(), rend = m.end();
737 if (!r->info(info_flags::numeric))
738 numeric_flag = false;
739 exmap srl; // symbol replacement list
740 ex rtest = r->to_rational(srl);
741 if (!rtest.is_zero())
743 if (!rtest.info(info_flags::crational_polynomial) &&
744 rtest.info(info_flags::rational_function))
749 // Here is the heuristics in case this routine has to decide:
750 if (algo == determinant_algo::automatic) {
751 // Minor expansion is generally a good guess:
752 algo = determinant_algo::laplace;
753 // Does anybody know when a matrix is really sparse?
754 // Maybe <~row/2.236 nonzero elements average in a row?
755 if (row>3 && 5*sparse_count<=row*col)
756 algo = determinant_algo::bareiss;
757 // Purely numeric matrix can be handled by Gauss elimination.
758 // This overrides any prior decisions.
760 algo = determinant_algo::gauss;
763 // Trap the trivial case here, since some algorithms don't like it
765 // for consistency with non-trivial determinants...
767 return m[0].normal();
769 return m[0].expand();
772 // Compute the determinant
774 case determinant_algo::gauss: {
777 int sign = tmp.gauss_elimination(true);
778 for (unsigned d=0; d<row; ++d)
779 det *= tmp.m[d*col+d];
781 return (sign*det).normal();
783 return (sign*det).normal().expand();
785 case determinant_algo::bareiss: {
788 sign = tmp.fraction_free_elimination(true);
790 return (sign*tmp.m[row*col-1]).normal();
792 return (sign*tmp.m[row*col-1]).expand();
794 case determinant_algo::divfree: {
797 sign = tmp.division_free_elimination(true);
800 ex det = tmp.m[row*col-1];
801 // factor out accumulated bogus slag
802 for (unsigned d=0; d<row-2; ++d)
803 for (unsigned j=0; j<row-d-2; ++j)
804 det = (det/tmp.m[d*col+d]).normal();
807 case determinant_algo::laplace:
809 // This is the minor expansion scheme. We always develop such
810 // that the smallest minors (i.e, the trivial 1x1 ones) are on the
811 // rightmost column. For this to be efficient, empirical tests
812 // have shown that the emptiest columns (i.e. the ones with most
813 // zeros) should be the ones on the right hand side -- although
814 // this might seem counter-intuitive (and in contradiction to some
815 // literature like the FORM manual). Please go ahead and test it
816 // if you don't believe me! Therefore we presort the columns of
818 typedef std::pair<unsigned,unsigned> uintpair;
819 std::vector<uintpair> c_zeros; // number of zeros in column
820 for (unsigned c=0; c<col; ++c) {
822 for (unsigned r=0; r<row; ++r)
823 if (m[r*col+c].is_zero())
825 c_zeros.push_back(uintpair(acc,c));
827 std::sort(c_zeros.begin(),c_zeros.end());
828 std::vector<unsigned> pre_sort;
829 for (std::vector<uintpair>::const_iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
830 pre_sort.push_back(i->second);
831 std::vector<unsigned> pre_sort_test(pre_sort); // permutation_sign() modifies the vector so we make a copy here
832 int sign = permutation_sign(pre_sort_test.begin(), pre_sort_test.end());
833 exvector result(row*col); // represents sorted matrix
835 for (std::vector<unsigned>::const_iterator i=pre_sort.begin();
838 for (unsigned r=0; r<row; ++r)
839 result[r*col+c] = m[r*col+(*i)];
843 return (sign*matrix(row,col,result).determinant_minor()).normal();
845 return sign*matrix(row,col,result).determinant_minor();
851 /** Trace of a matrix. The result is normalized if it is in some quotient
852 * field and expanded only otherwise. This implies that the trace of the
853 * symbolic 2x2 matrix [[a/(a-b),x],[y,b/(b-a)]] is recognized to be unity.
855 * @return the sum of diagonal elements
856 * @exception logic_error (matrix not square) */
857 ex matrix::trace() const
860 throw (std::logic_error("matrix::trace(): matrix not square"));
863 for (unsigned r=0; r<col; ++r)
866 if (tr.info(info_flags::rational_function) &&
867 !tr.info(info_flags::crational_polynomial))
874 /** Characteristic Polynomial. Following mathematica notation the
875 * characteristic polynomial of a matrix M is defined as the determiant of
876 * (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
877 * as M. Note that some CASs define it with a sign inside the determinant
878 * which gives rise to an overall sign if the dimension is odd. This method
879 * returns the characteristic polynomial collected in powers of lambda as a
882 * @return characteristic polynomial as new expression
883 * @exception logic_error (matrix not square)
884 * @see matrix::determinant() */
885 ex matrix::charpoly(const ex & lambda) const
888 throw (std::logic_error("matrix::charpoly(): matrix not square"));
890 bool numeric_flag = true;
891 exvector::const_iterator r = m.begin(), rend = m.end();
892 while (r!=rend && numeric_flag==true) {
893 if (!r->info(info_flags::numeric))
894 numeric_flag = false;
898 // The pure numeric case is traditionally rather common. Hence, it is
899 // trapped and we use Leverrier's algorithm which goes as row^3 for
900 // every coefficient. The expensive part is the matrix multiplication.
905 ex poly = power(lambda, row) - c*power(lambda, row-1);
906 for (unsigned i=1; i<row; ++i) {
907 for (unsigned j=0; j<row; ++j)
910 c = B.trace() / ex(i+1);
911 poly -= c*power(lambda, row-i-1);
921 for (unsigned r=0; r<col; ++r)
922 M.m[r*col+r] -= lambda;
924 return M.determinant().collect(lambda);
929 /** Inverse of this matrix.
931 * @return the inverted matrix
932 * @exception logic_error (matrix not square)
933 * @exception runtime_error (singular matrix) */
934 matrix matrix::inverse() const
937 throw (std::logic_error("matrix::inverse(): matrix not square"));
939 // This routine actually doesn't do anything fancy at all. We compute the
940 // inverse of the matrix A by solving the system A * A^{-1} == Id.
942 // First populate the identity matrix supposed to become the right hand side.
943 matrix identity(row,col);
944 for (unsigned i=0; i<row; ++i)
945 identity(i,i) = _ex1;
947 // Populate a dummy matrix of variables, just because of compatibility with
948 // matrix::solve() which wants this (for compatibility with under-determined
949 // systems of equations).
950 matrix vars(row,col);
951 for (unsigned r=0; r<row; ++r)
952 for (unsigned c=0; c<col; ++c)
953 vars(r,c) = symbol();
957 sol = this->solve(vars,identity);
958 } catch (const std::runtime_error & e) {
959 if (e.what()==std::string("matrix::solve(): inconsistent linear system"))
960 throw (std::runtime_error("matrix::inverse(): singular matrix"));
968 /** Solve a linear system consisting of a m x n matrix and a m x p right hand
969 * side by applying an elimination scheme to the augmented matrix.
971 * @param vars n x p matrix, all elements must be symbols
972 * @param rhs m x p matrix
973 * @param algo selects the solving algorithm
974 * @return n x p solution matrix
975 * @exception logic_error (incompatible matrices)
976 * @exception invalid_argument (1st argument must be matrix of symbols)
977 * @exception runtime_error (inconsistent linear system)
979 matrix matrix::solve(const matrix & vars,
983 const unsigned m = this->rows();
984 const unsigned n = this->cols();
985 const unsigned p = rhs.cols();
988 if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
989 throw (std::logic_error("matrix::solve(): incompatible matrices"));
990 for (unsigned ro=0; ro<n; ++ro)
991 for (unsigned co=0; co<p; ++co)
992 if (!vars(ro,co).info(info_flags::symbol))
993 throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
995 // build the augmented matrix of *this with rhs attached to the right
997 for (unsigned r=0; r<m; ++r) {
998 for (unsigned c=0; c<n; ++c)
999 aug.m[r*(n+p)+c] = this->m[r*n+c];
1000 for (unsigned c=0; c<p; ++c)
1001 aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
1004 // Gather some statistical information about the augmented matrix:
1005 bool numeric_flag = true;
1006 exvector::const_iterator r = aug.m.begin(), rend = aug.m.end();
1007 while (r!=rend && numeric_flag==true) {
1008 if (!r->info(info_flags::numeric))
1009 numeric_flag = false;
1013 // Here is the heuristics in case this routine has to decide:
1014 if (algo == solve_algo::automatic) {
1015 // Bareiss (fraction-free) elimination is generally a good guess:
1016 algo = solve_algo::bareiss;
1017 // For m<3, Bareiss elimination is equivalent to division free
1018 // elimination but has more logistic overhead
1020 algo = solve_algo::divfree;
1021 // This overrides any prior decisions.
1023 algo = solve_algo::gauss;
1026 // Eliminate the augmented matrix:
1028 case solve_algo::gauss:
1029 aug.gauss_elimination();
1031 case solve_algo::divfree:
1032 aug.division_free_elimination();
1034 case solve_algo::bareiss:
1036 aug.fraction_free_elimination();
1039 // assemble the solution matrix:
1041 for (unsigned co=0; co<p; ++co) {
1042 unsigned last_assigned_sol = n+1;
1043 for (int r=m-1; r>=0; --r) {
1044 unsigned fnz = 1; // first non-zero in row
1045 while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
1048 // row consists only of zeros, corresponding rhs must be 0, too
1049 if (!aug.m[r*(n+p)+n+co].is_zero()) {
1050 throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
1053 // assign solutions for vars between fnz+1 and
1054 // last_assigned_sol-1: free parameters
1055 for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
1056 sol(c,co) = vars.m[c*p+co];
1057 ex e = aug.m[r*(n+p)+n+co];
1058 for (unsigned c=fnz; c<n; ++c)
1059 e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
1060 sol(fnz-1,co) = (e/(aug.m[r*(n+p)+(fnz-1)])).normal();
1061 last_assigned_sol = fnz;
1064 // assign solutions for vars between 1 and
1065 // last_assigned_sol-1: free parameters
1066 for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
1067 sol(ro,co) = vars(ro,co);
1074 /** Compute the rank of this matrix. */
1075 unsigned matrix::rank() const
1078 // Transform this matrix into upper echelon form and then count the
1079 // number of non-zero rows.
1081 GINAC_ASSERT(row*col==m.capacity());
1083 // Actually, any elimination scheme will do since we are only
1084 // interested in the echelon matrix' zeros.
1085 matrix to_eliminate = *this;
1086 to_eliminate.fraction_free_elimination();
1088 unsigned r = row*col; // index of last non-zero element
1090 if (!to_eliminate.m[r].is_zero())
1099 /** Recursive determinant for small matrices having at least one symbolic
1100 * entry. The basic algorithm, known as Laplace-expansion, is enhanced by
1101 * some bookkeeping to avoid calculation of the same submatrices ("minors")
1102 * more than once. According to W.M.Gentleman and S.C.Johnson this algorithm
1103 * is better than elimination schemes for matrices of sparse multivariate
1104 * polynomials and also for matrices of dense univariate polynomials if the
1105 * matrix' dimesion is larger than 7.
1107 * @return the determinant as a new expression (in expanded form)
1108 * @see matrix::determinant() */
1109 ex matrix::determinant_minor() const
1111 // for small matrices the algorithm does not make any sense:
1112 const unsigned n = this->cols();
1114 return m[0].expand();
1116 return (m[0]*m[3]-m[2]*m[1]).expand();
1118 return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
1119 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
1120 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
1122 // This algorithm can best be understood by looking at a naive
1123 // implementation of Laplace-expansion, like this one:
1125 // matrix minorM(this->rows()-1,this->cols()-1);
1126 // for (unsigned r1=0; r1<this->rows(); ++r1) {
1127 // // shortcut if element(r1,0) vanishes
1128 // if (m[r1*col].is_zero())
1130 // // assemble the minor matrix
1131 // for (unsigned r=0; r<minorM.rows(); ++r) {
1132 // for (unsigned c=0; c<minorM.cols(); ++c) {
1134 // minorM(r,c) = m[r*col+c+1];
1136 // minorM(r,c) = m[(r+1)*col+c+1];
1139 // // recurse down and care for sign:
1141 // det -= m[r1*col] * minorM.determinant_minor();
1143 // det += m[r1*col] * minorM.determinant_minor();
1145 // return det.expand();
1146 // What happens is that while proceeding down many of the minors are
1147 // computed more than once. In particular, there are binomial(n,k)
1148 // kxk minors and each one is computed factorial(n-k) times. Therefore
1149 // it is reasonable to store the results of the minors. We proceed from
1150 // right to left. At each column c we only need to retrieve the minors
1151 // calculated in step c-1. We therefore only have to store at most
1152 // 2*binomial(n,n/2) minors.
1154 // Unique flipper counter for partitioning into minors
1155 std::vector<unsigned> Pkey;
1157 // key for minor determinant (a subpartition of Pkey)
1158 std::vector<unsigned> Mkey;
1160 // we store our subminors in maps, keys being the rows they arise from
1161 typedef std::map<std::vector<unsigned>,class ex> Rmap;
1162 typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
1166 // initialize A with last column:
1167 for (unsigned r=0; r<n; ++r) {
1168 Pkey.erase(Pkey.begin(),Pkey.end());
1170 A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
1172 // proceed from right to left through matrix
1173 for (int c=n-2; c>=0; --c) {
1174 Pkey.erase(Pkey.begin(),Pkey.end()); // don't change capacity
1175 Mkey.erase(Mkey.begin(),Mkey.end());
1176 for (unsigned i=0; i<n-c; ++i)
1178 unsigned fc = 0; // controls logic for our strange flipper counter
1181 for (unsigned r=0; r<n-c; ++r) {
1182 // maybe there is nothing to do?
1183 if (m[Pkey[r]*n+c].is_zero())
1185 // create the sorted key for all possible minors
1186 Mkey.erase(Mkey.begin(),Mkey.end());
1187 for (unsigned i=0; i<n-c; ++i)
1189 Mkey.push_back(Pkey[i]);
1190 // Fetch the minors and compute the new determinant
1192 det -= m[Pkey[r]*n+c]*A[Mkey];
1194 det += m[Pkey[r]*n+c]*A[Mkey];
1196 // prevent build-up of deep nesting of expressions saves time:
1198 // store the new determinant at its place in B:
1200 B.insert(Rmap_value(Pkey,det));
1201 // increment our strange flipper counter
1202 for (fc=n-c; fc>0; --fc) {
1204 if (Pkey[fc-1]<fc+c)
1208 for (unsigned j=fc; j<n-c; ++j)
1209 Pkey[j] = Pkey[j-1]+1;
1211 // next column, so change the role of A and B:
1220 /** Perform the steps of an ordinary Gaussian elimination to bring the m x n
1221 * matrix into an upper echelon form. The algorithm is ok for matrices
1222 * with numeric coefficients but quite unsuited for symbolic matrices.
1224 * @param det may be set to true to save a lot of space if one is only
1225 * interested in the diagonal elements (i.e. for calculating determinants).
1226 * The others are set to zero in this case.
1227 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1228 * number of rows was swapped and 0 if the matrix is singular. */
1229 int matrix::gauss_elimination(const bool det)
1231 ensure_if_modifiable();
1232 const unsigned m = this->rows();
1233 const unsigned n = this->cols();
1234 GINAC_ASSERT(!det || n==m);
1238 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1239 int indx = pivot(r0, c0, true);
1243 return 0; // leaves *this in a messy state
1248 for (unsigned r2=r0+1; r2<m; ++r2) {
1249 if (!this->m[r2*n+c0].is_zero()) {
1250 // yes, there is something to do in this row
1251 ex piv = this->m[r2*n+c0] / this->m[r0*n+c0];
1252 for (unsigned c=c0+1; c<n; ++c) {
1253 this->m[r2*n+c] -= piv * this->m[r0*n+c];
1254 if (!this->m[r2*n+c].info(info_flags::numeric))
1255 this->m[r2*n+c] = this->m[r2*n+c].normal();
1258 // fill up left hand side with zeros
1259 for (unsigned c=r0; c<=c0; ++c)
1260 this->m[r2*n+c] = _ex0;
1263 // save space by deleting no longer needed elements
1264 for (unsigned c=r0+1; c<n; ++c)
1265 this->m[r0*n+c] = _ex0;
1270 // clear remaining rows
1271 for (unsigned r=r0+1; r<m; ++r) {
1272 for (unsigned c=0; c<n; ++c)
1273 this->m[r*n+c] = _ex0;
1280 /** Perform the steps of division free elimination to bring the m x n matrix
1281 * into an upper echelon form.
1283 * @param det may be set to true to save a lot of space if one is only
1284 * interested in the diagonal elements (i.e. for calculating determinants).
1285 * The others are set to zero in this case.
1286 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1287 * number of rows was swapped and 0 if the matrix is singular. */
1288 int matrix::division_free_elimination(const bool det)
1290 ensure_if_modifiable();
1291 const unsigned m = this->rows();
1292 const unsigned n = this->cols();
1293 GINAC_ASSERT(!det || n==m);
1297 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1298 int indx = pivot(r0, c0, true);
1302 return 0; // leaves *this in a messy state
1307 for (unsigned r2=r0+1; r2<m; ++r2) {
1308 for (unsigned c=c0+1; c<n; ++c)
1309 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();
1310 // fill up left hand side with zeros
1311 for (unsigned c=r0; c<=c0; ++c)
1312 this->m[r2*n+c] = _ex0;
1315 // save space by deleting no longer needed elements
1316 for (unsigned c=r0+1; c<n; ++c)
1317 this->m[r0*n+c] = _ex0;
1322 // clear remaining rows
1323 for (unsigned r=r0+1; r<m; ++r) {
1324 for (unsigned c=0; c<n; ++c)
1325 this->m[r*n+c] = _ex0;
1332 /** Perform the steps of Bareiss' one-step fraction free elimination to bring
1333 * the matrix into an upper echelon form. Fraction free elimination means
1334 * that divide is used straightforwardly, without computing GCDs first. This
1335 * is possible, since we know the divisor at each step.
1337 * @param det may be set to true to save a lot of space if one is only
1338 * interested in the last element (i.e. for calculating determinants). The
1339 * others are set to zero in this case.
1340 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1341 * number of rows was swapped and 0 if the matrix is singular. */
1342 int matrix::fraction_free_elimination(const bool det)
1345 // (single-step fraction free elimination scheme, already known to Jordan)
1347 // Usual division-free elimination sets m[0](r,c) = m(r,c) and then sets
1348 // m[k+1](r,c) = m[k](k,k) * m[k](r,c) - m[k](r,k) * m[k](k,c).
1350 // Bareiss (fraction-free) elimination in addition divides that element
1351 // by m[k-1](k-1,k-1) for k>1, where it can be shown by means of the
1352 // Sylvester identity that this really divides m[k+1](r,c).
1354 // We also allow rational functions where the original prove still holds.
1355 // However, we must care for numerator and denominator separately and
1356 // "manually" work in the integral domains because of subtle cancellations
1357 // (see below). This blows up the bookkeeping a bit and the formula has
1358 // to be modified to expand like this (N{x} stands for numerator of x,
1359 // D{x} for denominator of x):
1360 // 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)}
1361 // -N{m[k](r,k)}*N{m[k](k,c)}*D{m[k](k,k)}*D{m[k](r,c)}
1362 // 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)}
1363 // where for k>1 we now divide N{m[k+1](r,c)} by
1364 // N{m[k-1](k-1,k-1)}
1365 // and D{m[k+1](r,c)} by
1366 // D{m[k-1](k-1,k-1)}.
1368 ensure_if_modifiable();
1369 const unsigned m = this->rows();
1370 const unsigned n = this->cols();
1371 GINAC_ASSERT(!det || n==m);
1380 // We populate temporary matrices to subsequently operate on. There is
1381 // one holding numerators and another holding denominators of entries.
1382 // This is a must since the evaluator (or even earlier mul's constructor)
1383 // might cancel some trivial element which causes divide() to fail. The
1384 // elements are normalized first (yes, even though this algorithm doesn't
1385 // need GCDs) since the elements of *this might be unnormalized, which
1386 // makes things more complicated than they need to be.
1387 matrix tmp_n(*this);
1388 matrix tmp_d(m,n); // for denominators, if needed
1389 exmap srl; // symbol replacement list
1390 exvector::const_iterator cit = this->m.begin(), citend = this->m.end();
1391 exvector::iterator tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
1392 while (cit != citend) {
1393 ex nd = cit->normal().to_rational(srl).numer_denom();
1395 *tmp_n_it++ = nd.op(0);
1396 *tmp_d_it++ = nd.op(1);
1400 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1401 int indx = tmp_n.pivot(r0, c0, true);
1410 // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
1411 for (unsigned c=c0; c<n; ++c)
1412 tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
1414 for (unsigned r2=r0+1; r2<m; ++r2) {
1415 for (unsigned c=c0+1; c<n; ++c) {
1416 dividend_n = (tmp_n.m[r0*n+c0]*tmp_n.m[r2*n+c]*
1417 tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]
1418 -tmp_n.m[r2*n+c0]*tmp_n.m[r0*n+c]*
1419 tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
1420 dividend_d = (tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]*
1421 tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
1422 bool check = divide(dividend_n, divisor_n,
1423 tmp_n.m[r2*n+c], true);
1424 check &= divide(dividend_d, divisor_d,
1425 tmp_d.m[r2*n+c], true);
1426 GINAC_ASSERT(check);
1428 // fill up left hand side with zeros
1429 for (unsigned c=r0; c<=c0; ++c)
1430 tmp_n.m[r2*n+c] = _ex0;
1432 if (c0<n && r0<m-1) {
1433 // compute next iteration's divisor
1434 divisor_n = tmp_n.m[r0*n+c0].expand();
1435 divisor_d = tmp_d.m[r0*n+c0].expand();
1437 // save space by deleting no longer needed elements
1438 for (unsigned c=0; c<n; ++c) {
1439 tmp_n.m[r0*n+c] = _ex0;
1440 tmp_d.m[r0*n+c] = _ex1;
1447 // clear remaining rows
1448 for (unsigned r=r0+1; r<m; ++r) {
1449 for (unsigned c=0; c<n; ++c)
1450 tmp_n.m[r*n+c] = _ex0;
1453 // repopulate *this matrix:
1454 exvector::iterator it = this->m.begin(), itend = this->m.end();
1455 tmp_n_it = tmp_n.m.begin();
1456 tmp_d_it = tmp_d.m.begin();
1458 *it++ = ((*tmp_n_it++)/(*tmp_d_it++)).subs(srl, subs_options::no_pattern);
1464 /** Partial pivoting method for matrix elimination schemes.
1465 * Usual pivoting (symbolic==false) returns the index to the element with the
1466 * largest absolute value in column ro and swaps the current row with the one
1467 * where the element was found. With (symbolic==true) it does the same thing
1468 * with the first non-zero element.
1470 * @param ro is the row from where to begin
1471 * @param co is the column to be inspected
1472 * @param symbolic signal if we want the first non-zero element to be pivoted
1473 * (true) or the one with the largest absolute value (false).
1474 * @return 0 if no interchange occured, -1 if all are zero (usually signaling
1475 * a degeneracy) and positive integer k means that rows ro and k were swapped.
1477 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
1481 // search first non-zero element in column co beginning at row ro
1482 while ((k<row) && (this->m[k*col+co].expand().is_zero()))
1485 // search largest element in column co beginning at row ro
1486 GINAC_ASSERT(is_exactly_a<numeric>(this->m[k*col+co]));
1487 unsigned kmax = k+1;
1488 numeric mmax = abs_function::eval_numeric(ex_to<numeric>(m[kmax*col+co]));
1490 GINAC_ASSERT(is_exactly_a<numeric>(this->m[kmax*col+co]));
1491 numeric tmp = ex_to<numeric>(this->m[kmax*col+co]);
1492 if (abs_function::eval_numeric(tmp) > mmax) {
1498 if (!mmax.is_zero())
1502 // all elements in column co below row ro vanish
1505 // matrix needs no pivoting
1507 // matrix needs pivoting, so swap rows k and ro
1508 ensure_if_modifiable();
1509 for (unsigned c=0; c<col; ++c)
1510 this->m[k*col+c].swap(this->m[ro*col+c]);
1515 ex lst_to_matrix(const lst & l)
1517 lst::const_iterator itr, itc;
1519 // Find number of rows and columns
1520 size_t rows = l.nops(), cols = 0;
1521 for (itr = l.begin(); itr != l.end(); ++itr) {
1522 if (!is_a<lst>(*itr))
1523 throw (std::invalid_argument("lst_to_matrix: argument must be a list of lists"));
1524 if (itr->nops() > cols)
1528 // Allocate and fill matrix
1529 matrix &M = *new matrix(rows, cols);
1530 M.setflag(status_flags::dynallocated);
1533 for (itr = l.begin(), i = 0; itr != l.end(); ++itr, ++i) {
1535 for (itc = ex_to<lst>(*itr).begin(), j = 0; itc != ex_to<lst>(*itr).end(); ++itc, ++j)
1542 ex diag_matrix(const lst & l)
1544 lst::const_iterator it;
1545 size_t dim = l.nops();
1547 // Allocate and fill matrix
1548 matrix &M = *new matrix(dim, dim);
1549 M.setflag(status_flags::dynallocated);
1552 for (it = l.begin(), i = 0; it != l.end(); ++it, ++i)
1558 ex unit_matrix(unsigned r, unsigned c)
1560 matrix &Id = *new matrix(r, c);
1561 Id.setflag(status_flags::dynallocated);
1562 for (unsigned i=0; i<r && i<c; i++)
1568 ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const std::string & tex_base_name)
1570 matrix &M = *new matrix(r, c);
1571 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1573 bool long_format = (r > 10 || c > 10);
1574 bool single_row = (r == 1 || c == 1);
1576 for (unsigned i=0; i<r; i++) {
1577 for (unsigned j=0; j<c; j++) {
1578 std::ostringstream s1, s2;
1580 s2 << tex_base_name << "_{";
1591 s1 << '_' << i << '_' << j;
1592 s2 << i << ';' << j << "}";
1595 s2 << i << j << '}';
1598 M(i, j) = symbol(s1.str(), s2.str());
1605 ex reduced_matrix(const matrix& m, unsigned r, unsigned c)
1607 if (r+1>m.rows() || c+1>m.cols() || m.cols()<2 || m.rows()<2)
1608 throw std::runtime_error("minor_matrix(): index out of bounds");
1610 const unsigned rows = m.rows()-1;
1611 const unsigned cols = m.cols()-1;
1612 matrix &M = *new matrix(rows, cols);
1613 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1625 M(ro2,co2) = m(ro, co);
1636 ex sub_matrix(const matrix&m, unsigned r, unsigned nr, unsigned c, unsigned nc)
1638 if (r+nr>m.rows() || c+nc>m.cols())
1639 throw std::runtime_error("sub_matrix(): index out of bounds");
1641 matrix &M = *new matrix(nr, nc);
1642 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1644 for (unsigned ro=0; ro<nr; ++ro) {
1645 for (unsigned co=0; co<nc; ++co) {
1646 M(ro,co) = m(ro+r,co+c);
1653 } // namespace GiNaC