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