]> www.ginac.de Git - ginac.git/blobdiff - ginac/normal.cpp
Univariate Hensel lifting now uses upoly.
[ginac.git] / ginac / normal.cpp
index 8d9f462140dd42fdcdefa4e00e50eb602b85e049..09773d37489693a4854964f129002abb186b2c5b 100644 (file)
@@ -6,7 +6,7 @@
  *  computation, square-free factorization and rational function normalization. */
 
 /*
- *  GiNaC Copyright (C) 1999-2003 Johannes Gutenberg University Mainz, Germany
+ *  GiNaC Copyright (C) 1999-2008 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
@@ -20,7 +20,7 @@
  *
  *  You should have received a copy of the GNU General Public License
  *  along with this program; if not, write to the Free Software
- *  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
+ *  Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
  */
 
 #include <algorithm>
@@ -84,17 +84,17 @@ static struct _stat_print {
 #endif
 
 
-/** Return pointer to first symbol found in expression.  Due to GiNaC´s
+/** Return pointer to first symbol found in expression.  Due to GiNaC's
  *  internal ordering of terms, it may not be obvious which symbol this
  *  function returns for a given expression.
  *
  *  @param e  expression to search
- *  @param x  pointer to first symbol found (returned)
+ *  @param x  first symbol found (returned)
  *  @return "false" if no symbol was found, "true" otherwise */
-static bool get_first_symbol(const ex &e, const symbol *&x)
+static bool get_first_symbol(const ex &e, ex &x)
 {
        if (is_a<symbol>(e)) {
-               x = &ex_to<symbol>(e);
+               x = e;
                return true;
        } else if (is_exactly_a<add>(e) || is_exactly_a<mul>(e)) {
                for (size_t i=0; i<e.nops(); i++)
@@ -119,8 +119,8 @@ static bool get_first_symbol(const ex &e, const symbol *&x)
  *
  *  @see get_symbol_stats */
 struct sym_desc {
-       /** Pointer to symbol */
-       const symbol *sym;
+       /** Reference to symbol */
+       ex sym;
 
        /** Highest degree of symbol in polynomial "a" */
        int deg_a;
@@ -154,11 +154,11 @@ struct sym_desc {
 typedef std::vector<sym_desc> sym_desc_vec;
 
 // Add symbol the sym_desc_vec (used internally by get_symbol_stats())
-static void add_symbol(const symbol *s, sym_desc_vec &v)
+static void add_symbol(const ex &s, sym_desc_vec &v)
 {
        sym_desc_vec::const_iterator it = v.begin(), itend = v.end();
        while (it != itend) {
-               if (it->sym->compare(*s) == 0)  // If it's already in there, don't add it a second time
+               if (it->sym.is_equal(s))  // If it's already in there, don't add it a second time
                        return;
                ++it;
        }
@@ -171,7 +171,7 @@ static void add_symbol(const symbol *s, sym_desc_vec &v)
 static void collect_symbols(const ex &e, sym_desc_vec &v)
 {
        if (is_a<symbol>(e)) {
-               add_symbol(&ex_to<symbol>(e), v);
+               add_symbol(e, v);
        } else if (is_exactly_a<add>(e) || is_exactly_a<mul>(e)) {
                for (size_t i=0; i<e.nops(); i++)
                        collect_symbols(e.op(i), v);
@@ -198,14 +198,14 @@ static void get_symbol_stats(const ex &a, const ex &b, sym_desc_vec &v)
        collect_symbols(b.eval(), v);
        sym_desc_vec::iterator it = v.begin(), itend = v.end();
        while (it != itend) {
-               int deg_a = a.degree(*(it->sym));
-               int deg_b = b.degree(*(it->sym));
+               int deg_a = a.degree(it->sym);
+               int deg_b = b.degree(it->sym);
                it->deg_a = deg_a;
                it->deg_b = deg_b;
                it->max_deg = std::max(deg_a, deg_b);
-               it->max_lcnops = std::max(a.lcoeff(*(it->sym)).nops(), b.lcoeff(*(it->sym)).nops());
-               it->ldeg_a = a.ldegree(*(it->sym));
-               it->ldeg_b = b.ldegree(*(it->sym));
+               it->max_lcnops = std::max(a.lcoeff(it->sym).nops(), b.lcoeff(it->sym).nops());
+               it->ldeg_a = a.ldegree(it->sym);
+               it->ldeg_b = b.ldegree(it->sym);
                ++it;
        }
        std::sort(v.begin(), v.end());
@@ -214,8 +214,8 @@ static void get_symbol_stats(const ex &a, const ex &b, sym_desc_vec &v)
        std::clog << "Symbols:\n";
        it = v.begin(); itend = v.end();
        while (it != itend) {
-               std::clog << " " << *it->sym << ": deg_a=" << it->deg_a << ", deg_b=" << it->deg_b << ", ldeg_a=" << it->ldeg_a << ", ldeg_b=" << it->ldeg_b << ", max_deg=" << it->max_deg << ", max_lcnops=" << it->max_lcnops << endl;
-               std::clog << "  lcoeff_a=" << a.lcoeff(*(it->sym)) << ", lcoeff_b=" << b.lcoeff(*(it->sym)) << endl;
+               std::clog << " " << it->sym << ": deg_a=" << it->deg_a << ", deg_b=" << it->deg_b << ", ldeg_a=" << it->ldeg_a << ", ldeg_b=" << it->ldeg_b << ", max_deg=" << it->max_deg << ", max_lcnops=" << it->max_lcnops << endl;
+               std::clog << "  lcoeff_a=" << a.lcoeff(it->sym) << ", lcoeff_b=" << b.lcoeff(it->sym) << endl;
                ++it;
        }
 #endif
@@ -233,14 +233,14 @@ static numeric lcmcoeff(const ex &e, const numeric &l)
        if (e.info(info_flags::rational))
                return lcm(ex_to<numeric>(e).denom(), l);
        else if (is_exactly_a<add>(e)) {
-               numeric c = _num1;
+               numeric c = *_num1_p;
                for (size_t i=0; i<e.nops(); i++)
                        c = lcmcoeff(e.op(i), c);
                return lcm(c, l);
        } else if (is_exactly_a<mul>(e)) {
-               numeric c = _num1;
+               numeric c = *_num1_p;
                for (size_t i=0; i<e.nops(); i++)
-                       c *= lcmcoeff(e.op(i), _num1);
+                       c *= lcmcoeff(e.op(i), *_num1_p);
                return lcm(c, l);
        } else if (is_exactly_a<power>(e)) {
                if (is_a<symbol>(e.op(0)))
@@ -260,7 +260,7 @@ static numeric lcmcoeff(const ex &e, const numeric &l)
  *  @return LCM of denominators of coefficients */
 static numeric lcm_of_coefficients_denominators(const ex &e)
 {
-       return lcmcoeff(e, _num1);
+       return lcmcoeff(e, *_num1_p);
 }
 
 /** Bring polynomial from Q[X] to Z[X] by multiplying in the previously
@@ -273,9 +273,9 @@ static ex multiply_lcm(const ex &e, const numeric &lcm)
        if (is_exactly_a<mul>(e)) {
                size_t num = e.nops();
                exvector v; v.reserve(num + 1);
-               numeric lcm_accum = _num1;
+               numeric lcm_accum = *_num1_p;
                for (size_t i=0; i<num; i++) {
-                       numeric op_lcm = lcmcoeff(e.op(i), _num1);
+                       numeric op_lcm = lcmcoeff(e.op(i), *_num1_p);
                        v.push_back(multiply_lcm(e.op(i), op_lcm));
                        lcm_accum *= op_lcm;
                }
@@ -298,9 +298,10 @@ static ex multiply_lcm(const ex &e, const numeric &lcm)
 
 
 /** Compute the integer content (= GCD of all numeric coefficients) of an
- *  expanded polynomial.
+ *  expanded polynomial. For a polynomial with rational coefficients, this
+ *  returns g/l where g is the GCD of the coefficients' numerators and l
+ *  is the LCM of the coefficients' denominators.
  *
- *  @param e  expanded polynomial
  *  @return integer content */
 numeric ex::integer_content() const
 {
@@ -309,7 +310,7 @@ numeric ex::integer_content() const
 
 numeric basic::integer_content() const
 {
-       return _num1;
+       return *_num1_p;
 }
 
 numeric numeric::integer_content() const
@@ -321,16 +322,18 @@ numeric add::integer_content() const
 {
        epvector::const_iterator it = seq.begin();
        epvector::const_iterator itend = seq.end();
-       numeric c = _num0;
+       numeric c = *_num0_p, l = *_num1_p;
        while (it != itend) {
                GINAC_ASSERT(!is_exactly_a<numeric>(it->rest));
                GINAC_ASSERT(is_exactly_a<numeric>(it->coeff));
-               c = gcd(ex_to<numeric>(it->coeff), c);
+               c = gcd(ex_to<numeric>(it->coeff).numer(), c);
+               l = lcm(ex_to<numeric>(it->coeff).denom(), l);
                it++;
        }
        GINAC_ASSERT(is_exactly_a<numeric>(overall_coeff));
-       c = gcd(ex_to<numeric>(overall_coeff),c);
-       return c;
+       c = gcd(ex_to<numeric>(overall_coeff).numer(), c);
+       l = lcm(ex_to<numeric>(overall_coeff).denom(), l);
+       return c/l;
 }
 
 numeric mul::integer_content() const
@@ -361,7 +364,7 @@ numeric mul::integer_content() const
  *  @param check_args  check whether a and b are polynomials with rational
  *         coefficients (defaults to "true")
  *  @return quotient of a and b in Q[x] */
-ex quo(const ex &a, const ex &b, const symbol &x, bool check_args)
+ex quo(const ex &a, const ex &b, const ex &x, bool check_args)
 {
        if (b.is_zero())
                throw(std::overflow_error("quo: division by zero"));
@@ -411,7 +414,7 @@ ex quo(const ex &a, const ex &b, const symbol &x, bool check_args)
  *  @param check_args  check whether a and b are polynomials with rational
  *         coefficients (defaults to "true")
  *  @return remainder of a(x) and b(x) in Q[x] */
-ex rem(const ex &a, const ex &b, const symbol &x, bool check_args)
+ex rem(const ex &a, const ex &b, const ex &x, bool check_args)
 {
        if (b.is_zero())
                throw(std::overflow_error("rem: division by zero"));
@@ -460,7 +463,7 @@ ex rem(const ex &a, const ex &b, const symbol &x, bool check_args)
  *  @param a rational function in x
  *  @param x a is a function of x
  *  @return decomposed function. */
-ex decomp_rational(const ex &a, const symbol &x)
+ex decomp_rational(const ex &a, const ex &x)
 {
        ex nd = numer_denom(a);
        ex numer = nd.op(0), denom = nd.op(1);
@@ -480,7 +483,7 @@ ex decomp_rational(const ex &a, const symbol &x)
  *  @param check_args  check whether a and b are polynomials with rational
  *         coefficients (defaults to "true")
  *  @return pseudo-remainder of a(x) and b(x) in Q[x] */
-ex prem(const ex &a, const ex &b, const symbol &x, bool check_args)
+ex prem(const ex &a, const ex &b, const ex &x, bool check_args)
 {
        if (b.is_zero())
                throw(std::overflow_error("prem: division by zero"));
@@ -532,7 +535,7 @@ ex prem(const ex &a, const ex &b, const symbol &x, bool check_args)
  *  @param check_args  check whether a and b are polynomials with rational
  *         coefficients (defaults to "true")
  *  @return sparse pseudo-remainder of a(x) and b(x) in Q[x] */
-ex sprem(const ex &a, const ex &b, const symbol &x, bool check_args)
+ex sprem(const ex &a, const ex &b, const ex &x, bool check_args)
 {
        if (b.is_zero())
                throw(std::overflow_error("prem: division by zero"));
@@ -607,36 +610,102 @@ bool divide(const ex &a, const ex &b, ex &q, bool check_args)
                throw(std::invalid_argument("divide: arguments must be polynomials over the rationals"));
 
        // Find first symbol
-       const symbol *x;
+       ex x;
        if (!get_first_symbol(a, x) && !get_first_symbol(b, x))
                throw(std::invalid_argument("invalid expression in divide()"));
 
+       // Try to avoid expanding partially factored expressions.
+       if (is_exactly_a<mul>(b)) {
+       // Divide sequentially by each term
+               ex rem_new, rem_old = a;
+               for (size_t i=0; i < b.nops(); i++) {
+                       if (! divide(rem_old, b.op(i), rem_new, false))
+                               return false;
+                       rem_old = rem_new;
+               }
+               q = rem_new;
+               return true;
+       } else if (is_exactly_a<power>(b)) {
+               const ex& bb(b.op(0));
+               int exp_b = ex_to<numeric>(b.op(1)).to_int();
+               ex rem_new, rem_old = a;
+               for (int i=exp_b; i>0; i--) {
+                       if (! divide(rem_old, bb, rem_new, false))
+                               return false;
+                       rem_old = rem_new;
+               }
+               q = rem_new;
+               return true;
+       } 
+       
+       if (is_exactly_a<mul>(a)) {
+               // Divide sequentially each term. If some term in a is divisible 
+               // by b we are done... and if not, we can't really say anything.
+               size_t i;
+               ex rem_i;
+               bool divisible_p = false;
+               for (i=0; i < a.nops(); ++i) {
+                       if (divide(a.op(i), b, rem_i, false)) {
+                               divisible_p = true;
+                               break;
+                       }
+               }
+               if (divisible_p) {
+                       exvector resv;
+                       resv.reserve(a.nops());
+                       for (size_t j=0; j < a.nops(); j++) {
+                               if (j==i)
+                                       resv.push_back(rem_i);
+                               else
+                                       resv.push_back(a.op(j));
+                       }
+                       q = (new mul(resv))->setflag(status_flags::dynallocated);
+                       return true;
+               }
+       } else if (is_exactly_a<power>(a)) {
+               // The base itself might be divisible by b, in that case we don't
+               // need to expand a
+               const ex& ab(a.op(0));
+               int a_exp = ex_to<numeric>(a.op(1)).to_int();
+               ex rem_i;
+               if (divide(ab, b, rem_i, false)) {
+                       q = rem_i*power(ab, a_exp - 1);
+                       return true;
+               }
+               for (int i=2; i < a_exp; i++) {
+                       if (divide(power(ab, i), b, rem_i, false)) {
+                               q = rem_i*power(ab, a_exp - i);
+                               return true;
+                       }
+               } // ... so we *really* need to expand expression.
+       }
+       
        // Polynomial long division (recursive)
        ex r = a.expand();
        if (r.is_zero()) {
                q = _ex0;
                return true;
        }
-       int bdeg = b.degree(*x);
-       int rdeg = r.degree(*x);
-       ex blcoeff = b.expand().coeff(*x, bdeg);
+       int bdeg = b.degree(x);
+       int rdeg = r.degree(x);
+       ex blcoeff = b.expand().coeff(x, bdeg);
        bool blcoeff_is_numeric = is_exactly_a<numeric>(blcoeff);
        exvector v; v.reserve(std::max(rdeg - bdeg + 1, 0));
        while (rdeg >= bdeg) {
-               ex term, rcoeff = r.coeff(*x, rdeg);
+               ex term, rcoeff = r.coeff(x, rdeg);
                if (blcoeff_is_numeric)
                        term = rcoeff / blcoeff;
                else
                        if (!divide(rcoeff, blcoeff, term, false))
                                return false;
-               term *= power(*x, rdeg - bdeg);
+               term *= power(x, rdeg - bdeg);
                v.push_back(term);
                r -= (term * b).expand();
                if (r.is_zero()) {
                        q = (new add(v))->setflag(status_flags::dynallocated);
                        return true;
                }
-               rdeg = r.degree(*x);
+               rdeg = r.degree(x);
        }
        return false;
 }
@@ -665,7 +734,7 @@ typedef std::map<ex2, exbool, ex2_less> ex2_exbool_remember;
 /** Exact polynomial division of a(X) by b(X) in Z[X].
  *  This functions works like divide() but the input and output polynomials are
  *  in Z[X] instead of Q[X] (i.e. they have integer coefficients). Unlike
- *  divide(), it doesn´t check whether the input polynomials really are integer
+ *  divide(), it doesn't check whether the input polynomials really are integer
  *  polynomials, so be careful of what you pass in. Also, you have to run
  *  get_symbol_stats() over the input polynomials before calling this function
  *  and pass an iterator to the first element of the sym_desc vector. This
@@ -711,11 +780,36 @@ static bool divide_in_z(const ex &a, const ex &b, ex &q, sym_desc_vec::const_ite
        }
 #endif
 
+       if (is_exactly_a<power>(b)) {
+               const ex& bb(b.op(0));
+               ex qbar = a;
+               int exp_b = ex_to<numeric>(b.op(1)).to_int();
+               for (int i=exp_b; i>0; i--) {
+                       if (!divide_in_z(qbar, bb, q, var))
+                               return false;
+                       qbar = q;
+               }
+               return true;
+       }
+
+       if (is_exactly_a<mul>(b)) {
+               ex qbar = a;
+               for (const_iterator itrb = b.begin(); itrb != b.end(); ++itrb) {
+                       sym_desc_vec sym_stats;
+                       get_symbol_stats(a, *itrb, sym_stats);
+                       if (!divide_in_z(qbar, *itrb, q, sym_stats.begin()))
+                               return false;
+
+                       qbar = q;
+               }
+               return true;
+       }
+
        // Main symbol
-       const symbol *x = var->sym;
+       const ex &x = var->sym;
 
        // Compare degrees
-       int adeg = a.degree(*x), bdeg = b.degree(*x);
+       int adeg = a.degree(x), bdeg = b.degree(x);
        if (bdeg > adeg)
                return false;
 
@@ -727,24 +821,24 @@ static bool divide_in_z(const ex &a, const ex &b, ex &q, sym_desc_vec::const_ite
        // Compute values at evaluation points 0..adeg
        vector<numeric> alpha; alpha.reserve(adeg + 1);
        exvector u; u.reserve(adeg + 1);
-       numeric point = _num0;
+       numeric point = *_num0_p;
        ex c;
        for (i=0; i<=adeg; i++) {
-               ex bs = b.subs(*x == point, subs_options::no_pattern);
+               ex bs = b.subs(x == point, subs_options::no_pattern);
                while (bs.is_zero()) {
-                       point += _num1;
-                       bs = b.subs(*x == point, subs_options::no_pattern);
+                       point += *_num1_p;
+                       bs = b.subs(x == point, subs_options::no_pattern);
                }
-               if (!divide_in_z(a.subs(*x == point, subs_options::no_pattern), bs, c, var+1))
+               if (!divide_in_z(a.subs(x == point, subs_options::no_pattern), bs, c, var+1))
                        return false;
                alpha.push_back(point);
                u.push_back(c);
-               point += _num1;
+               point += *_num1_p;
        }
 
        // Compute inverses
        vector<numeric> rcp; rcp.reserve(adeg + 1);
-       rcp.push_back(_num0);
+       rcp.push_back(*_num0_p);
        for (k=1; k<=adeg; k++) {
                numeric product = alpha[k] - alpha[0];
                for (i=1; i<k; i++)
@@ -765,9 +859,9 @@ static bool divide_in_z(const ex &a, const ex &b, ex &q, sym_desc_vec::const_ite
        // Convert from Newton form to standard form
        c = v[adeg];
        for (k=adeg-1; k>=0; k--)
-               c = c * (*x - alpha[k]) + v[k];
+               c = c * (x - alpha[k]) + v[k];
 
-       if (c.degree(*x) == (adeg - bdeg)) {
+       if (c.degree(x) == (adeg - bdeg)) {
                q = c.expand();
                return true;
        } else
@@ -781,13 +875,13 @@ static bool divide_in_z(const ex &a, const ex &b, ex &q, sym_desc_vec::const_ite
                return true;
        int rdeg = adeg;
        ex eb = b.expand();
-       ex blcoeff = eb.coeff(*x, bdeg);
+       ex blcoeff = eb.coeff(x, bdeg);
        exvector v; v.reserve(std::max(rdeg - bdeg + 1, 0));
        while (rdeg >= bdeg) {
-               ex term, rcoeff = r.coeff(*x, rdeg);
+               ex term, rcoeff = r.coeff(x, rdeg);
                if (!divide_in_z(rcoeff, blcoeff, term, var+1))
                        break;
-               term = (term * power(*x, rdeg - bdeg)).expand();
+               term = (term * power(x, rdeg - bdeg)).expand();
                v.push_back(term);
                r -= (term * eb).expand();
                if (r.is_zero()) {
@@ -797,7 +891,7 @@ static bool divide_in_z(const ex &a, const ex &b, ex &q, sym_desc_vec::const_ite
 #endif
                        return true;
                }
-               rdeg = r.degree(*x);
+               rdeg = r.degree(x);
        }
 #if USE_REMEMBER
        dr_remember[ex2(a, b)] = exbool(q, false);
@@ -813,21 +907,21 @@ static bool divide_in_z(const ex &a, const ex &b, ex &q, sym_desc_vec::const_ite
  */
 
 /** Compute unit part (= sign of leading coefficient) of a multivariate
- *  polynomial in Z[x]. The product of unit part, content part, and primitive
+ *  polynomial in Q[x]. The product of unit part, content part, and primitive
  *  part is the polynomial itself.
  *
- *  @param x  variable in which to compute the unit part
+ *  @param x  main variable
  *  @return unit part
- *  @see ex::content, ex::primpart */
-ex ex::unit(const symbol &x) const
+ *  @see ex::content, ex::primpart, ex::unitcontprim */
+ex ex::unit(const ex &x) const
 {
        ex c = expand().lcoeff(x);
        if (is_exactly_a<numeric>(c))
-               return c < _ex0 ? _ex_1 : _ex1;
+               return c.info(info_flags::negative) ?_ex_1 : _ex1;
        else {
-               const symbol *y;
+               ex y;
                if (get_first_symbol(c, y))
-                       return c.unit(*y);
+                       return c.unit(y);
                else
                        throw(std::invalid_argument("invalid expression in unit()"));
        }
@@ -835,82 +929,72 @@ ex ex::unit(const symbol &x) const
 
 
 /** Compute content part (= unit normal GCD of all coefficients) of a
- *  multivariate polynomial in Z[x].  The product of unit part, content part,
+ *  multivariate polynomial in Q[x]. The product of unit part, content part,
  *  and primitive part is the polynomial itself.
  *
- *  @param x  variable in which to compute the content part
+ *  @param x  main variable
  *  @return content part
- *  @see ex::unit, ex::primpart */
-ex ex::content(const symbol &x) const
+ *  @see ex::unit, ex::primpart, ex::unitcontprim */
+ex ex::content(const ex &x) const
 {
-       if (is_zero())
-               return _ex0;
        if (is_exactly_a<numeric>(*this))
                return info(info_flags::negative) ? -*this : *this;
+
        ex e = expand();
        if (e.is_zero())
                return _ex0;
 
-       // First, try the integer content
+       // First, divide out the integer content (which we can calculate very efficiently).
+       // If the leading coefficient of the quotient is an integer, we are done.
        ex c = e.integer_content();
        ex r = e / c;
-       ex lcoeff = r.lcoeff(x);
+       int deg = r.degree(x);
+       ex lcoeff = r.coeff(x, deg);
        if (lcoeff.info(info_flags::integer))
                return c;
 
        // GCD of all coefficients
-       int deg = e.degree(x);
-       int ldeg = e.ldegree(x);
+       int ldeg = r.ldegree(x);
        if (deg == ldeg)
-               return e.lcoeff(x) / e.unit(x);
-       c = _ex0;
+               return lcoeff * c / lcoeff.unit(x);
+       ex cont = _ex0;
        for (int i=ldeg; i<=deg; i++)
-               c = gcd(e.coeff(x, i), c, NULL, NULL, false);
-       return c;
+               cont = gcd(r.coeff(x, i), cont, NULL, NULL, false);
+       return cont * c;
 }
 
 
-/** Compute primitive part of a multivariate polynomial in Z[x].
- *  The product of unit part, content part, and primitive part is the
- *  polynomial itself.
+/** Compute primitive part of a multivariate polynomial in Q[x]. The result
+ *  will be a unit-normal polynomial with a content part of 1. The product
+ *  of unit part, content part, and primitive part is the polynomial itself.
  *
- *  @param x  variable in which to compute the primitive part
+ *  @param x  main variable
  *  @return primitive part
- *  @see ex::unit, ex::content */
-ex ex::primpart(const symbol &x) const
+ *  @see ex::unit, ex::content, ex::unitcontprim */
+ex ex::primpart(const ex &x) const
 {
-       if (is_zero())
-               return _ex0;
-       if (is_exactly_a<numeric>(*this))
-               return _ex1;
-
-       ex c = content(x);
-       if (c.is_zero())
-               return _ex0;
-       ex u = unit(x);
-       if (is_exactly_a<numeric>(c))
-               return *this / (c * u);
-       else
-               return quo(*this, c * u, x, false);
+       // We need to compute the unit and content anyway, so call unitcontprim()
+       ex u, c, p;
+       unitcontprim(x, u, c, p);
+       return p;
 }
 
 
-/** Compute primitive part of a multivariate polynomial in Z[x] when the
+/** Compute primitive part of a multivariate polynomial in Q[x] when the
  *  content part is already known. This function is faster in computing the
  *  primitive part than the previous function.
  *
- *  @param x  variable in which to compute the primitive part
+ *  @param x  main variable
  *  @param c  previously computed content part
  *  @return primitive part */
-ex ex::primpart(const symbol &x, const ex &c) const
+ex ex::primpart(const ex &x, const ex &c) const
 {
-       if (is_zero())
-               return _ex0;
-       if (c.is_zero())
+       if (is_zero() || c.is_zero())
                return _ex0;
        if (is_exactly_a<numeric>(*this))
                return _ex1;
 
+       // Divide by unit and content to get primitive part
        ex u = unit(x);
        if (is_exactly_a<numeric>(c))
                return *this / (c * u);
@@ -919,6 +1003,61 @@ ex ex::primpart(const symbol &x, const ex &c) const
 }
 
 
+/** Compute unit part, content part, and primitive part of a multivariate
+ *  polynomial in Q[x]. The product of the three parts is the polynomial
+ *  itself.
+ *
+ *  @param x  main variable
+ *  @param u  unit part (returned)
+ *  @param c  content part (returned)
+ *  @param p  primitive part (returned)
+ *  @see ex::unit, ex::content, ex::primpart */
+void ex::unitcontprim(const ex &x, ex &u, ex &c, ex &p) const
+{
+       // Quick check for zero (avoid expanding)
+       if (is_zero()) {
+               u = _ex1;
+               c = p = _ex0;
+               return;
+       }
+
+       // Special case: input is a number
+       if (is_exactly_a<numeric>(*this)) {
+               if (info(info_flags::negative)) {
+                       u = _ex_1;
+                       c = abs(ex_to<numeric>(*this));
+               } else {
+                       u = _ex1;
+                       c = *this;
+               }
+               p = _ex1;
+               return;
+       }
+
+       // Expand input polynomial
+       ex e = expand();
+       if (e.is_zero()) {
+               u = _ex1;
+               c = p = _ex0;
+               return;
+       }
+
+       // Compute unit and content
+       u = unit(x);
+       c = content(x);
+
+       // Divide by unit and content to get primitive part
+       if (c.is_zero()) {
+               p = _ex0;
+               return;
+       }
+       if (is_exactly_a<numeric>(c))
+               p = *this / (c * u);
+       else
+               p = quo(e, c * u, x, false);
+}
+
+
 /*
  *  GCD of multivariate polynomials
  */
@@ -939,7 +1078,7 @@ static ex sr_gcd(const ex &a, const ex &b, sym_desc_vec::const_iterator var)
 #endif
 
        // The first symbol is our main variable
-       const symbol &x = *(var->sym);
+       const ex &x = var->sym;
 
        // Sort c and d so that c has higher degree
        ex c, d;
@@ -1003,7 +1142,6 @@ static ex sr_gcd(const ex &a, const ex &b, sym_desc_vec::const_iterator var)
 /** Return maximum (absolute value) coefficient of a polynomial.
  *  This function is used internally by heur_gcd().
  *
- *  @param e  expanded multivariate polynomial
  *  @return maximum coefficient
  *  @see heur_gcd */
 numeric ex::max_coefficient() const
@@ -1015,7 +1153,7 @@ numeric ex::max_coefficient() const
  *  @see heur_gcd */
 numeric basic::max_coefficient() const
 {
-       return _num1;
+       return *_num1_p;
 }
 
 numeric numeric::max_coefficient() const
@@ -1109,7 +1247,7 @@ ex mul::smod(const numeric &xi) const
 
 
 /** xi-adic polynomial interpolation */
-static ex interpolate(const ex &gamma, const numeric &xi, const symbol &x, int degree_hint = 1)
+static ex interpolate(const ex &gamma, const numeric &xi, const ex &x, int degree_hint = 1)
 {
        exvector g; g.reserve(degree_hint);
        ex e = gamma;
@@ -1130,17 +1268,19 @@ class gcdheu_failed {};
  *  polynomials and an iterator to the first element of the sym_desc vector
  *  passed in. This function is used internally by gcd().
  *
- *  @param a  first multivariate polynomial (expanded)
- *  @param b  second multivariate polynomial (expanded)
+ *  @param a  first integer multivariate polynomial (expanded)
+ *  @param b  second integer multivariate polynomial (expanded)
  *  @param ca  cofactor of polynomial a (returned), NULL to suppress
  *             calculation of cofactor
  *  @param cb  cofactor of polynomial b (returned), NULL to suppress
  *             calculation of cofactor
  *  @param var iterator to first element of vector of sym_desc structs
- *  @return the GCD as a new expression
+ *  @param res the GCD (returned)
+ *  @return true if GCD was computed, false otherwise.
  *  @see gcd
  *  @exception gcdheu_failed() */
-static ex heur_gcd(const ex &a, const ex &b, ex *ca, ex *cb, sym_desc_vec::const_iterator var)
+static bool heur_gcd_z(ex& res, const ex &a, const ex &b, ex *ca, ex *cb,
+                      sym_desc_vec::const_iterator var)
 {
 #if STATISTICS
        heur_gcd_called++;
@@ -1148,7 +1288,7 @@ static ex heur_gcd(const ex &a, const ex &b, ex *ca, ex *cb, sym_desc_vec::const
 
        // Algorithm only works for non-vanishing input polynomials
        if (a.is_zero() || b.is_zero())
-               return (new fail())->setflag(status_flags::dynallocated);
+               return false;
 
        // GCD of two numeric values -> CLN
        if (is_exactly_a<numeric>(a) && is_exactly_a<numeric>(b)) {
@@ -1157,11 +1297,12 @@ static ex heur_gcd(const ex &a, const ex &b, ex *ca, ex *cb, sym_desc_vec::const
                        *ca = ex_to<numeric>(a) / g;
                if (cb)
                        *cb = ex_to<numeric>(b) / g;
-               return g;
+               res = g;
+               return true;
        }
 
        // The first symbol is our main variable
-       const symbol &x = *(var->sym);
+       const ex &x = var->sym;
 
        // Remove integer content
        numeric gc = gcd(a.integer_content(), b.integer_content());
@@ -1175,9 +1316,9 @@ static ex heur_gcd(const ex &a, const ex &b, ex *ca, ex *cb, sym_desc_vec::const
        numeric mq = q.max_coefficient();
        numeric xi;
        if (mp > mq)
-               xi = mq * _num2 + _num2;
+               xi = mq * (*_num2_p) + (*_num2_p);
        else
-               xi = mp * _num2 + _num2;
+               xi = mp * (*_num2_p) + (*_num2_p);
 
        // 6 tries maximum
        for (int t=0; t<6; t++) {
@@ -1187,9 +1328,13 @@ static ex heur_gcd(const ex &a, const ex &b, ex *ca, ex *cb, sym_desc_vec::const
 
                // Apply evaluation homomorphism and calculate GCD
                ex cp, cq;
-               ex gamma = heur_gcd(p.subs(x == xi, subs_options::no_pattern), q.subs(x == xi, subs_options::no_pattern), &cp, &cq, var+1).expand();
-               if (!is_exactly_a<fail>(gamma)) {
-
+               ex gamma;
+               bool found = heur_gcd_z(gamma,
+                                       p.subs(x == xi, subs_options::no_pattern),
+                                       q.subs(x == xi, subs_options::no_pattern),
+                                       &cp, &cq, var+1);
+               if (found) {
+                       gamma = gamma.expand();
                        // Reconstruct polynomial from GCD of mapped polynomials
                        ex g = interpolate(gamma, xi, x, maxdeg);
 
@@ -1200,30 +1345,96 @@ static ex heur_gcd(const ex &a, const ex &b, ex *ca, ex *cb, sym_desc_vec::const
                        ex dummy;
                        if (divide_in_z(p, g, ca ? *ca : dummy, var) && divide_in_z(q, g, cb ? *cb : dummy, var)) {
                                g *= gc;
-                               ex lc = g.lcoeff(x);
-                               if (is_exactly_a<numeric>(lc) && ex_to<numeric>(lc).is_negative())
-                                       return -g;
-                               else
-                                       return g;
+                               res = g;
+                               return true;
                        }
                }
 
                // Next evaluation point
                xi = iquo(xi * isqrt(isqrt(xi)) * numeric(73794), numeric(27011));
        }
-       return (new fail())->setflag(status_flags::dynallocated);
+       return false;
 }
 
+/** Compute GCD of multivariate polynomials using the heuristic GCD algorithm.
+ *  get_symbol_stats() must have been called previously with the input
+ *  polynomials and an iterator to the first element of the sym_desc vector
+ *  passed in. This function is used internally by gcd().
+ *
+ *  @param a  first rational multivariate polynomial (expanded)
+ *  @param b  second rational multivariate polynomial (expanded)
+ *  @param ca  cofactor of polynomial a (returned), NULL to suppress
+ *             calculation of cofactor
+ *  @param cb  cofactor of polynomial b (returned), NULL to suppress
+ *             calculation of cofactor
+ *  @param var iterator to first element of vector of sym_desc structs
+ *  @param res the GCD (returned)
+ *  @return true if GCD was computed, false otherwise.
+ *  @see heur_gcd_z
+ *  @see gcd
+ */
+static bool heur_gcd(ex& res, const ex& a, const ex& b, ex *ca, ex *cb,
+                    sym_desc_vec::const_iterator var)
+{
+       if (a.info(info_flags::integer_polynomial) && 
+           b.info(info_flags::integer_polynomial)) {
+               try {
+                       return heur_gcd_z(res, a, b, ca, cb, var);
+               } catch (gcdheu_failed) {
+                       return false;
+               }
+       }
+
+       // convert polynomials to Z[X]
+       const numeric a_lcm = lcm_of_coefficients_denominators(a);
+       const numeric ab_lcm = lcmcoeff(b, a_lcm);
+
+       const ex ai = a*ab_lcm;
+       const ex bi = b*ab_lcm;
+       if (!ai.info(info_flags::integer_polynomial))
+               throw std::logic_error("heur_gcd: not an integer polynomial [1]");
+
+       if (!bi.info(info_flags::integer_polynomial))
+               throw std::logic_error("heur_gcd: not an integer polynomial [2]");
+
+       bool found = false;
+       try {
+               found = heur_gcd_z(res, ai, bi, ca, cb, var);
+       } catch (gcdheu_failed) {
+               return false;
+       }
+       
+       // GCD is not unique, it's defined up to a unit (i.e. invertible
+       // element). If the coefficient ring is a field, every its element is
+       // invertible, so one can multiply the polynomial GCD with any element
+       // of the coefficient field. We use this ambiguity to make cofactors
+       // integer polynomials.
+       if (found)
+               res /= ab_lcm;
+       return found;
+}
+
+
+// gcd helper to handle partially factored polynomials (to avoid expanding
+// large expressions). At least one of the arguments should be a power.
+static ex gcd_pf_pow(const ex& a, const ex& b, ex* ca, ex* cb);
+
+// gcd helper to handle partially factored polynomials (to avoid expanding
+// large expressions). At least one of the arguments should be a product.
+static ex gcd_pf_mul(const ex& a, const ex& b, ex* ca, ex* cb);
 
 /** Compute GCD (Greatest Common Divisor) of multivariate polynomials a(X)
- *  and b(X) in Z[X].
+ *  and b(X) in Z[X]. Optionally also compute the cofactors of a and b,
+ *  defined by a = ca * gcd(a, b) and b = cb * gcd(a, b).
  *
  *  @param a  first multivariate polynomial
  *  @param b  second multivariate polynomial
+ *  @param ca pointer to expression that will receive the cofactor of a, or NULL
+ *  @param cb pointer to expression that will receive the cofactor of b, or NULL
  *  @param check_args  check whether a and b are polynomials with rational
  *         coefficients (defaults to "true")
  *  @return the GCD as a new expression */
-ex gcd(const ex &a, const ex &b, ex *ca, ex *cb, bool check_args)
+ex gcd(const ex &a, const ex &b, ex *ca, ex *cb, bool check_args, unsigned options)
 {
 #if STATISTICS
        gcd_called++;
@@ -1254,90 +1465,14 @@ ex gcd(const ex &a, const ex &b, ex *ca, ex *cb, bool check_args)
        }
 
        // Partially factored cases (to avoid expanding large expressions)
-       if (is_exactly_a<mul>(a)) {
-               if (is_exactly_a<mul>(b) && b.nops() > a.nops())
-                       goto factored_b;
-factored_a:
-               size_t num = a.nops();
-               exvector g; g.reserve(num);
-               exvector acc_ca; acc_ca.reserve(num);
-               ex part_b = b;
-               for (size_t i=0; i<num; i++) {
-                       ex part_ca, part_cb;
-                       g.push_back(gcd(a.op(i), part_b, &part_ca, &part_cb, check_args));
-                       acc_ca.push_back(part_ca);
-                       part_b = part_cb;
-               }
-               if (ca)
-                       *ca = (new mul(acc_ca))->setflag(status_flags::dynallocated);
-               if (cb)
-                       *cb = part_b;
-               return (new mul(g))->setflag(status_flags::dynallocated);
-       } else if (is_exactly_a<mul>(b)) {
-               if (is_exactly_a<mul>(a) && a.nops() > b.nops())
-                       goto factored_a;
-factored_b:
-               size_t num = b.nops();
-               exvector g; g.reserve(num);
-               exvector acc_cb; acc_cb.reserve(num);
-               ex part_a = a;
-               for (size_t i=0; i<num; i++) {
-                       ex part_ca, part_cb;
-                       g.push_back(gcd(part_a, b.op(i), &part_ca, &part_cb, check_args));
-                       acc_cb.push_back(part_cb);
-                       part_a = part_ca;
-               }
-               if (ca)
-                       *ca = part_a;
-               if (cb)
-                       *cb = (new mul(acc_cb))->setflag(status_flags::dynallocated);
-               return (new mul(g))->setflag(status_flags::dynallocated);
-       }
-
+       if (!(options & gcd_options::no_part_factored)) {
+               if (is_exactly_a<mul>(a) || is_exactly_a<mul>(b))
+                       return gcd_pf_mul(a, b, ca, cb);
 #if FAST_COMPARE
-       // Input polynomials of the form poly^n are sometimes also trivial
-       if (is_exactly_a<power>(a)) {
-               ex p = a.op(0);
-               if (is_exactly_a<power>(b)) {
-                       if (p.is_equal(b.op(0))) {
-                               // a = p^n, b = p^m, gcd = p^min(n, m)
-                               ex exp_a = a.op(1), exp_b = b.op(1);
-                               if (exp_a < exp_b) {
-                                       if (ca)
-                                               *ca = _ex1;
-                                       if (cb)
-                                               *cb = power(p, exp_b - exp_a);
-                                       return power(p, exp_a);
-                               } else {
-                                       if (ca)
-                                               *ca = power(p, exp_a - exp_b);
-                                       if (cb)
-                                               *cb = _ex1;
-                                       return power(p, exp_b);
-                               }
-                       }
-               } else {
-                       if (p.is_equal(b)) {
-                               // a = p^n, b = p, gcd = p
-                               if (ca)
-                                       *ca = power(p, a.op(1) - 1);
-                               if (cb)
-                                       *cb = _ex1;
-                               return p;
-                       }
-               }
-       } else if (is_exactly_a<power>(b)) {
-               ex p = b.op(0);
-               if (p.is_equal(a)) {
-                       // a = p, b = p^n, gcd = p
-                       if (ca)
-                               *ca = _ex1;
-                       if (cb)
-                               *cb = power(p, b.op(1) - 1);
-                       return p;
-               }
-       }
+               if (is_exactly_a<power>(a) || is_exactly_a<power>(b))
+                       return gcd_pf_pow(a, b, ca, cb);
 #endif
+       }
 
        // Some trivial cases
        ex aex = a.expand(), bex = b.expand();
@@ -1372,13 +1507,73 @@ factored_b:
        }
 #endif
 
+       if (is_a<symbol>(aex)) {
+               if (! bex.subs(aex==_ex0, subs_options::no_pattern).is_zero()) {
+                       if (ca)
+                               *ca = a;
+                       if (cb)
+                               *cb = b;
+                       return _ex1;
+               }
+       }
+
+       if (is_a<symbol>(bex)) {
+               if (! aex.subs(bex==_ex0, subs_options::no_pattern).is_zero()) {
+                       if (ca)
+                               *ca = a;
+                       if (cb)
+                               *cb = b;
+                       return _ex1;
+               }
+       }
+
+       if (is_exactly_a<numeric>(aex)) {
+               numeric bcont = bex.integer_content();
+               numeric g = gcd(ex_to<numeric>(aex), bcont);
+               if (ca)
+                       *ca = ex_to<numeric>(aex)/g;
+               if (cb)
+                       *cb = bex/g;
+               return g;
+       }
+
+       if (is_exactly_a<numeric>(bex)) {
+               numeric acont = aex.integer_content();
+               numeric g = gcd(ex_to<numeric>(bex), acont);
+               if (ca)
+                       *ca = aex/g;
+               if (cb)
+                       *cb = ex_to<numeric>(bex)/g;
+               return g;
+       }
+
        // Gather symbol statistics
        sym_desc_vec sym_stats;
        get_symbol_stats(a, b, sym_stats);
 
-       // The symbol with least degree is our main variable
+       // The symbol with least degree which is contained in both polynomials
+       // is our main variable
+       sym_desc_vec::iterator vari = sym_stats.begin();
+       while ((vari != sym_stats.end()) && 
+              (((vari->ldeg_b == 0) && (vari->deg_b == 0)) ||
+               ((vari->ldeg_a == 0) && (vari->deg_a == 0))))
+               vari++;
+
+       // No common symbols at all, just return 1:
+       if (vari == sym_stats.end()) {
+               // N.B: keep cofactors factored
+               if (ca)
+                       *ca = a;
+               if (cb)
+                       *cb = b;
+               return _ex1;
+       }
+       // move symbols which contained only in one of the polynomials
+       // to the end:
+       rotate(sym_stats.begin(), vari, sym_stats.end());
+
        sym_desc_vec::const_iterator var = sym_stats.begin();
-       const symbol &x = *(var->sym);
+       const ex &x = var->sym;
 
        // Cancel trivial common factor
        int ldeg_a = var->ldeg_a;
@@ -1390,57 +1585,173 @@ factored_b:
        }
 
        // Try to eliminate variables
-       if (var->deg_a == 0) {
-               ex c = bex.content(x);
-               ex g = gcd(aex, c, ca, cb, false);
+       if (var->deg_a == 0 && var->deg_b != 0 ) {
+               ex bex_u, bex_c, bex_p;
+               bex.unitcontprim(x, bex_u, bex_c, bex_p);
+               ex g = gcd(aex, bex_c, ca, cb, false);
                if (cb)
-                       *cb *= bex.unit(x) * bex.primpart(x, c);
+                       *cb *= bex_u * bex_p;
                return g;
-       } else if (var->deg_b == 0) {
-               ex c = aex.content(x);
-               ex g = gcd(c, bex, ca, cb, false);
+       } else if (var->deg_b == 0 && var->deg_a != 0) {
+               ex aex_u, aex_c, aex_p;
+               aex.unitcontprim(x, aex_u, aex_c, aex_p);
+               ex g = gcd(aex_c, bex, ca, cb, false);
                if (ca)
-                       *ca *= aex.unit(x) * aex.primpart(x, c);
+                       *ca *= aex_u * aex_p;
                return g;
        }
 
        // Try heuristic algorithm first, fall back to PRS if that failed
        ex g;
-       try {
-               g = heur_gcd(aex, bex, ca, cb, var);
-       } catch (gcdheu_failed) {
-               g = fail();
-       }
-       if (is_exactly_a<fail>(g)) {
+       if (!(options & gcd_options::no_heur_gcd)) {
+               bool found = heur_gcd(g, aex, bex, ca, cb, var);
+               if (found) {
+                       // heur_gcd have already computed cofactors...
+                       if (g.is_equal(_ex1)) {
+                               // ... but we want to keep them factored if possible.
+                               if (ca)
+                                       *ca = a;
+                               if (cb)
+                                       *cb = b;
+                       }
+                       return g;
+               }
 #if STATISTICS
-               heur_gcd_failed++;
+               else {
+                       heur_gcd_failed++;
+               }
 #endif
-               g = sr_gcd(aex, bex, var);
-               if (g.is_equal(_ex1)) {
-                       // Keep cofactors factored if possible
+       }
+
+       g = sr_gcd(aex, bex, var);
+       if (g.is_equal(_ex1)) {
+               // Keep cofactors factored if possible
+               if (ca)
+                       *ca = a;
+               if (cb)
+                       *cb = b;
+       } else {
+               if (ca)
+                       divide(aex, g, *ca, false);
+               if (cb)
+                       divide(bex, g, *cb, false);
+       }
+       return g;
+}
+
+// gcd helper to handle partially factored polynomials (to avoid expanding
+// large expressions). Both arguments should be powers.
+static ex gcd_pf_pow_pow(const ex& a, const ex& b, ex* ca, ex* cb)
+{
+       ex p = a.op(0);
+       const ex& exp_a = a.op(1);
+       ex pb = b.op(0);
+       const ex& exp_b = b.op(1);
+
+       // a = p^n, b = p^m, gcd = p^min(n, m)
+       if (p.is_equal(pb)) {
+               if (exp_a < exp_b) {
                        if (ca)
-                               *ca = a;
+                               *ca = _ex1;
                        if (cb)
-                               *cb = b;
+                               *cb = power(p, exp_b - exp_a);
+                       return power(p, exp_a);
                } else {
                        if (ca)
-                               divide(aex, g, *ca, false);
+                               *ca = power(p, exp_a - exp_b);
                        if (cb)
-                               divide(bex, g, *cb, false);
+                               *cb = _ex1;
+                       return power(p, exp_b);
                }
-       } else {
-               if (g.is_equal(_ex1)) {
-                       // Keep cofactors factored if possible
+       }
+
+       ex p_co, pb_co;
+       ex p_gcd = gcd(p, pb, &p_co, &pb_co, false);
+       // a(x) = p(x)^n, b(x) = p_b(x)^m, gcd (p, p_b) = 1 ==> gcd(a,b) = 1
+       if (p_gcd.is_equal(_ex1)) {
                        if (ca)
                                *ca = a;
                        if (cb)
                                *cb = b;
-               }
+                       return _ex1;
+                       // XXX: do I need to check for p_gcd = -1?
        }
 
-       return g;
+       // there are common factors:
+       // a(x) = g(x)^n A(x)^n, b(x) = g(x)^m B(x)^m ==>
+       // gcd(a, b) = g(x)^n gcd(A(x)^n, g(x)^(n-m) B(x)^m
+       if (exp_a < exp_b) {
+               ex pg =  gcd(power(p_co, exp_a), power(p_gcd, exp_b-exp_a)*power(pb_co, exp_b), ca, cb, false);
+               return power(p_gcd, exp_a)*pg;
+       } else {
+               ex pg = gcd(power(p_gcd, exp_a - exp_b)*power(p_co, exp_a), power(pb_co, exp_b), ca, cb, false);
+               return power(p_gcd, exp_b)*pg;
+       }
 }
 
+static ex gcd_pf_pow(const ex& a, const ex& b, ex* ca, ex* cb)
+{
+       if (is_exactly_a<power>(a) && is_exactly_a<power>(b))
+               return gcd_pf_pow_pow(a, b, ca, cb);
+
+       if (is_exactly_a<power>(b) && (! is_exactly_a<power>(a)))
+               return gcd_pf_pow(b, a, cb, ca);
+
+       GINAC_ASSERT(is_exactly_a<power>(a));
+
+       ex p = a.op(0);
+       const ex& exp_a = a.op(1);
+       if (p.is_equal(b)) {
+               // a = p^n, b = p, gcd = p
+               if (ca)
+                       *ca = power(p, a.op(1) - 1);
+               if (cb)
+                       *cb = _ex1;
+               return p;
+       } 
+
+       ex p_co, bpart_co;
+       ex p_gcd = gcd(p, b, &p_co, &bpart_co, false);
+
+       // a(x) = p(x)^n, gcd(p, b) = 1 ==> gcd(a, b) = 1
+       if (p_gcd.is_equal(_ex1)) {
+               if (ca)
+                       *ca = a;
+               if (cb)
+                       *cb = b;
+               return _ex1;
+       }
+       // a(x) = g(x)^n A(x)^n, b(x) = g(x) B(x) ==> gcd(a, b) = g(x) gcd(g(x)^(n-1) A(x)^n, B(x))
+       ex rg = gcd(power(p_gcd, exp_a-1)*power(p_co, exp_a), bpart_co, ca, cb, false);
+       return p_gcd*rg;
+}
+
+static ex gcd_pf_mul(const ex& a, const ex& b, ex* ca, ex* cb)
+{
+       if (is_exactly_a<mul>(a) && is_exactly_a<mul>(b)
+                                && (b.nops() >  a.nops()))
+               return gcd_pf_mul(b, a, cb, ca);
+
+       if (is_exactly_a<mul>(b) && (!is_exactly_a<mul>(a)))
+               return gcd_pf_mul(b, a, cb, ca);
+
+       GINAC_ASSERT(is_exactly_a<mul>(a));
+       size_t num = a.nops();
+       exvector g; g.reserve(num);
+       exvector acc_ca; acc_ca.reserve(num);
+       ex part_b = b;
+       for (size_t i=0; i<num; i++) {
+               ex part_ca, part_cb;
+               g.push_back(gcd(a.op(i), part_b, &part_ca, &part_cb, false));
+               acc_ca.push_back(part_ca);
+               part_b = part_cb;
+       }
+       if (ca)
+               *ca = (new mul(acc_ca))->setflag(status_flags::dynallocated);
+       if (cb)
+               *cb = part_b;
+       return (new mul(g))->setflag(status_flags::dynallocated);
+}
 
 /** Compute LCM (Least Common Multiple) of multivariate polynomials in Z[X].
  *
@@ -1467,7 +1778,7 @@ ex lcm(const ex &a, const ex &b, bool check_args)
  */
 
 /** Compute square-free factorization of multivariate polynomial a(x) using
- *  Yun´s algorithm.  Used internally by sqrfree().
+ *  Yun's algorithm.  Used internally by sqrfree().
  *
  *  @param a  multivariate polynomial over Z[X], treated here as univariate
  *            polynomial in x.
@@ -1498,7 +1809,7 @@ static exvector sqrfree_yun(const ex &a, const symbol &x)
 /** Compute a square-free factorization of a multivariate polynomial in Q[X].
  *
  *  @param a  multivariate polynomial over Q[X]
- *  @param x  lst of variables to factor in, may be left empty for autodetection
+ *  @param l  lst of variables to factor in, may be left empty for autodetection
  *  @return   a square-free factorization of \p a.
  *
  * \note
@@ -1545,7 +1856,7 @@ ex sqrfree(const ex &a, const lst &l)
                get_symbol_stats(a, _ex0, sdv);
                sym_desc_vec::const_iterator it = sdv.begin(), itend = sdv.end();
                while (it != itend) {
-                       args.append(*it->sym);
+                       args.append(it->sym);
                        ++it;
                }
        } else {
@@ -1562,7 +1873,7 @@ ex sqrfree(const ex &a, const lst &l)
        const ex tmp = multiply_lcm(a,lcm);
 
        // find the factors
-       exvector factors = sqrfree_yun(tmp,x);
+       exvector factors = sqrfree_yun(tmp, x);
 
        // construct the next list of symbols with the first element popped
        lst newargs = args;
@@ -1702,23 +2013,23 @@ static ex replace_with_symbol(const ex & e, exmap & repl, exmap & rev_lookup)
 }
 
 /** Create a symbol for replacing the expression "e" (or return a previously
- *  assigned symbol). An expression of the form "symbol == expression" is added
- *  to repl_lst and the symbol is returned.
+ *  assigned symbol). The symbol and expression are appended to repl, and the
+ *  symbol is returned.
  *  @see basic::to_rational
  *  @see basic::to_polynomial */
-static ex replace_with_symbol(const ex & e, lst & repl_lst)
+static ex replace_with_symbol(const ex & e, exmap & repl)
 {
-       // Expression already in repl_lst? Then return the assigned symbol
-       for (lst::const_iterator it = repl_lst.begin(); it != repl_lst.end(); ++it)
-               if (it->op(1).is_equal(e))
-                       return it->op(0);
+       // Expression already replaced? Then return the assigned symbol
+       for (exmap::const_iterator it = repl.begin(); it != repl.end(); ++it)
+               if (it->second.is_equal(e))
+                       return it->first;
        
        // Otherwise create new symbol and add to list, taking care that the
-       // replacement expression doesn't itself contain symbols from the repl_lst,
+       // replacement expression doesn't itself contain symbols from repl,
        // because subs() is not recursive
        ex es = (new symbol)->setflag(status_flags::dynallocated);
-       ex e_replaced = e.subs(repl_lst, subs_options::no_pattern);
-       repl_lst.append(es == e_replaced);
+       ex e_replaced = e.subs(repl, subs_options::no_pattern);
+       repl.insert(std::make_pair(es, e_replaced));
        return es;
 }
 
@@ -1790,7 +2101,7 @@ static ex frac_cancel(const ex &n, const ex &d)
 {
        ex num = n;
        ex den = d;
-       numeric pre_factor = _num1;
+       numeric pre_factor = *_num1_p;
 
 //std::clog << "frac_cancel num = " << num << ", den = " << den << std::endl;
 
@@ -1827,10 +2138,10 @@ static ex frac_cancel(const ex &n, const ex &d)
                        den *= _ex_1;
                }
        } else {
-               const symbol *x;
+               ex x;
                if (get_first_symbol(den, x)) {
-                       GINAC_ASSERT(is_exactly_a<numeric>(den.unit(*x)));
-                       if (ex_to<numeric>(den.unit(*x)).is_negative()) {
+                       GINAC_ASSERT(is_exactly_a<numeric>(den.unit(x)));
+                       if (ex_to<numeric>(den.unit(x)).is_negative()) {
                                num *= _ex_1;
                                den *= _ex_1;
                        }
@@ -2103,46 +2414,80 @@ ex ex::numer_denom() const
  *  on non-rational functions by applying to_rational() on the arguments,
  *  calling the desired function and re-substituting the temporary symbols
  *  in the result. To make the last step possible, all temporary symbols and
- *  their associated expressions are collected in the list specified by the
- *  repl_lst parameter in the form {symbol == expression}, ready to be passed
- *  as an argument to ex::subs().
+ *  their associated expressions are collected in the map specified by the
+ *  repl parameter, ready to be passed as an argument to ex::subs().
  *
- *  @param repl_lst collects a list of all temporary symbols and their replacements
+ *  @param repl collects all temporary symbols and their replacements
  *  @return rationalized expression */
-ex ex::to_rational(lst &repl_lst) const
+ex ex::to_rational(exmap & repl) const
+{
+       return bp->to_rational(repl);
+}
+
+// GiNaC 1.1 compatibility function
+ex ex::to_rational(lst & repl_lst) const
 {
-       return bp->to_rational(repl_lst);
+       // Convert lst to exmap
+       exmap m;
+       for (lst::const_iterator it = repl_lst.begin(); it != repl_lst.end(); ++it)
+               m.insert(std::make_pair(it->op(0), it->op(1)));
+
+       ex ret = bp->to_rational(m);
+
+       // Convert exmap back to lst
+       repl_lst.remove_all();
+       for (exmap::const_iterator it = m.begin(); it != m.end(); ++it)
+               repl_lst.append(it->first == it->second);
+
+       return ret;
 }
 
-ex ex::to_polynomial(lst &repl_lst) const
+ex ex::to_polynomial(exmap & repl) const
 {
-       return bp->to_polynomial(repl_lst);
+       return bp->to_polynomial(repl);
 }
 
+// GiNaC 1.1 compatibility function
+ex ex::to_polynomial(lst & repl_lst) const
+{
+       // Convert lst to exmap
+       exmap m;
+       for (lst::const_iterator it = repl_lst.begin(); it != repl_lst.end(); ++it)
+               m.insert(std::make_pair(it->op(0), it->op(1)));
+
+       ex ret = bp->to_polynomial(m);
+
+       // Convert exmap back to lst
+       repl_lst.remove_all();
+       for (exmap::const_iterator it = m.begin(); it != m.end(); ++it)
+               repl_lst.append(it->first == it->second);
+
+       return ret;
+}
 
 /** Default implementation of ex::to_rational(). This replaces the object with
  *  a temporary symbol. */
-ex basic::to_rational(lst &repl_lst) const
+ex basic::to_rational(exmap & repl) const
 {
-       return replace_with_symbol(*this, repl_lst);
+       return replace_with_symbol(*this, repl);
 }
 
-ex basic::to_polynomial(lst &repl_lst) const
+ex basic::to_polynomial(exmap & repl) const
 {
-       return replace_with_symbol(*this, repl_lst);
+       return replace_with_symbol(*this, repl);
 }
 
 
 /** Implementation of ex::to_rational() for symbols. This returns the
  *  unmodified symbol. */
-ex symbol::to_rational(lst &repl_lst) const
+ex symbol::to_rational(exmap & repl) const
 {
        return *this;
 }
 
 /** Implementation of ex::to_polynomial() for symbols. This returns the
  *  unmodified symbol. */
-ex symbol::to_polynomial(lst &repl_lst) const
+ex symbol::to_polynomial(exmap & repl) const
 {
        return *this;
 }
@@ -2151,17 +2496,17 @@ ex symbol::to_polynomial(lst &repl_lst) const
 /** Implementation of ex::to_rational() for a numeric. It splits complex
  *  numbers into re+I*im and replaces I and non-rational real numbers with a
  *  temporary symbol. */
-ex numeric::to_rational(lst &repl_lst) const
+ex numeric::to_rational(exmap & repl) const
 {
        if (is_real()) {
                if (!is_rational())
-                       return replace_with_symbol(*this, repl_lst);
+                       return replace_with_symbol(*this, repl);
        } else { // complex
                numeric re = real();
                numeric im = imag();
-               ex re_ex = re.is_rational() ? re : replace_with_symbol(re, repl_lst);
-               ex im_ex = im.is_rational() ? im : replace_with_symbol(im, repl_lst);
-               return re_ex + im_ex * replace_with_symbol(I, repl_lst);
+               ex re_ex = re.is_rational() ? re : replace_with_symbol(re, repl);
+               ex im_ex = im.is_rational() ? im : replace_with_symbol(im, repl);
+               return re_ex + im_ex * replace_with_symbol(I, repl);
        }
        return *this;
 }
@@ -2169,17 +2514,17 @@ ex numeric::to_rational(lst &repl_lst) const
 /** Implementation of ex::to_polynomial() for a numeric. It splits complex
  *  numbers into re+I*im and replaces I and non-integer real numbers with a
  *  temporary symbol. */
-ex numeric::to_polynomial(lst &repl_lst) const
+ex numeric::to_polynomial(exmap & repl) const
 {
        if (is_real()) {
                if (!is_integer())
-                       return replace_with_symbol(*this, repl_lst);
+                       return replace_with_symbol(*this, repl);
        } else { // complex
                numeric re = real();
                numeric im = imag();
-               ex re_ex = re.is_integer() ? re : replace_with_symbol(re, repl_lst);
-               ex im_ex = im.is_integer() ? im : replace_with_symbol(im, repl_lst);
-               return re_ex + im_ex * replace_with_symbol(I, repl_lst);
+               ex re_ex = re.is_integer() ? re : replace_with_symbol(re, repl);
+               ex im_ex = im.is_integer() ? im : replace_with_symbol(im, repl);
+               return re_ex + im_ex * replace_with_symbol(I, repl);
        }
        return *this;
 }
@@ -2187,36 +2532,47 @@ ex numeric::to_polynomial(lst &repl_lst) const
 
 /** Implementation of ex::to_rational() for powers. It replaces non-integer
  *  powers by temporary symbols. */
-ex power::to_rational(lst &repl_lst) const
+ex power::to_rational(exmap & repl) const
 {
        if (exponent.info(info_flags::integer))
-               return power(basis.to_rational(repl_lst), exponent);
+               return power(basis.to_rational(repl), exponent);
        else
-               return replace_with_symbol(*this, repl_lst);
+               return replace_with_symbol(*this, repl);
 }
 
 /** Implementation of ex::to_polynomial() for powers. It replaces non-posint
  *  powers by temporary symbols. */
-ex power::to_polynomial(lst &repl_lst) const
+ex power::to_polynomial(exmap & repl) const
 {
        if (exponent.info(info_flags::posint))
-               return power(basis.to_rational(repl_lst), exponent);
+               return power(basis.to_rational(repl), exponent);
+       else if (exponent.info(info_flags::negint))
+       {
+               ex basis_pref = collect_common_factors(basis);
+               if (is_exactly_a<mul>(basis_pref) || is_exactly_a<power>(basis_pref)) {
+                       // (A*B)^n will be automagically transformed to A^n*B^n
+                       ex t = power(basis_pref, exponent);
+                       return t.to_polynomial(repl);
+               }
+               else
+                       return power(replace_with_symbol(power(basis, _ex_1), repl), -exponent);
+       } 
        else
-               return replace_with_symbol(*this, repl_lst);
+               return replace_with_symbol(*this, repl);
 }
 
 
 /** Implementation of ex::to_rational() for expairseqs. */
-ex expairseq::to_rational(lst &repl_lst) const
+ex expairseq::to_rational(exmap & repl) const
 {
        epvector s;
        s.reserve(seq.size());
        epvector::const_iterator i = seq.begin(), end = seq.end();
        while (i != end) {
-               s.push_back(split_ex_to_pair(recombine_pair_to_ex(*i).to_rational(repl_lst)));
+               s.push_back(split_ex_to_pair(recombine_pair_to_ex(*i).to_rational(repl)));
                ++i;
        }
-       ex oc = overall_coeff.to_rational(repl_lst);
+       ex oc = overall_coeff.to_rational(repl);
        if (oc.info(info_flags::numeric))
                return thisexpairseq(s, overall_coeff);
        else
@@ -2225,16 +2581,16 @@ ex expairseq::to_rational(lst &repl_lst) const
 }
 
 /** Implementation of ex::to_polynomial() for expairseqs. */
-ex expairseq::to_polynomial(lst &repl_lst) const
+ex expairseq::to_polynomial(exmap & repl) const
 {
        epvector s;
        s.reserve(seq.size());
        epvector::const_iterator i = seq.begin(), end = seq.end();
        while (i != end) {
-               s.push_back(split_ex_to_pair(recombine_pair_to_ex(*i).to_polynomial(repl_lst)));
+               s.push_back(split_ex_to_pair(recombine_pair_to_ex(*i).to_polynomial(repl)));
                ++i;
        }
-       ex oc = overall_coeff.to_polynomial(repl_lst);
+       ex oc = overall_coeff.to_polynomial(repl);
        if (oc.info(info_flags::numeric))
                return thisexpairseq(s, overall_coeff);
        else
@@ -2246,7 +2602,7 @@ ex expairseq::to_polynomial(lst &repl_lst) const
 /** Remove the common factor in the terms of a sum 'e' by calculating the GCD,
  *  and multiply it into the expression 'factor' (which needs to be initialized
  *  to 1, unless you're accumulating factors). */
-static ex find_common_factor(const ex & e, ex & factor, lst & repl)
+static ex find_common_factor(const ex & e, ex & factor, exmap & repl)
 {
        if (is_exactly_a<add>(e)) {
 
@@ -2258,7 +2614,7 @@ static ex find_common_factor(const ex & e, ex & factor, lst & repl)
                for (size_t i=0; i<num; i++) {
                        ex x = e.op(i).to_polynomial(repl);
 
-                       if (is_exactly_a<add>(x) || is_exactly_a<mul>(x)) {
+                       if (is_exactly_a<add>(x) || is_exactly_a<mul>(x) || is_a<power>(x)) {
                                ex f = 1;
                                x = find_common_factor(x, f, repl);
                                x *= f;
@@ -2317,8 +2673,16 @@ term_done:       ;
                return (new mul(v))->setflag(status_flags::dynallocated);
 
        } else if (is_exactly_a<power>(e)) {
-
-               return e.to_polynomial(repl);
+               const ex e_exp(e.op(1));
+               if (e_exp.info(info_flags::integer)) {
+                       ex eb = e.op(0).to_polynomial(repl);
+                       ex factor_local(_ex1);
+                       ex pre_res = find_common_factor(eb, factor_local, repl);
+                       factor *= power(factor_local, e_exp);
+                       return power(pre_res, e_exp);
+                       
+               } else
+                       return e.to_polynomial(repl);
 
        } else
                return e;
@@ -2329,9 +2693,9 @@ term_done:        ;
  *  'a*(b*x+b*y)' to 'a*b*(x+y)'. */
 ex collect_common_factors(const ex & e)
 {
-       if (is_exactly_a<add>(e) || is_exactly_a<mul>(e)) {
+       if (is_exactly_a<add>(e) || is_exactly_a<mul>(e) || is_exactly_a<power>(e)) {
 
-               lst repl;
+               exmap repl;
                ex factor = 1;
                ex r = find_common_factor(e, factor, repl);
                return factor.subs(repl, subs_options::no_pattern) * r.subs(repl, subs_options::no_pattern);
@@ -2341,4 +2705,37 @@ ex collect_common_factors(const ex & e)
 }
 
 
+/** Resultant of two expressions e1,e2 with respect to symbol s.
+ *  Method: Compute determinant of Sylvester matrix of e1,e2,s.  */
+ex resultant(const ex & e1, const ex & e2, const ex & s)
+{
+       const ex ee1 = e1.expand();
+       const ex ee2 = e2.expand();
+       if (!ee1.info(info_flags::polynomial) ||
+           !ee2.info(info_flags::polynomial))
+               throw(std::runtime_error("resultant(): arguments must be polynomials"));
+
+       const int h1 = ee1.degree(s);
+       const int l1 = ee1.ldegree(s);
+       const int h2 = ee2.degree(s);
+       const int l2 = ee2.ldegree(s);
+
+       const int msize = h1 + h2;
+       matrix m(msize, msize);
+
+       for (int l = h1; l >= l1; --l) {
+               const ex e = ee1.coeff(s, l);
+               for (int k = 0; k < h2; ++k)
+                       m(k, k+h1-l) = e;
+       }
+       for (int l = h2; l >= l2; --l) {
+               const ex e = ee2.coeff(s, l);
+               for (int k = 0; k < h1; ++k)
+                       m(k+h2, k+h2-l) = e;
+       }
+
+       return m.determinant();
+}
+
+
 } // namespace GiNaC