polyroots.py 36 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226
  1. """Algorithms for computing symbolic roots of polynomials. """
  2. import math
  3. from functools import reduce
  4. from sympy.core import S, I, pi
  5. from sympy.core.exprtools import factor_terms
  6. from sympy.core.function import _mexpand
  7. from sympy.core.logic import fuzzy_not
  8. from sympy.core.mul import expand_2arg, Mul
  9. from sympy.core.numbers import Rational, igcd, comp
  10. from sympy.core.power import Pow
  11. from sympy.core.relational import Eq
  12. from sympy.core.sorting import ordered
  13. from sympy.core.symbol import Dummy, Symbol, symbols
  14. from sympy.core.sympify import sympify
  15. from sympy.functions import exp, im, cos, acos, Piecewise
  16. from sympy.functions.elementary.miscellaneous import root, sqrt
  17. from sympy.ntheory import divisors, isprime, nextprime
  18. from sympy.polys.domains import EX
  19. from sympy.polys.polyerrors import (PolynomialError, GeneratorsNeeded,
  20. DomainError, UnsolvableFactorError)
  21. from sympy.polys.polyquinticconst import PolyQuintic
  22. from sympy.polys.polytools import Poly, cancel, factor, gcd_list, discriminant
  23. from sympy.polys.rationaltools import together
  24. from sympy.polys.specialpolys import cyclotomic_poly
  25. from sympy.utilities import public
  26. from sympy.utilities.misc import filldedent
  27. z = Symbol('z') # importing from abc cause O to be lost as clashing symbol
  28. def roots_linear(f):
  29. """Returns a list of roots of a linear polynomial."""
  30. r = -f.nth(0)/f.nth(1)
  31. dom = f.get_domain()
  32. if not dom.is_Numerical:
  33. if dom.is_Composite:
  34. r = factor(r)
  35. else:
  36. from sympy.simplify.simplify import simplify
  37. r = simplify(r)
  38. return [r]
  39. def roots_quadratic(f):
  40. """Returns a list of roots of a quadratic polynomial. If the domain is ZZ
  41. then the roots will be sorted with negatives coming before positives.
  42. The ordering will be the same for any numerical coefficients as long as
  43. the assumptions tested are correct, otherwise the ordering will not be
  44. sorted (but will be canonical).
  45. """
  46. a, b, c = f.all_coeffs()
  47. dom = f.get_domain()
  48. def _sqrt(d):
  49. # remove squares from square root since both will be represented
  50. # in the results; a similar thing is happening in roots() but
  51. # must be duplicated here because not all quadratics are binomials
  52. co = []
  53. other = []
  54. for di in Mul.make_args(d):
  55. if di.is_Pow and di.exp.is_Integer and di.exp % 2 == 0:
  56. co.append(Pow(di.base, di.exp//2))
  57. else:
  58. other.append(di)
  59. if co:
  60. d = Mul(*other)
  61. co = Mul(*co)
  62. return co*sqrt(d)
  63. return sqrt(d)
  64. def _simplify(expr):
  65. if dom.is_Composite:
  66. return factor(expr)
  67. else:
  68. from sympy.simplify.simplify import simplify
  69. return simplify(expr)
  70. if c is S.Zero:
  71. r0, r1 = S.Zero, -b/a
  72. if not dom.is_Numerical:
  73. r1 = _simplify(r1)
  74. elif r1.is_negative:
  75. r0, r1 = r1, r0
  76. elif b is S.Zero:
  77. r = -c/a
  78. if not dom.is_Numerical:
  79. r = _simplify(r)
  80. R = _sqrt(r)
  81. r0 = -R
  82. r1 = R
  83. else:
  84. d = b**2 - 4*a*c
  85. A = 2*a
  86. B = -b/A
  87. if not dom.is_Numerical:
  88. d = _simplify(d)
  89. B = _simplify(B)
  90. D = factor_terms(_sqrt(d)/A)
  91. r0 = B - D
  92. r1 = B + D
  93. if a.is_negative:
  94. r0, r1 = r1, r0
  95. elif not dom.is_Numerical:
  96. r0, r1 = [expand_2arg(i) for i in (r0, r1)]
  97. return [r0, r1]
  98. def roots_cubic(f, trig=False):
  99. """Returns a list of roots of a cubic polynomial.
  100. References
  101. ==========
  102. [1] https://en.wikipedia.org/wiki/Cubic_function, General formula for roots,
  103. (accessed November 17, 2014).
  104. """
  105. if trig:
  106. a, b, c, d = f.all_coeffs()
  107. p = (3*a*c - b**2)/(3*a**2)
  108. q = (2*b**3 - 9*a*b*c + 27*a**2*d)/(27*a**3)
  109. D = 18*a*b*c*d - 4*b**3*d + b**2*c**2 - 4*a*c**3 - 27*a**2*d**2
  110. if (D > 0) == True:
  111. rv = []
  112. for k in range(3):
  113. rv.append(2*sqrt(-p/3)*cos(acos(q/p*sqrt(-3/p)*Rational(3, 2))/3 - k*pi*Rational(2, 3)))
  114. return [i - b/3/a for i in rv]
  115. # a*x**3 + b*x**2 + c*x + d -> x**3 + a*x**2 + b*x + c
  116. _, a, b, c = f.monic().all_coeffs()
  117. if c is S.Zero:
  118. x1, x2 = roots([1, a, b], multiple=True)
  119. return [x1, S.Zero, x2]
  120. # x**3 + a*x**2 + b*x + c -> u**3 + p*u + q
  121. p = b - a**2/3
  122. q = c - a*b/3 + 2*a**3/27
  123. pon3 = p/3
  124. aon3 = a/3
  125. u1 = None
  126. if p is S.Zero:
  127. if q is S.Zero:
  128. return [-aon3]*3
  129. u1 = -root(q, 3) if q.is_positive else root(-q, 3)
  130. elif q is S.Zero:
  131. y1, y2 = roots([1, 0, p], multiple=True)
  132. return [tmp - aon3 for tmp in [y1, S.Zero, y2]]
  133. elif q.is_real and q.is_negative:
  134. u1 = -root(-q/2 + sqrt(q**2/4 + pon3**3), 3)
  135. coeff = I*sqrt(3)/2
  136. if u1 is None:
  137. u1 = S.One
  138. u2 = Rational(-1, 2) + coeff
  139. u3 = Rational(-1, 2) - coeff
  140. b, c, d = a, b, c # a, b, c, d = S.One, a, b, c
  141. D0 = b**2 - 3*c # b**2 - 3*a*c
  142. D1 = 2*b**3 - 9*b*c + 27*d # 2*b**3 - 9*a*b*c + 27*a**2*d
  143. C = root((D1 + sqrt(D1**2 - 4*D0**3))/2, 3)
  144. return [-(b + uk*C + D0/C/uk)/3 for uk in [u1, u2, u3]] # -(b + uk*C + D0/C/uk)/3/a
  145. u2 = u1*(Rational(-1, 2) + coeff)
  146. u3 = u1*(Rational(-1, 2) - coeff)
  147. if p is S.Zero:
  148. return [u1 - aon3, u2 - aon3, u3 - aon3]
  149. soln = [
  150. -u1 + pon3/u1 - aon3,
  151. -u2 + pon3/u2 - aon3,
  152. -u3 + pon3/u3 - aon3
  153. ]
  154. return soln
  155. def _roots_quartic_euler(p, q, r, a):
  156. """
  157. Descartes-Euler solution of the quartic equation
  158. Parameters
  159. ==========
  160. p, q, r: coefficients of ``x**4 + p*x**2 + q*x + r``
  161. a: shift of the roots
  162. Notes
  163. =====
  164. This is a helper function for ``roots_quartic``.
  165. Look for solutions of the form ::
  166. ``x1 = sqrt(R) - sqrt(A + B*sqrt(R))``
  167. ``x2 = -sqrt(R) - sqrt(A - B*sqrt(R))``
  168. ``x3 = -sqrt(R) + sqrt(A - B*sqrt(R))``
  169. ``x4 = sqrt(R) + sqrt(A + B*sqrt(R))``
  170. To satisfy the quartic equation one must have
  171. ``p = -2*(R + A); q = -4*B*R; r = (R - A)**2 - B**2*R``
  172. so that ``R`` must satisfy the Descartes-Euler resolvent equation
  173. ``64*R**3 + 32*p*R**2 + (4*p**2 - 16*r)*R - q**2 = 0``
  174. If the resolvent does not have a rational solution, return None;
  175. in that case it is likely that the Ferrari method gives a simpler
  176. solution.
  177. Examples
  178. ========
  179. >>> from sympy import S
  180. >>> from sympy.polys.polyroots import _roots_quartic_euler
  181. >>> p, q, r = -S(64)/5, -S(512)/125, -S(1024)/3125
  182. >>> _roots_quartic_euler(p, q, r, S(0))[0]
  183. -sqrt(32*sqrt(5)/125 + 16/5) + 4*sqrt(5)/5
  184. """
  185. # solve the resolvent equation
  186. x = Dummy('x')
  187. eq = 64*x**3 + 32*p*x**2 + (4*p**2 - 16*r)*x - q**2
  188. xsols = list(roots(Poly(eq, x), cubics=False).keys())
  189. xsols = [sol for sol in xsols if sol.is_rational and sol.is_nonzero]
  190. if not xsols:
  191. return None
  192. R = max(xsols)
  193. c1 = sqrt(R)
  194. B = -q*c1/(4*R)
  195. A = -R - p/2
  196. c2 = sqrt(A + B)
  197. c3 = sqrt(A - B)
  198. return [c1 - c2 - a, -c1 - c3 - a, -c1 + c3 - a, c1 + c2 - a]
  199. def roots_quartic(f):
  200. r"""
  201. Returns a list of roots of a quartic polynomial.
  202. There are many references for solving quartic expressions available [1-5].
  203. This reviewer has found that many of them require one to select from among
  204. 2 or more possible sets of solutions and that some solutions work when one
  205. is searching for real roots but do not work when searching for complex roots
  206. (though this is not always stated clearly). The following routine has been
  207. tested and found to be correct for 0, 2 or 4 complex roots.
  208. The quasisymmetric case solution [6] looks for quartics that have the form
  209. `x**4 + A*x**3 + B*x**2 + C*x + D = 0` where `(C/A)**2 = D`.
  210. Although no general solution that is always applicable for all
  211. coefficients is known to this reviewer, certain conditions are tested
  212. to determine the simplest 4 expressions that can be returned:
  213. 1) `f = c + a*(a**2/8 - b/2) == 0`
  214. 2) `g = d - a*(a*(3*a**2/256 - b/16) + c/4) = 0`
  215. 3) if `f != 0` and `g != 0` and `p = -d + a*c/4 - b**2/12` then
  216. a) `p == 0`
  217. b) `p != 0`
  218. Examples
  219. ========
  220. >>> from sympy import Poly
  221. >>> from sympy.polys.polyroots import roots_quartic
  222. >>> r = roots_quartic(Poly('x**4-6*x**3+17*x**2-26*x+20'))
  223. >>> # 4 complex roots: 1+-I*sqrt(3), 2+-I
  224. >>> sorted(str(tmp.evalf(n=2)) for tmp in r)
  225. ['1.0 + 1.7*I', '1.0 - 1.7*I', '2.0 + 1.0*I', '2.0 - 1.0*I']
  226. References
  227. ==========
  228. 1. http://mathforum.org/dr.math/faq/faq.cubic.equations.html
  229. 2. https://en.wikipedia.org/wiki/Quartic_function#Summary_of_Ferrari.27s_method
  230. 3. https://planetmath.org/encyclopedia/GaloisTheoreticDerivationOfTheQuarticFormula.html
  231. 4. https://people.bath.ac.uk/masjhd/JHD-CA.pdf
  232. 5. http://www.albmath.org/files/Math_5713.pdf
  233. 6. https://web.archive.org/web/20171002081448/http://www.statemaster.com/encyclopedia/Quartic-equation
  234. 7. https://eqworld.ipmnet.ru/en/solutions/ae/ae0108.pdf
  235. """
  236. _, a, b, c, d = f.monic().all_coeffs()
  237. if not d:
  238. return [S.Zero] + roots([1, a, b, c], multiple=True)
  239. elif (c/a)**2 == d:
  240. x, m = f.gen, c/a
  241. g = Poly(x**2 + a*x + b - 2*m, x)
  242. z1, z2 = roots_quadratic(g)
  243. h1 = Poly(x**2 - z1*x + m, x)
  244. h2 = Poly(x**2 - z2*x + m, x)
  245. r1 = roots_quadratic(h1)
  246. r2 = roots_quadratic(h2)
  247. return r1 + r2
  248. else:
  249. a2 = a**2
  250. e = b - 3*a2/8
  251. f = _mexpand(c + a*(a2/8 - b/2))
  252. aon4 = a/4
  253. g = _mexpand(d - aon4*(a*(3*a2/64 - b/4) + c))
  254. if f.is_zero:
  255. y1, y2 = [sqrt(tmp) for tmp in
  256. roots([1, e, g], multiple=True)]
  257. return [tmp - aon4 for tmp in [-y1, -y2, y1, y2]]
  258. if g.is_zero:
  259. y = [S.Zero] + roots([1, 0, e, f], multiple=True)
  260. return [tmp - aon4 for tmp in y]
  261. else:
  262. # Descartes-Euler method, see [7]
  263. sols = _roots_quartic_euler(e, f, g, aon4)
  264. if sols:
  265. return sols
  266. # Ferrari method, see [1, 2]
  267. p = -e**2/12 - g
  268. q = -e**3/108 + e*g/3 - f**2/8
  269. TH = Rational(1, 3)
  270. def _ans(y):
  271. w = sqrt(e + 2*y)
  272. arg1 = 3*e + 2*y
  273. arg2 = 2*f/w
  274. ans = []
  275. for s in [-1, 1]:
  276. root = sqrt(-(arg1 + s*arg2))
  277. for t in [-1, 1]:
  278. ans.append((s*w - t*root)/2 - aon4)
  279. return ans
  280. # whether a Piecewise is returned or not
  281. # depends on knowing p, so try to put
  282. # in a simple form
  283. p = _mexpand(p)
  284. # p == 0 case
  285. y1 = e*Rational(-5, 6) - q**TH
  286. if p.is_zero:
  287. return _ans(y1)
  288. # if p != 0 then u below is not 0
  289. root = sqrt(q**2/4 + p**3/27)
  290. r = -q/2 + root # or -q/2 - root
  291. u = r**TH # primary root of solve(x**3 - r, x)
  292. y2 = e*Rational(-5, 6) + u - p/u/3
  293. if fuzzy_not(p.is_zero):
  294. return _ans(y2)
  295. # sort it out once they know the values of the coefficients
  296. return [Piecewise((a1, Eq(p, 0)), (a2, True))
  297. for a1, a2 in zip(_ans(y1), _ans(y2))]
  298. def roots_binomial(f):
  299. """Returns a list of roots of a binomial polynomial. If the domain is ZZ
  300. then the roots will be sorted with negatives coming before positives.
  301. The ordering will be the same for any numerical coefficients as long as
  302. the assumptions tested are correct, otherwise the ordering will not be
  303. sorted (but will be canonical).
  304. """
  305. n = f.degree()
  306. a, b = f.nth(n), f.nth(0)
  307. base = -cancel(b/a)
  308. alpha = root(base, n)
  309. if alpha.is_number:
  310. alpha = alpha.expand(complex=True)
  311. # define some parameters that will allow us to order the roots.
  312. # If the domain is ZZ this is guaranteed to return roots sorted
  313. # with reals before non-real roots and non-real sorted according
  314. # to real part and imaginary part, e.g. -1, 1, -1 + I, 2 - I
  315. neg = base.is_negative
  316. even = n % 2 == 0
  317. if neg:
  318. if even == True and (base + 1).is_positive:
  319. big = True
  320. else:
  321. big = False
  322. # get the indices in the right order so the computed
  323. # roots will be sorted when the domain is ZZ
  324. ks = []
  325. imax = n//2
  326. if even:
  327. ks.append(imax)
  328. imax -= 1
  329. if not neg:
  330. ks.append(0)
  331. for i in range(imax, 0, -1):
  332. if neg:
  333. ks.extend([i, -i])
  334. else:
  335. ks.extend([-i, i])
  336. if neg:
  337. ks.append(0)
  338. if big:
  339. for i in range(0, len(ks), 2):
  340. pair = ks[i: i + 2]
  341. pair = list(reversed(pair))
  342. # compute the roots
  343. roots, d = [], 2*I*pi/n
  344. for k in ks:
  345. zeta = exp(k*d).expand(complex=True)
  346. roots.append((alpha*zeta).expand(power_base=False))
  347. return roots
  348. def _inv_totient_estimate(m):
  349. """
  350. Find ``(L, U)`` such that ``L <= phi^-1(m) <= U``.
  351. Examples
  352. ========
  353. >>> from sympy.polys.polyroots import _inv_totient_estimate
  354. >>> _inv_totient_estimate(192)
  355. (192, 840)
  356. >>> _inv_totient_estimate(400)
  357. (400, 1750)
  358. """
  359. primes = [ d + 1 for d in divisors(m) if isprime(d + 1) ]
  360. a, b = 1, 1
  361. for p in primes:
  362. a *= p
  363. b *= p - 1
  364. L = m
  365. U = int(math.ceil(m*(float(a)/b)))
  366. P = p = 2
  367. primes = []
  368. while P <= U:
  369. p = nextprime(p)
  370. primes.append(p)
  371. P *= p
  372. P //= p
  373. b = 1
  374. for p in primes[:-1]:
  375. b *= p - 1
  376. U = int(math.ceil(m*(float(P)/b)))
  377. return L, U
  378. def roots_cyclotomic(f, factor=False):
  379. """Compute roots of cyclotomic polynomials. """
  380. L, U = _inv_totient_estimate(f.degree())
  381. for n in range(L, U + 1):
  382. g = cyclotomic_poly(n, f.gen, polys=True)
  383. if f.expr == g.expr:
  384. break
  385. else: # pragma: no cover
  386. raise RuntimeError("failed to find index of a cyclotomic polynomial")
  387. roots = []
  388. if not factor:
  389. # get the indices in the right order so the computed
  390. # roots will be sorted
  391. h = n//2
  392. ks = [i for i in range(1, n + 1) if igcd(i, n) == 1]
  393. ks.sort(key=lambda x: (x, -1) if x <= h else (abs(x - n), 1))
  394. d = 2*I*pi/n
  395. for k in reversed(ks):
  396. roots.append(exp(k*d).expand(complex=True))
  397. else:
  398. g = Poly(f, extension=root(-1, n))
  399. for h, _ in ordered(g.factor_list()[1]):
  400. roots.append(-h.TC())
  401. return roots
  402. def roots_quintic(f):
  403. """
  404. Calculate exact roots of a solvable irreducible quintic with rational coefficients.
  405. Return an empty list if the quintic is reducible or not solvable.
  406. """
  407. result = []
  408. coeff_5, coeff_4, p_, q_, r_, s_ = f.all_coeffs()
  409. if not all(coeff.is_Rational for coeff in (coeff_5, coeff_4, p_, q_, r_, s_)):
  410. return result
  411. if coeff_5 != 1:
  412. f = Poly(f / coeff_5)
  413. _, coeff_4, p_, q_, r_, s_ = f.all_coeffs()
  414. # Cancel coeff_4 to form x^5 + px^3 + qx^2 + rx + s
  415. if coeff_4:
  416. p = p_ - 2*coeff_4*coeff_4/5
  417. q = q_ - 3*coeff_4*p_/5 + 4*coeff_4**3/25
  418. r = r_ - 2*coeff_4*q_/5 + 3*coeff_4**2*p_/25 - 3*coeff_4**4/125
  419. s = s_ - coeff_4*r_/5 + coeff_4**2*q_/25 - coeff_4**3*p_/125 + 4*coeff_4**5/3125
  420. x = f.gen
  421. f = Poly(x**5 + p*x**3 + q*x**2 + r*x + s)
  422. else:
  423. p, q, r, s = p_, q_, r_, s_
  424. quintic = PolyQuintic(f)
  425. # Eqn standardized. Algo for solving starts here
  426. if not f.is_irreducible:
  427. return result
  428. f20 = quintic.f20
  429. # Check if f20 has linear factors over domain Z
  430. if f20.is_irreducible:
  431. return result
  432. # Now, we know that f is solvable
  433. for _factor in f20.factor_list()[1]:
  434. if _factor[0].is_linear:
  435. theta = _factor[0].root(0)
  436. break
  437. d = discriminant(f)
  438. delta = sqrt(d)
  439. # zeta = a fifth root of unity
  440. zeta1, zeta2, zeta3, zeta4 = quintic.zeta
  441. T = quintic.T(theta, d)
  442. tol = S(1e-10)
  443. alpha = T[1] + T[2]*delta
  444. alpha_bar = T[1] - T[2]*delta
  445. beta = T[3] + T[4]*delta
  446. beta_bar = T[3] - T[4]*delta
  447. disc = alpha**2 - 4*beta
  448. disc_bar = alpha_bar**2 - 4*beta_bar
  449. l0 = quintic.l0(theta)
  450. Stwo = S(2)
  451. l1 = _quintic_simplify((-alpha + sqrt(disc)) / Stwo)
  452. l4 = _quintic_simplify((-alpha - sqrt(disc)) / Stwo)
  453. l2 = _quintic_simplify((-alpha_bar + sqrt(disc_bar)) / Stwo)
  454. l3 = _quintic_simplify((-alpha_bar - sqrt(disc_bar)) / Stwo)
  455. order = quintic.order(theta, d)
  456. test = (order*delta.n()) - ( (l1.n() - l4.n())*(l2.n() - l3.n()) )
  457. # Comparing floats
  458. if not comp(test, 0, tol):
  459. l2, l3 = l3, l2
  460. # Now we have correct order of l's
  461. R1 = l0 + l1*zeta1 + l2*zeta2 + l3*zeta3 + l4*zeta4
  462. R2 = l0 + l3*zeta1 + l1*zeta2 + l4*zeta3 + l2*zeta4
  463. R3 = l0 + l2*zeta1 + l4*zeta2 + l1*zeta3 + l3*zeta4
  464. R4 = l0 + l4*zeta1 + l3*zeta2 + l2*zeta3 + l1*zeta4
  465. Res = [None, [None]*5, [None]*5, [None]*5, [None]*5]
  466. Res_n = [None, [None]*5, [None]*5, [None]*5, [None]*5]
  467. # Simplifying improves performance a lot for exact expressions
  468. R1 = _quintic_simplify(R1)
  469. R2 = _quintic_simplify(R2)
  470. R3 = _quintic_simplify(R3)
  471. R4 = _quintic_simplify(R4)
  472. # hard-coded results for [factor(i) for i in _vsolve(x**5 - a - I*b, x)]
  473. x0 = z**(S(1)/5)
  474. x1 = sqrt(2)
  475. x2 = sqrt(5)
  476. x3 = sqrt(5 - x2)
  477. x4 = I*x2
  478. x5 = x4 + I
  479. x6 = I*x0/4
  480. x7 = x1*sqrt(x2 + 5)
  481. sol = [x0, -x6*(x1*x3 - x5), x6*(x1*x3 + x5), -x6*(x4 + x7 - I), x6*(-x4 + x7 + I)]
  482. R1 = R1.as_real_imag()
  483. R2 = R2.as_real_imag()
  484. R3 = R3.as_real_imag()
  485. R4 = R4.as_real_imag()
  486. for i, s in enumerate(sol):
  487. Res[1][i] = _quintic_simplify(s.xreplace({z: R1[0] + I*R1[1]}))
  488. Res[2][i] = _quintic_simplify(s.xreplace({z: R2[0] + I*R2[1]}))
  489. Res[3][i] = _quintic_simplify(s.xreplace({z: R3[0] + I*R3[1]}))
  490. Res[4][i] = _quintic_simplify(s.xreplace({z: R4[0] + I*R4[1]}))
  491. for i in range(1, 5):
  492. for j in range(5):
  493. Res_n[i][j] = Res[i][j].n()
  494. Res[i][j] = _quintic_simplify(Res[i][j])
  495. r1 = Res[1][0]
  496. r1_n = Res_n[1][0]
  497. for i in range(5):
  498. if comp(im(r1_n*Res_n[4][i]), 0, tol):
  499. r4 = Res[4][i]
  500. break
  501. # Now we have various Res values. Each will be a list of five
  502. # values. We have to pick one r value from those five for each Res
  503. u, v = quintic.uv(theta, d)
  504. testplus = (u + v*delta*sqrt(5)).n()
  505. testminus = (u - v*delta*sqrt(5)).n()
  506. # Evaluated numbers suffixed with _n
  507. # We will use evaluated numbers for calculation. Much faster.
  508. r4_n = r4.n()
  509. r2 = r3 = None
  510. for i in range(5):
  511. r2temp_n = Res_n[2][i]
  512. for j in range(5):
  513. # Again storing away the exact number and using
  514. # evaluated numbers in computations
  515. r3temp_n = Res_n[3][j]
  516. if (comp((r1_n*r2temp_n**2 + r4_n*r3temp_n**2 - testplus).n(), 0, tol) and
  517. comp((r3temp_n*r1_n**2 + r2temp_n*r4_n**2 - testminus).n(), 0, tol)):
  518. r2 = Res[2][i]
  519. r3 = Res[3][j]
  520. break
  521. if r2 is not None:
  522. break
  523. else:
  524. return [] # fall back to normal solve
  525. # Now, we have r's so we can get roots
  526. x1 = (r1 + r2 + r3 + r4)/5
  527. x2 = (r1*zeta4 + r2*zeta3 + r3*zeta2 + r4*zeta1)/5
  528. x3 = (r1*zeta3 + r2*zeta1 + r3*zeta4 + r4*zeta2)/5
  529. x4 = (r1*zeta2 + r2*zeta4 + r3*zeta1 + r4*zeta3)/5
  530. x5 = (r1*zeta1 + r2*zeta2 + r3*zeta3 + r4*zeta4)/5
  531. result = [x1, x2, x3, x4, x5]
  532. # Now check if solutions are distinct
  533. saw = set()
  534. for r in result:
  535. r = r.n(2)
  536. if r in saw:
  537. # Roots were identical. Abort, return []
  538. # and fall back to usual solve
  539. return []
  540. saw.add(r)
  541. # Restore to original equation where coeff_4 is nonzero
  542. if coeff_4:
  543. result = [x - coeff_4 / 5 for x in result]
  544. return result
  545. def _quintic_simplify(expr):
  546. from sympy.simplify.simplify import powsimp
  547. expr = powsimp(expr)
  548. expr = cancel(expr)
  549. return together(expr)
  550. def _integer_basis(poly):
  551. """Compute coefficient basis for a polynomial over integers.
  552. Returns the integer ``div`` such that substituting ``x = div*y``
  553. ``p(x) = m*q(y)`` where the coefficients of ``q`` are smaller
  554. than those of ``p``.
  555. For example ``x**5 + 512*x + 1024 = 0``
  556. with ``div = 4`` becomes ``y**5 + 2*y + 1 = 0``
  557. Returns the integer ``div`` or ``None`` if there is no possible scaling.
  558. Examples
  559. ========
  560. >>> from sympy.polys import Poly
  561. >>> from sympy.abc import x
  562. >>> from sympy.polys.polyroots import _integer_basis
  563. >>> p = Poly(x**5 + 512*x + 1024, x, domain='ZZ')
  564. >>> _integer_basis(p)
  565. 4
  566. """
  567. monoms, coeffs = list(zip(*poly.terms()))
  568. monoms, = list(zip(*monoms))
  569. coeffs = list(map(abs, coeffs))
  570. if coeffs[0] < coeffs[-1]:
  571. coeffs = list(reversed(coeffs))
  572. n = monoms[0]
  573. monoms = [n - i for i in reversed(monoms)]
  574. else:
  575. return None
  576. monoms = monoms[:-1]
  577. coeffs = coeffs[:-1]
  578. # Special case for two-term polynominals
  579. if len(monoms) == 1:
  580. r = Pow(coeffs[0], S.One/monoms[0])
  581. if r.is_Integer:
  582. return int(r)
  583. else:
  584. return None
  585. divs = reversed(divisors(gcd_list(coeffs))[1:])
  586. try:
  587. div = next(divs)
  588. except StopIteration:
  589. return None
  590. while True:
  591. for monom, coeff in zip(monoms, coeffs):
  592. if coeff % div**monom != 0:
  593. try:
  594. div = next(divs)
  595. except StopIteration:
  596. return None
  597. else:
  598. break
  599. else:
  600. return div
  601. def preprocess_roots(poly):
  602. """Try to get rid of symbolic coefficients from ``poly``. """
  603. coeff = S.One
  604. poly_func = poly.func
  605. try:
  606. _, poly = poly.clear_denoms(convert=True)
  607. except DomainError:
  608. return coeff, poly
  609. poly = poly.primitive()[1]
  610. poly = poly.retract()
  611. # TODO: This is fragile. Figure out how to make this independent of construct_domain().
  612. if poly.get_domain().is_Poly and all(c.is_term for c in poly.rep.coeffs()):
  613. poly = poly.inject()
  614. strips = list(zip(*poly.monoms()))
  615. gens = list(poly.gens[1:])
  616. base, strips = strips[0], strips[1:]
  617. for gen, strip in zip(list(gens), strips):
  618. reverse = False
  619. if strip[0] < strip[-1]:
  620. strip = reversed(strip)
  621. reverse = True
  622. ratio = None
  623. for a, b in zip(base, strip):
  624. if not a and not b:
  625. continue
  626. elif not a or not b:
  627. break
  628. elif b % a != 0:
  629. break
  630. else:
  631. _ratio = b // a
  632. if ratio is None:
  633. ratio = _ratio
  634. elif ratio != _ratio:
  635. break
  636. else:
  637. if reverse:
  638. ratio = -ratio
  639. poly = poly.eval(gen, 1)
  640. coeff *= gen**(-ratio)
  641. gens.remove(gen)
  642. if gens:
  643. poly = poly.eject(*gens)
  644. if poly.is_univariate and poly.get_domain().is_ZZ:
  645. basis = _integer_basis(poly)
  646. if basis is not None:
  647. n = poly.degree()
  648. def func(k, coeff):
  649. return coeff//basis**(n - k[0])
  650. poly = poly.termwise(func)
  651. coeff *= basis
  652. if not isinstance(poly, poly_func):
  653. poly = poly_func(poly)
  654. return coeff, poly
  655. @public
  656. def roots(f, *gens,
  657. auto=True,
  658. cubics=True,
  659. trig=False,
  660. quartics=True,
  661. quintics=False,
  662. multiple=False,
  663. filter=None,
  664. predicate=None,
  665. strict=False,
  666. **flags):
  667. """
  668. Computes symbolic roots of a univariate polynomial.
  669. Given a univariate polynomial f with symbolic coefficients (or
  670. a list of the polynomial's coefficients), returns a dictionary
  671. with its roots and their multiplicities.
  672. Only roots expressible via radicals will be returned. To get
  673. a complete set of roots use RootOf class or numerical methods
  674. instead. By default cubic and quartic formulas are used in
  675. the algorithm. To disable them because of unreadable output
  676. set ``cubics=False`` or ``quartics=False`` respectively. If cubic
  677. roots are real but are expressed in terms of complex numbers
  678. (casus irreducibilis [1]) the ``trig`` flag can be set to True to
  679. have the solutions returned in terms of cosine and inverse cosine
  680. functions.
  681. To get roots from a specific domain set the ``filter`` flag with
  682. one of the following specifiers: Z, Q, R, I, C. By default all
  683. roots are returned (this is equivalent to setting ``filter='C'``).
  684. By default a dictionary is returned giving a compact result in
  685. case of multiple roots. However to get a list containing all
  686. those roots set the ``multiple`` flag to True; the list will
  687. have identical roots appearing next to each other in the result.
  688. (For a given Poly, the all_roots method will give the roots in
  689. sorted numerical order.)
  690. If the ``strict`` flag is True, ``UnsolvableFactorError`` will be
  691. raised if the roots found are known to be incomplete (because
  692. some roots are not expressible in radicals).
  693. Examples
  694. ========
  695. >>> from sympy import Poly, roots, degree
  696. >>> from sympy.abc import x, y
  697. >>> roots(x**2 - 1, x)
  698. {-1: 1, 1: 1}
  699. >>> p = Poly(x**2-1, x)
  700. >>> roots(p)
  701. {-1: 1, 1: 1}
  702. >>> p = Poly(x**2-y, x, y)
  703. >>> roots(Poly(p, x))
  704. {-sqrt(y): 1, sqrt(y): 1}
  705. >>> roots(x**2 - y, x)
  706. {-sqrt(y): 1, sqrt(y): 1}
  707. >>> roots([1, 0, -1])
  708. {-1: 1, 1: 1}
  709. ``roots`` will only return roots expressible in radicals. If
  710. the given polynomial has some or all of its roots inexpressible in
  711. radicals, the result of ``roots`` will be incomplete or empty
  712. respectively.
  713. Example where result is incomplete:
  714. >>> roots((x-1)*(x**5-x+1), x)
  715. {1: 1}
  716. In this case, the polynomial has an unsolvable quintic factor
  717. whose roots cannot be expressed by radicals. The polynomial has a
  718. rational root (due to the factor `(x-1)`), which is returned since
  719. ``roots`` always finds all rational roots.
  720. Example where result is empty:
  721. >>> roots(x**7-3*x**2+1, x)
  722. {}
  723. Here, the polynomial has no roots expressible in radicals, so
  724. ``roots`` returns an empty dictionary.
  725. The result produced by ``roots`` is complete if and only if the
  726. sum of the multiplicity of each root is equal to the degree of
  727. the polynomial. If strict=True, UnsolvableFactorError will be
  728. raised if the result is incomplete.
  729. The result can be be checked for completeness as follows:
  730. >>> f = x**3-2*x**2+1
  731. >>> sum(roots(f, x).values()) == degree(f, x)
  732. True
  733. >>> f = (x-1)*(x**5-x+1)
  734. >>> sum(roots(f, x).values()) == degree(f, x)
  735. False
  736. References
  737. ==========
  738. .. [1] https://en.wikipedia.org/wiki/Cubic_equation#Trigonometric_and_hyperbolic_solutions
  739. """
  740. from sympy.polys.polytools import to_rational_coeffs
  741. flags = dict(flags)
  742. if isinstance(f, list):
  743. if gens:
  744. raise ValueError('redundant generators given')
  745. x = Dummy('x')
  746. poly, i = {}, len(f) - 1
  747. for coeff in f:
  748. poly[i], i = sympify(coeff), i - 1
  749. f = Poly(poly, x, field=True)
  750. else:
  751. try:
  752. F = Poly(f, *gens, **flags)
  753. if not isinstance(f, Poly) and not F.gen.is_Symbol:
  754. raise PolynomialError("generator must be a Symbol")
  755. f = F
  756. except GeneratorsNeeded:
  757. if multiple:
  758. return []
  759. else:
  760. return {}
  761. else:
  762. n = f.degree()
  763. if f.length() == 2 and n > 2:
  764. # check for foo**n in constant if dep is c*gen**m
  765. con, dep = f.as_expr().as_independent(*f.gens)
  766. fcon = -(-con).factor()
  767. if fcon != con:
  768. con = fcon
  769. bases = []
  770. for i in Mul.make_args(con):
  771. if i.is_Pow:
  772. b, e = i.as_base_exp()
  773. if e.is_Integer and b.is_Add:
  774. bases.append((b, Dummy(positive=True)))
  775. if bases:
  776. rv = roots(Poly((dep + con).xreplace(dict(bases)),
  777. *f.gens), *F.gens,
  778. auto=auto,
  779. cubics=cubics,
  780. trig=trig,
  781. quartics=quartics,
  782. quintics=quintics,
  783. multiple=multiple,
  784. filter=filter,
  785. predicate=predicate,
  786. **flags)
  787. return {factor_terms(k.xreplace(
  788. {v: k for k, v in bases})
  789. ): v for k, v in rv.items()}
  790. if f.is_multivariate:
  791. raise PolynomialError('multivariate polynomials are not supported')
  792. def _update_dict(result, zeros, currentroot, k):
  793. if currentroot == S.Zero:
  794. if S.Zero in zeros:
  795. zeros[S.Zero] += k
  796. else:
  797. zeros[S.Zero] = k
  798. if currentroot in result:
  799. result[currentroot] += k
  800. else:
  801. result[currentroot] = k
  802. def _try_decompose(f):
  803. """Find roots using functional decomposition. """
  804. factors, roots = f.decompose(), []
  805. for currentroot in _try_heuristics(factors[0]):
  806. roots.append(currentroot)
  807. for currentfactor in factors[1:]:
  808. previous, roots = list(roots), []
  809. for currentroot in previous:
  810. g = currentfactor - Poly(currentroot, f.gen)
  811. for currentroot in _try_heuristics(g):
  812. roots.append(currentroot)
  813. return roots
  814. def _try_heuristics(f):
  815. """Find roots using formulas and some tricks. """
  816. if f.is_ground:
  817. return []
  818. if f.is_monomial:
  819. return [S.Zero]*f.degree()
  820. if f.length() == 2:
  821. if f.degree() == 1:
  822. return list(map(cancel, roots_linear(f)))
  823. else:
  824. return roots_binomial(f)
  825. result = []
  826. for i in [-1, 1]:
  827. if not f.eval(i):
  828. f = f.quo(Poly(f.gen - i, f.gen))
  829. result.append(i)
  830. break
  831. n = f.degree()
  832. if n == 1:
  833. result += list(map(cancel, roots_linear(f)))
  834. elif n == 2:
  835. result += list(map(cancel, roots_quadratic(f)))
  836. elif f.is_cyclotomic:
  837. result += roots_cyclotomic(f)
  838. elif n == 3 and cubics:
  839. result += roots_cubic(f, trig=trig)
  840. elif n == 4 and quartics:
  841. result += roots_quartic(f)
  842. elif n == 5 and quintics:
  843. result += roots_quintic(f)
  844. return result
  845. # Convert the generators to symbols
  846. dumgens = symbols('x:%d' % len(f.gens), cls=Dummy)
  847. f = f.per(f.rep, dumgens)
  848. (k,), f = f.terms_gcd()
  849. if not k:
  850. zeros = {}
  851. else:
  852. zeros = {S.Zero: k}
  853. coeff, f = preprocess_roots(f)
  854. if auto and f.get_domain().is_Ring:
  855. f = f.to_field()
  856. # Use EX instead of ZZ_I or QQ_I
  857. if f.get_domain().is_QQ_I:
  858. f = f.per(f.rep.convert(EX))
  859. rescale_x = None
  860. translate_x = None
  861. result = {}
  862. if not f.is_ground:
  863. dom = f.get_domain()
  864. if not dom.is_Exact and dom.is_Numerical:
  865. for r in f.nroots():
  866. _update_dict(result, zeros, r, 1)
  867. elif f.degree() == 1:
  868. _update_dict(result, zeros, roots_linear(f)[0], 1)
  869. elif f.length() == 2:
  870. roots_fun = roots_quadratic if f.degree() == 2 else roots_binomial
  871. for r in roots_fun(f):
  872. _update_dict(result, zeros, r, 1)
  873. else:
  874. _, factors = Poly(f.as_expr()).factor_list()
  875. if len(factors) == 1 and f.degree() == 2:
  876. for r in roots_quadratic(f):
  877. _update_dict(result, zeros, r, 1)
  878. else:
  879. if len(factors) == 1 and factors[0][1] == 1:
  880. if f.get_domain().is_EX:
  881. res = to_rational_coeffs(f)
  882. if res:
  883. if res[0] is None:
  884. translate_x, f = res[2:]
  885. else:
  886. rescale_x, f = res[1], res[-1]
  887. result = roots(f)
  888. if not result:
  889. for currentroot in _try_decompose(f):
  890. _update_dict(result, zeros, currentroot, 1)
  891. else:
  892. for r in _try_heuristics(f):
  893. _update_dict(result, zeros, r, 1)
  894. else:
  895. for currentroot in _try_decompose(f):
  896. _update_dict(result, zeros, currentroot, 1)
  897. else:
  898. for currentfactor, k in factors:
  899. for r in _try_heuristics(Poly(currentfactor, f.gen, field=True)):
  900. _update_dict(result, zeros, r, k)
  901. if coeff is not S.One:
  902. _result, result, = result, {}
  903. for currentroot, k in _result.items():
  904. result[coeff*currentroot] = k
  905. if filter not in [None, 'C']:
  906. handlers = {
  907. 'Z': lambda r: r.is_Integer,
  908. 'Q': lambda r: r.is_Rational,
  909. 'R': lambda r: all(a.is_real for a in r.as_numer_denom()),
  910. 'I': lambda r: r.is_imaginary,
  911. }
  912. try:
  913. query = handlers[filter]
  914. except KeyError:
  915. raise ValueError("Invalid filter: %s" % filter)
  916. for zero in dict(result).keys():
  917. if not query(zero):
  918. del result[zero]
  919. if predicate is not None:
  920. for zero in dict(result).keys():
  921. if not predicate(zero):
  922. del result[zero]
  923. if rescale_x:
  924. result1 = {}
  925. for k, v in result.items():
  926. result1[k*rescale_x] = v
  927. result = result1
  928. if translate_x:
  929. result1 = {}
  930. for k, v in result.items():
  931. result1[k + translate_x] = v
  932. result = result1
  933. # adding zero roots after non-trivial roots have been translated
  934. result.update(zeros)
  935. if strict and sum(result.values()) < f.degree():
  936. raise UnsolvableFactorError(filldedent('''
  937. Strict mode: some factors cannot be solved in radicals, so
  938. a complete list of solutions cannot be returned. Call
  939. roots with strict=False to get solutions expressible in
  940. radicals (if there are any).
  941. '''))
  942. if not multiple:
  943. return result
  944. else:
  945. zeros = []
  946. for zero in ordered(result):
  947. zeros.extend([zero]*result[zero])
  948. return zeros
  949. def root_factors(f, *gens, filter=None, **args):
  950. """
  951. Returns all factors of a univariate polynomial.
  952. Examples
  953. ========
  954. >>> from sympy.abc import x, y
  955. >>> from sympy.polys.polyroots import root_factors
  956. >>> root_factors(x**2 - y, x)
  957. [x - sqrt(y), x + sqrt(y)]
  958. """
  959. args = dict(args)
  960. F = Poly(f, *gens, **args)
  961. if not F.is_Poly:
  962. return [f]
  963. if F.is_multivariate:
  964. raise ValueError('multivariate polynomials are not supported')
  965. x = F.gens[0]
  966. zeros = roots(F, filter=filter)
  967. if not zeros:
  968. factors = [F]
  969. else:
  970. factors, N = [], 0
  971. for r, n in ordered(zeros.items()):
  972. factors, N = factors + [Poly(x - r, x)]*n, N + n
  973. if N < F.degree():
  974. G = reduce(lambda p, q: p*q, factors)
  975. factors.append(F.quo(G))
  976. if not isinstance(f, Poly):
  977. factors = [ f.as_expr() for f in factors ]
  978. return factors