GiNaC  1.6.2
matrix.cpp
Go to the documentation of this file.
00001 
00005 /*
00006  *  GiNaC Copyright (C) 1999-2011 Johannes Gutenberg University Mainz, Germany
00007  *
00008  *  This program is free software; you can redistribute it and/or modify
00009  *  it under the terms of the GNU General Public License as published by
00010  *  the Free Software Foundation; either version 2 of the License, or
00011  *  (at your option) any later version.
00012  *
00013  *  This program is distributed in the hope that it will be useful,
00014  *  but WITHOUT ANY WARRANTY; without even the implied warranty of
00015  *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016  *  GNU General Public License for more details.
00017  *
00018  *  You should have received a copy of the GNU General Public License
00019  *  along with this program; if not, write to the Free Software
00020  *  Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
00021  */
00022 
00023 #include "matrix.h"
00024 #include "numeric.h"
00025 #include "lst.h"
00026 #include "idx.h"
00027 #include "indexed.h"
00028 #include "add.h"
00029 #include "power.h"
00030 #include "symbol.h"
00031 #include "operators.h"
00032 #include "normal.h"
00033 #include "archive.h"
00034 #include "utils.h"
00035 
00036 #include <algorithm>
00037 #include <iostream>
00038 #include <map>
00039 #include <sstream>
00040 #include <stdexcept>
00041 #include <string>
00042 
00043 namespace GiNaC {
00044 
00045 GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(matrix, basic,
00046   print_func<print_context>(&matrix::do_print).
00047   print_func<print_latex>(&matrix::do_print_latex).
00048   print_func<print_tree>(&matrix::do_print_tree).
00049   print_func<print_python_repr>(&matrix::do_print_python_repr))
00050 
00051 
00052 // default constructor
00054 
00056 matrix::matrix() : row(1), col(1), m(1, _ex0)
00057 {
00058     setflag(status_flags::not_shareable);
00059 }
00060 
00062 // other constructors
00064 
00065 // public
00066 
00071 matrix::matrix(unsigned r, unsigned c) : row(r), col(c), m(r*c, _ex0)
00072 {
00073     setflag(status_flags::not_shareable);
00074 }
00075 
00076 // protected
00077 
00079 matrix::matrix(unsigned r, unsigned c, const exvector & m2)
00080   : row(r), col(c), m(m2)
00081 {
00082     setflag(status_flags::not_shareable);
00083 }
00084 
00089 matrix::matrix(unsigned r, unsigned c, const lst & l)
00090   : row(r), col(c), m(r*c, _ex0)
00091 {
00092     setflag(status_flags::not_shareable);
00093 
00094     size_t i = 0;
00095     for (lst::const_iterator it = l.begin(); it != l.end(); ++it, ++i) {
00096         size_t x = i % c;
00097         size_t y = i / c;
00098         if (y >= r)
00099             break; // matrix smaller than list: throw away excessive elements
00100         m[y*c+x] = *it;
00101     }
00102 }
00103 
00105 // archiving
00107 
00108 void matrix::read_archive(const archive_node &n, lst &sym_lst)
00109 {
00110     inherited::read_archive(n, sym_lst);
00111 
00112     if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
00113         throw (std::runtime_error("unknown matrix dimensions in archive"));
00114     m.reserve(row * col);
00115     // XXX: default ctor inserts a zero element, we need to erase it here.
00116     m.pop_back();
00117     archive_node::archive_node_cit first = n.find_first("m");
00118     archive_node::archive_node_cit last = n.find_last("m");
00119     ++last;
00120     for (archive_node::archive_node_cit i=first; i != last; ++i) {
00121         ex e;
00122         n.find_ex_by_loc(i, e, sym_lst);
00123         m.push_back(e);
00124     }
00125 }
00126 GINAC_BIND_UNARCHIVER(matrix);
00127 
00128 void matrix::archive(archive_node &n) const
00129 {
00130     inherited::archive(n);
00131     n.add_unsigned("row", row);
00132     n.add_unsigned("col", col);
00133     exvector::const_iterator i = m.begin(), iend = m.end();
00134     while (i != iend) {
00135         n.add_ex("m", *i);
00136         ++i;
00137     }
00138 }
00139 
00141 // functions overriding virtual functions from base classes
00143 
00144 // public
00145 
00146 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
00147 {
00148     for (unsigned ro=0; ro<row; ++ro) {
00149         c.s << row_start;
00150         for (unsigned co=0; co<col; ++co) {
00151             m[ro*col+co].print(c);
00152             if (co < col-1)
00153                 c.s << col_sep;
00154             else
00155                 c.s << row_end;
00156         }
00157         if (ro < row-1)
00158             c.s << row_sep;
00159     }
00160 }
00161 
00162 void matrix::do_print(const print_context & c, unsigned level) const
00163 {
00164     c.s << "[";
00165     print_elements(c, "[", "]", ",", ",");
00166     c.s << "]";
00167 }
00168 
00169 void matrix::do_print_latex(const print_latex & c, unsigned level) const
00170 {
00171     c.s << "\\left(\\begin{array}{" << std::string(col,'c') << "}";
00172     print_elements(c, "", "", "\\\\", "&");
00173     c.s << "\\end{array}\\right)";
00174 }
00175 
00176 void matrix::do_print_python_repr(const print_python_repr & c, unsigned level) const
00177 {
00178     c.s << class_name() << '(';
00179     print_elements(c, "[", "]", ",", ",");
00180     c.s << ')';
00181 }
00182 
00184 size_t matrix::nops() const
00185 {
00186     return static_cast<size_t>(row) * static_cast<size_t>(col);
00187 }
00188 
00190 ex matrix::op(size_t i) const
00191 {
00192     GINAC_ASSERT(i<nops());
00193     
00194     return m[i];
00195 }
00196 
00198 ex & matrix::let_op(size_t i)
00199 {
00200     GINAC_ASSERT(i<nops());
00201     
00202     ensure_if_modifiable();
00203     return m[i];
00204 }
00205 
00207 ex matrix::eval(int level) const
00208 {
00209     // check if we have to do anything at all
00210     if ((level==1)&&(flags & status_flags::evaluated))
00211         return *this;
00212     
00213     // emergency break
00214     if (level == -max_recursion_level)
00215         throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
00216     
00217     // eval() entry by entry
00218     exvector m2(row*col);
00219     --level;
00220     for (unsigned r=0; r<row; ++r)
00221         for (unsigned c=0; c<col; ++c)
00222             m2[r*col+c] = m[r*col+c].eval(level);
00223     
00224     return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
00225                                                status_flags::evaluated);
00226 }
00227 
00228 ex matrix::subs(const exmap & mp, unsigned options) const
00229 {
00230     exvector m2(row * col);
00231     for (unsigned r=0; r<row; ++r)
00232         for (unsigned c=0; c<col; ++c)
00233             m2[r*col+c] = m[r*col+c].subs(mp, options);
00234 
00235     return matrix(row, col, m2).subs_one_level(mp, options);
00236 }
00237 
00239 ex matrix::conjugate() const
00240 {
00241     exvector * ev = 0;
00242     for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i) {
00243         ex x = i->conjugate();
00244         if (ev) {
00245             ev->push_back(x);
00246             continue;
00247         }
00248         if (are_ex_trivially_equal(x, *i)) {
00249             continue;
00250         }
00251         ev = new exvector;
00252         ev->reserve(m.size());
00253         for (exvector::const_iterator j=m.begin(); j!=i; ++j) {
00254             ev->push_back(*j);
00255         }
00256         ev->push_back(x);
00257     }
00258     if (ev) {
00259         ex result = matrix(row, col, *ev);
00260         delete ev;
00261         return result;
00262     }
00263     return *this;
00264 }
00265 
00266 ex matrix::real_part() const
00267 {
00268     exvector v;
00269     v.reserve(m.size());
00270     for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
00271         v.push_back(i->real_part());
00272     return matrix(row, col, v);
00273 }
00274 
00275 ex matrix::imag_part() const
00276 {
00277     exvector v;
00278     v.reserve(m.size());
00279     for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
00280         v.push_back(i->imag_part());
00281     return matrix(row, col, v);
00282 }
00283 
00284 // protected
00285 
00286 int matrix::compare_same_type(const basic & other) const
00287 {
00288     GINAC_ASSERT(is_exactly_a<matrix>(other));
00289     const matrix &o = static_cast<const matrix &>(other);
00290     
00291     // compare number of rows
00292     if (row != o.rows())
00293         return row < o.rows() ? -1 : 1;
00294     
00295     // compare number of columns
00296     if (col != o.cols())
00297         return col < o.cols() ? -1 : 1;
00298     
00299     // equal number of rows and columns, compare individual elements
00300     int cmpval;
00301     for (unsigned r=0; r<row; ++r) {
00302         for (unsigned c=0; c<col; ++c) {
00303             cmpval = ((*this)(r,c)).compare(o(r,c));
00304             if (cmpval!=0) return cmpval;
00305         }
00306     }
00307     // all elements are equal => matrices are equal;
00308     return 0;
00309 }
00310 
00311 bool matrix::match_same_type(const basic & other) const
00312 {
00313     GINAC_ASSERT(is_exactly_a<matrix>(other));
00314     const matrix & o = static_cast<const matrix &>(other);
00315     
00316     // The number of rows and columns must be the same. This is necessary to
00317     // prevent a 2x3 matrix from matching a 3x2 one.
00318     return row == o.rows() && col == o.cols();
00319 }
00320 
00322 ex matrix::eval_indexed(const basic & i) const
00323 {
00324     GINAC_ASSERT(is_a<indexed>(i));
00325     GINAC_ASSERT(is_a<matrix>(i.op(0)));
00326 
00327     bool all_indices_unsigned = static_cast<const indexed &>(i).all_index_values_are(info_flags::nonnegint);
00328 
00329     // Check indices
00330     if (i.nops() == 2) {
00331 
00332         // One index, must be one-dimensional vector
00333         if (row != 1 && col != 1)
00334             throw (std::runtime_error("matrix::eval_indexed(): vector must have exactly 1 index"));
00335 
00336         const idx & i1 = ex_to<idx>(i.op(1));
00337 
00338         if (col == 1) {
00339 
00340             // Column vector
00341             if (!i1.get_dim().is_equal(row))
00342                 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
00343 
00344             // Index numeric -> return vector element
00345             if (all_indices_unsigned) {
00346                 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
00347                 if (n1 >= row)
00348                     throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
00349                 return (*this)(n1, 0);
00350             }
00351 
00352         } else {
00353 
00354             // Row vector
00355             if (!i1.get_dim().is_equal(col))
00356                 throw (std::runtime_error("matrix::eval_indexed(): dimension of index must match number of vector elements"));
00357 
00358             // Index numeric -> return vector element
00359             if (all_indices_unsigned) {
00360                 unsigned n1 = ex_to<numeric>(i1.get_value()).to_int();
00361                 if (n1 >= col)
00362                     throw (std::runtime_error("matrix::eval_indexed(): value of index exceeds number of vector elements"));
00363                 return (*this)(0, n1);
00364             }
00365         }
00366 
00367     } else if (i.nops() == 3) {
00368 
00369         // Two indices
00370         const idx & i1 = ex_to<idx>(i.op(1));
00371         const idx & i2 = ex_to<idx>(i.op(2));
00372 
00373         if (!i1.get_dim().is_equal(row))
00374             throw (std::runtime_error("matrix::eval_indexed(): dimension of first index must match number of rows"));
00375         if (!i2.get_dim().is_equal(col))
00376             throw (std::runtime_error("matrix::eval_indexed(): dimension of second index must match number of columns"));
00377 
00378         // Pair of dummy indices -> compute trace
00379         if (is_dummy_pair(i1, i2))
00380             return trace();
00381 
00382         // Both indices numeric -> return matrix element
00383         if (all_indices_unsigned) {
00384             unsigned n1 = ex_to<numeric>(i1.get_value()).to_int(), n2 = ex_to<numeric>(i2.get_value()).to_int();
00385             if (n1 >= row)
00386                 throw (std::runtime_error("matrix::eval_indexed(): value of first index exceeds number of rows"));
00387             if (n2 >= col)
00388                 throw (std::runtime_error("matrix::eval_indexed(): value of second index exceeds number of columns"));
00389             return (*this)(n1, n2);
00390         }
00391 
00392     } else
00393         throw (std::runtime_error("matrix::eval_indexed(): matrix must have exactly 2 indices"));
00394 
00395     return i.hold();
00396 }
00397 
00399 ex matrix::add_indexed(const ex & self, const ex & other) const
00400 {
00401     GINAC_ASSERT(is_a<indexed>(self));
00402     GINAC_ASSERT(is_a<matrix>(self.op(0)));
00403     GINAC_ASSERT(is_a<indexed>(other));
00404     GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
00405 
00406     // Only add two matrices
00407     if (is_a<matrix>(other.op(0))) {
00408         GINAC_ASSERT(other.nops() == 2 || other.nops() == 3);
00409 
00410         const matrix &self_matrix = ex_to<matrix>(self.op(0));
00411         const matrix &other_matrix = ex_to<matrix>(other.op(0));
00412 
00413         if (self.nops() == 2 && other.nops() == 2) { // vector + vector
00414 
00415             if (self_matrix.row == other_matrix.row)
00416                 return indexed(self_matrix.add(other_matrix), self.op(1));
00417             else if (self_matrix.row == other_matrix.col)
00418                 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1));
00419 
00420         } else if (self.nops() == 3 && other.nops() == 3) { // matrix + matrix
00421 
00422             if (self.op(1).is_equal(other.op(1)) && self.op(2).is_equal(other.op(2)))
00423                 return indexed(self_matrix.add(other_matrix), self.op(1), self.op(2));
00424             else if (self.op(1).is_equal(other.op(2)) && self.op(2).is_equal(other.op(1)))
00425                 return indexed(self_matrix.add(other_matrix.transpose()), self.op(1), self.op(2));
00426 
00427         }
00428     }
00429 
00430     // Don't know what to do, return unevaluated sum
00431     return self + other;
00432 }
00433 
00435 ex matrix::scalar_mul_indexed(const ex & self, const numeric & other) const
00436 {
00437     GINAC_ASSERT(is_a<indexed>(self));
00438     GINAC_ASSERT(is_a<matrix>(self.op(0)));
00439     GINAC_ASSERT(self.nops() == 2 || self.nops() == 3);
00440 
00441     const matrix &self_matrix = ex_to<matrix>(self.op(0));
00442 
00443     if (self.nops() == 2)
00444         return indexed(self_matrix.mul(other), self.op(1));
00445     else // self.nops() == 3
00446         return indexed(self_matrix.mul(other), self.op(1), self.op(2));
00447 }
00448 
00450 bool matrix::contract_with(exvector::iterator self, exvector::iterator other, exvector & v) const
00451 {
00452     GINAC_ASSERT(is_a<indexed>(*self));
00453     GINAC_ASSERT(is_a<indexed>(*other));
00454     GINAC_ASSERT(self->nops() == 2 || self->nops() == 3);
00455     GINAC_ASSERT(is_a<matrix>(self->op(0)));
00456 
00457     // Only contract with other matrices
00458     if (!is_a<matrix>(other->op(0)))
00459         return false;
00460 
00461     GINAC_ASSERT(other->nops() == 2 || other->nops() == 3);
00462 
00463     const matrix &self_matrix = ex_to<matrix>(self->op(0));
00464     const matrix &other_matrix = ex_to<matrix>(other->op(0));
00465 
00466     if (self->nops() == 2) {
00467 
00468         if (other->nops() == 2) { // vector * vector (scalar product)
00469 
00470             if (self_matrix.col == 1) {
00471                 if (other_matrix.col == 1) {
00472                     // Column vector * column vector, transpose first vector
00473                     *self = self_matrix.transpose().mul(other_matrix)(0, 0);
00474                 } else {
00475                     // Column vector * row vector, swap factors
00476                     *self = other_matrix.mul(self_matrix)(0, 0);
00477                 }
00478             } else {
00479                 if (other_matrix.col == 1) {
00480                     // Row vector * column vector, perfect
00481                     *self = self_matrix.mul(other_matrix)(0, 0);
00482                 } else {
00483                     // Row vector * row vector, transpose second vector
00484                     *self = self_matrix.mul(other_matrix.transpose())(0, 0);
00485                 }
00486             }
00487             *other = _ex1;
00488             return true;
00489 
00490         } else { // vector * matrix
00491 
00492             // B_i * A_ij = (B*A)_j (B is row vector)
00493             if (is_dummy_pair(self->op(1), other->op(1))) {
00494                 if (self_matrix.row == 1)
00495                     *self = indexed(self_matrix.mul(other_matrix), other->op(2));
00496                 else
00497                     *self = indexed(self_matrix.transpose().mul(other_matrix), other->op(2));
00498                 *other = _ex1;
00499                 return true;
00500             }
00501 
00502             // B_j * A_ij = (A*B)_i (B is column vector)
00503             if (is_dummy_pair(self->op(1), other->op(2))) {
00504                 if (self_matrix.col == 1)
00505                     *self = indexed(other_matrix.mul(self_matrix), other->op(1));
00506                 else
00507                     *self = indexed(other_matrix.mul(self_matrix.transpose()), other->op(1));
00508                 *other = _ex1;
00509                 return true;
00510             }
00511         }
00512 
00513     } else if (other->nops() == 3) { // matrix * matrix
00514 
00515         // A_ij * B_jk = (A*B)_ik
00516         if (is_dummy_pair(self->op(2), other->op(1))) {
00517             *self = indexed(self_matrix.mul(other_matrix), self->op(1), other->op(2));
00518             *other = _ex1;
00519             return true;
00520         }
00521 
00522         // A_ij * B_kj = (A*Btrans)_ik
00523         if (is_dummy_pair(self->op(2), other->op(2))) {
00524             *self = indexed(self_matrix.mul(other_matrix.transpose()), self->op(1), other->op(1));
00525             *other = _ex1;
00526             return true;
00527         }
00528 
00529         // A_ji * B_jk = (Atrans*B)_ik
00530         if (is_dummy_pair(self->op(1), other->op(1))) {
00531             *self = indexed(self_matrix.transpose().mul(other_matrix), self->op(2), other->op(2));
00532             *other = _ex1;
00533             return true;
00534         }
00535 
00536         // A_ji * B_kj = (B*A)_ki
00537         if (is_dummy_pair(self->op(1), other->op(2))) {
00538             *self = indexed(other_matrix.mul(self_matrix), other->op(1), self->op(2));
00539             *other = _ex1;
00540             return true;
00541         }
00542     }
00543 
00544     return false;
00545 }
00546 
00547 
00549 // non-virtual functions in this class
00551 
00552 // public
00553 
00557 matrix matrix::add(const matrix & other) const
00558 {
00559     if (col != other.col || row != other.row)
00560         throw std::logic_error("matrix::add(): incompatible matrices");
00561     
00562     exvector sum(this->m);
00563     exvector::iterator i = sum.begin(), end = sum.end();
00564     exvector::const_iterator ci = other.m.begin();
00565     while (i != end)
00566         *i++ += *ci++;
00567     
00568     return matrix(row,col,sum);
00569 }
00570 
00571 
00575 matrix matrix::sub(const matrix & other) const
00576 {
00577     if (col != other.col || row != other.row)
00578         throw std::logic_error("matrix::sub(): incompatible matrices");
00579     
00580     exvector dif(this->m);
00581     exvector::iterator i = dif.begin(), end = dif.end();
00582     exvector::const_iterator ci = other.m.begin();
00583     while (i != end)
00584         *i++ -= *ci++;
00585     
00586     return matrix(row,col,dif);
00587 }
00588 
00589 
00593 matrix matrix::mul(const matrix & other) const
00594 {
00595     if (this->cols() != other.rows())
00596         throw std::logic_error("matrix::mul(): incompatible matrices");
00597     
00598     exvector prod(this->rows()*other.cols());
00599     
00600     for (unsigned r1=0; r1<this->rows(); ++r1) {
00601         for (unsigned c=0; c<this->cols(); ++c) {
00602             // Quick test: can we shortcut?
00603             if (m[r1*col+c].is_zero())
00604                 continue;
00605             for (unsigned r2=0; r2<other.cols(); ++r2)
00606                 prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]);
00607         }
00608     }
00609     return matrix(row, other.col, prod);
00610 }
00611 
00612 
00614 matrix matrix::mul(const numeric & other) const
00615 {
00616     exvector prod(row * col);
00617 
00618     for (unsigned r=0; r<row; ++r)
00619         for (unsigned c=0; c<col; ++c)
00620             prod[r*col+c] = m[r*col+c] * other;
00621 
00622     return matrix(row, col, prod);
00623 }
00624 
00625 
00627 matrix matrix::mul_scalar(const ex & other) const
00628 {
00629     if (other.return_type() != return_types::commutative)
00630         throw std::runtime_error("matrix::mul_scalar(): non-commutative scalar");
00631 
00632     exvector prod(row * col);
00633 
00634     for (unsigned r=0; r<row; ++r)
00635         for (unsigned c=0; c<col; ++c)
00636             prod[r*col+c] = m[r*col+c] * other;
00637 
00638     return matrix(row, col, prod);
00639 }
00640 
00641 
00643 matrix matrix::pow(const ex & expn) const
00644 {
00645     if (col!=row)
00646         throw (std::logic_error("matrix::pow(): matrix not square"));
00647     
00648     if (is_exactly_a<numeric>(expn)) {
00649         // Integer cases are computed by successive multiplication, using the
00650         // obvious shortcut of storing temporaries, like A^4 == (A*A)*(A*A).
00651         if (expn.info(info_flags::integer)) {
00652             numeric b = ex_to<numeric>(expn);
00653             matrix A(row,col);
00654             if (expn.info(info_flags::negative)) {
00655                 b *= -1;
00656                 A = this->inverse();
00657             } else {
00658                 A = *this;
00659             }
00660             matrix C(row,col);
00661             for (unsigned r=0; r<row; ++r)
00662                 C(r,r) = _ex1;
00663             if (b.is_zero())
00664                 return C;
00665             // This loop computes the representation of b in base 2 from right
00666             // to left and multiplies the factors whenever needed.  Note
00667             // that this is not entirely optimal but close to optimal and
00668             // "better" algorithms are much harder to implement.  (See Knuth,
00669             // TAoCP2, section "Evaluation of Powers" for a good discussion.)
00670             while (b!=*_num1_p) {
00671                 if (b.is_odd()) {
00672                     C = C.mul(A);
00673                     --b;
00674                 }
00675                 b /= *_num2_p;  // still integer.
00676                 A = A.mul(A);
00677             }
00678             return A.mul(C);
00679         }
00680     }
00681     throw (std::runtime_error("matrix::pow(): don't know how to handle exponent"));
00682 }
00683 
00684 
00690 const ex & matrix::operator() (unsigned ro, unsigned co) const
00691 {
00692     if (ro>=row || co>=col)
00693         throw (std::range_error("matrix::operator(): index out of range"));
00694 
00695     return m[ro*col+co];
00696 }
00697 
00698 
00704 ex & matrix::operator() (unsigned ro, unsigned co)
00705 {
00706     if (ro>=row || co>=col)
00707         throw (std::range_error("matrix::operator(): index out of range"));
00708 
00709     ensure_if_modifiable();
00710     return m[ro*col+co];
00711 }
00712 
00713 
00716 matrix matrix::transpose() const
00717 {
00718     exvector trans(this->cols()*this->rows());
00719     
00720     for (unsigned r=0; r<this->cols(); ++r)
00721         for (unsigned c=0; c<this->rows(); ++c)
00722             trans[r*this->rows()+c] = m[c*this->cols()+r];
00723     
00724     return matrix(this->cols(),this->rows(),trans);
00725 }
00726 
00741 ex matrix::determinant(unsigned algo) const
00742 {
00743     if (row!=col)
00744         throw (std::logic_error("matrix::determinant(): matrix not square"));
00745     GINAC_ASSERT(row*col==m.capacity());
00746     
00747     // Gather some statistical information about this matrix:
00748     bool numeric_flag = true;
00749     bool normal_flag = false;
00750     unsigned sparse_count = 0;  // counts non-zero elements
00751     exvector::const_iterator r = m.begin(), rend = m.end();
00752     while (r != rend) {
00753         if (!r->info(info_flags::numeric))
00754             numeric_flag = false;
00755         exmap srl;  // symbol replacement list
00756         ex rtest = r->to_rational(srl);
00757         if (!rtest.is_zero())
00758             ++sparse_count;
00759         if (!rtest.info(info_flags::crational_polynomial) &&
00760              rtest.info(info_flags::rational_function))
00761             normal_flag = true;
00762         ++r;
00763     }
00764     
00765     // Here is the heuristics in case this routine has to decide:
00766     if (algo == determinant_algo::automatic) {
00767         // Minor expansion is generally a good guess:
00768         algo = determinant_algo::laplace;
00769         // Does anybody know when a matrix is really sparse?
00770         // Maybe <~row/2.236 nonzero elements average in a row?
00771         if (row>3 && 5*sparse_count<=row*col)
00772             algo = determinant_algo::bareiss;
00773         // Purely numeric matrix can be handled by Gauss elimination.
00774         // This overrides any prior decisions.
00775         if (numeric_flag)
00776             algo = determinant_algo::gauss;
00777     }
00778     
00779     // Trap the trivial case here, since some algorithms don't like it
00780     if (this->row==1) {
00781         // for consistency with non-trivial determinants...
00782         if (normal_flag)
00783             return m[0].normal();
00784         else
00785             return m[0].expand();
00786     }
00787 
00788     // Compute the determinant
00789     switch(algo) {
00790         case determinant_algo::gauss: {
00791             ex det = 1;
00792             matrix tmp(*this);
00793             int sign = tmp.gauss_elimination(true);
00794             for (unsigned d=0; d<row; ++d)
00795                 det *= tmp.m[d*col+d];
00796             if (normal_flag)
00797                 return (sign*det).normal();
00798             else
00799                 return (sign*det).normal().expand();
00800         }
00801         case determinant_algo::bareiss: {
00802             matrix tmp(*this);
00803             int sign;
00804             sign = tmp.fraction_free_elimination(true);
00805             if (normal_flag)
00806                 return (sign*tmp.m[row*col-1]).normal();
00807             else
00808                 return (sign*tmp.m[row*col-1]).expand();
00809         }
00810         case determinant_algo::divfree: {
00811             matrix tmp(*this);
00812             int sign;
00813             sign = tmp.division_free_elimination(true);
00814             if (sign==0)
00815                 return _ex0;
00816             ex det = tmp.m[row*col-1];
00817             // factor out accumulated bogus slag
00818             for (unsigned d=0; d<row-2; ++d)
00819                 for (unsigned j=0; j<row-d-2; ++j)
00820                     det = (det/tmp.m[d*col+d]).normal();
00821             return (sign*det);
00822         }
00823         case determinant_algo::laplace:
00824         default: {
00825             // This is the minor expansion scheme.  We always develop such
00826             // that the smallest minors (i.e, the trivial 1x1 ones) are on the
00827             // rightmost column.  For this to be efficient, empirical tests
00828             // have shown that the emptiest columns (i.e. the ones with most
00829             // zeros) should be the ones on the right hand side -- although
00830             // this might seem counter-intuitive (and in contradiction to some
00831             // literature like the FORM manual).  Please go ahead and test it
00832             // if you don't believe me!  Therefore we presort the columns of
00833             // the matrix:
00834             typedef std::pair<unsigned,unsigned> uintpair;
00835             std::vector<uintpair> c_zeros;  // number of zeros in column
00836             for (unsigned c=0; c<col; ++c) {
00837                 unsigned acc = 0;
00838                 for (unsigned r=0; r<row; ++r)
00839                     if (m[r*col+c].is_zero())
00840                         ++acc;
00841                 c_zeros.push_back(uintpair(acc,c));
00842             }
00843             std::sort(c_zeros.begin(),c_zeros.end());
00844             std::vector<unsigned> pre_sort;
00845             for (std::vector<uintpair>::const_iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
00846                 pre_sort.push_back(i->second);
00847             std::vector<unsigned> pre_sort_test(pre_sort); // permutation_sign() modifies the vector so we make a copy here
00848             int sign = permutation_sign(pre_sort_test.begin(), pre_sort_test.end());
00849             exvector result(row*col);  // represents sorted matrix
00850             unsigned c = 0;
00851             for (std::vector<unsigned>::const_iterator i=pre_sort.begin();
00852                  i!=pre_sort.end();
00853                  ++i,++c) {
00854                 for (unsigned r=0; r<row; ++r)
00855                     result[r*col+c] = m[r*col+(*i)];
00856             }
00857             
00858             if (normal_flag)
00859                 return (sign*matrix(row,col,result).determinant_minor()).normal();
00860             else
00861                 return sign*matrix(row,col,result).determinant_minor();
00862         }
00863     }
00864 }
00865 
00866 
00873 ex matrix::trace() const
00874 {
00875     if (row != col)
00876         throw (std::logic_error("matrix::trace(): matrix not square"));
00877     
00878     ex tr;
00879     for (unsigned r=0; r<col; ++r)
00880         tr += m[r*col+r];
00881     
00882     if (tr.info(info_flags::rational_function) &&
00883        !tr.info(info_flags::crational_polynomial))
00884         return tr.normal();
00885     else
00886         return tr.expand();
00887 }
00888 
00889 
00901 ex matrix::charpoly(const ex & lambda) const
00902 {
00903     if (row != col)
00904         throw (std::logic_error("matrix::charpoly(): matrix not square"));
00905     
00906     bool numeric_flag = true;
00907     exvector::const_iterator r = m.begin(), rend = m.end();
00908     while (r!=rend && numeric_flag==true) {
00909         if (!r->info(info_flags::numeric))
00910             numeric_flag = false;
00911         ++r;
00912     }
00913     
00914     // The pure numeric case is traditionally rather common.  Hence, it is
00915     // trapped and we use Leverrier's algorithm which goes as row^3 for
00916     // every coefficient.  The expensive part is the matrix multiplication.
00917     if (numeric_flag) {
00918 
00919         matrix B(*this);
00920         ex c = B.trace();
00921         ex poly = power(lambda, row) - c*power(lambda, row-1);
00922         for (unsigned i=1; i<row; ++i) {
00923             for (unsigned j=0; j<row; ++j)
00924                 B.m[j*col+j] -= c;
00925             B = this->mul(B);
00926             c = B.trace() / ex(i+1);
00927             poly -= c*power(lambda, row-i-1);
00928         }
00929         if (row%2)
00930             return -poly;
00931         else
00932             return poly;
00933 
00934     } else {
00935     
00936         matrix M(*this);
00937         for (unsigned r=0; r<col; ++r)
00938             M.m[r*col+r] -= lambda;
00939     
00940         return M.determinant().collect(lambda);
00941     }
00942 }
00943 
00944 
00950 matrix matrix::inverse() const
00951 {
00952     if (row != col)
00953         throw (std::logic_error("matrix::inverse(): matrix not square"));
00954     
00955     // This routine actually doesn't do anything fancy at all.  We compute the
00956     // inverse of the matrix A by solving the system A * A^{-1} == Id.
00957     
00958     // First populate the identity matrix supposed to become the right hand side.
00959     matrix identity(row,col);
00960     for (unsigned i=0; i<row; ++i)
00961         identity(i,i) = _ex1;
00962     
00963     // Populate a dummy matrix of variables, just because of compatibility with
00964     // matrix::solve() which wants this (for compatibility with under-determined
00965     // systems of equations).
00966     matrix vars(row,col);
00967     for (unsigned r=0; r<row; ++r)
00968         for (unsigned c=0; c<col; ++c)
00969             vars(r,c) = symbol();
00970     
00971     matrix sol(row,col);
00972     try {
00973         sol = this->solve(vars,identity);
00974     } catch (const std::runtime_error & e) {
00975         if (e.what()==std::string("matrix::solve(): inconsistent linear system"))
00976             throw (std::runtime_error("matrix::inverse(): singular matrix"));
00977         else
00978             throw;
00979     }
00980     return sol;
00981 }
00982 
00983 
00995 matrix matrix::solve(const matrix & vars,
00996                      const matrix & rhs,
00997                      unsigned algo) const
00998 {
00999     const unsigned m = this->rows();
01000     const unsigned n = this->cols();
01001     const unsigned p = rhs.cols();
01002     
01003     // syntax checks    
01004     if ((rhs.rows() != m) || (vars.rows() != n) || (vars.col != p))
01005         throw (std::logic_error("matrix::solve(): incompatible matrices"));
01006     for (unsigned ro=0; ro<n; ++ro)
01007         for (unsigned co=0; co<p; ++co)
01008             if (!vars(ro,co).info(info_flags::symbol))
01009                 throw (std::invalid_argument("matrix::solve(): 1st argument must be matrix of symbols"));
01010     
01011     // build the augmented matrix of *this with rhs attached to the right
01012     matrix aug(m,n+p);
01013     for (unsigned r=0; r<m; ++r) {
01014         for (unsigned c=0; c<n; ++c)
01015             aug.m[r*(n+p)+c] = this->m[r*n+c];
01016         for (unsigned c=0; c<p; ++c)
01017             aug.m[r*(n+p)+c+n] = rhs.m[r*p+c];
01018     }
01019     
01020     // Gather some statistical information about the augmented matrix:
01021     bool numeric_flag = true;
01022     exvector::const_iterator r = aug.m.begin(), rend = aug.m.end();
01023     while (r!=rend && numeric_flag==true) {
01024         if (!r->info(info_flags::numeric))
01025             numeric_flag = false;
01026         ++r;
01027     }
01028     
01029     // Here is the heuristics in case this routine has to decide:
01030     if (algo == solve_algo::automatic) {
01031         // Bareiss (fraction-free) elimination is generally a good guess:
01032         algo = solve_algo::bareiss;
01033         // For m<3, Bareiss elimination is equivalent to division free
01034         // elimination but has more logistic overhead
01035         if (m<3)
01036             algo = solve_algo::divfree;
01037         // This overrides any prior decisions.
01038         if (numeric_flag)
01039             algo = solve_algo::gauss;
01040     }
01041     
01042     // Eliminate the augmented matrix:
01043     switch(algo) {
01044         case solve_algo::gauss:
01045             aug.gauss_elimination();
01046             break;
01047         case solve_algo::divfree:
01048             aug.division_free_elimination();
01049             break;
01050         case solve_algo::bareiss:
01051         default:
01052             aug.fraction_free_elimination();
01053     }
01054     
01055     // assemble the solution matrix:
01056     matrix sol(n,p);
01057     for (unsigned co=0; co<p; ++co) {
01058         unsigned last_assigned_sol = n+1;
01059         for (int r=m-1; r>=0; --r) {
01060             unsigned fnz = 1;    // first non-zero in row
01061             while ((fnz<=n) && (aug.m[r*(n+p)+(fnz-1)].is_zero()))
01062                 ++fnz;
01063             if (fnz>n) {
01064                 // row consists only of zeros, corresponding rhs must be 0, too
01065                 if (!aug.m[r*(n+p)+n+co].is_zero()) {
01066                     throw (std::runtime_error("matrix::solve(): inconsistent linear system"));
01067                 }
01068             } else {
01069                 // assign solutions for vars between fnz+1 and
01070                 // last_assigned_sol-1: free parameters
01071                 for (unsigned c=fnz; c<last_assigned_sol-1; ++c)
01072                     sol(c,co) = vars.m[c*p+co];
01073                 ex e = aug.m[r*(n+p)+n+co];
01074                 for (unsigned c=fnz; c<n; ++c)
01075                     e -= aug.m[r*(n+p)+c]*sol.m[c*p+co];
01076                 sol(fnz-1,co) = (e/(aug.m[r*(n+p)+(fnz-1)])).normal();
01077                 last_assigned_sol = fnz;
01078             }
01079         }
01080         // assign solutions for vars between 1 and
01081         // last_assigned_sol-1: free parameters
01082         for (unsigned ro=0; ro<last_assigned_sol-1; ++ro)
01083             sol(ro,co) = vars(ro,co);
01084     }
01085     
01086     return sol;
01087 }
01088 
01089 
01091 unsigned matrix::rank() const
01092 {
01093     // Method:
01094     // Transform this matrix into upper echelon form and then count the
01095     // number of non-zero rows.
01096 
01097     GINAC_ASSERT(row*col==m.capacity());
01098 
01099     // Actually, any elimination scheme will do since we are only
01100     // interested in the echelon matrix' zeros.
01101     matrix to_eliminate = *this;
01102     to_eliminate.fraction_free_elimination();
01103 
01104     unsigned r = row*col;  // index of last non-zero element
01105     while (r--) {
01106         if (!to_eliminate.m[r].is_zero())
01107             return 1+r/col;
01108     }
01109     return 0;
01110 }
01111 
01112 
01113 // protected
01114 
01125 ex matrix::determinant_minor() const
01126 {
01127     // for small matrices the algorithm does not make any sense:
01128     const unsigned n = this->cols();
01129     if (n==1)
01130         return m[0].expand();
01131     if (n==2)
01132         return (m[0]*m[3]-m[2]*m[1]).expand();
01133     if (n==3)
01134         return (m[0]*m[4]*m[8]-m[0]*m[5]*m[7]-
01135                 m[1]*m[3]*m[8]+m[2]*m[3]*m[7]+
01136                 m[1]*m[5]*m[6]-m[2]*m[4]*m[6]).expand();
01137     
01138     // This algorithm can best be understood by looking at a naive
01139     // implementation of Laplace-expansion, like this one:
01140     // ex det;
01141     // matrix minorM(this->rows()-1,this->cols()-1);
01142     // for (unsigned r1=0; r1<this->rows(); ++r1) {
01143     //     // shortcut if element(r1,0) vanishes
01144     //     if (m[r1*col].is_zero())
01145     //         continue;
01146     //     // assemble the minor matrix
01147     //     for (unsigned r=0; r<minorM.rows(); ++r) {
01148     //         for (unsigned c=0; c<minorM.cols(); ++c) {
01149     //             if (r<r1)
01150     //                 minorM(r,c) = m[r*col+c+1];
01151     //             else
01152     //                 minorM(r,c) = m[(r+1)*col+c+1];
01153     //         }
01154     //     }
01155     //     // recurse down and care for sign:
01156     //     if (r1%2)
01157     //         det -= m[r1*col] * minorM.determinant_minor();
01158     //     else
01159     //         det += m[r1*col] * minorM.determinant_minor();
01160     // }
01161     // return det.expand();
01162     // What happens is that while proceeding down many of the minors are
01163     // computed more than once.  In particular, there are binomial(n,k)
01164     // kxk minors and each one is computed factorial(n-k) times.  Therefore
01165     // it is reasonable to store the results of the minors.  We proceed from
01166     // right to left.  At each column c we only need to retrieve the minors
01167     // calculated in step c-1.  We therefore only have to store at most 
01168     // 2*binomial(n,n/2) minors.
01169     
01170     // Unique flipper counter for partitioning into minors
01171     std::vector<unsigned> Pkey;
01172     Pkey.reserve(n);
01173     // key for minor determinant (a subpartition of Pkey)
01174     std::vector<unsigned> Mkey;
01175     Mkey.reserve(n-1);
01176     // we store our subminors in maps, keys being the rows they arise from
01177     typedef std::map<std::vector<unsigned>,class ex> Rmap;
01178     typedef std::map<std::vector<unsigned>,class ex>::value_type Rmap_value;
01179     Rmap A;
01180     Rmap B;
01181     ex det;
01182     // initialize A with last column:
01183     for (unsigned r=0; r<n; ++r) {
01184         Pkey.erase(Pkey.begin(),Pkey.end());
01185         Pkey.push_back(r);
01186         A.insert(Rmap_value(Pkey,m[n*(r+1)-1]));
01187     }
01188     // proceed from right to left through matrix
01189     for (int c=n-2; c>=0; --c) {
01190         Pkey.erase(Pkey.begin(),Pkey.end());  // don't change capacity
01191         Mkey.erase(Mkey.begin(),Mkey.end());
01192         for (unsigned i=0; i<n-c; ++i)
01193             Pkey.push_back(i);
01194         unsigned fc = 0;  // controls logic for our strange flipper counter
01195         do {
01196             det = _ex0;
01197             for (unsigned r=0; r<n-c; ++r) {
01198                 // maybe there is nothing to do?
01199                 if (m[Pkey[r]*n+c].is_zero())
01200                     continue;
01201                 // create the sorted key for all possible minors
01202                 Mkey.erase(Mkey.begin(),Mkey.end());
01203                 for (unsigned i=0; i<n-c; ++i)
01204                     if (i!=r)
01205                         Mkey.push_back(Pkey[i]);
01206                 // Fetch the minors and compute the new determinant
01207                 if (r%2)
01208                     det -= m[Pkey[r]*n+c]*A[Mkey];
01209                 else
01210                     det += m[Pkey[r]*n+c]*A[Mkey];
01211             }
01212             // prevent build-up of deep nesting of expressions saves time:
01213             det = det.expand();
01214             // store the new determinant at its place in B:
01215             if (!det.is_zero())
01216                 B.insert(Rmap_value(Pkey,det));
01217             // increment our strange flipper counter
01218             for (fc=n-c; fc>0; --fc) {
01219                 ++Pkey[fc-1];
01220                 if (Pkey[fc-1]<fc+c)
01221                     break;
01222             }
01223             if (fc<n-c && fc>0)
01224                 for (unsigned j=fc; j<n-c; ++j)
01225                     Pkey[j] = Pkey[j-1]+1;
01226         } while(fc);
01227         // next column, so change the role of A and B:
01228         A.swap(B);
01229         B.clear();
01230     }
01231     
01232     return det;
01233 }
01234 
01235 
01245 int matrix::gauss_elimination(const bool det)
01246 {
01247     ensure_if_modifiable();
01248     const unsigned m = this->rows();
01249     const unsigned n = this->cols();
01250     GINAC_ASSERT(!det || n==m);
01251     int sign = 1;
01252     
01253     unsigned r0 = 0;
01254     for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
01255         int indx = pivot(r0, c0, true);
01256         if (indx == -1) {
01257             sign = 0;
01258             if (det)
01259                 return 0;  // leaves *this in a messy state
01260         }
01261         if (indx>=0) {
01262             if (indx > 0)
01263                 sign = -sign;
01264             for (unsigned r2=r0+1; r2<m; ++r2) {
01265                 if (!this->m[r2*n+c0].is_zero()) {
01266                     // yes, there is something to do in this row
01267                     ex piv = this->m[r2*n+c0] / this->m[r0*n+c0];
01268                     for (unsigned c=c0+1; c<n; ++c) {
01269                         this->m[r2*n+c] -= piv * this->m[r0*n+c];
01270                         if (!this->m[r2*n+c].info(info_flags::numeric))
01271                             this->m[r2*n+c] = this->m[r2*n+c].normal();
01272                     }
01273                 }
01274                 // fill up left hand side with zeros
01275                 for (unsigned c=r0; c<=c0; ++c)
01276                     this->m[r2*n+c] = _ex0;
01277             }
01278             if (det) {
01279                 // save space by deleting no longer needed elements
01280                 for (unsigned c=r0+1; c<n; ++c)
01281                     this->m[r0*n+c] = _ex0;
01282             }
01283             ++r0;
01284         }
01285     }
01286     // clear remaining rows
01287     for (unsigned r=r0+1; r<m; ++r) {
01288         for (unsigned c=0; c<n; ++c)
01289             this->m[r*n+c] = _ex0;
01290     }
01291 
01292     return sign;
01293 }
01294 
01295 
01304 int matrix::division_free_elimination(const bool det)
01305 {
01306     ensure_if_modifiable();
01307     const unsigned m = this->rows();
01308     const unsigned n = this->cols();
01309     GINAC_ASSERT(!det || n==m);
01310     int sign = 1;
01311     
01312     unsigned r0 = 0;
01313     for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
01314         int indx = pivot(r0, c0, true);
01315         if (indx==-1) {
01316             sign = 0;
01317             if (det)
01318                 return 0;  // leaves *this in a messy state
01319         }
01320         if (indx>=0) {
01321             if (indx>0)
01322                 sign = -sign;
01323             for (unsigned r2=r0+1; r2<m; ++r2) {
01324                 for (unsigned c=c0+1; c<n; ++c)
01325                     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();
01326                 // fill up left hand side with zeros
01327                 for (unsigned c=r0; c<=c0; ++c)
01328                     this->m[r2*n+c] = _ex0;
01329             }
01330             if (det) {
01331                 // save space by deleting no longer needed elements
01332                 for (unsigned c=r0+1; c<n; ++c)
01333                     this->m[r0*n+c] = _ex0;
01334             }
01335             ++r0;
01336         }
01337     }
01338     // clear remaining rows
01339     for (unsigned r=r0+1; r<m; ++r) {
01340         for (unsigned c=0; c<n; ++c)
01341             this->m[r*n+c] = _ex0;
01342     }
01343 
01344     return sign;
01345 }
01346 
01347 
01358 int matrix::fraction_free_elimination(const bool det)
01359 {
01360     // Method:
01361     // (single-step fraction free elimination scheme, already known to Jordan)
01362     //
01363     // Usual division-free elimination sets m[0](r,c) = m(r,c) and then sets
01364     //     m[k+1](r,c) = m[k](k,k) * m[k](r,c) - m[k](r,k) * m[k](k,c).
01365     //
01366     // Bareiss (fraction-free) elimination in addition divides that element
01367     // by m[k-1](k-1,k-1) for k>1, where it can be shown by means of the
01368     // Sylvester identity that this really divides m[k+1](r,c).
01369     //
01370     // We also allow rational functions where the original prove still holds.
01371     // However, we must care for numerator and denominator separately and
01372     // "manually" work in the integral domains because of subtle cancellations
01373     // (see below).  This blows up the bookkeeping a bit and the formula has
01374     // to be modified to expand like this (N{x} stands for numerator of x,
01375     // D{x} for denominator of x):
01376     //     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)}
01377     //                     -N{m[k](r,k)}*N{m[k](k,c)}*D{m[k](k,k)}*D{m[k](r,c)}
01378     //     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)}
01379     // where for k>1 we now divide N{m[k+1](r,c)} by
01380     //     N{m[k-1](k-1,k-1)}
01381     // and D{m[k+1](r,c)} by
01382     //     D{m[k-1](k-1,k-1)}.
01383     
01384     ensure_if_modifiable();
01385     const unsigned m = this->rows();
01386     const unsigned n = this->cols();
01387     GINAC_ASSERT(!det || n==m);
01388     int sign = 1;
01389     if (m==1)
01390         return 1;
01391     ex divisor_n = 1;
01392     ex divisor_d = 1;
01393     ex dividend_n;
01394     ex dividend_d;
01395     
01396     // We populate temporary matrices to subsequently operate on.  There is
01397     // one holding numerators and another holding denominators of entries.
01398     // This is a must since the evaluator (or even earlier mul's constructor)
01399     // might cancel some trivial element which causes divide() to fail.  The
01400     // elements are normalized first (yes, even though this algorithm doesn't
01401     // need GCDs) since the elements of *this might be unnormalized, which
01402     // makes things more complicated than they need to be.
01403     matrix tmp_n(*this);
01404     matrix tmp_d(m,n);  // for denominators, if needed
01405     exmap srl;  // symbol replacement list
01406     exvector::const_iterator cit = this->m.begin(), citend = this->m.end();
01407     exvector::iterator tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
01408     while (cit != citend) {
01409         ex nd = cit->normal().to_rational(srl).numer_denom();
01410         ++cit;
01411         *tmp_n_it++ = nd.op(0);
01412         *tmp_d_it++ = nd.op(1);
01413     }
01414     
01415     unsigned r0 = 0;
01416     for (unsigned c0=0; c0<n && r0<m-1; ++c0) {
01417         // When trying to find a pivot, we should try a bit harder than expand().
01418         // Searching the first non-zero element in-place here instead of calling
01419         // pivot() allows us to do no more substitutions and back-substitutions
01420         // than are actually necessary.
01421         unsigned indx = r0;
01422         while ((indx<m) &&
01423                (tmp_n[indx*n+c0].subs(srl, subs_options::no_pattern).expand().is_zero()))
01424             ++indx;
01425         if (indx==m) {
01426             // all elements in column c0 below row r0 vanish
01427             sign = 0;
01428             if (det)
01429                 return 0;
01430         } else {
01431             if (indx>r0) {
01432                 // Matrix needs pivoting, swap rows r0 and indx of tmp_n and tmp_d.
01433                 sign = -sign;
01434                 for (unsigned c=c0; c<n; ++c) {
01435                     tmp_n.m[n*indx+c].swap(tmp_n.m[n*r0+c]);
01436                     tmp_d.m[n*indx+c].swap(tmp_d.m[n*r0+c]);
01437                 }
01438             }
01439             for (unsigned r2=r0+1; r2<m; ++r2) {
01440                 for (unsigned c=c0+1; c<n; ++c) {
01441                     dividend_n = (tmp_n.m[r0*n+c0]*tmp_n.m[r2*n+c]*
01442                                   tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]
01443                                  -tmp_n.m[r2*n+c0]*tmp_n.m[r0*n+c]*
01444                                   tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
01445                     dividend_d = (tmp_d.m[r2*n+c0]*tmp_d.m[r0*n+c]*
01446                                   tmp_d.m[r0*n+c0]*tmp_d.m[r2*n+c]).expand();
01447                     bool check = divide(dividend_n, divisor_n,
01448                                         tmp_n.m[r2*n+c], true);
01449                     check &= divide(dividend_d, divisor_d,
01450                                     tmp_d.m[r2*n+c], true);
01451                     GINAC_ASSERT(check);
01452                 }
01453                 // fill up left hand side with zeros
01454                 for (unsigned c=r0; c<=c0; ++c)
01455                     tmp_n.m[r2*n+c] = _ex0;
01456             }
01457             if (c0<n && r0<m-1) {
01458                 // compute next iteration's divisor
01459                 divisor_n = tmp_n.m[r0*n+c0].expand();
01460                 divisor_d = tmp_d.m[r0*n+c0].expand();
01461                 if (det) {
01462                     // save space by deleting no longer needed elements
01463                     for (unsigned c=0; c<n; ++c) {
01464                         tmp_n.m[r0*n+c] = _ex0;
01465                         tmp_d.m[r0*n+c] = _ex1;
01466                     }
01467                 }
01468             }
01469             ++r0;
01470         }
01471     }
01472     // clear remaining rows
01473     for (unsigned r=r0+1; r<m; ++r) {
01474         for (unsigned c=0; c<n; ++c)
01475             tmp_n.m[r*n+c] = _ex0;
01476     }
01477 
01478     // repopulate *this matrix:
01479     exvector::iterator it = this->m.begin(), itend = this->m.end();
01480     tmp_n_it = tmp_n.m.begin();
01481     tmp_d_it = tmp_d.m.begin();
01482     while (it != itend)
01483         *it++ = ((*tmp_n_it++)/(*tmp_d_it++)).subs(srl, subs_options::no_pattern);
01484     
01485     return sign;
01486 }
01487 
01488 
01502 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
01503 {
01504     unsigned k = ro;
01505     if (symbolic) {
01506         // search first non-zero element in column co beginning at row ro
01507         while ((k<row) && (this->m[k*col+co].expand().is_zero()))
01508             ++k;
01509     } else {
01510         // search largest element in column co beginning at row ro
01511         GINAC_ASSERT(is_exactly_a<numeric>(this->m[k*col+co]));
01512         unsigned kmax = k+1;
01513         numeric mmax = abs(ex_to<numeric>(m[kmax*col+co]));
01514         while (kmax<row) {
01515             GINAC_ASSERT(is_exactly_a<numeric>(this->m[kmax*col+co]));
01516             numeric tmp = ex_to<numeric>(this->m[kmax*col+co]);
01517             if (abs(tmp) > mmax) {
01518                 mmax = tmp;
01519                 k = kmax;
01520             }
01521             ++kmax;
01522         }
01523         if (!mmax.is_zero())
01524             k = kmax;
01525     }
01526     if (k==row)
01527         // all elements in column co below row ro vanish
01528         return -1;
01529     if (k==ro)
01530         // matrix needs no pivoting
01531         return 0;
01532     // matrix needs pivoting, so swap rows k and ro
01533     ensure_if_modifiable();
01534     for (unsigned c=0; c<col; ++c)
01535         this->m[k*col+c].swap(this->m[ro*col+c]);
01536     
01537     return k;
01538 }
01539 
01542 bool matrix::is_zero_matrix() const
01543 {
01544     for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i) 
01545         if(!(i->is_zero()))
01546             return false;
01547     return true;
01548 }
01549 
01550 ex lst_to_matrix(const lst & l)
01551 {
01552     lst::const_iterator itr, itc;
01553 
01554     // Find number of rows and columns
01555     size_t rows = l.nops(), cols = 0;
01556     for (itr = l.begin(); itr != l.end(); ++itr) {
01557         if (!is_a<lst>(*itr))
01558             throw (std::invalid_argument("lst_to_matrix: argument must be a list of lists"));
01559         if (itr->nops() > cols)
01560             cols = itr->nops();
01561     }
01562 
01563     // Allocate and fill matrix
01564     matrix &M = *new matrix(rows, cols);
01565     M.setflag(status_flags::dynallocated);
01566 
01567     unsigned i;
01568     for (itr = l.begin(), i = 0; itr != l.end(); ++itr, ++i) {
01569         unsigned j;
01570         for (itc = ex_to<lst>(*itr).begin(), j = 0; itc != ex_to<lst>(*itr).end(); ++itc, ++j)
01571             M(i, j) = *itc;
01572     }
01573 
01574     return M;
01575 }
01576 
01577 ex diag_matrix(const lst & l)
01578 {
01579     lst::const_iterator it;
01580     size_t dim = l.nops();
01581 
01582     // Allocate and fill matrix
01583     matrix &M = *new matrix(dim, dim);
01584     M.setflag(status_flags::dynallocated);
01585 
01586     unsigned i;
01587     for (it = l.begin(), i = 0; it != l.end(); ++it, ++i)
01588         M(i, i) = *it;
01589 
01590     return M;
01591 }
01592 
01593 ex unit_matrix(unsigned r, unsigned c)
01594 {
01595     matrix &Id = *new matrix(r, c);
01596     Id.setflag(status_flags::dynallocated);
01597     for (unsigned i=0; i<r && i<c; i++)
01598         Id(i,i) = _ex1;
01599 
01600     return Id;
01601 }
01602 
01603 ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const std::string & tex_base_name)
01604 {
01605     matrix &M = *new matrix(r, c);
01606     M.setflag(status_flags::dynallocated | status_flags::evaluated);
01607 
01608     bool long_format = (r > 10 || c > 10);
01609     bool single_row = (r == 1 || c == 1);
01610 
01611     for (unsigned i=0; i<r; i++) {
01612         for (unsigned j=0; j<c; j++) {
01613             std::ostringstream s1, s2;
01614             s1 << base_name;
01615             s2 << tex_base_name << "_{";
01616             if (single_row) {
01617                 if (c == 1) {
01618                     s1 << i;
01619                     s2 << i << '}';
01620                 } else {
01621                     s1 << j;
01622                     s2 << j << '}';
01623                 }
01624             } else {
01625                 if (long_format) {
01626                     s1 << '_' << i << '_' << j;
01627                     s2 << i << ';' << j << "}";
01628                 } else {
01629                     s1 << i << j;
01630                     s2 << i << j << '}';
01631                 }
01632             }
01633             M(i, j) = symbol(s1.str(), s2.str());
01634         }
01635     }
01636 
01637     return M;
01638 }
01639 
01640 ex reduced_matrix(const matrix& m, unsigned r, unsigned c)
01641 {
01642     if (r+1>m.rows() || c+1>m.cols() || m.cols()<2 || m.rows()<2)
01643         throw std::runtime_error("minor_matrix(): index out of bounds");
01644 
01645     const unsigned rows = m.rows()-1;
01646     const unsigned cols = m.cols()-1;
01647     matrix &M = *new matrix(rows, cols);
01648     M.setflag(status_flags::dynallocated | status_flags::evaluated);
01649 
01650     unsigned ro = 0;
01651     unsigned ro2 = 0;
01652     while (ro2<rows) {
01653         if (ro==r)
01654             ++ro;
01655         unsigned co = 0;
01656         unsigned co2 = 0;
01657         while (co2<cols) {
01658             if (co==c)
01659                 ++co;
01660             M(ro2,co2) = m(ro, co);
01661             ++co;
01662             ++co2;
01663         }
01664         ++ro;
01665         ++ro2;
01666     }
01667 
01668     return M;
01669 }
01670 
01671 ex sub_matrix(const matrix&m, unsigned r, unsigned nr, unsigned c, unsigned nc)
01672 {
01673     if (r+nr>m.rows() || c+nc>m.cols())
01674         throw std::runtime_error("sub_matrix(): index out of bounds");
01675 
01676     matrix &M = *new matrix(nr, nc);
01677     M.setflag(status_flags::dynallocated | status_flags::evaluated);
01678 
01679     for (unsigned ro=0; ro<nr; ++ro) {
01680         for (unsigned co=0; co<nc; ++co) {
01681             M(ro,co) = m(ro+r,co+c);
01682         }
01683     }
01684 
01685     return M;
01686 }
01687 
01688 } // namespace GiNaC

This page is part of the GiNaC developer's reference. It was generated automatically by doxygen. For an introduction, see the tutorial.