]> www.ginac.de Git - ginac.git/blobdiff - ginac/matrix.cpp
Make .eval() evaluate top-level only.
[ginac.git] / ginac / matrix.cpp
index d108bd4c04be6bb315a3966c3f4a716d77a45504..c842885f8103c0a2cc6530b6fc65b6b203197c56 100644 (file)
@@ -3,7 +3,7 @@
  *  Implementation of symbolic matrices */
 
 /*
- *  GiNaC Copyright (C) 1999-2008 Johannes Gutenberg University Mainz, Germany
+ *  GiNaC Copyright (C) 1999-2015 Johannes Gutenberg University Mainz, Germany
  *
  *  This program is free software; you can redistribute it and/or modify
  *  it under the terms of the GNU General Public License as published by
  *  Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
  */
 
-#include <string>
-#include <iostream>
-#include <sstream>
-#include <algorithm>
-#include <map>
-#include <stdexcept>
-
 #include "matrix.h"
 #include "numeric.h"
 #include "lst.h"
 #include "archive.h"
 #include "utils.h"
 
+#include <algorithm>
+#include <iostream>
+#include <map>
+#include <sstream>
+#include <stdexcept>
+#include <string>
+
 namespace GiNaC {
 
 GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(matrix, basic,
@@ -53,7 +53,7 @@ GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(matrix, basic,
 //////////
 
 /** Default ctor.  Initializes to 1 x 1-dimensional zero-matrix. */
-matrix::matrix() : inherited(&matrix::tinfo_static), row(1), col(1), m(1, _ex0)
+matrix::matrix() : row(1), col(1), m(1, _ex0)
 {
        setflag(status_flags::not_shareable);
 }
@@ -68,17 +68,7 @@ matrix::matrix() : inherited(&matrix::tinfo_static), row(1), col(1), m(1, _ex0)
  *
  *  @param r number of rows
  *  @param c number of cols */
-matrix::matrix(unsigned r, unsigned c)
-  : inherited(&matrix::tinfo_static), row(r), col(c), m(r*c, _ex0)
-{
-       setflag(status_flags::not_shareable);
-}
-
-// protected
-
-/** Ctor from representation, for internal use only. */
-matrix::matrix(unsigned r, unsigned c, const exvector & m2)
-  : inherited(&matrix::tinfo_static), row(r), col(c), m(m2)
+matrix::matrix(unsigned r, unsigned c) : row(r), col(c), m(r*c, _ex0)
 {
        setflag(status_flags::not_shareable);
 }
@@ -88,55 +78,89 @@ matrix::matrix(unsigned r, unsigned c, const exvector & m2)
  *  If the list has more elements than the matrix, the excessive elements are
  *  thrown away. */
 matrix::matrix(unsigned r, unsigned c, const lst & l)
-  : inherited(&matrix::tinfo_static), row(r), col(c), m(r*c, _ex0)
+  : row(r), col(c), m(r*c, _ex0)
 {
        setflag(status_flags::not_shareable);
 
        size_t i = 0;
-       for (lst::const_iterator it = l.begin(); it != l.end(); ++it, ++i) {
+       for (auto & it : l) {
                size_t x = i % c;
                size_t y = i / c;
                if (y >= r)
                        break; // matrix smaller than list: throw away excessive elements
-               m[y*c+x] = *it;
+               m[y*c+x] = it;
+               ++i;
+       }
+}
+
+/** Construct a matrix from an 2 dimensional initializer list.
+ *  Throws an exception if some row has a different length than all the others.
+ */
+matrix::matrix(std::initializer_list<std::initializer_list<ex>> l)
+  : row(l.size()), col(l.begin()->size())
+{
+       setflag(status_flags::not_shareable);
+
+       m.reserve(row*col);
+       for (const auto & r : l) {
+               unsigned c = 0;
+               for (const auto & e : r) {
+                       m.push_back(e);
+                       ++c;
+               }
+               if (c != col)
+                       throw std::invalid_argument("matrix::matrix{{}}: wrong dimension");
        }
 }
 
+// protected
+
+/** Ctor from representation, for internal use only. */
+matrix::matrix(unsigned r, unsigned c, const exvector & m2)
+  : row(r), col(c), m(m2)
+{
+       setflag(status_flags::not_shareable);
+}
+matrix::matrix(unsigned r, unsigned c, exvector && m2)
+  : row(r), col(c), m(std::move(m2))
+{
+       setflag(status_flags::not_shareable);
+}
+
 //////////
 // archiving
 //////////
 
-matrix::matrix(const archive_node &n, lst &sym_lst) : inherited(n, sym_lst)
+void matrix::read_archive(const archive_node &n, lst &sym_lst)
 {
-       setflag(status_flags::not_shareable);
+       inherited::read_archive(n, sym_lst);
 
        if (!(n.find_unsigned("row", row)) || !(n.find_unsigned("col", col)))
                throw (std::runtime_error("unknown matrix dimensions in archive"));
        m.reserve(row * col);
-       archive_node::archive_node_cit first = n.find_first("m");
-       archive_node::archive_node_cit last = n.find_last("m");
+       // XXX: default ctor inserts a zero element, we need to erase it here.
+       m.pop_back();
+       auto first = n.find_first("m");
+       auto last = n.find_last("m");
        ++last;
-       for (archive_node::archive_node_cit i=first; i<last; ++i) {
+       for (auto i=first; i != last; ++i) {
                ex e;
                n.find_ex_by_loc(i, e, sym_lst);
                m.push_back(e);
        }
 }
+GINAC_BIND_UNARCHIVER(matrix);
 
 void matrix::archive(archive_node &n) const
 {
        inherited::archive(n);
        n.add_unsigned("row", row);
        n.add_unsigned("col", col);
-       exvector::const_iterator i = m.begin(), iend = m.end();
-       while (i != iend) {
-               n.add_ex("m", *i);
-               ++i;
+       for (auto & i : m) {
+               n.add_ex("m", i);
        }
 }
 
-DEFAULT_UNARCHIVE(matrix)
-
 //////////
 // functions overriding virtual functions from base classes
 //////////
@@ -203,28 +227,6 @@ ex & matrix::let_op(size_t i)
        return m[i];
 }
 
-/** Evaluate matrix entry by entry. */
-ex matrix::eval(int level) const
-{
-       // check if we have to do anything at all
-       if ((level==1)&&(flags & status_flags::evaluated))
-               return *this;
-       
-       // emergency break
-       if (level == -max_recursion_level)
-               throw (std::runtime_error("matrix::eval(): recursion limit exceeded"));
-       
-       // eval() entry by entry
-       exvector m2(row*col);
-       --level;
-       for (unsigned r=0; r<row; ++r)
-               for (unsigned c=0; c<col; ++c)
-                       m2[r*col+c] = m[r*col+c].eval(level);
-       
-       return (new matrix(row, col, m2))->setflag(status_flags::dynallocated |
-                                                  status_flags::evaluated);
-}
-
 ex matrix::subs(const exmap & mp, unsigned options) const
 {
        exvector m2(row * col);
@@ -232,14 +234,14 @@ ex matrix::subs(const exmap & mp, unsigned options) const
                for (unsigned c=0; c<col; ++c)
                        m2[r*col+c] = m[r*col+c].subs(mp, options);
 
-       return matrix(row, col, m2).subs_one_level(mp, options);
+       return matrix(row, col, std::move(m2)).subs_one_level(mp, options);
 }
 
 /** Complex conjugate every matrix entry. */
 ex matrix::conjugate() const
 {
-       exvector * ev = 0;
-       for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i) {
+       std::unique_ptr<exvector> ev(nullptr);
+       for (auto i=m.begin(); i!=m.end(); ++i) {
                ex x = i->conjugate();
                if (ev) {
                        ev->push_back(x);
@@ -248,17 +250,15 @@ ex matrix::conjugate() const
                if (are_ex_trivially_equal(x, *i)) {
                        continue;
                }
-               ev = new exvector;
+               ev.reset(new exvector);
                ev->reserve(m.size());
-               for (exvector::const_iterator j=m.begin(); j!=i; ++j) {
+               for (auto j=m.begin(); j!=i; ++j) {
                        ev->push_back(*j);
                }
                ev->push_back(x);
        }
        if (ev) {
-               ex result = matrix(row, col, *ev);
-               delete ev;
-               return result;
+               return matrix(row, col, std::move(*ev));
        }
        return *this;
 }
@@ -267,18 +267,18 @@ ex matrix::real_part() const
 {
        exvector v;
        v.reserve(m.size());
-       for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
-               v.push_back(i->real_part());
-       return matrix(row, col, v);
+       for (auto & i : m)
+               v.push_back(i.real_part());
+       return matrix(row, col, std::move(v));
 }
 
 ex matrix::imag_part() const
 {
        exvector v;
        v.reserve(m.size());
-       for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i)
-               v.push_back(i->imag_part());
-       return matrix(row, col, v);
+       for (auto & i : m)
+               v.push_back(i.imag_part());
+       return matrix(row, col, std::move(v));
 }
 
 // protected
@@ -560,12 +560,11 @@ matrix matrix::add(const matrix & other) const
                throw std::logic_error("matrix::add(): incompatible matrices");
        
        exvector sum(this->m);
-       exvector::iterator i = sum.begin(), end = sum.end();
-       exvector::const_iterator ci = other.m.begin();
-       while (i != end)
-               *i++ += *ci++;
+       auto ci = other.m.begin();
+       for (auto & i : sum)
+               i += *ci++;
        
-       return matrix(row,col,sum);
+       return matrix(row, col, std::move(sum));
 }
 
 
@@ -578,12 +577,11 @@ matrix matrix::sub(const matrix & other) const
                throw std::logic_error("matrix::sub(): incompatible matrices");
        
        exvector dif(this->m);
-       exvector::iterator i = dif.begin(), end = dif.end();
-       exvector::const_iterator ci = other.m.begin();
-       while (i != end)
-               *i++ -= *ci++;
+       auto ci = other.m.begin();
+       for (auto & i : dif)
+               i -= *ci++;
        
-       return matrix(row,col,dif);
+       return matrix(row, col, std::move(dif));
 }
 
 
@@ -606,7 +604,7 @@ matrix matrix::mul(const matrix & other) const
                                prod[r1*other.col+r2] += (m[r1*col+c] * other.m[c*other.col+r2]);
                }
        }
-       return matrix(row, other.col, prod);
+       return matrix(row, other.col, std::move(prod));
 }
 
 
@@ -619,7 +617,7 @@ matrix matrix::mul(const numeric & other) const
                for (unsigned c=0; c<col; ++c)
                        prod[r*col+c] = m[r*col+c] * other;
 
-       return matrix(row, col, prod);
+       return matrix(row, col, std::move(prod));
 }
 
 
@@ -635,7 +633,7 @@ matrix matrix::mul_scalar(const ex & other) const
                for (unsigned c=0; c<col; ++c)
                        prod[r*col+c] = m[r*col+c] * other;
 
-       return matrix(row, col, prod);
+       return matrix(row, col, std::move(prod));
 }
 
 
@@ -721,7 +719,7 @@ matrix matrix::transpose() const
                for (unsigned c=0; c<this->rows(); ++c)
                        trans[r*this->rows()+c] = m[c*this->cols()+r];
        
-       return matrix(this->cols(),this->rows(),trans);
+       return matrix(this->cols(), this->rows(), std::move(trans));
 }
 
 /** Determinant of square matrix.  This routine doesn't actually calculate the
@@ -748,18 +746,16 @@ ex matrix::determinant(unsigned algo) const
        bool numeric_flag = true;
        bool normal_flag = false;
        unsigned sparse_count = 0;  // counts non-zero elements
-       exvector::const_iterator r = m.begin(), rend = m.end();
-       while (r != rend) {
-               if (!r->info(info_flags::numeric))
+       for (auto r : m) {
+               if (!r.info(info_flags::numeric))
                        numeric_flag = false;
                exmap srl;  // symbol replacement list
-               ex rtest = r->to_rational(srl);
+               ex rtest = r.to_rational(srl);
                if (!rtest.is_zero())
                        ++sparse_count;
                if (!rtest.info(info_flags::crational_polynomial) &&
-                        rtest.info(info_flags::rational_function))
+                    rtest.info(info_flags::rational_function))
                        normal_flag = true;
-               ++r;
        }
        
        // Here is the heuristics in case this routine has to decide:
@@ -842,23 +838,22 @@ ex matrix::determinant(unsigned algo) const
                        }
                        std::sort(c_zeros.begin(),c_zeros.end());
                        std::vector<unsigned> pre_sort;
-                       for (std::vector<uintpair>::const_iterator i=c_zeros.begin(); i!=c_zeros.end(); ++i)
-                               pre_sort.push_back(i->second);
+                       for (auto & i : c_zeros)
+                               pre_sort.push_back(i.second);
                        std::vector<unsigned> pre_sort_test(pre_sort); // permutation_sign() modifies the vector so we make a copy here
                        int sign = permutation_sign(pre_sort_test.begin(), pre_sort_test.end());
                        exvector result(row*col);  // represents sorted matrix
                        unsigned c = 0;
-                       for (std::vector<unsigned>::const_iterator i=pre_sort.begin();
-                                i!=pre_sort.end();
-                                ++i,++c) {
+                       for (auto & it : pre_sort) {
                                for (unsigned r=0; r<row; ++r)
-                                       result[r*col+c] = m[r*col+(*i)];
+                                       result[r*col+c] = m[r*col+it];
+                               ++c;
                        }
                        
                        if (normal_flag)
-                               return (sign*matrix(row,col,result).determinant_minor()).normal();
+                               return (sign*matrix(row, col, std::move(result)).determinant_minor()).normal();
                        else
-                               return sign*matrix(row,col,result).determinant_minor();
+                               return sign*matrix(row, col, std::move(result)).determinant_minor();
                }
        }
 }
@@ -888,7 +883,7 @@ ex matrix::trace() const
 
 
 /** Characteristic Polynomial.  Following mathematica notation the
- *  characteristic polynomial of a matrix M is defined as the determiant of
+ *  characteristic polynomial of a matrix M is defined as the determinant of
  *  (M - lambda * 1) where 1 stands for the unit matrix of the same dimension
  *  as M.  Note that some CASs define it with a sign inside the determinant
  *  which gives rise to an overall sign if the dimension is odd.  This method
@@ -904,11 +899,11 @@ ex matrix::charpoly(const ex & lambda) const
                throw (std::logic_error("matrix::charpoly(): matrix not square"));
        
        bool numeric_flag = true;
-       exvector::const_iterator r = m.begin(), rend = m.end();
-       while (r!=rend && numeric_flag==true) {
-               if (!r->info(info_flags::numeric))
+       for (auto & r : m) {
+               if (!r.info(info_flags::numeric)) {
                        numeric_flag = false;
-               ++r;
+                       break;
+               }
        }
        
        // The pure numeric case is traditionally rather common.  Hence, it is
@@ -1019,11 +1014,11 @@ matrix matrix::solve(const matrix & vars,
        
        // Gather some statistical information about the augmented matrix:
        bool numeric_flag = true;
-       exvector::const_iterator r = aug.m.begin(), rend = aug.m.end();
-       while (r!=rend && numeric_flag==true) {
-               if (!r->info(info_flags::numeric))
+       for (auto & r : aug.m) {
+               if (!r.info(info_flags::numeric)) {
                        numeric_flag = false;
-               ++r;
+                       break;
+               }
        }
        
        // Here is the heuristics in case this routine has to decide:
@@ -1118,7 +1113,7 @@ unsigned matrix::rank() const
  *  more than once.  According to W.M.Gentleman and S.C.Johnson this algorithm
  *  is better than elimination schemes for matrices of sparse multivariate
  *  polynomials and also for matrices of dense univariate polynomials if the
- *  matrix' dimesion is larger than 7.
+ *  matrix' dimension is larger than 7.
  *
  *  @return the determinant as a new expression (in expanded form)
  *  @see matrix::determinant() */
@@ -1224,9 +1219,8 @@ ex matrix::determinant_minor() const
                                for (unsigned j=fc; j<n-c; ++j)
                                        Pkey[j] = Pkey[j-1]+1;
                } while(fc);
-               // next column, so change the role of A and B:
-               A.swap(B);
-               B.clear();
+               // next column, clear B and change the role of A and B:
+               A = std::move(B);
        }
        
        return det;
@@ -1403,11 +1397,9 @@ int matrix::fraction_free_elimination(const bool det)
        matrix tmp_n(*this);
        matrix tmp_d(m,n);  // for denominators, if needed
        exmap srl;  // symbol replacement list
-       exvector::const_iterator cit = this->m.begin(), citend = this->m.end();
-       exvector::iterator tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
-       while (cit != citend) {
-               ex nd = cit->normal().to_rational(srl).numer_denom();
-               ++cit;
+       auto tmp_n_it = tmp_n.m.begin(), tmp_d_it = tmp_d.m.begin();
+       for (auto & it : this->m) {
+               ex nd = it.normal().to_rational(srl).numer_denom();
                *tmp_n_it++ = nd.op(0);
                *tmp_d_it++ = nd.op(1);
        }
@@ -1476,11 +1468,10 @@ int matrix::fraction_free_elimination(const bool det)
        }
 
        // repopulate *this matrix:
-       exvector::iterator it = this->m.begin(), itend = this->m.end();
        tmp_n_it = tmp_n.m.begin();
        tmp_d_it = tmp_d.m.begin();
-       while (it != itend)
-               *it++ = ((*tmp_n_it++)/(*tmp_d_it++)).subs(srl, subs_options::no_pattern);
+       for (auto & it : this->m)
+               it = ((*tmp_n_it++)/(*tmp_d_it++)).subs(srl, subs_options::no_pattern);
        
        return sign;
 }
@@ -1496,7 +1487,7 @@ int matrix::fraction_free_elimination(const bool det)
  *  @param co is the column to be inspected
  *  @param symbolic signal if we want the first non-zero element to be pivoted
  *  (true) or the one with the largest absolute value (false).
- *  @return 0 if no interchange occured, -1 if all are zero (usually signaling
+ *  @return 0 if no interchange occurred, -1 if all are zero (usually signaling
  *  a degeneracy) and positive integer k means that rows ro and k were swapped.
  */
 int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
@@ -1541,34 +1532,34 @@ int matrix::pivot(unsigned ro, unsigned co, bool symbolic)
  */
 bool matrix::is_zero_matrix() const
 {
-       for (exvector::const_iterator i=m.begin(); i!=m.end(); ++i) 
-               if(!(i->is_zero()))
+       for (auto & i : m)
+               if (!i.is_zero())
                        return false;
        return true;
 }
 
 ex lst_to_matrix(const lst & l)
 {
-       lst::const_iterator itr, itc;
-
        // Find number of rows and columns
        size_t rows = l.nops(), cols = 0;
-       for (itr = l.begin(); itr != l.end(); ++itr) {
-               if (!is_a<lst>(*itr))
+       for (auto & itr : l) {
+               if (!is_a<lst>(itr))
                        throw (std::invalid_argument("lst_to_matrix: argument must be a list of lists"));
-               if (itr->nops() > cols)
-                       cols = itr->nops();
+               if (itr.nops() > cols)
+                       cols = itr.nops();
        }
 
        // Allocate and fill matrix
-       matrix &M = *new matrix(rows, cols);
-       M.setflag(status_flags::dynallocated);
-
-       unsigned i;
-       for (itr = l.begin(), i = 0; itr != l.end(); ++itr, ++i) {
-               unsigned j;
-               for (itc = ex_to<lst>(*itr).begin(), j = 0; itc != ex_to<lst>(*itr).end(); ++itc, ++j)
-                       M(i, j) = *itc;
+       matrix & M = dynallocate<matrix>(rows, cols);
+
+       unsigned i = 0;
+       for (auto & itr : l) {
+               unsigned j = 0;
+               for (auto & itc : ex_to<lst>(itr)) {
+                       M(i, j) = itc;
+                       ++j;
+               }
+               ++i;
        }
 
        return M;
@@ -1576,24 +1567,40 @@ ex lst_to_matrix(const lst & l)
 
 ex diag_matrix(const lst & l)
 {
-       lst::const_iterator it;
        size_t dim = l.nops();
 
        // Allocate and fill matrix
-       matrix &M = *new matrix(dim, dim);
-       M.setflag(status_flags::dynallocated);
+       matrix & M = dynallocate<matrix>(dim, dim);
+
+       unsigned i = 0;
+       for (auto & it : l) {
+               M(i, i) = it;
+               ++i;
+       }
+
+       return M;
+}
+
+ex diag_matrix(std::initializer_list<ex> l)
+{
+       size_t dim = l.size();
+
+       // Allocate and fill matrix
+       matrix & M = dynallocate<matrix>(dim, dim);
 
-       unsigned i;
-       for (it = l.begin(), i = 0; it != l.end(); ++it, ++i)
-               M(i, i) = *it;
+       unsigned i = 0;
+       for (auto & it : l) {
+               M(i, i) = it;
+               ++i;
+       }
 
        return M;
 }
 
 ex unit_matrix(unsigned r, unsigned c)
 {
-       matrix &Id = *new matrix(r, c);
-       Id.setflag(status_flags::dynallocated);
+       matrix & Id = dynallocate<matrix>(r, c);
+       Id.setflag(status_flags::evaluated);
        for (unsigned i=0; i<r && i<c; i++)
                Id(i,i) = _ex1;
 
@@ -1602,8 +1609,8 @@ ex unit_matrix(unsigned r, unsigned c)
 
 ex symbolic_matrix(unsigned r, unsigned c, const std::string & base_name, const std::string & tex_base_name)
 {
-       matrix &M = *new matrix(r, c);
-       M.setflag(status_flags::dynallocated | status_flags::evaluated);
+       matrix & M = dynallocate<matrix>(r, c);
+       M.setflag(status_flags::evaluated);
 
        bool long_format = (r > 10 || c > 10);
        bool single_row = (r == 1 || c == 1);
@@ -1644,8 +1651,8 @@ ex reduced_matrix(const matrix& m, unsigned r, unsigned c)
 
        const unsigned rows = m.rows()-1;
        const unsigned cols = m.cols()-1;
-       matrix &M = *new matrix(rows, cols);
-       M.setflag(status_flags::dynallocated | status_flags::evaluated);
+       matrix & M = dynallocate<matrix>(rows, cols);
+       M.setflag(status_flags::evaluated);
 
        unsigned ro = 0;
        unsigned ro2 = 0;
@@ -1673,8 +1680,8 @@ ex sub_matrix(const matrix&m, unsigned r, unsigned nr, unsigned c, unsigned nc)
        if (r+nr>m.rows() || c+nc>m.cols())
                throw std::runtime_error("sub_matrix(): index out of bounds");
 
-       matrix &M = *new matrix(nr, nc);
-       M.setflag(status_flags::dynallocated | status_flags::evaluated);
+       matrix & M = dynallocate<matrix>(nr, nc);
+       M.setflag(status_flags::evaluated);
 
        for (unsigned ro=0; ro<nr; ++ro) {
                for (unsigned co=0; co<nc; ++co) {