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