]> www.ginac.de Git - ginac.git/blob - ginac/normal.cpp
gcd(): allow user to override (some of) heuristics.
[ginac.git] / ginac / normal.cpp
1 /** @file normal.cpp
2  *
3  *  This file implements several functions that work on univariate and
4  *  multivariate polynomials and rational functions.
5  *  These functions include polynomial quotient and remainder, GCD and LCM
6  *  computation, square-free factorization and rational function normalization. */
7
8 /*
9  *  GiNaC Copyright (C) 1999-2008 Johannes Gutenberg University Mainz, Germany
10  *
11  *  This program is free software; you can redistribute it and/or modify
12  *  it under the terms of the GNU General Public License as published by
13  *  the Free Software Foundation; either version 2 of the License, or
14  *  (at your option) any later version.
15  *
16  *  This program is distributed in the hope that it will be useful,
17  *  but WITHOUT ANY WARRANTY; without even the implied warranty of
18  *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
19  *  GNU General Public License for more details.
20  *
21  *  You should have received a copy of the GNU General Public License
22  *  along with this program; if not, write to the Free Software
23  *  Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
24  */
25
26 #include <algorithm>
27 #include <map>
28
29 #include "normal.h"
30 #include "basic.h"
31 #include "ex.h"
32 #include "add.h"
33 #include "constant.h"
34 #include "expairseq.h"
35 #include "fail.h"
36 #include "inifcns.h"
37 #include "lst.h"
38 #include "mul.h"
39 #include "numeric.h"
40 #include "power.h"
41 #include "relational.h"
42 #include "operators.h"
43 #include "matrix.h"
44 #include "pseries.h"
45 #include "symbol.h"
46 #include "utils.h"
47
48 namespace GiNaC {
49
50 // If comparing expressions (ex::compare()) is fast, you can set this to 1.
51 // Some routines like quo(), rem() and gcd() will then return a quick answer
52 // when they are called with two identical arguments.
53 #define FAST_COMPARE 1
54
55 // Set this if you want divide_in_z() to use remembering
56 #define USE_REMEMBER 0
57
58 // Set this if you want divide_in_z() to use trial division followed by
59 // polynomial interpolation (always slower except for completely dense
60 // polynomials)
61 #define USE_TRIAL_DIVISION 0
62
63 // Set this to enable some statistical output for the GCD routines
64 #define STATISTICS 0
65
66
67 #if STATISTICS
68 // Statistics variables
69 static int gcd_called = 0;
70 static int sr_gcd_called = 0;
71 static int heur_gcd_called = 0;
72 static int heur_gcd_failed = 0;
73
74 // Print statistics at end of program
75 static struct _stat_print {
76         _stat_print() {}
77         ~_stat_print() {
78                 std::cout << "gcd() called " << gcd_called << " times\n";
79                 std::cout << "sr_gcd() called " << sr_gcd_called << " times\n";
80                 std::cout << "heur_gcd() called " << heur_gcd_called << " times\n";
81                 std::cout << "heur_gcd() failed " << heur_gcd_failed << " times\n";
82         }
83 } stat_print;
84 #endif
85
86
87 /** Return pointer to first symbol found in expression.  Due to GiNaC's
88  *  internal ordering of terms, it may not be obvious which symbol this
89  *  function returns for a given expression.
90  *
91  *  @param e  expression to search
92  *  @param x  first symbol found (returned)
93  *  @return "false" if no symbol was found, "true" otherwise */
94 static bool get_first_symbol(const ex &e, ex &x)
95 {
96         if (is_a<symbol>(e)) {
97                 x = e;
98                 return true;
99         } else if (is_exactly_a<add>(e) || is_exactly_a<mul>(e)) {
100                 for (size_t i=0; i<e.nops(); i++)
101                         if (get_first_symbol(e.op(i), x))
102                                 return true;
103         } else if (is_exactly_a<power>(e)) {
104                 if (get_first_symbol(e.op(0), x))
105                         return true;
106         }
107         return false;
108 }
109
110
111 /*
112  *  Statistical information about symbols in polynomials
113  */
114
115 /** This structure holds information about the highest and lowest degrees
116  *  in which a symbol appears in two multivariate polynomials "a" and "b".
117  *  A vector of these structures with information about all symbols in
118  *  two polynomials can be created with the function get_symbol_stats().
119  *
120  *  @see get_symbol_stats */
121 struct sym_desc {
122         /** Reference to symbol */
123         ex sym;
124
125         /** Highest degree of symbol in polynomial "a" */
126         int deg_a;
127
128         /** Highest degree of symbol in polynomial "b" */
129         int deg_b;
130
131         /** Lowest degree of symbol in polynomial "a" */
132         int ldeg_a;
133
134         /** Lowest degree of symbol in polynomial "b" */
135         int ldeg_b;
136
137         /** Maximum of deg_a and deg_b (Used for sorting) */
138         int max_deg;
139
140         /** Maximum number of terms of leading coefficient of symbol in both polynomials */
141         size_t max_lcnops;
142
143         /** Commparison operator for sorting */
144         bool operator<(const sym_desc &x) const
145         {
146                 if (max_deg == x.max_deg)
147                         return max_lcnops < x.max_lcnops;
148                 else
149                         return max_deg < x.max_deg;
150         }
151 };
152
153 // Vector of sym_desc structures
154 typedef std::vector<sym_desc> sym_desc_vec;
155
156 // Add symbol the sym_desc_vec (used internally by get_symbol_stats())
157 static void add_symbol(const ex &s, sym_desc_vec &v)
158 {
159         sym_desc_vec::const_iterator it = v.begin(), itend = v.end();
160         while (it != itend) {
161                 if (it->sym.is_equal(s))  // If it's already in there, don't add it a second time
162                         return;
163                 ++it;
164         }
165         sym_desc d;
166         d.sym = s;
167         v.push_back(d);
168 }
169
170 // Collect all symbols of an expression (used internally by get_symbol_stats())
171 static void collect_symbols(const ex &e, sym_desc_vec &v)
172 {
173         if (is_a<symbol>(e)) {
174                 add_symbol(e, v);
175         } else if (is_exactly_a<add>(e) || is_exactly_a<mul>(e)) {
176                 for (size_t i=0; i<e.nops(); i++)
177                         collect_symbols(e.op(i), v);
178         } else if (is_exactly_a<power>(e)) {
179                 collect_symbols(e.op(0), v);
180         }
181 }
182
183 /** Collect statistical information about symbols in polynomials.
184  *  This function fills in a vector of "sym_desc" structs which contain
185  *  information about the highest and lowest degrees of all symbols that
186  *  appear in two polynomials. The vector is then sorted by minimum
187  *  degree (lowest to highest). The information gathered by this
188  *  function is used by the GCD routines to identify trivial factors
189  *  and to determine which variable to choose as the main variable
190  *  for GCD computation.
191  *
192  *  @param a  first multivariate polynomial
193  *  @param b  second multivariate polynomial
194  *  @param v  vector of sym_desc structs (filled in) */
195 static void get_symbol_stats(const ex &a, const ex &b, sym_desc_vec &v)
196 {
197         collect_symbols(a.eval(), v);   // eval() to expand assigned symbols
198         collect_symbols(b.eval(), v);
199         sym_desc_vec::iterator it = v.begin(), itend = v.end();
200         while (it != itend) {
201                 int deg_a = a.degree(it->sym);
202                 int deg_b = b.degree(it->sym);
203                 it->deg_a = deg_a;
204                 it->deg_b = deg_b;
205                 it->max_deg = std::max(deg_a, deg_b);
206                 it->max_lcnops = std::max(a.lcoeff(it->sym).nops(), b.lcoeff(it->sym).nops());
207                 it->ldeg_a = a.ldegree(it->sym);
208                 it->ldeg_b = b.ldegree(it->sym);
209                 ++it;
210         }
211         std::sort(v.begin(), v.end());
212
213 #if 0
214         std::clog << "Symbols:\n";
215         it = v.begin(); itend = v.end();
216         while (it != itend) {
217                 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;
218                 std::clog << "  lcoeff_a=" << a.lcoeff(it->sym) << ", lcoeff_b=" << b.lcoeff(it->sym) << endl;
219                 ++it;
220         }
221 #endif
222 }
223
224
225 /*
226  *  Computation of LCM of denominators of coefficients of a polynomial
227  */
228
229 // Compute LCM of denominators of coefficients by going through the
230 // expression recursively (used internally by lcm_of_coefficients_denominators())
231 static numeric lcmcoeff(const ex &e, const numeric &l)
232 {
233         if (e.info(info_flags::rational))
234                 return lcm(ex_to<numeric>(e).denom(), l);
235         else if (is_exactly_a<add>(e)) {
236                 numeric c = *_num1_p;
237                 for (size_t i=0; i<e.nops(); i++)
238                         c = lcmcoeff(e.op(i), c);
239                 return lcm(c, l);
240         } else if (is_exactly_a<mul>(e)) {
241                 numeric c = *_num1_p;
242                 for (size_t i=0; i<e.nops(); i++)
243                         c *= lcmcoeff(e.op(i), *_num1_p);
244                 return lcm(c, l);
245         } else if (is_exactly_a<power>(e)) {
246                 if (is_a<symbol>(e.op(0)))
247                         return l;
248                 else
249                         return pow(lcmcoeff(e.op(0), l), ex_to<numeric>(e.op(1)));
250         }
251         return l;
252 }
253
254 /** Compute LCM of denominators of coefficients of a polynomial.
255  *  Given a polynomial with rational coefficients, this function computes
256  *  the LCM of the denominators of all coefficients. This can be used
257  *  to bring a polynomial from Q[X] to Z[X].
258  *
259  *  @param e  multivariate polynomial (need not be expanded)
260  *  @return LCM of denominators of coefficients */
261 static numeric lcm_of_coefficients_denominators(const ex &e)
262 {
263         return lcmcoeff(e, *_num1_p);
264 }
265
266 /** Bring polynomial from Q[X] to Z[X] by multiplying in the previously
267  *  determined LCM of the coefficient's denominators.
268  *
269  *  @param e  multivariate polynomial (need not be expanded)
270  *  @param lcm  LCM to multiply in */
271 static ex multiply_lcm(const ex &e, const numeric &lcm)
272 {
273         if (is_exactly_a<mul>(e)) {
274                 size_t num = e.nops();
275                 exvector v; v.reserve(num + 1);
276                 numeric lcm_accum = *_num1_p;
277                 for (size_t i=0; i<num; i++) {
278                         numeric op_lcm = lcmcoeff(e.op(i), *_num1_p);
279                         v.push_back(multiply_lcm(e.op(i), op_lcm));
280                         lcm_accum *= op_lcm;
281                 }
282                 v.push_back(lcm / lcm_accum);
283                 return (new mul(v))->setflag(status_flags::dynallocated);
284         } else if (is_exactly_a<add>(e)) {
285                 size_t num = e.nops();
286                 exvector v; v.reserve(num);
287                 for (size_t i=0; i<num; i++)
288                         v.push_back(multiply_lcm(e.op(i), lcm));
289                 return (new add(v))->setflag(status_flags::dynallocated);
290         } else if (is_exactly_a<power>(e)) {
291                 if (is_a<symbol>(e.op(0)))
292                         return e * lcm;
293                 else
294                         return pow(multiply_lcm(e.op(0), lcm.power(ex_to<numeric>(e.op(1)).inverse())), e.op(1));
295         } else
296                 return e * lcm;
297 }
298
299
300 /** Compute the integer content (= GCD of all numeric coefficients) of an
301  *  expanded polynomial. For a polynomial with rational coefficients, this
302  *  returns g/l where g is the GCD of the coefficients' numerators and l
303  *  is the LCM of the coefficients' denominators.
304  *
305  *  @return integer content */
306 numeric ex::integer_content() const
307 {
308         return bp->integer_content();
309 }
310
311 numeric basic::integer_content() const
312 {
313         return *_num1_p;
314 }
315
316 numeric numeric::integer_content() const
317 {
318         return abs(*this);
319 }
320
321 numeric add::integer_content() const
322 {
323         epvector::const_iterator it = seq.begin();
324         epvector::const_iterator itend = seq.end();
325         numeric c = *_num0_p, l = *_num1_p;
326         while (it != itend) {
327                 GINAC_ASSERT(!is_exactly_a<numeric>(it->rest));
328                 GINAC_ASSERT(is_exactly_a<numeric>(it->coeff));
329                 c = gcd(ex_to<numeric>(it->coeff).numer(), c);
330                 l = lcm(ex_to<numeric>(it->coeff).denom(), l);
331                 it++;
332         }
333         GINAC_ASSERT(is_exactly_a<numeric>(overall_coeff));
334         c = gcd(ex_to<numeric>(overall_coeff).numer(), c);
335         l = lcm(ex_to<numeric>(overall_coeff).denom(), l);
336         return c/l;
337 }
338
339 numeric mul::integer_content() const
340 {
341 #ifdef DO_GINAC_ASSERT
342         epvector::const_iterator it = seq.begin();
343         epvector::const_iterator itend = seq.end();
344         while (it != itend) {
345                 GINAC_ASSERT(!is_exactly_a<numeric>(recombine_pair_to_ex(*it)));
346                 ++it;
347         }
348 #endif // def DO_GINAC_ASSERT
349         GINAC_ASSERT(is_exactly_a<numeric>(overall_coeff));
350         return abs(ex_to<numeric>(overall_coeff));
351 }
352
353
354 /*
355  *  Polynomial quotients and remainders
356  */
357
358 /** Quotient q(x) of polynomials a(x) and b(x) in Q[x].
359  *  It satisfies a(x)=b(x)*q(x)+r(x).
360  *
361  *  @param a  first polynomial in x (dividend)
362  *  @param b  second polynomial in x (divisor)
363  *  @param x  a and b are polynomials in x
364  *  @param check_args  check whether a and b are polynomials with rational
365  *         coefficients (defaults to "true")
366  *  @return quotient of a and b in Q[x] */
367 ex quo(const ex &a, const ex &b, const ex &x, bool check_args)
368 {
369         if (b.is_zero())
370                 throw(std::overflow_error("quo: division by zero"));
371         if (is_exactly_a<numeric>(a) && is_exactly_a<numeric>(b))
372                 return a / b;
373 #if FAST_COMPARE
374         if (a.is_equal(b))
375                 return _ex1;
376 #endif
377         if (check_args && (!a.info(info_flags::rational_polynomial) || !b.info(info_flags::rational_polynomial)))
378                 throw(std::invalid_argument("quo: arguments must be polynomials over the rationals"));
379
380         // Polynomial long division
381         ex r = a.expand();
382         if (r.is_zero())
383                 return r;
384         int bdeg = b.degree(x);
385         int rdeg = r.degree(x);
386         ex blcoeff = b.expand().coeff(x, bdeg);
387         bool blcoeff_is_numeric = is_exactly_a<numeric>(blcoeff);
388         exvector v; v.reserve(std::max(rdeg - bdeg + 1, 0));
389         while (rdeg >= bdeg) {
390                 ex term, rcoeff = r.coeff(x, rdeg);
391                 if (blcoeff_is_numeric)
392                         term = rcoeff / blcoeff;
393                 else {
394                         if (!divide(rcoeff, blcoeff, term, false))
395                                 return (new fail())->setflag(status_flags::dynallocated);
396                 }
397                 term *= power(x, rdeg - bdeg);
398                 v.push_back(term);
399                 r -= (term * b).expand();
400                 if (r.is_zero())
401                         break;
402                 rdeg = r.degree(x);
403         }
404         return (new add(v))->setflag(status_flags::dynallocated);
405 }
406
407
408 /** Remainder r(x) of polynomials a(x) and b(x) in Q[x].
409  *  It satisfies a(x)=b(x)*q(x)+r(x).
410  *
411  *  @param a  first polynomial in x (dividend)
412  *  @param b  second polynomial in x (divisor)
413  *  @param x  a and b are polynomials in x
414  *  @param check_args  check whether a and b are polynomials with rational
415  *         coefficients (defaults to "true")
416  *  @return remainder of a(x) and b(x) in Q[x] */
417 ex rem(const ex &a, const ex &b, const ex &x, bool check_args)
418 {
419         if (b.is_zero())
420                 throw(std::overflow_error("rem: division by zero"));
421         if (is_exactly_a<numeric>(a)) {
422                 if  (is_exactly_a<numeric>(b))
423                         return _ex0;
424                 else
425                         return a;
426         }
427 #if FAST_COMPARE
428         if (a.is_equal(b))
429                 return _ex0;
430 #endif
431         if (check_args && (!a.info(info_flags::rational_polynomial) || !b.info(info_flags::rational_polynomial)))
432                 throw(std::invalid_argument("rem: arguments must be polynomials over the rationals"));
433
434         // Polynomial long division
435         ex r = a.expand();
436         if (r.is_zero())
437                 return r;
438         int bdeg = b.degree(x);
439         int rdeg = r.degree(x);
440         ex blcoeff = b.expand().coeff(x, bdeg);
441         bool blcoeff_is_numeric = is_exactly_a<numeric>(blcoeff);
442         while (rdeg >= bdeg) {
443                 ex term, rcoeff = r.coeff(x, rdeg);
444                 if (blcoeff_is_numeric)
445                         term = rcoeff / blcoeff;
446                 else {
447                         if (!divide(rcoeff, blcoeff, term, false))
448                                 return (new fail())->setflag(status_flags::dynallocated);
449                 }
450                 term *= power(x, rdeg - bdeg);
451                 r -= (term * b).expand();
452                 if (r.is_zero())
453                         break;
454                 rdeg = r.degree(x);
455         }
456         return r;
457 }
458
459
460 /** Decompose rational function a(x)=N(x)/D(x) into P(x)+n(x)/D(x)
461  *  with degree(n, x) < degree(D, x).
462  *
463  *  @param a rational function in x
464  *  @param x a is a function of x
465  *  @return decomposed function. */
466 ex decomp_rational(const ex &a, const ex &x)
467 {
468         ex nd = numer_denom(a);
469         ex numer = nd.op(0), denom = nd.op(1);
470         ex q = quo(numer, denom, x);
471         if (is_exactly_a<fail>(q))
472                 return a;
473         else
474                 return q + rem(numer, denom, x) / denom;
475 }
476
477
478 /** Pseudo-remainder of polynomials a(x) and b(x) in Q[x].
479  *
480  *  @param a  first polynomial in x (dividend)
481  *  @param b  second polynomial in x (divisor)
482  *  @param x  a and b are polynomials in x
483  *  @param check_args  check whether a and b are polynomials with rational
484  *         coefficients (defaults to "true")
485  *  @return pseudo-remainder of a(x) and b(x) in Q[x] */
486 ex prem(const ex &a, const ex &b, const ex &x, bool check_args)
487 {
488         if (b.is_zero())
489                 throw(std::overflow_error("prem: division by zero"));
490         if (is_exactly_a<numeric>(a)) {
491                 if (is_exactly_a<numeric>(b))
492                         return _ex0;
493                 else
494                         return b;
495         }
496         if (check_args && (!a.info(info_flags::rational_polynomial) || !b.info(info_flags::rational_polynomial)))
497                 throw(std::invalid_argument("prem: arguments must be polynomials over the rationals"));
498
499         // Polynomial long division
500         ex r = a.expand();
501         ex eb = b.expand();
502         int rdeg = r.degree(x);
503         int bdeg = eb.degree(x);
504         ex blcoeff;
505         if (bdeg <= rdeg) {
506                 blcoeff = eb.coeff(x, bdeg);
507                 if (bdeg == 0)
508                         eb = _ex0;
509                 else
510                         eb -= blcoeff * power(x, bdeg);
511         } else
512                 blcoeff = _ex1;
513
514         int delta = rdeg - bdeg + 1, i = 0;
515         while (rdeg >= bdeg && !r.is_zero()) {
516                 ex rlcoeff = r.coeff(x, rdeg);
517                 ex term = (power(x, rdeg - bdeg) * eb * rlcoeff).expand();
518                 if (rdeg == 0)
519                         r = _ex0;
520                 else
521                         r -= rlcoeff * power(x, rdeg);
522                 r = (blcoeff * r).expand() - term;
523                 rdeg = r.degree(x);
524                 i++;
525         }
526         return power(blcoeff, delta - i) * r;
527 }
528
529
530 /** Sparse pseudo-remainder of polynomials a(x) and b(x) in Q[x].
531  *
532  *  @param a  first polynomial in x (dividend)
533  *  @param b  second polynomial in x (divisor)
534  *  @param x  a and b are polynomials in x
535  *  @param check_args  check whether a and b are polynomials with rational
536  *         coefficients (defaults to "true")
537  *  @return sparse pseudo-remainder of a(x) and b(x) in Q[x] */
538 ex sprem(const ex &a, const ex &b, const ex &x, bool check_args)
539 {
540         if (b.is_zero())
541                 throw(std::overflow_error("prem: division by zero"));
542         if (is_exactly_a<numeric>(a)) {
543                 if (is_exactly_a<numeric>(b))
544                         return _ex0;
545                 else
546                         return b;
547         }
548         if (check_args && (!a.info(info_flags::rational_polynomial) || !b.info(info_flags::rational_polynomial)))
549                 throw(std::invalid_argument("prem: arguments must be polynomials over the rationals"));
550
551         // Polynomial long division
552         ex r = a.expand();
553         ex eb = b.expand();
554         int rdeg = r.degree(x);
555         int bdeg = eb.degree(x);
556         ex blcoeff;
557         if (bdeg <= rdeg) {
558                 blcoeff = eb.coeff(x, bdeg);
559                 if (bdeg == 0)
560                         eb = _ex0;
561                 else
562                         eb -= blcoeff * power(x, bdeg);
563         } else
564                 blcoeff = _ex1;
565
566         while (rdeg >= bdeg && !r.is_zero()) {
567                 ex rlcoeff = r.coeff(x, rdeg);
568                 ex term = (power(x, rdeg - bdeg) * eb * rlcoeff).expand();
569                 if (rdeg == 0)
570                         r = _ex0;
571                 else
572                         r -= rlcoeff * power(x, rdeg);
573                 r = (blcoeff * r).expand() - term;
574                 rdeg = r.degree(x);
575         }
576         return r;
577 }
578
579
580 /** Exact polynomial division of a(X) by b(X) in Q[X].
581  *  
582  *  @param a  first multivariate polynomial (dividend)
583  *  @param b  second multivariate polynomial (divisor)
584  *  @param q  quotient (returned)
585  *  @param check_args  check whether a and b are polynomials with rational
586  *         coefficients (defaults to "true")
587  *  @return "true" when exact division succeeds (quotient returned in q),
588  *          "false" otherwise (q left untouched) */
589 bool divide(const ex &a, const ex &b, ex &q, bool check_args)
590 {
591         if (b.is_zero())
592                 throw(std::overflow_error("divide: division by zero"));
593         if (a.is_zero()) {
594                 q = _ex0;
595                 return true;
596         }
597         if (is_exactly_a<numeric>(b)) {
598                 q = a / b;
599                 return true;
600         } else if (is_exactly_a<numeric>(a))
601                 return false;
602 #if FAST_COMPARE
603         if (a.is_equal(b)) {
604                 q = _ex1;
605                 return true;
606         }
607 #endif
608         if (check_args && (!a.info(info_flags::rational_polynomial) ||
609                            !b.info(info_flags::rational_polynomial)))
610                 throw(std::invalid_argument("divide: arguments must be polynomials over the rationals"));
611
612         // Find first symbol
613         ex x;
614         if (!get_first_symbol(a, x) && !get_first_symbol(b, x))
615                 throw(std::invalid_argument("invalid expression in divide()"));
616
617         // Try to avoid expanding partially factored expressions.
618         if (is_exactly_a<mul>(b)) {
619         // Divide sequentially by each term
620                 ex rem_new, rem_old = a;
621                 for (size_t i=0; i < b.nops(); i++) {
622                         if (! divide(rem_old, b.op(i), rem_new, false))
623                                 return false;
624                         rem_old = rem_new;
625                 }
626                 q = rem_new;
627                 return true;
628         } else if (is_exactly_a<power>(b)) {
629                 const ex& bb(b.op(0));
630                 int exp_b = ex_to<numeric>(b.op(1)).to_int();
631                 ex rem_new, rem_old = a;
632                 for (int i=exp_b; i>0; i--) {
633                         if (! divide(rem_old, bb, rem_new, false))
634                                 return false;
635                         rem_old = rem_new;
636                 }
637                 q = rem_new;
638                 return true;
639         } 
640         
641         if (is_exactly_a<mul>(a)) {
642                 // Divide sequentially each term. If some term in a is divisible 
643                 // by b we are done... and if not, we can't really say anything.
644                 size_t i;
645                 ex rem_i;
646                 bool divisible_p = false;
647                 for (i=0; i < a.nops(); ++i) {
648                         if (divide(a.op(i), b, rem_i, false)) {
649                                 divisible_p = true;
650                                 break;
651                         }
652                 }
653                 if (divisible_p) {
654                         exvector resv;
655                         resv.reserve(a.nops());
656                         for (size_t j=0; j < a.nops(); j++) {
657                                 if (j==i)
658                                         resv.push_back(rem_i);
659                                 else
660                                         resv.push_back(a.op(j));
661                         }
662                         q = (new mul(resv))->setflag(status_flags::dynallocated);
663                         return true;
664                 }
665         } else if (is_exactly_a<power>(a)) {
666                 // The base itself might be divisible by b, in that case we don't
667                 // need to expand a
668                 const ex& ab(a.op(0));
669                 int a_exp = ex_to<numeric>(a.op(1)).to_int();
670                 ex rem_i;
671                 if (divide(ab, b, rem_i, false)) {
672                         q = rem_i*power(ab, a_exp - 1);
673                         return true;
674                 }
675                 for (int i=2; i < a_exp; i++) {
676                         if (divide(power(ab, i), b, rem_i, false)) {
677                                 q = rem_i*power(ab, a_exp - i);
678                                 return true;
679                         }
680                 } // ... so we *really* need to expand expression.
681         }
682         
683         // Polynomial long division (recursive)
684         ex r = a.expand();
685         if (r.is_zero()) {
686                 q = _ex0;
687                 return true;
688         }
689         int bdeg = b.degree(x);
690         int rdeg = r.degree(x);
691         ex blcoeff = b.expand().coeff(x, bdeg);
692         bool blcoeff_is_numeric = is_exactly_a<numeric>(blcoeff);
693         exvector v; v.reserve(std::max(rdeg - bdeg + 1, 0));
694         while (rdeg >= bdeg) {
695                 ex term, rcoeff = r.coeff(x, rdeg);
696                 if (blcoeff_is_numeric)
697                         term = rcoeff / blcoeff;
698                 else
699                         if (!divide(rcoeff, blcoeff, term, false))
700                                 return false;
701                 term *= power(x, rdeg - bdeg);
702                 v.push_back(term);
703                 r -= (term * b).expand();
704                 if (r.is_zero()) {
705                         q = (new add(v))->setflag(status_flags::dynallocated);
706                         return true;
707                 }
708                 rdeg = r.degree(x);
709         }
710         return false;
711 }
712
713
714 #if USE_REMEMBER
715 /*
716  *  Remembering
717  */
718
719 typedef std::pair<ex, ex> ex2;
720 typedef std::pair<ex, bool> exbool;
721
722 struct ex2_less {
723         bool operator() (const ex2 &p, const ex2 &q) const 
724         {
725                 int cmp = p.first.compare(q.first);
726                 return ((cmp<0) || (!(cmp>0) && p.second.compare(q.second)<0));
727         }
728 };
729
730 typedef std::map<ex2, exbool, ex2_less> ex2_exbool_remember;
731 #endif
732
733
734 /** Exact polynomial division of a(X) by b(X) in Z[X].
735  *  This functions works like divide() but the input and output polynomials are
736  *  in Z[X] instead of Q[X] (i.e. they have integer coefficients). Unlike
737  *  divide(), it doesn't check whether the input polynomials really are integer
738  *  polynomials, so be careful of what you pass in. Also, you have to run
739  *  get_symbol_stats() over the input polynomials before calling this function
740  *  and pass an iterator to the first element of the sym_desc vector. This
741  *  function is used internally by the heur_gcd().
742  *  
743  *  @param a  first multivariate polynomial (dividend)
744  *  @param b  second multivariate polynomial (divisor)
745  *  @param q  quotient (returned)
746  *  @param var  iterator to first element of vector of sym_desc structs
747  *  @return "true" when exact division succeeds (the quotient is returned in
748  *          q), "false" otherwise.
749  *  @see get_symbol_stats, heur_gcd */
750 static bool divide_in_z(const ex &a, const ex &b, ex &q, sym_desc_vec::const_iterator var)
751 {
752         q = _ex0;
753         if (b.is_zero())
754                 throw(std::overflow_error("divide_in_z: division by zero"));
755         if (b.is_equal(_ex1)) {
756                 q = a;
757                 return true;
758         }
759         if (is_exactly_a<numeric>(a)) {
760                 if (is_exactly_a<numeric>(b)) {
761                         q = a / b;
762                         return q.info(info_flags::integer);
763                 } else
764                         return false;
765         }
766 #if FAST_COMPARE
767         if (a.is_equal(b)) {
768                 q = _ex1;
769                 return true;
770         }
771 #endif
772
773 #if USE_REMEMBER
774         // Remembering
775         static ex2_exbool_remember dr_remember;
776         ex2_exbool_remember::const_iterator remembered = dr_remember.find(ex2(a, b));
777         if (remembered != dr_remember.end()) {
778                 q = remembered->second.first;
779                 return remembered->second.second;
780         }
781 #endif
782
783         if (is_exactly_a<power>(b)) {
784                 const ex& bb(b.op(0));
785                 ex qbar = a;
786                 int exp_b = ex_to<numeric>(b.op(1)).to_int();
787                 for (int i=exp_b; i>0; i--) {
788                         if (!divide_in_z(qbar, bb, q, var))
789                                 return false;
790                         qbar = q;
791                 }
792                 return true;
793         }
794
795         if (is_exactly_a<mul>(b)) {
796                 ex qbar = a;
797                 for (const_iterator itrb = b.begin(); itrb != b.end(); ++itrb) {
798                         sym_desc_vec sym_stats;
799                         get_symbol_stats(a, *itrb, sym_stats);
800                         if (!divide_in_z(qbar, *itrb, q, sym_stats.begin()))
801                                 return false;
802
803                         qbar = q;
804                 }
805                 return true;
806         }
807
808         // Main symbol
809         const ex &x = var->sym;
810
811         // Compare degrees
812         int adeg = a.degree(x), bdeg = b.degree(x);
813         if (bdeg > adeg)
814                 return false;
815
816 #if USE_TRIAL_DIVISION
817
818         // Trial division with polynomial interpolation
819         int i, k;
820
821         // Compute values at evaluation points 0..adeg
822         vector<numeric> alpha; alpha.reserve(adeg + 1);
823         exvector u; u.reserve(adeg + 1);
824         numeric point = *_num0_p;
825         ex c;
826         for (i=0; i<=adeg; i++) {
827                 ex bs = b.subs(x == point, subs_options::no_pattern);
828                 while (bs.is_zero()) {
829                         point += *_num1_p;
830                         bs = b.subs(x == point, subs_options::no_pattern);
831                 }
832                 if (!divide_in_z(a.subs(x == point, subs_options::no_pattern), bs, c, var+1))
833                         return false;
834                 alpha.push_back(point);
835                 u.push_back(c);
836                 point += *_num1_p;
837         }
838
839         // Compute inverses
840         vector<numeric> rcp; rcp.reserve(adeg + 1);
841         rcp.push_back(*_num0_p);
842         for (k=1; k<=adeg; k++) {
843                 numeric product = alpha[k] - alpha[0];
844                 for (i=1; i<k; i++)
845                         product *= alpha[k] - alpha[i];
846                 rcp.push_back(product.inverse());
847         }
848
849         // Compute Newton coefficients
850         exvector v; v.reserve(adeg + 1);
851         v.push_back(u[0]);
852         for (k=1; k<=adeg; k++) {
853                 ex temp = v[k - 1];
854                 for (i=k-2; i>=0; i--)
855                         temp = temp * (alpha[k] - alpha[i]) + v[i];
856                 v.push_back((u[k] - temp) * rcp[k]);
857         }
858
859         // Convert from Newton form to standard form
860         c = v[adeg];
861         for (k=adeg-1; k>=0; k--)
862                 c = c * (x - alpha[k]) + v[k];
863
864         if (c.degree(x) == (adeg - bdeg)) {
865                 q = c.expand();
866                 return true;
867         } else
868                 return false;
869
870 #else
871
872         // Polynomial long division (recursive)
873         ex r = a.expand();
874         if (r.is_zero())
875                 return true;
876         int rdeg = adeg;
877         ex eb = b.expand();
878         ex blcoeff = eb.coeff(x, bdeg);
879         exvector v; v.reserve(std::max(rdeg - bdeg + 1, 0));
880         while (rdeg >= bdeg) {
881                 ex term, rcoeff = r.coeff(x, rdeg);
882                 if (!divide_in_z(rcoeff, blcoeff, term, var+1))
883                         break;
884                 term = (term * power(x, rdeg - bdeg)).expand();
885                 v.push_back(term);
886                 r -= (term * eb).expand();
887                 if (r.is_zero()) {
888                         q = (new add(v))->setflag(status_flags::dynallocated);
889 #if USE_REMEMBER
890                         dr_remember[ex2(a, b)] = exbool(q, true);
891 #endif
892                         return true;
893                 }
894                 rdeg = r.degree(x);
895         }
896 #if USE_REMEMBER
897         dr_remember[ex2(a, b)] = exbool(q, false);
898 #endif
899         return false;
900
901 #endif
902 }
903
904
905 /*
906  *  Separation of unit part, content part and primitive part of polynomials
907  */
908
909 /** Compute unit part (= sign of leading coefficient) of a multivariate
910  *  polynomial in Q[x]. The product of unit part, content part, and primitive
911  *  part is the polynomial itself.
912  *
913  *  @param x  main variable
914  *  @return unit part
915  *  @see ex::content, ex::primpart, ex::unitcontprim */
916 ex ex::unit(const ex &x) const
917 {
918         ex c = expand().lcoeff(x);
919         if (is_exactly_a<numeric>(c))
920                 return c.info(info_flags::negative) ?_ex_1 : _ex1;
921         else {
922                 ex y;
923                 if (get_first_symbol(c, y))
924                         return c.unit(y);
925                 else
926                         throw(std::invalid_argument("invalid expression in unit()"));
927         }
928 }
929
930
931 /** Compute content part (= unit normal GCD of all coefficients) of a
932  *  multivariate polynomial in Q[x]. The product of unit part, content part,
933  *  and primitive part is the polynomial itself.
934  *
935  *  @param x  main variable
936  *  @return content part
937  *  @see ex::unit, ex::primpart, ex::unitcontprim */
938 ex ex::content(const ex &x) const
939 {
940         if (is_exactly_a<numeric>(*this))
941                 return info(info_flags::negative) ? -*this : *this;
942
943         ex e = expand();
944         if (e.is_zero())
945                 return _ex0;
946
947         // First, divide out the integer content (which we can calculate very efficiently).
948         // If the leading coefficient of the quotient is an integer, we are done.
949         ex c = e.integer_content();
950         ex r = e / c;
951         int deg = r.degree(x);
952         ex lcoeff = r.coeff(x, deg);
953         if (lcoeff.info(info_flags::integer))
954                 return c;
955
956         // GCD of all coefficients
957         int ldeg = r.ldegree(x);
958         if (deg == ldeg)
959                 return lcoeff * c / lcoeff.unit(x);
960         ex cont = _ex0;
961         for (int i=ldeg; i<=deg; i++)
962                 cont = gcd(r.coeff(x, i), cont, NULL, NULL, false);
963         return cont * c;
964 }
965
966
967 /** Compute primitive part of a multivariate polynomial in Q[x]. The result
968  *  will be a unit-normal polynomial with a content part of 1. The product
969  *  of unit part, content part, and primitive part is the polynomial itself.
970  *
971  *  @param x  main variable
972  *  @return primitive part
973  *  @see ex::unit, ex::content, ex::unitcontprim */
974 ex ex::primpart(const ex &x) const
975 {
976         // We need to compute the unit and content anyway, so call unitcontprim()
977         ex u, c, p;
978         unitcontprim(x, u, c, p);
979         return p;
980 }
981
982
983 /** Compute primitive part of a multivariate polynomial in Q[x] when the
984  *  content part is already known. This function is faster in computing the
985  *  primitive part than the previous function.
986  *
987  *  @param x  main variable
988  *  @param c  previously computed content part
989  *  @return primitive part */
990 ex ex::primpart(const ex &x, const ex &c) const
991 {
992         if (is_zero() || c.is_zero())
993                 return _ex0;
994         if (is_exactly_a<numeric>(*this))
995                 return _ex1;
996
997         // Divide by unit and content to get primitive part
998         ex u = unit(x);
999         if (is_exactly_a<numeric>(c))
1000                 return *this / (c * u);
1001         else
1002                 return quo(*this, c * u, x, false);
1003 }
1004
1005
1006 /** Compute unit part, content part, and primitive part of a multivariate
1007  *  polynomial in Q[x]. The product of the three parts is the polynomial
1008  *  itself.
1009  *
1010  *  @param x  main variable
1011  *  @param u  unit part (returned)
1012  *  @param c  content part (returned)
1013  *  @param p  primitive part (returned)
1014  *  @see ex::unit, ex::content, ex::primpart */
1015 void ex::unitcontprim(const ex &x, ex &u, ex &c, ex &p) const
1016 {
1017         // Quick check for zero (avoid expanding)
1018         if (is_zero()) {
1019                 u = _ex1;
1020                 c = p = _ex0;
1021                 return;
1022         }
1023
1024         // Special case: input is a number
1025         if (is_exactly_a<numeric>(*this)) {
1026                 if (info(info_flags::negative)) {
1027                         u = _ex_1;
1028                         c = abs(ex_to<numeric>(*this));
1029                 } else {
1030                         u = _ex1;
1031                         c = *this;
1032                 }
1033                 p = _ex1;
1034                 return;
1035         }
1036
1037         // Expand input polynomial
1038         ex e = expand();
1039         if (e.is_zero()) {
1040                 u = _ex1;
1041                 c = p = _ex0;
1042                 return;
1043         }
1044
1045         // Compute unit and content
1046         u = unit(x);
1047         c = content(x);
1048
1049         // Divide by unit and content to get primitive part
1050         if (c.is_zero()) {
1051                 p = _ex0;
1052                 return;
1053         }
1054         if (is_exactly_a<numeric>(c))
1055                 p = *this / (c * u);
1056         else
1057                 p = quo(e, c * u, x, false);
1058 }
1059
1060
1061 /*
1062  *  GCD of multivariate polynomials
1063  */
1064
1065 /** Compute GCD of multivariate polynomials using the subresultant PRS
1066  *  algorithm. This function is used internally by gcd().
1067  *
1068  *  @param a   first multivariate polynomial
1069  *  @param b   second multivariate polynomial
1070  *  @param var iterator to first element of vector of sym_desc structs
1071  *  @return the GCD as a new expression
1072  *  @see gcd */
1073
1074 static ex sr_gcd(const ex &a, const ex &b, sym_desc_vec::const_iterator var)
1075 {
1076 #if STATISTICS
1077         sr_gcd_called++;
1078 #endif
1079
1080         // The first symbol is our main variable
1081         const ex &x = var->sym;
1082
1083         // Sort c and d so that c has higher degree
1084         ex c, d;
1085         int adeg = a.degree(x), bdeg = b.degree(x);
1086         int cdeg, ddeg;
1087         if (adeg >= bdeg) {
1088                 c = a;
1089                 d = b;
1090                 cdeg = adeg;
1091                 ddeg = bdeg;
1092         } else {
1093                 c = b;
1094                 d = a;
1095                 cdeg = bdeg;
1096                 ddeg = adeg;
1097         }
1098
1099         // Remove content from c and d, to be attached to GCD later
1100         ex cont_c = c.content(x);
1101         ex cont_d = d.content(x);
1102         ex gamma = gcd(cont_c, cont_d, NULL, NULL, false);
1103         if (ddeg == 0)
1104                 return gamma;
1105         c = c.primpart(x, cont_c);
1106         d = d.primpart(x, cont_d);
1107
1108         // First element of subresultant sequence
1109         ex r = _ex0, ri = _ex1, psi = _ex1;
1110         int delta = cdeg - ddeg;
1111
1112         for (;;) {
1113
1114                 // Calculate polynomial pseudo-remainder
1115                 r = prem(c, d, x, false);
1116                 if (r.is_zero())
1117                         return gamma * d.primpart(x);
1118
1119                 c = d;
1120                 cdeg = ddeg;
1121                 if (!divide_in_z(r, ri * pow(psi, delta), d, var))
1122                         throw(std::runtime_error("invalid expression in sr_gcd(), division failed"));
1123                 ddeg = d.degree(x);
1124                 if (ddeg == 0) {
1125                         if (is_exactly_a<numeric>(r))
1126                                 return gamma;
1127                         else
1128                                 return gamma * r.primpart(x);
1129                 }
1130
1131                 // Next element of subresultant sequence
1132                 ri = c.expand().lcoeff(x);
1133                 if (delta == 1)
1134                         psi = ri;
1135                 else if (delta)
1136                         divide_in_z(pow(ri, delta), pow(psi, delta-1), psi, var+1);
1137                 delta = cdeg - ddeg;
1138         }
1139 }
1140
1141
1142 /** Return maximum (absolute value) coefficient of a polynomial.
1143  *  This function is used internally by heur_gcd().
1144  *
1145  *  @return maximum coefficient
1146  *  @see heur_gcd */
1147 numeric ex::max_coefficient() const
1148 {
1149         return bp->max_coefficient();
1150 }
1151
1152 /** Implementation ex::max_coefficient().
1153  *  @see heur_gcd */
1154 numeric basic::max_coefficient() const
1155 {
1156         return *_num1_p;
1157 }
1158
1159 numeric numeric::max_coefficient() const
1160 {
1161         return abs(*this);
1162 }
1163
1164 numeric add::max_coefficient() const
1165 {
1166         epvector::const_iterator it = seq.begin();
1167         epvector::const_iterator itend = seq.end();
1168         GINAC_ASSERT(is_exactly_a<numeric>(overall_coeff));
1169         numeric cur_max = abs(ex_to<numeric>(overall_coeff));
1170         while (it != itend) {
1171                 numeric a;
1172                 GINAC_ASSERT(!is_exactly_a<numeric>(it->rest));
1173                 a = abs(ex_to<numeric>(it->coeff));
1174                 if (a > cur_max)
1175                         cur_max = a;
1176                 it++;
1177         }
1178         return cur_max;
1179 }
1180
1181 numeric mul::max_coefficient() const
1182 {
1183 #ifdef DO_GINAC_ASSERT
1184         epvector::const_iterator it = seq.begin();
1185         epvector::const_iterator itend = seq.end();
1186         while (it != itend) {
1187                 GINAC_ASSERT(!is_exactly_a<numeric>(recombine_pair_to_ex(*it)));
1188                 it++;
1189         }
1190 #endif // def DO_GINAC_ASSERT
1191         GINAC_ASSERT(is_exactly_a<numeric>(overall_coeff));
1192         return abs(ex_to<numeric>(overall_coeff));
1193 }
1194
1195
1196 /** Apply symmetric modular homomorphism to an expanded multivariate
1197  *  polynomial.  This function is usually used internally by heur_gcd().
1198  *
1199  *  @param xi  modulus
1200  *  @return mapped polynomial
1201  *  @see heur_gcd */
1202 ex basic::smod(const numeric &xi) const
1203 {
1204         return *this;
1205 }
1206
1207 ex numeric::smod(const numeric &xi) const
1208 {
1209         return GiNaC::smod(*this, xi);
1210 }
1211
1212 ex add::smod(const numeric &xi) const
1213 {
1214         epvector newseq;
1215         newseq.reserve(seq.size()+1);
1216         epvector::const_iterator it = seq.begin();
1217         epvector::const_iterator itend = seq.end();
1218         while (it != itend) {
1219                 GINAC_ASSERT(!is_exactly_a<numeric>(it->rest));
1220                 numeric coeff = GiNaC::smod(ex_to<numeric>(it->coeff), xi);
1221                 if (!coeff.is_zero())
1222                         newseq.push_back(expair(it->rest, coeff));
1223                 it++;
1224         }
1225         GINAC_ASSERT(is_exactly_a<numeric>(overall_coeff));
1226         numeric coeff = GiNaC::smod(ex_to<numeric>(overall_coeff), xi);
1227         return (new add(newseq,coeff))->setflag(status_flags::dynallocated);
1228 }
1229
1230 ex mul::smod(const numeric &xi) const
1231 {
1232 #ifdef DO_GINAC_ASSERT
1233         epvector::const_iterator it = seq.begin();
1234         epvector::const_iterator itend = seq.end();
1235         while (it != itend) {
1236                 GINAC_ASSERT(!is_exactly_a<numeric>(recombine_pair_to_ex(*it)));
1237                 it++;
1238         }
1239 #endif // def DO_GINAC_ASSERT
1240         mul * mulcopyp = new mul(*this);
1241         GINAC_ASSERT(is_exactly_a<numeric>(overall_coeff));
1242         mulcopyp->overall_coeff = GiNaC::smod(ex_to<numeric>(overall_coeff),xi);
1243         mulcopyp->clearflag(status_flags::evaluated);
1244         mulcopyp->clearflag(status_flags::hash_calculated);
1245         return mulcopyp->setflag(status_flags::dynallocated);
1246 }
1247
1248
1249 /** xi-adic polynomial interpolation */
1250 static ex interpolate(const ex &gamma, const numeric &xi, const ex &x, int degree_hint = 1)
1251 {
1252         exvector g; g.reserve(degree_hint);
1253         ex e = gamma;
1254         numeric rxi = xi.inverse();
1255         for (int i=0; !e.is_zero(); i++) {
1256                 ex gi = e.smod(xi);
1257                 g.push_back(gi * power(x, i));
1258                 e = (e - gi) * rxi;
1259         }
1260         return (new add(g))->setflag(status_flags::dynallocated);
1261 }
1262
1263 /** Exception thrown by heur_gcd() to signal failure. */
1264 class gcdheu_failed {};
1265
1266 /** Compute GCD of multivariate polynomials using the heuristic GCD algorithm.
1267  *  get_symbol_stats() must have been called previously with the input
1268  *  polynomials and an iterator to the first element of the sym_desc vector
1269  *  passed in. This function is used internally by gcd().
1270  *
1271  *  @param a  first integer multivariate polynomial (expanded)
1272  *  @param b  second integer multivariate polynomial (expanded)
1273  *  @param ca  cofactor of polynomial a (returned), NULL to suppress
1274  *             calculation of cofactor
1275  *  @param cb  cofactor of polynomial b (returned), NULL to suppress
1276  *             calculation of cofactor
1277  *  @param var iterator to first element of vector of sym_desc structs
1278  *  @param res the GCD (returned)
1279  *  @return true if GCD was computed, false otherwise.
1280  *  @see gcd
1281  *  @exception gcdheu_failed() */
1282 static bool heur_gcd_z(ex& res, const ex &a, const ex &b, ex *ca, ex *cb,
1283                        sym_desc_vec::const_iterator var)
1284 {
1285 #if STATISTICS
1286         heur_gcd_called++;
1287 #endif
1288
1289         // Algorithm only works for non-vanishing input polynomials
1290         if (a.is_zero() || b.is_zero())
1291                 return false;
1292
1293         // GCD of two numeric values -> CLN
1294         if (is_exactly_a<numeric>(a) && is_exactly_a<numeric>(b)) {
1295                 numeric g = gcd(ex_to<numeric>(a), ex_to<numeric>(b));
1296                 if (ca)
1297                         *ca = ex_to<numeric>(a) / g;
1298                 if (cb)
1299                         *cb = ex_to<numeric>(b) / g;
1300                 res = g;
1301                 return true;
1302         }
1303
1304         // The first symbol is our main variable
1305         const ex &x = var->sym;
1306
1307         // Remove integer content
1308         numeric gc = gcd(a.integer_content(), b.integer_content());
1309         numeric rgc = gc.inverse();
1310         ex p = a * rgc;
1311         ex q = b * rgc;
1312         int maxdeg =  std::max(p.degree(x), q.degree(x));
1313         
1314         // Find evaluation point
1315         numeric mp = p.max_coefficient();
1316         numeric mq = q.max_coefficient();
1317         numeric xi;
1318         if (mp > mq)
1319                 xi = mq * (*_num2_p) + (*_num2_p);
1320         else
1321                 xi = mp * (*_num2_p) + (*_num2_p);
1322
1323         // 6 tries maximum
1324         for (int t=0; t<6; t++) {
1325                 if (xi.int_length() * maxdeg > 100000) {
1326                         throw gcdheu_failed();
1327                 }
1328
1329                 // Apply evaluation homomorphism and calculate GCD
1330                 ex cp, cq;
1331                 ex gamma;
1332                 bool found = heur_gcd_z(gamma,
1333                                         p.subs(x == xi, subs_options::no_pattern),
1334                                         q.subs(x == xi, subs_options::no_pattern),
1335                                         &cp, &cq, var+1);
1336                 if (found) {
1337                         gamma = gamma.expand();
1338                         // Reconstruct polynomial from GCD of mapped polynomials
1339                         ex g = interpolate(gamma, xi, x, maxdeg);
1340
1341                         // Remove integer content
1342                         g /= g.integer_content();
1343
1344                         // If the calculated polynomial divides both p and q, this is the GCD
1345                         ex dummy;
1346                         if (divide_in_z(p, g, ca ? *ca : dummy, var) && divide_in_z(q, g, cb ? *cb : dummy, var)) {
1347                                 g *= gc;
1348                                 res = g;
1349                                 return true;
1350                         }
1351                 }
1352
1353                 // Next evaluation point
1354                 xi = iquo(xi * isqrt(isqrt(xi)) * numeric(73794), numeric(27011));
1355         }
1356         return false;
1357 }
1358
1359 /** Compute GCD of multivariate polynomials using the heuristic GCD algorithm.
1360  *  get_symbol_stats() must have been called previously with the input
1361  *  polynomials and an iterator to the first element of the sym_desc vector
1362  *  passed in. This function is used internally by gcd().
1363  *
1364  *  @param a  first rational multivariate polynomial (expanded)
1365  *  @param b  second rational multivariate polynomial (expanded)
1366  *  @param ca  cofactor of polynomial a (returned), NULL to suppress
1367  *             calculation of cofactor
1368  *  @param cb  cofactor of polynomial b (returned), NULL to suppress
1369  *             calculation of cofactor
1370  *  @param var iterator to first element of vector of sym_desc structs
1371  *  @param res the GCD (returned)
1372  *  @return true if GCD was computed, false otherwise.
1373  *  @see heur_gcd_z
1374  *  @see gcd
1375  */
1376 static bool heur_gcd(ex& res, const ex& a, const ex& b, ex *ca, ex *cb,
1377                      sym_desc_vec::const_iterator var)
1378 {
1379         if (a.info(info_flags::integer_polynomial) && 
1380             b.info(info_flags::integer_polynomial)) {
1381                 try {
1382                         return heur_gcd_z(res, a, b, ca, cb, var);
1383                 } catch (gcdheu_failed) {
1384                         return false;
1385                 }
1386         }
1387
1388         // convert polynomials to Z[X]
1389         const numeric a_lcm = lcm_of_coefficients_denominators(a);
1390         const numeric ab_lcm = lcmcoeff(b, a_lcm);
1391
1392         const ex ai = a*ab_lcm;
1393         const ex bi = b*ab_lcm;
1394         if (!ai.info(info_flags::integer_polynomial))
1395                 throw std::logic_error("heur_gcd: not an integer polynomial [1]");
1396
1397         if (!bi.info(info_flags::integer_polynomial))
1398                 throw std::logic_error("heur_gcd: not an integer polynomial [2]");
1399
1400         bool found = false;
1401         try {
1402                 found = heur_gcd_z(res, ai, bi, ca, cb, var);
1403         } catch (gcdheu_failed) {
1404                 return false;
1405         }
1406         
1407         // GCD is not unique, it's defined up to a unit (i.e. invertible
1408         // element). If the coefficient ring is a field, every its element is
1409         // invertible, so one can multiply the polynomial GCD with any element
1410         // of the coefficient field. We use this ambiguity to make cofactors
1411         // integer polynomials.
1412         if (found)
1413                 res /= ab_lcm;
1414         return found;
1415 }
1416
1417
1418 // gcd helper to handle partially factored polynomials (to avoid expanding
1419 // large expressions). At least one of the arguments should be a power.
1420 static ex gcd_pf_pow(const ex& a, const ex& b, ex* ca, ex* cb, bool check_args);
1421
1422 // gcd helper to handle partially factored polynomials (to avoid expanding
1423 // large expressions). At least one of the arguments should be a product.
1424 static ex gcd_pf_mul(const ex& a, const ex& b, ex* ca, ex* cb, bool check_args);
1425
1426 /** Compute GCD (Greatest Common Divisor) of multivariate polynomials a(X)
1427  *  and b(X) in Z[X]. Optionally also compute the cofactors of a and b,
1428  *  defined by a = ca * gcd(a, b) and b = cb * gcd(a, b).
1429  *
1430  *  @param a  first multivariate polynomial
1431  *  @param b  second multivariate polynomial
1432  *  @param ca pointer to expression that will receive the cofactor of a, or NULL
1433  *  @param cb pointer to expression that will receive the cofactor of b, or NULL
1434  *  @param check_args  check whether a and b are polynomials with rational
1435  *         coefficients (defaults to "true")
1436  *  @return the GCD as a new expression */
1437 ex gcd(const ex &a, const ex &b, ex *ca, ex *cb, bool check_args, unsigned options)
1438 {
1439 #if STATISTICS
1440         gcd_called++;
1441 #endif
1442
1443         // GCD of numerics -> CLN
1444         if (is_exactly_a<numeric>(a) && is_exactly_a<numeric>(b)) {
1445                 numeric g = gcd(ex_to<numeric>(a), ex_to<numeric>(b));
1446                 if (ca || cb) {
1447                         if (g.is_zero()) {
1448                                 if (ca)
1449                                         *ca = _ex0;
1450                                 if (cb)
1451                                         *cb = _ex0;
1452                         } else {
1453                                 if (ca)
1454                                         *ca = ex_to<numeric>(a) / g;
1455                                 if (cb)
1456                                         *cb = ex_to<numeric>(b) / g;
1457                         }
1458                 }
1459                 return g;
1460         }
1461
1462         // Check arguments
1463         if (check_args && (!a.info(info_flags::rational_polynomial) || !b.info(info_flags::rational_polynomial))) {
1464                 throw(std::invalid_argument("gcd: arguments must be polynomials over the rationals"));
1465         }
1466
1467         // Partially factored cases (to avoid expanding large expressions)
1468         if (!(options & gcd_options::no_part_factored)) {
1469                 if (is_exactly_a<mul>(a) || is_exactly_a<mul>(b))
1470                         return gcd_pf_mul(a, b, ca, cb, check_args);
1471 #if FAST_COMPARE
1472                 if (is_exactly_a<power>(a) || is_exactly_a<power>(b))
1473                         return gcd_pf_pow(a, b, ca, cb, check_args);
1474 #endif
1475         }
1476
1477         // Some trivial cases
1478         ex aex = a.expand(), bex = b.expand();
1479         if (aex.is_zero()) {
1480                 if (ca)
1481                         *ca = _ex0;
1482                 if (cb)
1483                         *cb = _ex1;
1484                 return b;
1485         }
1486         if (bex.is_zero()) {
1487                 if (ca)
1488                         *ca = _ex1;
1489                 if (cb)
1490                         *cb = _ex0;
1491                 return a;
1492         }
1493         if (aex.is_equal(_ex1) || bex.is_equal(_ex1)) {
1494                 if (ca)
1495                         *ca = a;
1496                 if (cb)
1497                         *cb = b;
1498                 return _ex1;
1499         }
1500 #if FAST_COMPARE
1501         if (a.is_equal(b)) {
1502                 if (ca)
1503                         *ca = _ex1;
1504                 if (cb)
1505                         *cb = _ex1;
1506                 return a;
1507         }
1508 #endif
1509
1510         if (is_a<symbol>(aex)) {
1511                 if (! bex.subs(aex==_ex0, subs_options::no_pattern).is_zero()) {
1512                         if (ca)
1513                                 *ca = a;
1514                         if (cb)
1515                                 *cb = b;
1516                         return _ex1;
1517                 }
1518         }
1519
1520         if (is_a<symbol>(bex)) {
1521                 if (! aex.subs(bex==_ex0, subs_options::no_pattern).is_zero()) {
1522                         if (ca)
1523                                 *ca = a;
1524                         if (cb)
1525                                 *cb = b;
1526                         return _ex1;
1527                 }
1528         }
1529
1530         if (is_exactly_a<numeric>(aex)) {
1531                 numeric bcont = bex.integer_content();
1532                 numeric g = gcd(ex_to<numeric>(aex), bcont);
1533                 if (ca)
1534                         *ca = ex_to<numeric>(aex)/g;
1535                 if (cb)
1536                         *cb = bex/g;
1537                 return g;
1538         }
1539
1540         if (is_exactly_a<numeric>(bex)) {
1541                 numeric acont = aex.integer_content();
1542                 numeric g = gcd(ex_to<numeric>(bex), acont);
1543                 if (ca)
1544                         *ca = aex/g;
1545                 if (cb)
1546                         *cb = ex_to<numeric>(bex)/g;
1547                 return g;
1548         }
1549
1550         // Gather symbol statistics
1551         sym_desc_vec sym_stats;
1552         get_symbol_stats(a, b, sym_stats);
1553
1554         // The symbol with least degree which is contained in both polynomials
1555         // is our main variable
1556         sym_desc_vec::iterator vari = sym_stats.begin();
1557         while ((vari != sym_stats.end()) && 
1558                (((vari->ldeg_b == 0) && (vari->deg_b == 0)) ||
1559                 ((vari->ldeg_a == 0) && (vari->deg_a == 0))))
1560                 vari++;
1561
1562         // No common symbols at all, just return 1:
1563         if (vari == sym_stats.end()) {
1564                 // N.B: keep cofactors factored
1565                 if (ca)
1566                         *ca = a;
1567                 if (cb)
1568                         *cb = b;
1569                 return _ex1;
1570         }
1571         // move symbols which contained only in one of the polynomials
1572         // to the end:
1573         rotate(sym_stats.begin(), vari, sym_stats.end());
1574
1575         sym_desc_vec::const_iterator var = sym_stats.begin();
1576         const ex &x = var->sym;
1577
1578         // Cancel trivial common factor
1579         int ldeg_a = var->ldeg_a;
1580         int ldeg_b = var->ldeg_b;
1581         int min_ldeg = std::min(ldeg_a,ldeg_b);
1582         if (min_ldeg > 0) {
1583                 ex common = power(x, min_ldeg);
1584                 return gcd((aex / common).expand(), (bex / common).expand(), ca, cb, false) * common;
1585         }
1586
1587         // Try to eliminate variables
1588         if (var->deg_a == 0 && var->deg_b != 0 ) {
1589                 ex bex_u, bex_c, bex_p;
1590                 bex.unitcontprim(x, bex_u, bex_c, bex_p);
1591                 ex g = gcd(aex, bex_c, ca, cb, false);
1592                 if (cb)
1593                         *cb *= bex_u * bex_p;
1594                 return g;
1595         } else if (var->deg_b == 0 && var->deg_a != 0) {
1596                 ex aex_u, aex_c, aex_p;
1597                 aex.unitcontprim(x, aex_u, aex_c, aex_p);
1598                 ex g = gcd(aex_c, bex, ca, cb, false);
1599                 if (ca)
1600                         *ca *= aex_u * aex_p;
1601                 return g;
1602         }
1603
1604         // Try heuristic algorithm first, fall back to PRS if that failed
1605         ex g;
1606         if (!(options & gcd_options::no_heur_gcd)) {
1607                 bool found = heur_gcd(g, aex, bex, ca, cb, var);
1608                 if (found) {
1609                         // heur_gcd have already computed cofactors...
1610                         if (g.is_equal(_ex1)) {
1611                                 // ... but we want to keep them factored if possible.
1612                                 if (ca)
1613                                         *ca = a;
1614                                 if (cb)
1615                                         *cb = b;
1616                         }
1617                         return g;
1618                 }
1619 #if STATISTICS
1620                 else {
1621                         heur_gcd_failed++;
1622                 }
1623 #endif
1624         }
1625
1626         g = sr_gcd(aex, bex, var);
1627         if (g.is_equal(_ex1)) {
1628                 // Keep cofactors factored if possible
1629                 if (ca)
1630                         *ca = a;
1631                 if (cb)
1632                         *cb = b;
1633         } else {
1634                 if (ca)
1635                         divide(aex, g, *ca, false);
1636                 if (cb)
1637                         divide(bex, g, *cb, false);
1638         }
1639         return g;
1640 }
1641
1642 static ex gcd_pf_pow(const ex& a, const ex& b, ex* ca, ex* cb, bool check_args)
1643 {
1644         if (is_exactly_a<power>(a)) {
1645                 ex p = a.op(0);
1646                 const ex& exp_a = a.op(1);
1647                 if (is_exactly_a<power>(b)) {
1648                         ex pb = b.op(0);
1649                         const ex& exp_b = b.op(1);
1650                         if (p.is_equal(pb)) {
1651                                 // a = p^n, b = p^m, gcd = p^min(n, m)
1652                                 if (exp_a < exp_b) {
1653                                         if (ca)
1654                                                 *ca = _ex1;
1655                                         if (cb)
1656                                                 *cb = power(p, exp_b - exp_a);
1657                                         return power(p, exp_a);
1658                                 } else {
1659                                         if (ca)
1660                                                 *ca = power(p, exp_a - exp_b);
1661                                         if (cb)
1662                                                 *cb = _ex1;
1663                                         return power(p, exp_b);
1664                                 }
1665                         } else {
1666                                 ex p_co, pb_co;
1667                                 ex p_gcd = gcd(p, pb, &p_co, &pb_co, check_args);
1668                                 if (p_gcd.is_equal(_ex1)) {
1669                                         // a(x) = p(x)^n, b(x) = p_b(x)^m, gcd (p, p_b) = 1 ==>
1670                                         // gcd(a,b) = 1
1671                                         if (ca)
1672                                                 *ca = a;
1673                                         if (cb)
1674                                                 *cb = b;
1675                                         return _ex1;
1676                                         // XXX: do I need to check for p_gcd = -1?
1677                                 } else {
1678                                         // there are common factors:
1679                                         // a(x) = g(x)^n A(x)^n, b(x) = g(x)^m B(x)^m ==>
1680                                         // gcd(a, b) = g(x)^n gcd(A(x)^n, g(x)^(n-m) B(x)^m
1681                                         if (exp_a < exp_b) {
1682                                                 return power(p_gcd, exp_a)*
1683                                                         gcd(power(p_co, exp_a), power(p_gcd, exp_b-exp_a)*power(pb_co, exp_b), ca, cb, false);
1684                                         } else {
1685                                                 return power(p_gcd, exp_b)*
1686                                                         gcd(power(p_gcd, exp_a - exp_b)*power(p_co, exp_a), power(pb_co, exp_b), ca, cb, false);
1687                                         }
1688                                 } // p_gcd.is_equal(_ex1)
1689                         } // p.is_equal(pb)
1690
1691                 } else {
1692                         if (p.is_equal(b)) {
1693                                 // a = p^n, b = p, gcd = p
1694                                 if (ca)
1695                                         *ca = power(p, a.op(1) - 1);
1696                                 if (cb)
1697                                         *cb = _ex1;
1698                                 return p;
1699                         } 
1700
1701                         ex p_co, bpart_co;
1702                         ex p_gcd = gcd(p, b, &p_co, &bpart_co, false);
1703
1704                         if (p_gcd.is_equal(_ex1)) {
1705                                 // a(x) = p(x)^n, gcd(p, b) = 1 ==> gcd(a, b) = 1
1706                                 if (ca)
1707                                         *ca = a;
1708                                 if (cb)
1709                                         *cb = b;
1710                                 return _ex1;
1711                         } else {
1712                                 // 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))
1713                                 return p_gcd*gcd(power(p_gcd, exp_a-1)*power(p_co, exp_a), bpart_co, ca, cb, false);
1714                         }
1715                 } // is_exactly_a<power>(b)
1716
1717         } else if (is_exactly_a<power>(b)) {
1718                 ex p = b.op(0);
1719                 if (p.is_equal(a)) {
1720                         // a = p, b = p^n, gcd = p
1721                         if (ca)
1722                                 *ca = _ex1;
1723                         if (cb)
1724                                 *cb = power(p, b.op(1) - 1);
1725                         return p;
1726                 }
1727
1728                 ex p_co, apart_co;
1729                 const ex& exp_b(b.op(1));
1730                 ex p_gcd = gcd(a, p, &apart_co, &p_co, false);
1731                 if (p_gcd.is_equal(_ex1)) {
1732                         // b=p(x)^n, gcd(a, p) = 1 ==> gcd(a, b) == 1
1733                         if (ca)
1734                                 *ca = a;
1735                         if (cb)
1736                                 *cb = b;
1737                         return _ex1;
1738                 } else {
1739                         // there are common factors:
1740                         // a(x) = g(x) A(x), b(x) = g(x)^n B(x)^n ==> gcd = g(x) gcd(g(x)^(n-1) A(x)^n, B(x))
1741
1742                         return p_gcd*gcd(apart_co, power(p_gcd, exp_b-1)*power(p_co, exp_b), ca, cb, false);
1743                 } // p_gcd.is_equal(_ex1)
1744         }
1745 }
1746
1747 static ex gcd_pf_mul(const ex& a, const ex& b, ex* ca, ex* cb, bool check_args)
1748 {
1749         if (is_exactly_a<mul>(a)) {
1750                 if (is_exactly_a<mul>(b) && b.nops() > a.nops())
1751                         goto factored_b;
1752 factored_a:
1753                 size_t num = a.nops();
1754                 exvector g; g.reserve(num);
1755                 exvector acc_ca; acc_ca.reserve(num);
1756                 ex part_b = b;
1757                 for (size_t i=0; i<num; i++) {
1758                         ex part_ca, part_cb;
1759                         g.push_back(gcd(a.op(i), part_b, &part_ca, &part_cb, check_args));
1760                         acc_ca.push_back(part_ca);
1761                         part_b = part_cb;
1762                 }
1763                 if (ca)
1764                         *ca = (new mul(acc_ca))->setflag(status_flags::dynallocated);
1765                 if (cb)
1766                         *cb = part_b;
1767                 return (new mul(g))->setflag(status_flags::dynallocated);
1768         } else if (is_exactly_a<mul>(b)) {
1769                 if (is_exactly_a<mul>(a) && a.nops() > b.nops())
1770                         goto factored_a;
1771 factored_b:
1772                 size_t num = b.nops();
1773                 exvector g; g.reserve(num);
1774                 exvector acc_cb; acc_cb.reserve(num);
1775                 ex part_a = a;
1776                 for (size_t i=0; i<num; i++) {
1777                         ex part_ca, part_cb;
1778                         g.push_back(gcd(part_a, b.op(i), &part_ca, &part_cb, check_args));
1779                         acc_cb.push_back(part_cb);
1780                         part_a = part_ca;
1781                 }
1782                 if (ca)
1783                         *ca = part_a;
1784                 if (cb)
1785                         *cb = (new mul(acc_cb))->setflag(status_flags::dynallocated);
1786                 return (new mul(g))->setflag(status_flags::dynallocated);
1787         }
1788 }
1789
1790 /** Compute LCM (Least Common Multiple) of multivariate polynomials in Z[X].
1791  *
1792  *  @param a  first multivariate polynomial
1793  *  @param b  second multivariate polynomial
1794  *  @param check_args  check whether a and b are polynomials with rational
1795  *         coefficients (defaults to "true")
1796  *  @return the LCM as a new expression */
1797 ex lcm(const ex &a, const ex &b, bool check_args)
1798 {
1799         if (is_exactly_a<numeric>(a) && is_exactly_a<numeric>(b))
1800                 return lcm(ex_to<numeric>(a), ex_to<numeric>(b));
1801         if (check_args && (!a.info(info_flags::rational_polynomial) || !b.info(info_flags::rational_polynomial)))
1802                 throw(std::invalid_argument("lcm: arguments must be polynomials over the rationals"));
1803         
1804         ex ca, cb;
1805         ex g = gcd(a, b, &ca, &cb, false);
1806         return ca * cb * g;
1807 }
1808
1809
1810 /*
1811  *  Square-free factorization
1812  */
1813
1814 /** Compute square-free factorization of multivariate polynomial a(x) using
1815  *  Yun's algorithm.  Used internally by sqrfree().
1816  *
1817  *  @param a  multivariate polynomial over Z[X], treated here as univariate
1818  *            polynomial in x.
1819  *  @param x  variable to factor in
1820  *  @return   vector of factors sorted in ascending degree */
1821 static exvector sqrfree_yun(const ex &a, const symbol &x)
1822 {
1823         exvector res;
1824         ex w = a;
1825         ex z = w.diff(x);
1826         ex g = gcd(w, z);
1827         if (g.is_equal(_ex1)) {
1828                 res.push_back(a);
1829                 return res;
1830         }
1831         ex y;
1832         do {
1833                 w = quo(w, g, x);
1834                 y = quo(z, g, x);
1835                 z = y - w.diff(x);
1836                 g = gcd(w, z);
1837                 res.push_back(g);
1838         } while (!z.is_zero());
1839         return res;
1840 }
1841
1842
1843 /** Compute a square-free factorization of a multivariate polynomial in Q[X].
1844  *
1845  *  @param a  multivariate polynomial over Q[X]
1846  *  @param l  lst of variables to factor in, may be left empty for autodetection
1847  *  @return   a square-free factorization of \p a.
1848  *
1849  * \note
1850  * A polynomial \f$p(X) \in C[X]\f$ is said <EM>square-free</EM>
1851  * if, whenever any two polynomials \f$q(X)\f$ and \f$r(X)\f$
1852  * are such that
1853  * \f[
1854  *     p(X) = q(X)^2 r(X),
1855  * \f]
1856  * we have \f$q(X) \in C\f$.
1857  * This means that \f$p(X)\f$ has no repeated factors, apart
1858  * eventually from constants.
1859  * Given a polynomial \f$p(X) \in C[X]\f$, we say that the
1860  * decomposition
1861  * \f[
1862  *   p(X) = b \cdot p_1(X)^{a_1} \cdot p_2(X)^{a_2} \cdots p_r(X)^{a_r}
1863  * \f]
1864  * is a <EM>square-free factorization</EM> of \f$p(X)\f$ if the
1865  * following conditions hold:
1866  * -#  \f$b \in C\f$ and \f$b \neq 0\f$;
1867  * -#  \f$a_i\f$ is a positive integer for \f$i = 1, \ldots, r\f$;
1868  * -#  the degree of the polynomial \f$p_i\f$ is strictly positive
1869  *     for \f$i = 1, \ldots, r\f$;
1870  * -#  the polynomial \f$\Pi_{i=1}^r p_i(X)\f$ is square-free.
1871  *
1872  * Square-free factorizations need not be unique.  For example, if
1873  * \f$a_i\f$ is even, we could change the polynomial \f$p_i(X)\f$
1874  * into \f$-p_i(X)\f$.
1875  * Observe also that the factors \f$p_i(X)\f$ need not be irreducible
1876  * polynomials.
1877  */
1878 ex sqrfree(const ex &a, const lst &l)
1879 {
1880         if (is_exactly_a<numeric>(a) ||     // algorithm does not trap a==0
1881             is_a<symbol>(a))        // shortcut
1882                 return a;
1883
1884         // If no lst of variables to factorize in was specified we have to
1885         // invent one now.  Maybe one can optimize here by reversing the order
1886         // or so, I don't know.
1887         lst args;
1888         if (l.nops()==0) {
1889                 sym_desc_vec sdv;
1890                 get_symbol_stats(a, _ex0, sdv);
1891                 sym_desc_vec::const_iterator it = sdv.begin(), itend = sdv.end();
1892                 while (it != itend) {
1893                         args.append(it->sym);
1894                         ++it;
1895                 }
1896         } else {
1897                 args = l;
1898         }
1899
1900         // Find the symbol to factor in at this stage
1901         if (!is_a<symbol>(args.op(0)))
1902                 throw (std::runtime_error("sqrfree(): invalid factorization variable"));
1903         const symbol &x = ex_to<symbol>(args.op(0));
1904
1905         // convert the argument from something in Q[X] to something in Z[X]
1906         const numeric lcm = lcm_of_coefficients_denominators(a);
1907         const ex tmp = multiply_lcm(a,lcm);
1908
1909         // find the factors
1910         exvector factors = sqrfree_yun(tmp, x);
1911
1912         // construct the next list of symbols with the first element popped
1913         lst newargs = args;
1914         newargs.remove_first();
1915
1916         // recurse down the factors in remaining variables
1917         if (newargs.nops()>0) {
1918                 exvector::iterator i = factors.begin();
1919                 while (i != factors.end()) {
1920                         *i = sqrfree(*i, newargs);
1921                         ++i;
1922                 }
1923         }
1924
1925         // Done with recursion, now construct the final result
1926         ex result = _ex1;
1927         exvector::const_iterator it = factors.begin(), itend = factors.end();
1928         for (int p = 1; it!=itend; ++it, ++p)
1929                 result *= power(*it, p);
1930
1931         // Yun's algorithm does not account for constant factors.  (For univariate
1932         // polynomials it works only in the monic case.)  We can correct this by
1933         // inserting what has been lost back into the result.  For completeness
1934         // we'll also have to recurse down that factor in the remaining variables.
1935         if (newargs.nops()>0)
1936                 result *= sqrfree(quo(tmp, result, x), newargs);
1937         else
1938                 result *= quo(tmp, result, x);
1939
1940         // Put in the reational overall factor again and return
1941         return result * lcm.inverse();
1942 }
1943
1944
1945 /** Compute square-free partial fraction decomposition of rational function
1946  *  a(x).
1947  *
1948  *  @param a rational function over Z[x], treated as univariate polynomial
1949  *           in x
1950  *  @param x variable to factor in
1951  *  @return decomposed rational function */
1952 ex sqrfree_parfrac(const ex & a, const symbol & x)
1953 {
1954         // Find numerator and denominator
1955         ex nd = numer_denom(a);
1956         ex numer = nd.op(0), denom = nd.op(1);
1957 //clog << "numer = " << numer << ", denom = " << denom << endl;
1958
1959         // Convert N(x)/D(x) -> Q(x) + R(x)/D(x), so degree(R) < degree(D)
1960         ex red_poly = quo(numer, denom, x), red_numer = rem(numer, denom, x).expand();
1961 //clog << "red_poly = " << red_poly << ", red_numer = " << red_numer << endl;
1962
1963         // Factorize denominator and compute cofactors
1964         exvector yun = sqrfree_yun(denom, x);
1965 //clog << "yun factors: " << exprseq(yun) << endl;
1966         size_t num_yun = yun.size();
1967         exvector factor; factor.reserve(num_yun);
1968         exvector cofac; cofac.reserve(num_yun);
1969         for (size_t i=0; i<num_yun; i++) {
1970                 if (!yun[i].is_equal(_ex1)) {
1971                         for (size_t j=0; j<=i; j++) {
1972                                 factor.push_back(pow(yun[i], j+1));
1973                                 ex prod = _ex1;
1974                                 for (size_t k=0; k<num_yun; k++) {
1975                                         if (k == i)
1976                                                 prod *= pow(yun[k], i-j);
1977                                         else
1978                                                 prod *= pow(yun[k], k+1);
1979                                 }
1980                                 cofac.push_back(prod.expand());
1981                         }
1982                 }
1983         }
1984         size_t num_factors = factor.size();
1985 //clog << "factors  : " << exprseq(factor) << endl;
1986 //clog << "cofactors: " << exprseq(cofac) << endl;
1987
1988         // Construct coefficient matrix for decomposition
1989         int max_denom_deg = denom.degree(x);
1990         matrix sys(max_denom_deg + 1, num_factors);
1991         matrix rhs(max_denom_deg + 1, 1);
1992         for (int i=0; i<=max_denom_deg; i++) {
1993                 for (size_t j=0; j<num_factors; j++)
1994                         sys(i, j) = cofac[j].coeff(x, i);
1995                 rhs(i, 0) = red_numer.coeff(x, i);
1996         }
1997 //clog << "coeffs: " << sys << endl;
1998 //clog << "rhs   : " << rhs << endl;
1999
2000         // Solve resulting linear system
2001         matrix vars(num_factors, 1);
2002         for (size_t i=0; i<num_factors; i++)
2003                 vars(i, 0) = symbol();
2004         matrix sol = sys.solve(vars, rhs);
2005
2006         // Sum up decomposed fractions
2007         ex sum = 0;
2008         for (size_t i=0; i<num_factors; i++)
2009                 sum += sol(i, 0) / factor[i];
2010
2011         return red_poly + sum;
2012 }
2013
2014
2015 /*
2016  *  Normal form of rational functions
2017  */
2018
2019 /*
2020  *  Note: The internal normal() functions (= basic::normal() and overloaded
2021  *  functions) all return lists of the form {numerator, denominator}. This
2022  *  is to get around mul::eval()'s automatic expansion of numeric coefficients.
2023  *  E.g. (a+b)/3 is automatically converted to a/3+b/3 but we want to keep
2024  *  the information that (a+b) is the numerator and 3 is the denominator.
2025  */
2026
2027
2028 /** Create a symbol for replacing the expression "e" (or return a previously
2029  *  assigned symbol). The symbol and expression are appended to repl, for
2030  *  a later application of subs().
2031  *  @see ex::normal */
2032 static ex replace_with_symbol(const ex & e, exmap & repl, exmap & rev_lookup)
2033 {
2034         // Expression already replaced? Then return the assigned symbol
2035         exmap::const_iterator it = rev_lookup.find(e);
2036         if (it != rev_lookup.end())
2037                 return it->second;
2038         
2039         // Otherwise create new symbol and add to list, taking care that the
2040         // replacement expression doesn't itself contain symbols from repl,
2041         // because subs() is not recursive
2042         ex es = (new symbol)->setflag(status_flags::dynallocated);
2043         ex e_replaced = e.subs(repl, subs_options::no_pattern);
2044         repl.insert(std::make_pair(es, e_replaced));
2045         rev_lookup.insert(std::make_pair(e_replaced, es));
2046         return es;
2047 }
2048
2049 /** Create a symbol for replacing the expression "e" (or return a previously
2050  *  assigned symbol). The symbol and expression are appended to repl, and the
2051  *  symbol is returned.
2052  *  @see basic::to_rational
2053  *  @see basic::to_polynomial */
2054 static ex replace_with_symbol(const ex & e, exmap & repl)
2055 {
2056         // Expression already replaced? Then return the assigned symbol
2057         for (exmap::const_iterator it = repl.begin(); it != repl.end(); ++it)
2058                 if (it->second.is_equal(e))
2059                         return it->first;
2060         
2061         // Otherwise create new symbol and add to list, taking care that the
2062         // replacement expression doesn't itself contain symbols from repl,
2063         // because subs() is not recursive
2064         ex es = (new symbol)->setflag(status_flags::dynallocated);
2065         ex e_replaced = e.subs(repl, subs_options::no_pattern);
2066         repl.insert(std::make_pair(es, e_replaced));
2067         return es;
2068 }
2069
2070
2071 /** Function object to be applied by basic::normal(). */
2072 struct normal_map_function : public map_function {
2073         int level;
2074         normal_map_function(int l) : level(l) {}
2075         ex operator()(const ex & e) { return normal(e, level); }
2076 };
2077
2078 /** Default implementation of ex::normal(). It normalizes the children and
2079  *  replaces the object with a temporary symbol.
2080  *  @see ex::normal */
2081 ex basic::normal(exmap & repl, exmap & rev_lookup, int level) const
2082 {
2083         if (nops() == 0)
2084                 return (new lst(replace_with_symbol(*this, repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2085         else {
2086                 if (level == 1)
2087                         return (new lst(replace_with_symbol(*this, repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2088                 else if (level == -max_recursion_level)
2089                         throw(std::runtime_error("max recursion level reached"));
2090                 else {
2091                         normal_map_function map_normal(level - 1);
2092                         return (new lst(replace_with_symbol(map(map_normal), repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2093                 }
2094         }
2095 }
2096
2097
2098 /** Implementation of ex::normal() for symbols. This returns the unmodified symbol.
2099  *  @see ex::normal */
2100 ex symbol::normal(exmap & repl, exmap & rev_lookup, int level) const
2101 {
2102         return (new lst(*this, _ex1))->setflag(status_flags::dynallocated);
2103 }
2104
2105
2106 /** Implementation of ex::normal() for a numeric. It splits complex numbers
2107  *  into re+I*im and replaces I and non-rational real numbers with a temporary
2108  *  symbol.
2109  *  @see ex::normal */
2110 ex numeric::normal(exmap & repl, exmap & rev_lookup, int level) const
2111 {
2112         numeric num = numer();
2113         ex numex = num;
2114
2115         if (num.is_real()) {
2116                 if (!num.is_integer())
2117                         numex = replace_with_symbol(numex, repl, rev_lookup);
2118         } else { // complex
2119                 numeric re = num.real(), im = num.imag();
2120                 ex re_ex = re.is_rational() ? re : replace_with_symbol(re, repl, rev_lookup);
2121                 ex im_ex = im.is_rational() ? im : replace_with_symbol(im, repl, rev_lookup);
2122                 numex = re_ex + im_ex * replace_with_symbol(I, repl, rev_lookup);
2123         }
2124
2125         // Denominator is always a real integer (see numeric::denom())
2126         return (new lst(numex, denom()))->setflag(status_flags::dynallocated);
2127 }
2128
2129
2130 /** Fraction cancellation.
2131  *  @param n  numerator
2132  *  @param d  denominator
2133  *  @return cancelled fraction {n, d} as a list */
2134 static ex frac_cancel(const ex &n, const ex &d)
2135 {
2136         ex num = n;
2137         ex den = d;
2138         numeric pre_factor = *_num1_p;
2139
2140 //std::clog << "frac_cancel num = " << num << ", den = " << den << std::endl;
2141
2142         // Handle trivial case where denominator is 1
2143         if (den.is_equal(_ex1))
2144                 return (new lst(num, den))->setflag(status_flags::dynallocated);
2145
2146         // Handle special cases where numerator or denominator is 0
2147         if (num.is_zero())
2148                 return (new lst(num, _ex1))->setflag(status_flags::dynallocated);
2149         if (den.expand().is_zero())
2150                 throw(std::overflow_error("frac_cancel: division by zero in frac_cancel"));
2151
2152         // Bring numerator and denominator to Z[X] by multiplying with
2153         // LCM of all coefficients' denominators
2154         numeric num_lcm = lcm_of_coefficients_denominators(num);
2155         numeric den_lcm = lcm_of_coefficients_denominators(den);
2156         num = multiply_lcm(num, num_lcm);
2157         den = multiply_lcm(den, den_lcm);
2158         pre_factor = den_lcm / num_lcm;
2159
2160         // Cancel GCD from numerator and denominator
2161         ex cnum, cden;
2162         if (gcd(num, den, &cnum, &cden, false) != _ex1) {
2163                 num = cnum;
2164                 den = cden;
2165         }
2166
2167         // Make denominator unit normal (i.e. coefficient of first symbol
2168         // as defined by get_first_symbol() is made positive)
2169         if (is_exactly_a<numeric>(den)) {
2170                 if (ex_to<numeric>(den).is_negative()) {
2171                         num *= _ex_1;
2172                         den *= _ex_1;
2173                 }
2174         } else {
2175                 ex x;
2176                 if (get_first_symbol(den, x)) {
2177                         GINAC_ASSERT(is_exactly_a<numeric>(den.unit(x)));
2178                         if (ex_to<numeric>(den.unit(x)).is_negative()) {
2179                                 num *= _ex_1;
2180                                 den *= _ex_1;
2181                         }
2182                 }
2183         }
2184
2185         // Return result as list
2186 //std::clog << " returns num = " << num << ", den = " << den << ", pre_factor = " << pre_factor << std::endl;
2187         return (new lst(num * pre_factor.numer(), den * pre_factor.denom()))->setflag(status_flags::dynallocated);
2188 }
2189
2190
2191 /** Implementation of ex::normal() for a sum. It expands terms and performs
2192  *  fractional addition.
2193  *  @see ex::normal */
2194 ex add::normal(exmap & repl, exmap & rev_lookup, int level) const
2195 {
2196         if (level == 1)
2197                 return (new lst(replace_with_symbol(*this, repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2198         else if (level == -max_recursion_level)
2199                 throw(std::runtime_error("max recursion level reached"));
2200
2201         // Normalize children and split each one into numerator and denominator
2202         exvector nums, dens;
2203         nums.reserve(seq.size()+1);
2204         dens.reserve(seq.size()+1);
2205         epvector::const_iterator it = seq.begin(), itend = seq.end();
2206         while (it != itend) {
2207                 ex n = ex_to<basic>(recombine_pair_to_ex(*it)).normal(repl, rev_lookup, level-1);
2208                 nums.push_back(n.op(0));
2209                 dens.push_back(n.op(1));
2210                 it++;
2211         }
2212         ex n = ex_to<numeric>(overall_coeff).normal(repl, rev_lookup, level-1);
2213         nums.push_back(n.op(0));
2214         dens.push_back(n.op(1));
2215         GINAC_ASSERT(nums.size() == dens.size());
2216
2217         // Now, nums is a vector of all numerators and dens is a vector of
2218         // all denominators
2219 //std::clog << "add::normal uses " << nums.size() << " summands:\n";
2220
2221         // Add fractions sequentially
2222         exvector::const_iterator num_it = nums.begin(), num_itend = nums.end();
2223         exvector::const_iterator den_it = dens.begin(), den_itend = dens.end();
2224 //std::clog << " num = " << *num_it << ", den = " << *den_it << std::endl;
2225         ex num = *num_it++, den = *den_it++;
2226         while (num_it != num_itend) {
2227 //std::clog << " num = " << *num_it << ", den = " << *den_it << std::endl;
2228                 ex next_num = *num_it++, next_den = *den_it++;
2229
2230                 // Trivially add sequences of fractions with identical denominators
2231                 while ((den_it != den_itend) && next_den.is_equal(*den_it)) {
2232                         next_num += *num_it;
2233                         num_it++; den_it++;
2234                 }
2235
2236                 // Additiion of two fractions, taking advantage of the fact that
2237                 // the heuristic GCD algorithm computes the cofactors at no extra cost
2238                 ex co_den1, co_den2;
2239                 ex g = gcd(den, next_den, &co_den1, &co_den2, false);
2240                 num = ((num * co_den2) + (next_num * co_den1)).expand();
2241                 den *= co_den2;         // this is the lcm(den, next_den)
2242         }
2243 //std::clog << " common denominator = " << den << std::endl;
2244
2245         // Cancel common factors from num/den
2246         return frac_cancel(num, den);
2247 }
2248
2249
2250 /** Implementation of ex::normal() for a product. It cancels common factors
2251  *  from fractions.
2252  *  @see ex::normal() */
2253 ex mul::normal(exmap & repl, exmap & rev_lookup, int level) const
2254 {
2255         if (level == 1)
2256                 return (new lst(replace_with_symbol(*this, repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2257         else if (level == -max_recursion_level)
2258                 throw(std::runtime_error("max recursion level reached"));
2259
2260         // Normalize children, separate into numerator and denominator
2261         exvector num; num.reserve(seq.size());
2262         exvector den; den.reserve(seq.size());
2263         ex n;
2264         epvector::const_iterator it = seq.begin(), itend = seq.end();
2265         while (it != itend) {
2266                 n = ex_to<basic>(recombine_pair_to_ex(*it)).normal(repl, rev_lookup, level-1);
2267                 num.push_back(n.op(0));
2268                 den.push_back(n.op(1));
2269                 it++;
2270         }
2271         n = ex_to<numeric>(overall_coeff).normal(repl, rev_lookup, level-1);
2272         num.push_back(n.op(0));
2273         den.push_back(n.op(1));
2274
2275         // Perform fraction cancellation
2276         return frac_cancel((new mul(num))->setflag(status_flags::dynallocated),
2277                            (new mul(den))->setflag(status_flags::dynallocated));
2278 }
2279
2280
2281 /** Implementation of ex::normal([B) for powers. It normalizes the basis,
2282  *  distributes integer exponents to numerator and denominator, and replaces
2283  *  non-integer powers by temporary symbols.
2284  *  @see ex::normal */
2285 ex power::normal(exmap & repl, exmap & rev_lookup, int level) const
2286 {
2287         if (level == 1)
2288                 return (new lst(replace_with_symbol(*this, repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2289         else if (level == -max_recursion_level)
2290                 throw(std::runtime_error("max recursion level reached"));
2291
2292         // Normalize basis and exponent (exponent gets reassembled)
2293         ex n_basis = ex_to<basic>(basis).normal(repl, rev_lookup, level-1);
2294         ex n_exponent = ex_to<basic>(exponent).normal(repl, rev_lookup, level-1);
2295         n_exponent = n_exponent.op(0) / n_exponent.op(1);
2296
2297         if (n_exponent.info(info_flags::integer)) {
2298
2299                 if (n_exponent.info(info_flags::positive)) {
2300
2301                         // (a/b)^n -> {a^n, b^n}
2302                         return (new lst(power(n_basis.op(0), n_exponent), power(n_basis.op(1), n_exponent)))->setflag(status_flags::dynallocated);
2303
2304                 } else if (n_exponent.info(info_flags::negative)) {
2305
2306                         // (a/b)^-n -> {b^n, a^n}
2307                         return (new lst(power(n_basis.op(1), -n_exponent), power(n_basis.op(0), -n_exponent)))->setflag(status_flags::dynallocated);
2308                 }
2309
2310         } else {
2311
2312                 if (n_exponent.info(info_flags::positive)) {
2313
2314                         // (a/b)^x -> {sym((a/b)^x), 1}
2315                         return (new lst(replace_with_symbol(power(n_basis.op(0) / n_basis.op(1), n_exponent), repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2316
2317                 } else if (n_exponent.info(info_flags::negative)) {
2318
2319                         if (n_basis.op(1).is_equal(_ex1)) {
2320
2321                                 // a^-x -> {1, sym(a^x)}
2322                                 return (new lst(_ex1, replace_with_symbol(power(n_basis.op(0), -n_exponent), repl, rev_lookup)))->setflag(status_flags::dynallocated);
2323
2324                         } else {
2325
2326                                 // (a/b)^-x -> {sym((b/a)^x), 1}
2327                                 return (new lst(replace_with_symbol(power(n_basis.op(1) / n_basis.op(0), -n_exponent), repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2328                         }
2329                 }
2330         }
2331
2332         // (a/b)^x -> {sym((a/b)^x, 1}
2333         return (new lst(replace_with_symbol(power(n_basis.op(0) / n_basis.op(1), n_exponent), repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2334 }
2335
2336
2337 /** Implementation of ex::normal() for pseries. It normalizes each coefficient
2338  *  and replaces the series by a temporary symbol.
2339  *  @see ex::normal */
2340 ex pseries::normal(exmap & repl, exmap & rev_lookup, int level) const
2341 {
2342         epvector newseq;
2343         epvector::const_iterator i = seq.begin(), end = seq.end();
2344         while (i != end) {
2345                 ex restexp = i->rest.normal();
2346                 if (!restexp.is_zero())
2347                         newseq.push_back(expair(restexp, i->coeff));
2348                 ++i;
2349         }
2350         ex n = pseries(relational(var,point), newseq);
2351         return (new lst(replace_with_symbol(n, repl, rev_lookup), _ex1))->setflag(status_flags::dynallocated);
2352 }
2353
2354
2355 /** Normalization of rational functions.
2356  *  This function converts an expression to its normal form
2357  *  "numerator/denominator", where numerator and denominator are (relatively
2358  *  prime) polynomials. Any subexpressions which are not rational functions
2359  *  (like non-rational numbers, non-integer powers or functions like sin(),
2360  *  cos() etc.) are replaced by temporary symbols which are re-substituted by
2361  *  the (normalized) subexpressions before normal() returns (this way, any
2362  *  expression can be treated as a rational function). normal() is applied
2363  *  recursively to arguments of functions etc.
2364  *
2365  *  @param level maximum depth of recursion
2366  *  @return normalized expression */
2367 ex ex::normal(int level) const
2368 {
2369         exmap repl, rev_lookup;
2370
2371         ex e = bp->normal(repl, rev_lookup, level);
2372         GINAC_ASSERT(is_a<lst>(e));
2373
2374         // Re-insert replaced symbols
2375         if (!repl.empty())
2376                 e = e.subs(repl, subs_options::no_pattern);
2377
2378         // Convert {numerator, denominator} form back to fraction
2379         return e.op(0) / e.op(1);
2380 }
2381
2382 /** Get numerator of an expression. If the expression is not of the normal
2383  *  form "numerator/denominator", it is first converted to this form and
2384  *  then the numerator is returned.
2385  *
2386  *  @see ex::normal
2387  *  @return numerator */
2388 ex ex::numer() const
2389 {
2390         exmap repl, rev_lookup;
2391
2392         ex e = bp->normal(repl, rev_lookup, 0);
2393         GINAC_ASSERT(is_a<lst>(e));
2394
2395         // Re-insert replaced symbols
2396         if (repl.empty())
2397                 return e.op(0);
2398         else
2399                 return e.op(0).subs(repl, subs_options::no_pattern);
2400 }
2401
2402 /** Get denominator of an expression. If the expression is not of the normal
2403  *  form "numerator/denominator", it is first converted to this form and
2404  *  then the denominator is returned.
2405  *
2406  *  @see ex::normal
2407  *  @return denominator */
2408 ex ex::denom() const
2409 {
2410         exmap repl, rev_lookup;
2411
2412         ex e = bp->normal(repl, rev_lookup, 0);
2413         GINAC_ASSERT(is_a<lst>(e));
2414
2415         // Re-insert replaced symbols
2416         if (repl.empty())
2417                 return e.op(1);
2418         else
2419                 return e.op(1).subs(repl, subs_options::no_pattern);
2420 }
2421
2422 /** Get numerator and denominator of an expression. If the expresison is not
2423  *  of the normal form "numerator/denominator", it is first converted to this
2424  *  form and then a list [numerator, denominator] is returned.
2425  *
2426  *  @see ex::normal
2427  *  @return a list [numerator, denominator] */
2428 ex ex::numer_denom() const
2429 {
2430         exmap repl, rev_lookup;
2431
2432         ex e = bp->normal(repl, rev_lookup, 0);
2433         GINAC_ASSERT(is_a<lst>(e));
2434
2435         // Re-insert replaced symbols
2436         if (repl.empty())
2437                 return e;
2438         else
2439                 return e.subs(repl, subs_options::no_pattern);
2440 }
2441
2442
2443 /** Rationalization of non-rational functions.
2444  *  This function converts a general expression to a rational function
2445  *  by replacing all non-rational subexpressions (like non-rational numbers,
2446  *  non-integer powers or functions like sin(), cos() etc.) to temporary
2447  *  symbols. This makes it possible to use functions like gcd() and divide()
2448  *  on non-rational functions by applying to_rational() on the arguments,
2449  *  calling the desired function and re-substituting the temporary symbols
2450  *  in the result. To make the last step possible, all temporary symbols and
2451  *  their associated expressions are collected in the map specified by the
2452  *  repl parameter, ready to be passed as an argument to ex::subs().
2453  *
2454  *  @param repl collects all temporary symbols and their replacements
2455  *  @return rationalized expression */
2456 ex ex::to_rational(exmap & repl) const
2457 {
2458         return bp->to_rational(repl);
2459 }
2460
2461 // GiNaC 1.1 compatibility function
2462 ex ex::to_rational(lst & repl_lst) const
2463 {
2464         // Convert lst to exmap
2465         exmap m;
2466         for (lst::const_iterator it = repl_lst.begin(); it != repl_lst.end(); ++it)
2467                 m.insert(std::make_pair(it->op(0), it->op(1)));
2468
2469         ex ret = bp->to_rational(m);
2470
2471         // Convert exmap back to lst
2472         repl_lst.remove_all();
2473         for (exmap::const_iterator it = m.begin(); it != m.end(); ++it)
2474                 repl_lst.append(it->first == it->second);
2475
2476         return ret;
2477 }
2478
2479 ex ex::to_polynomial(exmap & repl) const
2480 {
2481         return bp->to_polynomial(repl);
2482 }
2483
2484 // GiNaC 1.1 compatibility function
2485 ex ex::to_polynomial(lst & repl_lst) const
2486 {
2487         // Convert lst to exmap
2488         exmap m;
2489         for (lst::const_iterator it = repl_lst.begin(); it != repl_lst.end(); ++it)
2490                 m.insert(std::make_pair(it->op(0), it->op(1)));
2491
2492         ex ret = bp->to_polynomial(m);
2493
2494         // Convert exmap back to lst
2495         repl_lst.remove_all();
2496         for (exmap::const_iterator it = m.begin(); it != m.end(); ++it)
2497                 repl_lst.append(it->first == it->second);
2498
2499         return ret;
2500 }
2501
2502 /** Default implementation of ex::to_rational(). This replaces the object with
2503  *  a temporary symbol. */
2504 ex basic::to_rational(exmap & repl) const
2505 {
2506         return replace_with_symbol(*this, repl);
2507 }
2508
2509 ex basic::to_polynomial(exmap & repl) const
2510 {
2511         return replace_with_symbol(*this, repl);
2512 }
2513
2514
2515 /** Implementation of ex::to_rational() for symbols. This returns the
2516  *  unmodified symbol. */
2517 ex symbol::to_rational(exmap & repl) const
2518 {
2519         return *this;
2520 }
2521
2522 /** Implementation of ex::to_polynomial() for symbols. This returns the
2523  *  unmodified symbol. */
2524 ex symbol::to_polynomial(exmap & repl) const
2525 {
2526         return *this;
2527 }
2528
2529
2530 /** Implementation of ex::to_rational() for a numeric. It splits complex
2531  *  numbers into re+I*im and replaces I and non-rational real numbers with a
2532  *  temporary symbol. */
2533 ex numeric::to_rational(exmap & repl) const
2534 {
2535         if (is_real()) {
2536                 if (!is_rational())
2537                         return replace_with_symbol(*this, repl);
2538         } else { // complex
2539                 numeric re = real();
2540                 numeric im = imag();
2541                 ex re_ex = re.is_rational() ? re : replace_with_symbol(re, repl);
2542                 ex im_ex = im.is_rational() ? im : replace_with_symbol(im, repl);
2543                 return re_ex + im_ex * replace_with_symbol(I, repl);
2544         }
2545         return *this;
2546 }
2547
2548 /** Implementation of ex::to_polynomial() for a numeric. It splits complex
2549  *  numbers into re+I*im and replaces I and non-integer real numbers with a
2550  *  temporary symbol. */
2551 ex numeric::to_polynomial(exmap & repl) const
2552 {
2553         if (is_real()) {
2554                 if (!is_integer())
2555                         return replace_with_symbol(*this, repl);
2556         } else { // complex
2557                 numeric re = real();
2558                 numeric im = imag();
2559                 ex re_ex = re.is_integer() ? re : replace_with_symbol(re, repl);
2560                 ex im_ex = im.is_integer() ? im : replace_with_symbol(im, repl);
2561                 return re_ex + im_ex * replace_with_symbol(I, repl);
2562         }
2563         return *this;
2564 }
2565
2566
2567 /** Implementation of ex::to_rational() for powers. It replaces non-integer
2568  *  powers by temporary symbols. */
2569 ex power::to_rational(exmap & repl) const
2570 {
2571         if (exponent.info(info_flags::integer))
2572                 return power(basis.to_rational(repl), exponent);
2573         else
2574                 return replace_with_symbol(*this, repl);
2575 }
2576
2577 /** Implementation of ex::to_polynomial() for powers. It replaces non-posint
2578  *  powers by temporary symbols. */
2579 ex power::to_polynomial(exmap & repl) const
2580 {
2581         if (exponent.info(info_flags::posint))
2582                 return power(basis.to_rational(repl), exponent);
2583         else if (exponent.info(info_flags::negint))
2584         {
2585                 ex basis_pref = collect_common_factors(basis);
2586                 if (is_exactly_a<mul>(basis_pref) || is_exactly_a<power>(basis_pref)) {
2587                         // (A*B)^n will be automagically transformed to A^n*B^n
2588                         ex t = power(basis_pref, exponent);
2589                         return t.to_polynomial(repl);
2590                 }
2591                 else
2592                         return power(replace_with_symbol(power(basis, _ex_1), repl), -exponent);
2593         } 
2594         else
2595                 return replace_with_symbol(*this, repl);
2596 }
2597
2598
2599 /** Implementation of ex::to_rational() for expairseqs. */
2600 ex expairseq::to_rational(exmap & repl) const
2601 {
2602         epvector s;
2603         s.reserve(seq.size());
2604         epvector::const_iterator i = seq.begin(), end = seq.end();
2605         while (i != end) {
2606                 s.push_back(split_ex_to_pair(recombine_pair_to_ex(*i).to_rational(repl)));
2607                 ++i;
2608         }
2609         ex oc = overall_coeff.to_rational(repl);
2610         if (oc.info(info_flags::numeric))
2611                 return thisexpairseq(s, overall_coeff);
2612         else
2613                 s.push_back(combine_ex_with_coeff_to_pair(oc, _ex1));
2614         return thisexpairseq(s, default_overall_coeff());
2615 }
2616
2617 /** Implementation of ex::to_polynomial() for expairseqs. */
2618 ex expairseq::to_polynomial(exmap & repl) const
2619 {
2620         epvector s;
2621         s.reserve(seq.size());
2622         epvector::const_iterator i = seq.begin(), end = seq.end();
2623         while (i != end) {
2624                 s.push_back(split_ex_to_pair(recombine_pair_to_ex(*i).to_polynomial(repl)));
2625                 ++i;
2626         }
2627         ex oc = overall_coeff.to_polynomial(repl);
2628         if (oc.info(info_flags::numeric))
2629                 return thisexpairseq(s, overall_coeff);
2630         else
2631                 s.push_back(combine_ex_with_coeff_to_pair(oc, _ex1));
2632         return thisexpairseq(s, default_overall_coeff());
2633 }
2634
2635
2636 /** Remove the common factor in the terms of a sum 'e' by calculating the GCD,
2637  *  and multiply it into the expression 'factor' (which needs to be initialized
2638  *  to 1, unless you're accumulating factors). */
2639 static ex find_common_factor(const ex & e, ex & factor, exmap & repl)
2640 {
2641         if (is_exactly_a<add>(e)) {
2642
2643                 size_t num = e.nops();
2644                 exvector terms; terms.reserve(num);
2645                 ex gc;
2646
2647                 // Find the common GCD
2648                 for (size_t i=0; i<num; i++) {
2649                         ex x = e.op(i).to_polynomial(repl);
2650
2651                         if (is_exactly_a<add>(x) || is_exactly_a<mul>(x) || is_a<power>(x)) {
2652                                 ex f = 1;
2653                                 x = find_common_factor(x, f, repl);
2654                                 x *= f;
2655                         }
2656
2657                         if (i == 0)
2658                                 gc = x;
2659                         else
2660                                 gc = gcd(gc, x);
2661
2662                         terms.push_back(x);
2663                 }
2664
2665                 if (gc.is_equal(_ex1))
2666                         return e;
2667
2668                 // The GCD is the factor we pull out
2669                 factor *= gc;
2670
2671                 // Now divide all terms by the GCD
2672                 for (size_t i=0; i<num; i++) {
2673                         ex x;
2674
2675                         // Try to avoid divide() because it expands the polynomial
2676                         ex &t = terms[i];
2677                         if (is_exactly_a<mul>(t)) {
2678                                 for (size_t j=0; j<t.nops(); j++) {
2679                                         if (t.op(j).is_equal(gc)) {
2680                                                 exvector v; v.reserve(t.nops());
2681                                                 for (size_t k=0; k<t.nops(); k++) {
2682                                                         if (k == j)
2683                                                                 v.push_back(_ex1);
2684                                                         else
2685                                                                 v.push_back(t.op(k));
2686                                                 }
2687                                                 t = (new mul(v))->setflag(status_flags::dynallocated);
2688                                                 goto term_done;
2689                                         }
2690                                 }
2691                         }
2692
2693                         divide(t, gc, x);
2694                         t = x;
2695 term_done:      ;
2696                 }
2697                 return (new add(terms))->setflag(status_flags::dynallocated);
2698
2699         } else if (is_exactly_a<mul>(e)) {
2700
2701                 size_t num = e.nops();
2702                 exvector v; v.reserve(num);
2703
2704                 for (size_t i=0; i<num; i++)
2705                         v.push_back(find_common_factor(e.op(i), factor, repl));
2706
2707                 return (new mul(v))->setflag(status_flags::dynallocated);
2708
2709         } else if (is_exactly_a<power>(e)) {
2710                 const ex e_exp(e.op(1));
2711                 if (e_exp.info(info_flags::integer)) {
2712                         ex eb = e.op(0).to_polynomial(repl);
2713                         ex factor_local(_ex1);
2714                         ex pre_res = find_common_factor(eb, factor_local, repl);
2715                         factor *= power(factor_local, e_exp);
2716                         return power(pre_res, e_exp);
2717                         
2718                 } else
2719                         return e.to_polynomial(repl);
2720
2721         } else
2722                 return e;
2723 }
2724
2725
2726 /** Collect common factors in sums. This converts expressions like
2727  *  'a*(b*x+b*y)' to 'a*b*(x+y)'. */
2728 ex collect_common_factors(const ex & e)
2729 {
2730         if (is_exactly_a<add>(e) || is_exactly_a<mul>(e) || is_exactly_a<power>(e)) {
2731
2732                 exmap repl;
2733                 ex factor = 1;
2734                 ex r = find_common_factor(e, factor, repl);
2735                 return factor.subs(repl, subs_options::no_pattern) * r.subs(repl, subs_options::no_pattern);
2736
2737         } else
2738                 return e;
2739 }
2740
2741
2742 /** Resultant of two expressions e1,e2 with respect to symbol s.
2743  *  Method: Compute determinant of Sylvester matrix of e1,e2,s.  */
2744 ex resultant(const ex & e1, const ex & e2, const ex & s)
2745 {
2746         const ex ee1 = e1.expand();
2747         const ex ee2 = e2.expand();
2748         if (!ee1.info(info_flags::polynomial) ||
2749             !ee2.info(info_flags::polynomial))
2750                 throw(std::runtime_error("resultant(): arguments must be polynomials"));
2751
2752         const int h1 = ee1.degree(s);
2753         const int l1 = ee1.ldegree(s);
2754         const int h2 = ee2.degree(s);
2755         const int l2 = ee2.ldegree(s);
2756
2757         const int msize = h1 + h2;
2758         matrix m(msize, msize);
2759
2760         for (int l = h1; l >= l1; --l) {
2761                 const ex e = ee1.coeff(s, l);
2762                 for (int k = 0; k < h2; ++k)
2763                         m(k, k+h1-l) = e;
2764         }
2765         for (int l = h2; l >= l2; --l) {
2766                 const ex e = ee2.coeff(s, l);
2767                 for (int k = 0; k < h1; ++k)
2768                         m(k+h2, k+h2-l) = e;
2769         }
2770
2771         return m.determinant();
2772 }
2773
2774
2775 } // namespace GiNaC