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);
267 int matrix::compare_same_type(const basic & other) const
269 GINAC_ASSERT(is_exactly_a<matrix>(other));
270 const matrix &o = static_cast<const matrix &>(other);
272 // compare number of rows
274 return row < o.rows() ? -1 : 1;
276 // compare number of columns
278 return col < o.cols() ? -1 : 1;
280 // equal number of rows and columns, compare individual elements
282 for (unsigned r=0; r<row; ++r) {
283 for (unsigned c=0; c<col; ++c) {
284 cmpval = ((*this)(r,c)).compare(o(r,c));
285 if (cmpval!=0) return cmpval;
288 // all elements are equal => matrices are equal;
292 bool matrix::match_same_type(const basic & other) const
294 GINAC_ASSERT(is_exactly_a<matrix>(other));
295 const matrix & o = static_cast<const matrix &>(other);
297 // The number of rows and columns must be the same. This is necessary to
298 // prevent a 2x3 matrix from matching a 3x2 one.
299 return row == o.rows() && col == o.cols();
302 /** Automatic symbolic evaluation of an indexed matrix. */
303 ex matrix::eval_indexed(const basic & i) const
305 GINAC_ASSERT(is_a<indexed>(i));
306 GINAC_ASSERT(is_a<matrix>(i.op(0)));
308 bool all_indices_unsigned = static_cast<const indexed &>(i).all_index_values_are(info_flags::nonnegint);
313 // One index, must be one-dimensional vector
314 if (row != 1 && col != 1)
315 throw (std::runtime_error("matrix::eval_indexed(): vector must have exactly 1 index"));
317 const idx & i1 = ex_to<idx>(i.op(1));
322 if (!i1.get_dim().is_equal(row))
323 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
325 // Index numeric -> return vector element
326 if (all_indices_unsigned) {
327 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
329 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
330 return (*this)(n1, 0);
336 if (!i1.get_dim().is_equal(col))
337 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
339 // Index numeric -> return vector element
340 if (all_indices_unsigned) {
341 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
343 throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
344 return (*this)(0, n1);
348 } else if (i.nops() == 3) {
351 const idx & i1 = ex_to<idx>(i.op(1));
352 const idx & i2 = ex_to<idx>(i.op(2));
354 if (!i1.get_dim().is_equal(row))
355 throw (std::runtime_error("matrix::eval_indexed(): dimension of first index must match number of rows"));
356 if (!i2.get_dim().is_equal(col))
357 throw (std::runtime_error("matrix::eval_indexed(): dimension of second index must match number of columns"));
359 // Pair of dummy indices -> compute trace
360 if (is_dummy_pair(i1, i2))
363 // Both indices numeric -> return matrix element
364 if (all_indices_unsigned) {
365 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int(), n2 = ex_to<numeric>(i2.get_value()).to_int();
367 throw (std::runtime_error("matrix::eval_indexed(): value of first index exceeds number of rows"));
369 throw (std::runtime_error("matrix::eval_indexed(): value of second index exceeds number of columns"));
370 return (*this)(n1, n2);
374 throw (std::runtime_error("matrix::eval_indexed(): matrix must have exactly 2 indices"));
379 /** Sum of two indexed matrices. */
380 ex matrix::add_indexed(const ex & self, const ex & other) const
382 GINAC_ASSERT(is_a<indexed>(self));
383 GINAC_ASSERT(is_a<matrix>(self.op(0)));
384 GINAC_ASSERT(is_a<indexed>(other));
385 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
387 // Only add two matrices
388 if (is_a<matrix>(other.op(0))) {
389 GINAC_ASSERT(other.nops() == 2 || other.nops() == 3);
391 const matrix &self_matrix = ex_to<matrix>(self.op(0));
392 const matrix &other_matrix = ex_to<matrix>(other.op(0));
394 if (self.nops() == 2 && other.nops() == 2) { // vector + vector
396 if (self_matrix.row == other_matrix.row)
397 return indexed(self_matrix.add(other_matrix), self.op(1));
398 else if (self_matrix.row == other_matrix.col)
399 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1));
401 } else if (self.nops() == 3 && other.nops() == 3) { // matrix + matrix
403 if (self.op(1).is_equal(other.op(1)) && self.op(2).is_equal(other.op(2)))
404 return indexed(self_matrix.add(other_matrix), self.op(1), self.op(2));
405 else if (self.op(1).is_equal(other.op(2)) && self.op(2).is_equal(other.op(1)))
406 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1), self.op(2));
411 // Don't know what to do, return unevaluated sum
415 /** Product of an indexed matrix with a number. */
416 ex matrix::scalar_mul_indexed(const ex & self, const numeric & other) const
418 GINAC_ASSERT(is_a<indexed>(self));
419 GINAC_ASSERT(is_a<matrix>(self.op(0)));
420 GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
422 const matrix &self_matrix = ex_to<matrix>(self.op(0));
424 if (self.nops() == 2)
425 return indexed(self_matrix.mul(other), self.op(1));
426 else // self.nops() == 3
427 return indexed(self_matrix.mul(other), self.op(1), self.op(2));
430 /** Contraction of an indexed matrix with something else. */
431 bool matrix::contract_with(exvector::iterator self, exvector::iterator other, exvector & v) const
433 GINAC_ASSERT(is_a<indexed>(*self));
434 GINAC_ASSERT(is_a<indexed>(*other));
435 GINAC_ASSERT(self->nops() == 2 || self->nops() == 3);
436 GINAC_ASSERT(is_a<matrix>(self->op(0)));
438 // Only contract with other matrices
439 if (!is_a<matrix>(other->op(0)))
442 GINAC_ASSERT(other->nops() == 2 || other->nops() == 3);
444 const matrix &self_matrix = ex_to<matrix>(self->op(0));
445 const matrix &other_matrix = ex_to<matrix>(other->op(0));
447 if (self->nops() == 2) {
449 if (other->nops() == 2) { // vector * vector (scalar product)
451 if (self_matrix.col == 1) {
452 if (other_matrix.col == 1) {
453 // Column vector * column vector, transpose first vector
454 *self = self_matrix.transpose().mul(other_matrix)(0, 0);
456 // Column vector * row vector, swap factors
457 *self = other_matrix.mul(self_matrix)(0, 0);
460 if (other_matrix.col == 1) {
461 // Row vector * column vector, perfect
462 *self = self_matrix.mul(other_matrix)(0, 0);
464 // Row vector * row vector, transpose second vector
465 *self = self_matrix.mul(other_matrix.transpose())(0, 0);
471 } else { // vector * matrix
473 // B_i * A_ij = (B*A)_j (B is row vector)
474 if (is_dummy_pair(self->op(1), other->op(1))) {
475 if (self_matrix.row == 1)
476 *self = indexed(self_matrix.mul(other_matrix), other->op(2));
478 *self = indexed(self_matrix.transpose().mul(other_matrix), other->op(2));
483 // B_j * A_ij = (A*B)_i (B is column vector)
484 if (is_dummy_pair(self->op(1), other->op(2))) {
485 if (self_matrix.col == 1)
486 *self = indexed(other_matrix.mul(self_matrix), other->op(1));
488 *self = indexed(other_matrix.mul(self_matrix.transpose()), other->op(1));
494 } else if (other->nops() == 3) { // matrix * matrix
496 // A_ij * B_jk = (A*B)_ik
497 if (is_dummy_pair(self->op(2), other->op(1))) {
498 *self = indexed(self_matrix.mul(other_matrix), self->op(1), other->op(2));
503 // A_ij * B_kj = (A*Btrans)_ik
504 if (is_dummy_pair(self->op(2), other->op(2))) {
505 *self = indexed(self_matrix.mul(other_matrix.transpose()), self->op(1), other->op(1));
510 // A_ji * B_jk = (Atrans*B)_ik
511 if (is_dummy_pair(self->op(1), other->op(1))) {
512 *self = indexed(self_matrix.transpose().mul(other_matrix), self->op(2), other->op(2));
517 // A_ji * B_kj = (B*A)_ki
518 if (is_dummy_pair(self->op(1), other->op(2))) {
519 *self = indexed(other_matrix.mul(self_matrix), other->op(1), self->op(2));
530 // non-virtual functions in this class
537 * @exception logic_error (incompatible matrices) */
538 matrix matrix::add(const matrix & other) const
540 if (col != other.col || row != other.row)
541 throw std::logic_error("matrix::add(): incompatible matrices");
543 exvector sum(this->m);
544 exvector::iterator i = sum.begin(), end = sum.end();
545 exvector::const_iterator ci = other.m.begin();
549 return matrix(row,col,sum);
553 /** Difference of matrices.
555 * @exception logic_error (incompatible matrices) */
556 matrix matrix::sub(const matrix & other) const
558 if (col != other.col || row != other.row)
559 throw std::logic_error("matrix::sub(): incompatible matrices");
561 exvector dif(this->m);
562 exvector::iterator i = dif.begin(), end = dif.end();
563 exvector::const_iterator ci = other.m.begin();
567 return matrix(row,col,dif);
571 /** Product of matrices.
573 * @exception logic_error (incompatible matrices) */
574 matrix matrix::mul(const matrix & other) const
576 if (this->cols() != other.rows())
577 throw std::logic_error("matrix::mul(): incompatible matrices");
579 exvector prod(this->rows()*other.cols());
581 for (unsigned r1=0; r1<this->rows(); ++r1) {
582 for (unsigned c=0; c<this->cols(); ++c) {
583 // Quick test: can we shortcut?
584 if (m[r1*col+c].is_zero())
586 for (unsigned r2=0; r2<other.cols(); ++r2)
587 prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]);
590 return matrix(row, other.col, prod);
594 /** Product of matrix and scalar. */
595 matrix matrix::mul(const numeric & other) const
597 exvector prod(row * col);
599 for (unsigned r=0; r<row; ++r)
600 for (unsigned c=0; c<col; ++c)
601 prod[r*col+c] = m[r*col+c] * other;
603 return matrix(row, col, prod);
607 /** Product of matrix and scalar expression. */
608 matrix matrix::mul_scalar(const ex & other) const
610 if (other.return_type() != return_types::commutative)
611 throw std::runtime_error("matrix::mul_scalar(): non-commutative scalar");
613 exvector prod(row * col);
615 for (unsigned r=0; r<row; ++r)
616 for (unsigned c=0; c<col; ++c)
617 prod[r*col+c] = m[r*col+c] * other;
619 return matrix(row, col, prod);
623 /** Power of a matrix. Currently handles integer exponents only. */
624 matrix matrix::pow(const ex & expn) const
627 throw (std::logic_error("matrix::pow(): matrix not square"));
629 if (is_exactly_a<numeric>(expn)) {
630 // Integer cases are computed by successive multiplication, using the
631 // obvious shortcut of storing temporaries, like A^4 == (A*A)*(A*A).
632 if (expn.info(info_flags::integer)) {
633 numeric b = ex_to<numeric>(expn);
635 if (expn.info(info_flags::negative)) {
642 for (unsigned r=0; r<row; ++r)
646 // This loop computes the representation of b in base 2 from right
647 // to left and multiplies the factors whenever needed. Note
648 // that this is not entirely optimal but close to optimal and
649 // "better" algorithms are much harder to implement. (See Knuth,
650 // TAoCP2, section "Evaluation of Powers" for a good discussion.)
651 while (b!=*_num1_p) {
656 b /= *_num2_p; // still integer.
662 throw (std::runtime_error("matrix::pow(): don't know how to handle exponent"));
666 /** operator() to access elements for reading.
668 * @param ro row of element
669 * @param co column of element
670 * @exception range_error (index out of range) */
671 const ex & matrix::operator() (unsigned ro, unsigned co) const
673 if (ro>=row || co>=col)
674 throw (std::range_error("matrix::operator(): index out of range"));
680 /** operator() to access elements for writing.
682 * @param ro row of element
683 * @param co column of element
684 * @exception range_error (index out of range) */
685 ex & matrix::operator() (unsigned ro, unsigned co)
687 if (ro>=row || co>=col)
688 throw (std::range_error("matrix::operator(): index out of range"));
690 ensure_if_modifiable();
695 /** Transposed of an m x n matrix, producing a new n x m matrix object that
696 * represents the transposed. */
697 matrix matrix::transpose() const
699 exvector trans(this->cols()*this->rows());
701 for (unsigned r=0; r<this->cols(); ++r)
702 for (unsigned c=0; c<this->rows(); ++c)
703 trans[r*this->rows()+c] = m[c*this->cols()+r];
705 return matrix(this->cols(),this->rows(),trans);
708 /** Determinant of square matrix. This routine doesn't actually calculate the
709 * determinant, it only implements some heuristics about which algorithm to
710 * run. If all the elements of the matrix are elements of an integral domain
711 * the determinant is also in that integral domain and the result is expanded
712 * only. If one or more elements are from a quotient field the determinant is
713 * usually also in that quotient field and the result is normalized before it
714 * is returned. This implies that the determinant of the symbolic 2x2 matrix
715 * [[a/(a-b),1],[b/(a-b),1]] is returned as unity. (In this respect, it
716 * behaves like MapleV and unlike Mathematica.)
718 * @param algo allows to chose an algorithm
719 * @return the determinant as a new expression
720 * @exception logic_error (matrix not square)
721 * @see determinant_algo */
722 ex matrix::determinant(unsigned algo) const
725 throw (std::logic_error("matrix::determinant(): matrix not square"));
726 GINAC_ASSERT(row*col==m.capacity());
728 // Gather some statistical information about this matrix:
729 bool numeric_flag = true;
730 bool normal_flag = false;
731 unsigned sparse_count = 0; // counts non-zero elements
732 exvector::const_iterator r = m.begin(), rend = m.end();
734 if (!r->info(info_flags::numeric))
735 numeric_flag = false;
736 exmap srl; // symbol replacement list
737 ex rtest = r->to_rational(srl);
738 if (!rtest.is_zero())
740 if (!rtest.info(info_flags::crational_polynomial) &&
741 rtest.info(info_flags::rational_function))
746 // Here is the heuristics in case this routine has to decide:
747 if (algo == determinant_algo::automatic) {
748 // Minor expansion is generally a good guess:
749 algo = determinant_algo::laplace;
750 // Does anybody know when a matrix is really sparse?
751 // Maybe <~row/2.236 nonzero elements average in a row?
752 if (row>3 && 5*sparse_count<=row*col)
753 algo = determinant_algo::bareiss;
754 // Purely numeric matrix can be handled by Gauss elimination.
755 // This overrides any prior decisions.
757 algo = determinant_algo::gauss;
760 // Trap the trivial case here, since some algorithms don't like it
762 // for consistency with non-trivial determinants...
764 return m[0].normal();
766 return m[0].expand();
769 // Compute the determinant
771 case determinant_algo::gauss: {
774 int sign = tmp.gauss_elimination(true);
775 for (unsigned d=0; d<row; ++d)
776 det *= tmp.m[d*col+d];
778 return (sign*det).normal();
780 return (sign*det).normal().expand();
782 case determinant_algo::bareiss: {
785 sign = tmp.fraction_free_elimination(true);
787 return (sign*tmp.m[row*col-1]).normal();
789 return (sign*tmp.m[row*col-1]).expand();
791 case determinant_algo::divfree: {
794 sign = tmp.division_free_elimination(true);
797 ex det = tmp.m[row*col-1];
798 // factor out accumulated bogus slag
799 for (unsigned d=0; d<row-2; ++d)
800 for (unsigned j=0; j<row-d-2; ++j)
801 det = (det/tmp.m[d*col+d]).normal();
804 case determinant_algo::laplace:
806 // This is the minor expansion scheme. We always develop such
807 // that the smallest minors (i.e, the trivial 1x1 ones) are on the
808 // rightmost column. For this to be efficient, empirical tests
809 // have shown that the emptiest columns (i.e. the ones with most
810 // zeros) should be the ones on the right hand side -- although
811 // this might seem counter-intuitive (and in contradiction to some
812 // literature like the FORM manual). Please go ahead and test it
813 // if you don't believe me! Therefore we presort the columns of
815 typedef std::pair<unsigned,unsigned> uintpair;
816 std::vector<uintpair> c_zeros; // number of zeros in column
817 for (unsigned c=0; c<col; ++c) {
819 for (unsigned r=0; r<row; ++r)
820 if (m[r*col+c].is_zero())
822 c_zeros.push_back(uintpair(acc,c));
824 std::sort(c_zeros.begin(),c_zeros.end());
825 std::vector<unsigned> pre_sort;
826 for (std::vector<uintpair>::const_iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
827 pre_sort.push_back(i->second);
828 std::vector<unsigned> pre_sort_test(pre_sort); // permutation_sign() modifies the vector so we make a copy here
829 int sign = permutation_sign(pre_sort_test.begin(), pre_sort_test.end());
830 exvector result(row*col); // represents sorted matrix
832 for (std::vector<unsigned>::const_iterator i=pre_sort.begin();
835 for (unsigned r=0; r<row; ++r)
836 result[r*col+c] = m[r*col+(*i)];
840 return (sign*matrix(row,col,result).determinant_minor()).normal();
842 return sign*matrix(row,col,result).determinant_minor();
848 /** Trace of a matrix. The result is normalized if it is in some quotient
849 * field and expanded only otherwise. This implies that the trace of the
850 * symbolic 2x2 matrix [[a/(a-b),x],[y,b/(b-a)]] is recognized to be unity.
852 * @return the sum of diagonal elements
853 * @exception logic_error (matrix not square) */
854 ex matrix::trace() const
857 throw (std::logic_error("matrix::trace(): matrix not square"));
860 for (unsigned r=0; r<col; ++r)
863 if (tr.info(info_flags::rational_function) &&
864 !tr.info(info_flags::crational_polynomial))
871 /** Characteristic Polynomial. Following mathematica notation the
872 * characteristic polynomial of a matrix M is defined as the determiant of
873 * (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
874 * as M. Note that some CASs define it with a sign inside the determinant
875 * which gives rise to an overall sign if the dimension is odd. This method
876 * returns the characteristic polynomial collected in powers of lambda as a
879 * @return characteristic polynomial as new expression
880 * @exception logic_error (matrix not square)
881 * @see matrix::determinant() */
882 ex matrix::charpoly(const ex & lambda) const
885 throw (std::logic_error("matrix::charpoly(): matrix not square"));
887 bool numeric_flag = true;
888 exvector::const_iterator r = m.begin(), rend = m.end();
889 while (r!=rend && numeric_flag==true) {
890 if (!r->info(info_flags::numeric))
891 numeric_flag = false;
895 // The pure numeric case is traditionally rather common. Hence, it is
896 // trapped and we use Leverrier's algorithm which goes as row^3 for
897 // every coefficient. The expensive part is the matrix multiplication.
902 ex poly = power(lambda, row) - c*power(lambda, row-1);
903 for (unsigned i=1; i<row; ++i) {
904 for (unsigned j=0; j<row; ++j)
907 c = B.trace() / ex(i+1);
908 poly -= c*power(lambda, row-i-1);
918 for (unsigned r=0; r<col; ++r)
919 M.m[r*col+r] -= lambda;
921 return M.determinant().collect(lambda);
926 /** Inverse of this matrix.
928 * @return the inverted matrix
929 * @exception logic_error (matrix not square)
930 * @exception runtime_error (singular matrix) */
931 matrix matrix::inverse() const
934 throw (std::logic_error("matrix::inverse(): matrix not square"));
936 // This routine actually doesn't do anything fancy at all. We compute the
937 // inverse of the matrix A by solving the system A * A^{-1} == Id.
939 // First populate the identity matrix supposed to become the right hand side.
940 matrix identity(row,col);
941 for (unsigned i=0; i<row; ++i)
942 identity(i,i) = _ex1;
944 // Populate a dummy matrix of variables, just because of compatibility with
945 // matrix::solve() which wants this (for compatibility with under-determined
946 // systems of equations).
947 matrix vars(row,col);
948 for (unsigned r=0; r<row; ++r)
949 for (unsigned c=0; c<col; ++c)
950 vars(r,c) = symbol();
954 sol = this->solve(vars,identity);
955 } catch (const std::runtime_error & e) {
956 if (e.what()==std::string("matrix::solve(): inconsistent linear system"))
957 throw (std::runtime_error("matrix::inverse(): singular matrix"));
965 /** Solve a linear system consisting of a m x n matrix and a m x p right hand
966 * side by applying an elimination scheme to the augmented matrix.
968 * @param vars n x p matrix, all elements must be symbols
969 * @param rhs m x p matrix
970 * @param algo selects the solving algorithm
971 * @return n x p solution matrix
972 * @exception logic_error (incompatible matrices)
973 * @exception invalid_argument (1st argument must be matrix of symbols)
974 * @exception runtime_error (inconsistent linear system)
976 matrix matrix::solve(const matrix & vars,
980 const unsigned m = this->rows();
981 const unsigned n = this->cols();
982 const unsigned p = rhs.cols();
985 if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
986 throw (std::logic_error("matrix::solve(): incompatible matrices"));
987 for (unsigned ro=0; ro<n; ++ro)
988 for (unsigned co=0; co<p; ++co)
989 if (!vars(ro,co).info(info_flags::symbol))
990 throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
992 // build the augmented matrix of *this with rhs attached to the right
994 for (unsigned r=0; r<m; ++r) {
995 for (unsigned c=0; c<n; ++c)
996 aug.m[r*(n+p)+c] = this->m[r*n+c];
997 for (unsigned c=0; c<p; ++c)
998 aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
1001 // Gather some statistical information about the augmented matrix:
1002 bool numeric_flag = true;
1003 exvector::const_iterator r = aug.m.begin(), rend = aug.m.end();
1004 while (r!=rend && numeric_flag==true) {
1005 if (!r->info(info_flags::numeric))
1006 numeric_flag = false;
1010 // Here is the heuristics in case this routine has to decide:
1011 if (algo == solve_algo::automatic) {
1012 // Bareiss (fraction-free) elimination is generally a good guess:
1013 algo = solve_algo::bareiss;
1014 // For m<3, Bareiss elimination is equivalent to division free
1015 // elimination but has more logistic overhead
1017 algo = solve_algo::divfree;
1018 // This overrides any prior decisions.
1020 algo = solve_algo::gauss;
1023 // Eliminate the augmented matrix:
1025 case solve_algo::gauss:
1026 aug.gauss_elimination();
1028 case solve_algo::divfree:
1029 aug.division_free_elimination();
1031 case solve_algo::bareiss:
1033 aug.fraction_free_elimination();
1036 // assemble the solution matrix:
1038 for (unsigned co=0; co<p; ++co) {
1039 unsigned last_assigned_sol = n+1;
1040 for (int r=m-1; r>=0; --r) {
1041 unsigned fnz = 1; // first non-zero in row
1042 while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
1045 // row consists only of zeros, corresponding rhs must be 0, too
1046 if (!aug.m[r*(n+p)+n+co].is_zero()) {
1047 throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
1050 // assign solutions for vars between fnz+1 and
1051 // last_assigned_sol-1: free parameters
1052 for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
1053 sol(c,co) = vars.m[c*p+co];
1054 ex e = aug.m[r*(n+p)+n+co];
1055 for (unsigned c=fnz; c<n; ++c)
1056 e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
1057 sol(fnz-1,co) = (e/(aug.m[r*(n+p)+(fnz-1)])).normal();
1058 last_assigned_sol = fnz;
1061 // assign solutions for vars between 1 and
1062 // last_assigned_sol-1: free parameters
1063 for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
1064 sol(ro,co) = vars(ro,co);
1071 /** Compute the rank of this matrix. */
1072 unsigned matrix::rank() const
1075 // Transform this matrix into upper echelon form and then count the
1076 // number of non-zero rows.
1078 GINAC_ASSERT(row*col==m.capacity());
1080 // Actually, any elimination scheme will do since we are only
1081 // interested in the echelon matrix' zeros.
1082 matrix to_eliminate = *this;
1083 to_eliminate.fraction_free_elimination();
1085 unsigned r = row*col; // index of last non-zero element
1087 if (!to_eliminate.m[r].is_zero())
1096 /** Recursive determinant for small matrices having at least one symbolic
1097 * entry. The basic algorithm, known as Laplace-expansion, is enhanced by
1098 * some bookkeeping to avoid calculation of the same submatrices ("minors")
1099 * more than once. According to W.M.Gentleman and S.C.Johnson this algorithm
1100 * is better than elimination schemes for matrices of sparse multivariate
1101 * polynomials and also for matrices of dense univariate polynomials if the
1102 * matrix' dimesion is larger than 7.
1104 * @return the determinant as a new expression (in expanded form)
1105 * @see matrix::determinant() */
1106 ex matrix::determinant_minor() const
1108 // for small matrices the algorithm does not make any sense:
1109 const unsigned n = this->cols();
1111 return m[0].expand();
1113 return (m[0]*m[3]-m[2]*m[1]).expand();
1115 return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
1116 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
1117 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
1119 // This algorithm can best be understood by looking at a naive
1120 // implementation of Laplace-expansion, like this one:
1122 // matrix minorM(this->rows()-1,this->cols()-1);
1123 // for (unsigned r1=0; r1<this->rows(); ++r1) {
1124 // // shortcut if element(r1,0) vanishes
1125 // if (m[r1*col].is_zero())
1127 // // assemble the minor matrix
1128 // for (unsigned r=0; r<minorM.rows(); ++r) {
1129 // for (unsigned c=0; c<minorM.cols(); ++c) {
1131 // minorM(r,c) = m[r*col+c+1];
1133 // minorM(r,c) = m[(r+1)*col+c+1];
1136 // // recurse down and care for sign:
1138 // det -= m[r1*col] * minorM.determinant_minor();
1140 // det += m[r1*col] * minorM.determinant_minor();
1142 // return det.expand();
1143 // What happens is that while proceeding down many of the minors are
1144 // computed more than once. In particular, there are binomial(n,k)
1145 // kxk minors and each one is computed factorial(n-k) times. Therefore
1146 // it is reasonable to store the results of the minors. We proceed from
1147 // right to left. At each column c we only need to retrieve the minors
1148 // calculated in step c-1. We therefore only have to store at most
1149 // 2*binomial(n,n/2) minors.
1151 // Unique flipper counter for partitioning into minors
1152 std::vector<unsigned> Pkey;
1154 // key for minor determinant (a subpartition of Pkey)
1155 std::vector<unsigned> Mkey;
1157 // we store our subminors in maps, keys being the rows they arise from
1158 typedef std::map<std::vector<unsigned>,class ex> Rmap;
1159 typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
1163 // initialize A with last column:
1164 for (unsigned r=0; r<n; ++r) {
1165 Pkey.erase(Pkey.begin(),Pkey.end());
1167 A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
1169 // proceed from right to left through matrix
1170 for (int c=n-2; c>=0; --c) {
1171 Pkey.erase(Pkey.begin(),Pkey.end()); // don't change capacity
1172 Mkey.erase(Mkey.begin(),Mkey.end());
1173 for (unsigned i=0; i<n-c; ++i)
1175 unsigned fc = 0; // controls logic for our strange flipper counter
1178 for (unsigned r=0; r<n-c; ++r) {
1179 // maybe there is nothing to do?
1180 if (m[Pkey[r]*n+c].is_zero())
1182 // create the sorted key for all possible minors
1183 Mkey.erase(Mkey.begin(),Mkey.end());
1184 for (unsigned i=0; i<n-c; ++i)
1186 Mkey.push_back(Pkey[i]);
1187 // Fetch the minors and compute the new determinant
1189 det -= m[Pkey[r]*n+c]*A[Mkey];
1191 det += m[Pkey[r]*n+c]*A[Mkey];
1193 // prevent build-up of deep nesting of expressions saves time:
1195 // store the new determinant at its place in B:
1197 B.insert(Rmap_value(Pkey,det));
1198 // increment our strange flipper counter
1199 for (fc=n-c; fc>0; --fc) {
1201 if (Pkey[fc-1]<fc+c)
1205 for (unsigned j=fc; j<n-c; ++j)
1206 Pkey[j] = Pkey[j-1]+1;
1208 // next column, so change the role of A and B:
1217 /** Perform the steps of an ordinary Gaussian elimination to bring the m x n
1218 * matrix into an upper echelon form. The algorithm is ok for matrices
1219 * with numeric coefficients but quite unsuited for symbolic matrices.
1221 * @param det may be set to true to save a lot of space if one is only
1222 * interested in the diagonal elements (i.e. for calculating determinants).
1223 * The others are set to zero in this case.
1224 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1225 * number of rows was swapped and 0 if the matrix is singular. */
1226 int matrix::gauss_elimination(const bool det)
1228 ensure_if_modifiable();
1229 const unsigned m = this->rows();
1230 const unsigned n = this->cols();
1231 GINAC_ASSERT(!det || n==m);
1235 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1236 int indx = pivot(r0, c0, true);
1240 return 0; // leaves *this in a messy state
1245 for (unsigned r2=r0+1; r2<m; ++r2) {
1246 if (!this->m[r2*n+c0].is_zero()) {
1247 // yes, there is something to do in this row
1248 ex piv = this->m[r2*n+c0] / this->m[r0*n+c0];
1249 for (unsigned c=c0+1; c<n; ++c) {
1250 this->m[r2*n+c] -= piv * this->m[r0*n+c];
1251 if (!this->m[r2*n+c].info(info_flags::numeric))
1252 this->m[r2*n+c] = this->m[r2*n+c].normal();
1255 // fill up left hand side with zeros
1256 for (unsigned c=r0; c<=c0; ++c)
1257 this->m[r2*n+c] = _ex0;
1260 // save space by deleting no longer needed elements
1261 for (unsigned c=r0+1; c<n; ++c)
1262 this->m[r0*n+c] = _ex0;
1267 // clear remaining rows
1268 for (unsigned r=r0+1; r<m; ++r) {
1269 for (unsigned c=0; c<n; ++c)
1270 this->m[r*n+c] = _ex0;
1277 /** Perform the steps of division free elimination to bring the m x n matrix
1278 * into an upper echelon form.
1280 * @param det may be set to true to save a lot of space if one is only
1281 * interested in the diagonal elements (i.e. for calculating determinants).
1282 * The others are set to zero in this case.
1283 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1284 * number of rows was swapped and 0 if the matrix is singular. */
1285 int matrix::division_free_elimination(const bool det)
1287 ensure_if_modifiable();
1288 const unsigned m = this->rows();
1289 const unsigned n = this->cols();
1290 GINAC_ASSERT(!det || n==m);
1294 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1295 int indx = pivot(r0, c0, true);
1299 return 0; // leaves *this in a messy state
1304 for (unsigned r2=r0+1; r2<m; ++r2) {
1305 for (unsigned c=c0+1; c<n; ++c)
1306 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();
1307 // fill up left hand side with zeros
1308 for (unsigned c=r0; c<=c0; ++c)
1309 this->m[r2*n+c] = _ex0;
1312 // save space by deleting no longer needed elements
1313 for (unsigned c=r0+1; c<n; ++c)
1314 this->m[r0*n+c] = _ex0;
1319 // clear remaining rows
1320 for (unsigned r=r0+1; r<m; ++r) {
1321 for (unsigned c=0; c<n; ++c)
1322 this->m[r*n+c] = _ex0;
1329 /** Perform the steps of Bareiss' one-step fraction free elimination to bring
1330 * the matrix into an upper echelon form. Fraction free elimination means
1331 * that divide is used straightforwardly, without computing GCDs first. This
1332 * is possible, since we know the divisor at each step.
1334 * @param det may be set to true to save a lot of space if one is only
1335 * interested in the last element (i.e. for calculating determinants). The
1336 * others are set to zero in this case.
1337 * @return sign is 1 if an even number of rows was swapped, -1 if an odd
1338 * number of rows was swapped and 0 if the matrix is singular. */
1339 int matrix::fraction_free_elimination(const bool det)
1342 // (single-step fraction free elimination scheme, already known to Jordan)
1344 // Usual division-free elimination sets m[0](r,c) = m(r,c) and then sets
1345 // m[k+1](r,c) = m[k](k,k) * m[k](r,c) - m[k](r,k) * m[k](k,c).
1347 // Bareiss (fraction-free) elimination in addition divides that element
1348 // by m[k-1](k-1,k-1) for k>1, where it can be shown by means of the
1349 // Sylvester identity that this really divides m[k+1](r,c).
1351 // We also allow rational functions where the original prove still holds.
1352 // However, we must care for numerator and denominator separately and
1353 // "manually" work in the integral domains because of subtle cancellations
1354 // (see below). This blows up the bookkeeping a bit and the formula has
1355 // to be modified to expand like this (N{x} stands for numerator of x,
1356 // D{x} for denominator of x):
1357 // 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)}
1358 // -N{m[k](r,k)}*N{m[k](k,c)}*D{m[k](k,k)}*D{m[k](r,c)}
1359 // 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)}
1360 // where for k>1 we now divide N{m[k+1](r,c)} by
1361 // N{m[k-1](k-1,k-1)}
1362 // and D{m[k+1](r,c)} by
1363 // D{m[k-1](k-1,k-1)}.
1365 ensure_if_modifiable();
1366 const unsigned m = this->rows();
1367 const unsigned n = this->cols();
1368 GINAC_ASSERT(!det || n==m);
1377 // We populate temporary matrices to subsequently operate on. There is
1378 // one holding numerators and another holding denominators of entries.
1379 // This is a must since the evaluator (or even earlier mul's constructor)
1380 // might cancel some trivial element which causes divide() to fail. The
1381 // elements are normalized first (yes, even though this algorithm doesn't
1382 // need GCDs) since the elements of *this might be unnormalized, which
1383 // makes things more complicated than they need to be.
1384 matrix tmp_n(*this);
1385 matrix tmp_d(m,n); // for denominators, if needed
1386 exmap srl; // symbol replacement list
1387 exvector::const_iterator cit = this->m.begin(), citend = this->m.end();
1388 exvector::iterator tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
1389 while (cit != citend) {
1390 ex nd = cit->normal().to_rational(srl).numer_denom();
1392 *tmp_n_it++ = nd.op(0);
1393 *tmp_d_it++ = nd.op(1);
1397 for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
1398 int indx = tmp_n.pivot(r0, c0, true);
1407 // tmp_n's rows r0 and indx were swapped, do the same in tmp_d:
1408 for (unsigned c=c0; c<n; ++c)
1409 tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
1411 for (unsigned r2=r0+1; r2<m; ++r2) {
1412 for (unsigned c=c0+1; c<n; ++c) {
1413 dividend_n = (tmp_n.m[r0*n+c0]*tmp_n.m[r2*n+c]*
1414 tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]
1415 -tmp_n.m[r2*n+c0]*tmp_n.m[r0*n+c]*
1416 tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
1417 dividend_d = (tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]*
1418 tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
1419 bool check = divide(dividend_n, divisor_n,
1420 tmp_n.m[r2*n+c], true);
1421 check &= divide(dividend_d, divisor_d,
1422 tmp_d.m[r2*n+c], true);
1423 GINAC_ASSERT(check);
1425 // fill up left hand side with zeros
1426 for (unsigned c=r0; c<=c0; ++c)
1427 tmp_n.m[r2*n+c] = _ex0;
1429 if (c0<n && r0<m-1) {
1430 // compute next iteration's divisor
1431 divisor_n = tmp_n.m[r0*n+c0].expand();
1432 divisor_d = tmp_d.m[r0*n+c0].expand();
1434 // save space by deleting no longer needed elements
1435 for (unsigned c=0; c<n; ++c) {
1436 tmp_n.m[r0*n+c] = _ex0;
1437 tmp_d.m[r0*n+c] = _ex1;
1444 // clear remaining rows
1445 for (unsigned r=r0+1; r<m; ++r) {
1446 for (unsigned c=0; c<n; ++c)
1447 tmp_n.m[r*n+c] = _ex0;
1450 // repopulate *this matrix:
1451 exvector::iterator it = this->m.begin(), itend = this->m.end();
1452 tmp_n_it = tmp_n.m.begin();
1453 tmp_d_it = tmp_d.m.begin();
1455 *it++ = ((*tmp_n_it++)/(*tmp_d_it++)).subs(srl, subs_options::no_pattern);
1461 /** Partial pivoting method for matrix elimination schemes.
1462 * Usual pivoting (symbolic==false) returns the index to the element with the
1463 * largest absolute value in column ro and swaps the current row with the one
1464 * where the element was found. With (symbolic==true) it does the same thing
1465 * with the first non-zero element.
1467 * @param ro is the row from where to begin
1468 * @param co is the column to be inspected
1469 * @param symbolic signal if we want the first non-zero element to be pivoted
1470 * (true) or the one with the largest absolute value (false).
1471 * @return 0 if no interchange occured, -1 if all are zero (usually signaling
1472 * a degeneracy) and positive integer k means that rows ro and k were swapped.
1474 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
1478 // search first non-zero element in column co beginning at row ro
1479 while ((k<row) && (this->m[k*col+co].expand().is_zero()))
1482 // search largest element in column co beginning at row ro
1483 GINAC_ASSERT(is_exactly_a<numeric>(this->m[k*col+co]));
1484 unsigned kmax = k+1;
1485 numeric mmax = abs(ex_to<numeric>(m[kmax*col+co]));
1487 GINAC_ASSERT(is_exactly_a<numeric>(this->m[kmax*col+co]));
1488 numeric tmp = ex_to<numeric>(this->m[kmax*col+co]);
1489 if (abs(tmp) > mmax) {
1495 if (!mmax.is_zero())
1499 // all elements in column co below row ro vanish
1502 // matrix needs no pivoting
1504 // matrix needs pivoting, so swap rows k and ro
1505 ensure_if_modifiable();
1506 for (unsigned c=0; c<col; ++c)
1507 this->m[k*col+c].swap(this->m[ro*col+c]);
1512 ex lst_to_matrix(const lst & l)
1514 lst::const_iterator itr, itc;
1516 // Find number of rows and columns
1517 size_t rows = l.nops(), cols = 0;
1518 for (itr = l.begin(); itr != l.end(); ++itr) {
1519 if (!is_a<lst>(*itr))
1520 throw (std::invalid_argument("lst_to_matrix: argument must be a list of lists"));
1521 if (itr->nops() > cols)
1525 // Allocate and fill matrix
1526 matrix &M = *new matrix(rows, cols);
1527 M.setflag(status_flags::dynallocated);
1530 for (itr = l.begin(), i = 0; itr != l.end(); ++itr, ++i) {
1532 for (itc = ex_to<lst>(*itr).begin(), j = 0; itc != ex_to<lst>(*itr).end(); ++itc, ++j)
1539 ex diag_matrix(const lst & l)
1541 lst::const_iterator it;
1542 size_t dim = l.nops();
1544 // Allocate and fill matrix
1545 matrix &M = *new matrix(dim, dim);
1546 M.setflag(status_flags::dynallocated);
1549 for (it = l.begin(), i = 0; it != l.end(); ++it, ++i)
1555 ex unit_matrix(unsigned r, unsigned c)
1557 matrix &Id = *new matrix(r, c);
1558 Id.setflag(status_flags::dynallocated);
1559 for (unsigned i=0; i<r && i<c; i++)
1565 ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const std::string & tex_base_name)
1567 matrix &M = *new matrix(r, c);
1568 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1570 bool long_format = (r > 10 || c > 10);
1571 bool single_row = (r == 1 || c == 1);
1573 for (unsigned i=0; i<r; i++) {
1574 for (unsigned j=0; j<c; j++) {
1575 std::ostringstream s1, s2;
1577 s2 << tex_base_name << "_{";
1588 s1 << '_' << i << '_' << j;
1589 s2 << i << ';' << j << "}";
1592 s2 << i << j << '}';
1595 M(i, j) = symbol(s1.str(), s2.str());
1602 ex reduced_matrix(const matrix& m, unsigned r, unsigned c)
1604 if (r+1>m.rows() || c+1>m.cols() || m.cols()<2 || m.rows()<2)
1605 throw std::runtime_error("minor_matrix(): index out of bounds");
1607 const unsigned rows = m.rows()-1;
1608 const unsigned cols = m.cols()-1;
1609 matrix &M = *new matrix(rows, cols);
1610 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1622 M(ro2,co2) = m(ro, co);
1633 ex sub_matrix(const matrix&m, unsigned r, unsigned nr, unsigned c, unsigned nc)
1635 if (r+nr>m.rows() || c+nc>m.cols())
1636 throw std::runtime_error("sub_matrix(): index out of bounds");
1638 matrix &M = *new matrix(nr, nc);
1639 M.setflag(status_flags::dynallocated | status_flags::evaluated);
1641 for (unsigned ro=0; ro<nr; ++ro) {
1642 for (unsigned co=0; co<nc; ++co) {
1643 M(ro,co) = m(ro+r,co+c);
1650 } // namespace GiNaC