euclidtools.py 41 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789179017911792179317941795179617971798179918001801180218031804180518061807180818091810181118121813181418151816181718181819182018211822182318241825182618271828182918301831183218331834183518361837183818391840184118421843184418451846184718481849185018511852185318541855185618571858185918601861186218631864186518661867186818691870187118721873187418751876187718781879188018811882188318841885188618871888188918901891189218931894189518961897189818991900190119021903190419051906190719081909191019111912
  1. """Euclidean algorithms, GCDs, LCMs and polynomial remainder sequences. """
  2. from sympy.polys.densearith import (
  3. dup_sub_mul,
  4. dup_neg, dmp_neg,
  5. dmp_add,
  6. dmp_sub,
  7. dup_mul, dmp_mul,
  8. dmp_pow,
  9. dup_div, dmp_div,
  10. dup_rem,
  11. dup_quo, dmp_quo,
  12. dup_prem, dmp_prem,
  13. dup_mul_ground, dmp_mul_ground,
  14. dmp_mul_term,
  15. dup_quo_ground, dmp_quo_ground,
  16. dup_max_norm, dmp_max_norm)
  17. from sympy.polys.densebasic import (
  18. dup_strip, dmp_raise,
  19. dmp_zero, dmp_one, dmp_ground,
  20. dmp_one_p, dmp_zero_p,
  21. dmp_zeros,
  22. dup_degree, dmp_degree, dmp_degree_in,
  23. dup_LC, dmp_LC, dmp_ground_LC,
  24. dmp_multi_deflate, dmp_inflate,
  25. dup_convert, dmp_convert,
  26. dmp_apply_pairs)
  27. from sympy.polys.densetools import (
  28. dup_clear_denoms, dmp_clear_denoms,
  29. dup_diff, dmp_diff,
  30. dup_eval, dmp_eval, dmp_eval_in,
  31. dup_trunc, dmp_ground_trunc,
  32. dup_monic, dmp_ground_monic,
  33. dup_primitive, dmp_ground_primitive,
  34. dup_extract, dmp_ground_extract)
  35. from sympy.polys.galoistools import (
  36. gf_int, gf_crt)
  37. from sympy.polys.polyconfig import query
  38. from sympy.polys.polyerrors import (
  39. MultivariatePolynomialError,
  40. HeuristicGCDFailed,
  41. HomomorphismFailed,
  42. NotInvertible,
  43. DomainError)
  44. def dup_half_gcdex(f, g, K):
  45. """
  46. Half extended Euclidean algorithm in `F[x]`.
  47. Returns ``(s, h)`` such that ``h = gcd(f, g)`` and ``s*f = h (mod g)``.
  48. Examples
  49. ========
  50. >>> from sympy.polys import ring, QQ
  51. >>> R, x = ring("x", QQ)
  52. >>> f = x**4 - 2*x**3 - 6*x**2 + 12*x + 15
  53. >>> g = x**3 + x**2 - 4*x - 4
  54. >>> R.dup_half_gcdex(f, g)
  55. (-1/5*x + 3/5, x + 1)
  56. """
  57. if not K.is_Field:
  58. raise DomainError("Cannot compute half extended GCD over %s" % K)
  59. a, b = [K.one], []
  60. while g:
  61. q, r = dup_div(f, g, K)
  62. f, g = g, r
  63. a, b = b, dup_sub_mul(a, q, b, K)
  64. a = dup_quo_ground(a, dup_LC(f, K), K)
  65. f = dup_monic(f, K)
  66. return a, f
  67. def dmp_half_gcdex(f, g, u, K):
  68. """
  69. Half extended Euclidean algorithm in `F[X]`.
  70. Examples
  71. ========
  72. >>> from sympy.polys import ring, ZZ
  73. >>> R, x,y = ring("x,y", ZZ)
  74. """
  75. if not u:
  76. return dup_half_gcdex(f, g, K)
  77. else:
  78. raise MultivariatePolynomialError(f, g)
  79. def dup_gcdex(f, g, K):
  80. """
  81. Extended Euclidean algorithm in `F[x]`.
  82. Returns ``(s, t, h)`` such that ``h = gcd(f, g)`` and ``s*f + t*g = h``.
  83. Examples
  84. ========
  85. >>> from sympy.polys import ring, QQ
  86. >>> R, x = ring("x", QQ)
  87. >>> f = x**4 - 2*x**3 - 6*x**2 + 12*x + 15
  88. >>> g = x**3 + x**2 - 4*x - 4
  89. >>> R.dup_gcdex(f, g)
  90. (-1/5*x + 3/5, 1/5*x**2 - 6/5*x + 2, x + 1)
  91. """
  92. s, h = dup_half_gcdex(f, g, K)
  93. F = dup_sub_mul(h, s, f, K)
  94. t = dup_quo(F, g, K)
  95. return s, t, h
  96. def dmp_gcdex(f, g, u, K):
  97. """
  98. Extended Euclidean algorithm in `F[X]`.
  99. Examples
  100. ========
  101. >>> from sympy.polys import ring, ZZ
  102. >>> R, x,y = ring("x,y", ZZ)
  103. """
  104. if not u:
  105. return dup_gcdex(f, g, K)
  106. else:
  107. raise MultivariatePolynomialError(f, g)
  108. def dup_invert(f, g, K):
  109. """
  110. Compute multiplicative inverse of `f` modulo `g` in `F[x]`.
  111. Examples
  112. ========
  113. >>> from sympy.polys import ring, QQ
  114. >>> R, x = ring("x", QQ)
  115. >>> f = x**2 - 1
  116. >>> g = 2*x - 1
  117. >>> h = x - 1
  118. >>> R.dup_invert(f, g)
  119. -4/3
  120. >>> R.dup_invert(f, h)
  121. Traceback (most recent call last):
  122. ...
  123. NotInvertible: zero divisor
  124. """
  125. s, h = dup_half_gcdex(f, g, K)
  126. if h == [K.one]:
  127. return dup_rem(s, g, K)
  128. else:
  129. raise NotInvertible("zero divisor")
  130. def dmp_invert(f, g, u, K):
  131. """
  132. Compute multiplicative inverse of `f` modulo `g` in `F[X]`.
  133. Examples
  134. ========
  135. >>> from sympy.polys import ring, QQ
  136. >>> R, x = ring("x", QQ)
  137. """
  138. if not u:
  139. return dup_invert(f, g, K)
  140. else:
  141. raise MultivariatePolynomialError(f, g)
  142. def dup_euclidean_prs(f, g, K):
  143. """
  144. Euclidean polynomial remainder sequence (PRS) in `K[x]`.
  145. Examples
  146. ========
  147. >>> from sympy.polys import ring, QQ
  148. >>> R, x = ring("x", QQ)
  149. >>> f = x**8 + x**6 - 3*x**4 - 3*x**3 + 8*x**2 + 2*x - 5
  150. >>> g = 3*x**6 + 5*x**4 - 4*x**2 - 9*x + 21
  151. >>> prs = R.dup_euclidean_prs(f, g)
  152. >>> prs[0]
  153. x**8 + x**6 - 3*x**4 - 3*x**3 + 8*x**2 + 2*x - 5
  154. >>> prs[1]
  155. 3*x**6 + 5*x**4 - 4*x**2 - 9*x + 21
  156. >>> prs[2]
  157. -5/9*x**4 + 1/9*x**2 - 1/3
  158. >>> prs[3]
  159. -117/25*x**2 - 9*x + 441/25
  160. >>> prs[4]
  161. 233150/19773*x - 102500/6591
  162. >>> prs[5]
  163. -1288744821/543589225
  164. """
  165. prs = [f, g]
  166. h = dup_rem(f, g, K)
  167. while h:
  168. prs.append(h)
  169. f, g = g, h
  170. h = dup_rem(f, g, K)
  171. return prs
  172. def dmp_euclidean_prs(f, g, u, K):
  173. """
  174. Euclidean polynomial remainder sequence (PRS) in `K[X]`.
  175. Examples
  176. ========
  177. >>> from sympy.polys import ring, ZZ
  178. >>> R, x,y = ring("x,y", ZZ)
  179. """
  180. if not u:
  181. return dup_euclidean_prs(f, g, K)
  182. else:
  183. raise MultivariatePolynomialError(f, g)
  184. def dup_primitive_prs(f, g, K):
  185. """
  186. Primitive polynomial remainder sequence (PRS) in `K[x]`.
  187. Examples
  188. ========
  189. >>> from sympy.polys import ring, ZZ
  190. >>> R, x = ring("x", ZZ)
  191. >>> f = x**8 + x**6 - 3*x**4 - 3*x**3 + 8*x**2 + 2*x - 5
  192. >>> g = 3*x**6 + 5*x**4 - 4*x**2 - 9*x + 21
  193. >>> prs = R.dup_primitive_prs(f, g)
  194. >>> prs[0]
  195. x**8 + x**6 - 3*x**4 - 3*x**3 + 8*x**2 + 2*x - 5
  196. >>> prs[1]
  197. 3*x**6 + 5*x**4 - 4*x**2 - 9*x + 21
  198. >>> prs[2]
  199. -5*x**4 + x**2 - 3
  200. >>> prs[3]
  201. 13*x**2 + 25*x - 49
  202. >>> prs[4]
  203. 4663*x - 6150
  204. >>> prs[5]
  205. 1
  206. """
  207. prs = [f, g]
  208. _, h = dup_primitive(dup_prem(f, g, K), K)
  209. while h:
  210. prs.append(h)
  211. f, g = g, h
  212. _, h = dup_primitive(dup_prem(f, g, K), K)
  213. return prs
  214. def dmp_primitive_prs(f, g, u, K):
  215. """
  216. Primitive polynomial remainder sequence (PRS) in `K[X]`.
  217. Examples
  218. ========
  219. >>> from sympy.polys import ring, ZZ
  220. >>> R, x,y = ring("x,y", ZZ)
  221. """
  222. if not u:
  223. return dup_primitive_prs(f, g, K)
  224. else:
  225. raise MultivariatePolynomialError(f, g)
  226. def dup_inner_subresultants(f, g, K):
  227. """
  228. Subresultant PRS algorithm in `K[x]`.
  229. Computes the subresultant polynomial remainder sequence (PRS)
  230. and the non-zero scalar subresultants of `f` and `g`.
  231. By [1] Thm. 3, these are the constants '-c' (- to optimize
  232. computation of sign).
  233. The first subdeterminant is set to 1 by convention to match
  234. the polynomial and the scalar subdeterminants.
  235. If 'deg(f) < deg(g)', the subresultants of '(g,f)' are computed.
  236. Examples
  237. ========
  238. >>> from sympy.polys import ring, ZZ
  239. >>> R, x = ring("x", ZZ)
  240. >>> R.dup_inner_subresultants(x**2 + 1, x**2 - 1)
  241. ([x**2 + 1, x**2 - 1, -2], [1, 1, 4])
  242. References
  243. ==========
  244. .. [1] W.S. Brown, The Subresultant PRS Algorithm.
  245. ACM Transaction of Mathematical Software 4 (1978) 237-249
  246. """
  247. n = dup_degree(f)
  248. m = dup_degree(g)
  249. if n < m:
  250. f, g = g, f
  251. n, m = m, n
  252. if not f:
  253. return [], []
  254. if not g:
  255. return [f], [K.one]
  256. R = [f, g]
  257. d = n - m
  258. b = (-K.one)**(d + 1)
  259. h = dup_prem(f, g, K)
  260. h = dup_mul_ground(h, b, K)
  261. lc = dup_LC(g, K)
  262. c = lc**d
  263. # Conventional first scalar subdeterminant is 1
  264. S = [K.one, c]
  265. c = -c
  266. while h:
  267. k = dup_degree(h)
  268. R.append(h)
  269. f, g, m, d = g, h, k, m - k
  270. b = -lc * c**d
  271. h = dup_prem(f, g, K)
  272. h = dup_quo_ground(h, b, K)
  273. lc = dup_LC(g, K)
  274. if d > 1: # abnormal case
  275. q = c**(d - 1)
  276. c = K.quo((-lc)**d, q)
  277. else:
  278. c = -lc
  279. S.append(-c)
  280. return R, S
  281. def dup_subresultants(f, g, K):
  282. """
  283. Computes subresultant PRS of two polynomials in `K[x]`.
  284. Examples
  285. ========
  286. >>> from sympy.polys import ring, ZZ
  287. >>> R, x = ring("x", ZZ)
  288. >>> R.dup_subresultants(x**2 + 1, x**2 - 1)
  289. [x**2 + 1, x**2 - 1, -2]
  290. """
  291. return dup_inner_subresultants(f, g, K)[0]
  292. def dup_prs_resultant(f, g, K):
  293. """
  294. Resultant algorithm in `K[x]` using subresultant PRS.
  295. Examples
  296. ========
  297. >>> from sympy.polys import ring, ZZ
  298. >>> R, x = ring("x", ZZ)
  299. >>> R.dup_prs_resultant(x**2 + 1, x**2 - 1)
  300. (4, [x**2 + 1, x**2 - 1, -2])
  301. """
  302. if not f or not g:
  303. return (K.zero, [])
  304. R, S = dup_inner_subresultants(f, g, K)
  305. if dup_degree(R[-1]) > 0:
  306. return (K.zero, R)
  307. return S[-1], R
  308. def dup_resultant(f, g, K, includePRS=False):
  309. """
  310. Computes resultant of two polynomials in `K[x]`.
  311. Examples
  312. ========
  313. >>> from sympy.polys import ring, ZZ
  314. >>> R, x = ring("x", ZZ)
  315. >>> R.dup_resultant(x**2 + 1, x**2 - 1)
  316. 4
  317. """
  318. if includePRS:
  319. return dup_prs_resultant(f, g, K)
  320. return dup_prs_resultant(f, g, K)[0]
  321. def dmp_inner_subresultants(f, g, u, K):
  322. """
  323. Subresultant PRS algorithm in `K[X]`.
  324. Examples
  325. ========
  326. >>> from sympy.polys import ring, ZZ
  327. >>> R, x,y = ring("x,y", ZZ)
  328. >>> f = 3*x**2*y - y**3 - 4
  329. >>> g = x**2 + x*y**3 - 9
  330. >>> a = 3*x*y**4 + y**3 - 27*y + 4
  331. >>> b = -3*y**10 - 12*y**7 + y**6 - 54*y**4 + 8*y**3 + 729*y**2 - 216*y + 16
  332. >>> prs = [f, g, a, b]
  333. >>> sres = [[1], [1], [3, 0, 0, 0, 0], [-3, 0, 0, -12, 1, 0, -54, 8, 729, -216, 16]]
  334. >>> R.dmp_inner_subresultants(f, g) == (prs, sres)
  335. True
  336. """
  337. if not u:
  338. return dup_inner_subresultants(f, g, K)
  339. n = dmp_degree(f, u)
  340. m = dmp_degree(g, u)
  341. if n < m:
  342. f, g = g, f
  343. n, m = m, n
  344. if dmp_zero_p(f, u):
  345. return [], []
  346. v = u - 1
  347. if dmp_zero_p(g, u):
  348. return [f], [dmp_ground(K.one, v)]
  349. R = [f, g]
  350. d = n - m
  351. b = dmp_pow(dmp_ground(-K.one, v), d + 1, v, K)
  352. h = dmp_prem(f, g, u, K)
  353. h = dmp_mul_term(h, b, 0, u, K)
  354. lc = dmp_LC(g, K)
  355. c = dmp_pow(lc, d, v, K)
  356. S = [dmp_ground(K.one, v), c]
  357. c = dmp_neg(c, v, K)
  358. while not dmp_zero_p(h, u):
  359. k = dmp_degree(h, u)
  360. R.append(h)
  361. f, g, m, d = g, h, k, m - k
  362. b = dmp_mul(dmp_neg(lc, v, K),
  363. dmp_pow(c, d, v, K), v, K)
  364. h = dmp_prem(f, g, u, K)
  365. h = [ dmp_quo(ch, b, v, K) for ch in h ]
  366. lc = dmp_LC(g, K)
  367. if d > 1:
  368. p = dmp_pow(dmp_neg(lc, v, K), d, v, K)
  369. q = dmp_pow(c, d - 1, v, K)
  370. c = dmp_quo(p, q, v, K)
  371. else:
  372. c = dmp_neg(lc, v, K)
  373. S.append(dmp_neg(c, v, K))
  374. return R, S
  375. def dmp_subresultants(f, g, u, K):
  376. """
  377. Computes subresultant PRS of two polynomials in `K[X]`.
  378. Examples
  379. ========
  380. >>> from sympy.polys import ring, ZZ
  381. >>> R, x,y = ring("x,y", ZZ)
  382. >>> f = 3*x**2*y - y**3 - 4
  383. >>> g = x**2 + x*y**3 - 9
  384. >>> a = 3*x*y**4 + y**3 - 27*y + 4
  385. >>> b = -3*y**10 - 12*y**7 + y**6 - 54*y**4 + 8*y**3 + 729*y**2 - 216*y + 16
  386. >>> R.dmp_subresultants(f, g) == [f, g, a, b]
  387. True
  388. """
  389. return dmp_inner_subresultants(f, g, u, K)[0]
  390. def dmp_prs_resultant(f, g, u, K):
  391. """
  392. Resultant algorithm in `K[X]` using subresultant PRS.
  393. Examples
  394. ========
  395. >>> from sympy.polys import ring, ZZ
  396. >>> R, x,y = ring("x,y", ZZ)
  397. >>> f = 3*x**2*y - y**3 - 4
  398. >>> g = x**2 + x*y**3 - 9
  399. >>> a = 3*x*y**4 + y**3 - 27*y + 4
  400. >>> b = -3*y**10 - 12*y**7 + y**6 - 54*y**4 + 8*y**3 + 729*y**2 - 216*y + 16
  401. >>> res, prs = R.dmp_prs_resultant(f, g)
  402. >>> res == b # resultant has n-1 variables
  403. False
  404. >>> res == b.drop(x)
  405. True
  406. >>> prs == [f, g, a, b]
  407. True
  408. """
  409. if not u:
  410. return dup_prs_resultant(f, g, K)
  411. if dmp_zero_p(f, u) or dmp_zero_p(g, u):
  412. return (dmp_zero(u - 1), [])
  413. R, S = dmp_inner_subresultants(f, g, u, K)
  414. if dmp_degree(R[-1], u) > 0:
  415. return (dmp_zero(u - 1), R)
  416. return S[-1], R
  417. def dmp_zz_modular_resultant(f, g, p, u, K):
  418. """
  419. Compute resultant of `f` and `g` modulo a prime `p`.
  420. Examples
  421. ========
  422. >>> from sympy.polys import ring, ZZ
  423. >>> R, x,y = ring("x,y", ZZ)
  424. >>> f = x + y + 2
  425. >>> g = 2*x*y + x + 3
  426. >>> R.dmp_zz_modular_resultant(f, g, 5)
  427. -2*y**2 + 1
  428. """
  429. if not u:
  430. return gf_int(dup_prs_resultant(f, g, K)[0] % p, p)
  431. v = u - 1
  432. n = dmp_degree(f, u)
  433. m = dmp_degree(g, u)
  434. N = dmp_degree_in(f, 1, u)
  435. M = dmp_degree_in(g, 1, u)
  436. B = n*M + m*N
  437. D, a = [K.one], -K.one
  438. r = dmp_zero(v)
  439. while dup_degree(D) <= B:
  440. while True:
  441. a += K.one
  442. if a == p:
  443. raise HomomorphismFailed('no luck')
  444. F = dmp_eval_in(f, gf_int(a, p), 1, u, K)
  445. if dmp_degree(F, v) == n:
  446. G = dmp_eval_in(g, gf_int(a, p), 1, u, K)
  447. if dmp_degree(G, v) == m:
  448. break
  449. R = dmp_zz_modular_resultant(F, G, p, v, K)
  450. e = dmp_eval(r, a, v, K)
  451. if not v:
  452. R = dup_strip([R])
  453. e = dup_strip([e])
  454. else:
  455. R = [R]
  456. e = [e]
  457. d = K.invert(dup_eval(D, a, K), p)
  458. d = dup_mul_ground(D, d, K)
  459. d = dmp_raise(d, v, 0, K)
  460. c = dmp_mul(d, dmp_sub(R, e, v, K), v, K)
  461. r = dmp_add(r, c, v, K)
  462. r = dmp_ground_trunc(r, p, v, K)
  463. D = dup_mul(D, [K.one, -a], K)
  464. D = dup_trunc(D, p, K)
  465. return r
  466. def _collins_crt(r, R, P, p, K):
  467. """Wrapper of CRT for Collins's resultant algorithm. """
  468. return gf_int(gf_crt([r, R], [P, p], K), P*p)
  469. def dmp_zz_collins_resultant(f, g, u, K):
  470. """
  471. Collins's modular resultant algorithm in `Z[X]`.
  472. Examples
  473. ========
  474. >>> from sympy.polys import ring, ZZ
  475. >>> R, x,y = ring("x,y", ZZ)
  476. >>> f = x + y + 2
  477. >>> g = 2*x*y + x + 3
  478. >>> R.dmp_zz_collins_resultant(f, g)
  479. -2*y**2 - 5*y + 1
  480. """
  481. n = dmp_degree(f, u)
  482. m = dmp_degree(g, u)
  483. if n < 0 or m < 0:
  484. return dmp_zero(u - 1)
  485. A = dmp_max_norm(f, u, K)
  486. B = dmp_max_norm(g, u, K)
  487. a = dmp_ground_LC(f, u, K)
  488. b = dmp_ground_LC(g, u, K)
  489. v = u - 1
  490. B = K(2)*K.factorial(K(n + m))*A**m*B**n
  491. r, p, P = dmp_zero(v), K.one, K.one
  492. from sympy.ntheory import nextprime
  493. while P <= B:
  494. p = K(nextprime(p))
  495. while not (a % p) or not (b % p):
  496. p = K(nextprime(p))
  497. F = dmp_ground_trunc(f, p, u, K)
  498. G = dmp_ground_trunc(g, p, u, K)
  499. try:
  500. R = dmp_zz_modular_resultant(F, G, p, u, K)
  501. except HomomorphismFailed:
  502. continue
  503. if K.is_one(P):
  504. r = R
  505. else:
  506. r = dmp_apply_pairs(r, R, _collins_crt, (P, p, K), v, K)
  507. P *= p
  508. return r
  509. def dmp_qq_collins_resultant(f, g, u, K0):
  510. """
  511. Collins's modular resultant algorithm in `Q[X]`.
  512. Examples
  513. ========
  514. >>> from sympy.polys import ring, QQ
  515. >>> R, x,y = ring("x,y", QQ)
  516. >>> f = QQ(1,2)*x + y + QQ(2,3)
  517. >>> g = 2*x*y + x + 3
  518. >>> R.dmp_qq_collins_resultant(f, g)
  519. -2*y**2 - 7/3*y + 5/6
  520. """
  521. n = dmp_degree(f, u)
  522. m = dmp_degree(g, u)
  523. if n < 0 or m < 0:
  524. return dmp_zero(u - 1)
  525. K1 = K0.get_ring()
  526. cf, f = dmp_clear_denoms(f, u, K0, K1)
  527. cg, g = dmp_clear_denoms(g, u, K0, K1)
  528. f = dmp_convert(f, u, K0, K1)
  529. g = dmp_convert(g, u, K0, K1)
  530. r = dmp_zz_collins_resultant(f, g, u, K1)
  531. r = dmp_convert(r, u - 1, K1, K0)
  532. c = K0.convert(cf**m * cg**n, K1)
  533. return dmp_quo_ground(r, c, u - 1, K0)
  534. def dmp_resultant(f, g, u, K, includePRS=False):
  535. """
  536. Computes resultant of two polynomials in `K[X]`.
  537. Examples
  538. ========
  539. >>> from sympy.polys import ring, ZZ
  540. >>> R, x,y = ring("x,y", ZZ)
  541. >>> f = 3*x**2*y - y**3 - 4
  542. >>> g = x**2 + x*y**3 - 9
  543. >>> R.dmp_resultant(f, g)
  544. -3*y**10 - 12*y**7 + y**6 - 54*y**4 + 8*y**3 + 729*y**2 - 216*y + 16
  545. """
  546. if not u:
  547. return dup_resultant(f, g, K, includePRS=includePRS)
  548. if includePRS:
  549. return dmp_prs_resultant(f, g, u, K)
  550. if K.is_Field:
  551. if K.is_QQ and query('USE_COLLINS_RESULTANT'):
  552. return dmp_qq_collins_resultant(f, g, u, K)
  553. else:
  554. if K.is_ZZ and query('USE_COLLINS_RESULTANT'):
  555. return dmp_zz_collins_resultant(f, g, u, K)
  556. return dmp_prs_resultant(f, g, u, K)[0]
  557. def dup_discriminant(f, K):
  558. """
  559. Computes discriminant of a polynomial in `K[x]`.
  560. Examples
  561. ========
  562. >>> from sympy.polys import ring, ZZ
  563. >>> R, x = ring("x", ZZ)
  564. >>> R.dup_discriminant(x**2 + 2*x + 3)
  565. -8
  566. """
  567. d = dup_degree(f)
  568. if d <= 0:
  569. return K.zero
  570. else:
  571. s = (-1)**((d*(d - 1)) // 2)
  572. c = dup_LC(f, K)
  573. r = dup_resultant(f, dup_diff(f, 1, K), K)
  574. return K.quo(r, c*K(s))
  575. def dmp_discriminant(f, u, K):
  576. """
  577. Computes discriminant of a polynomial in `K[X]`.
  578. Examples
  579. ========
  580. >>> from sympy.polys import ring, ZZ
  581. >>> R, x,y,z,t = ring("x,y,z,t", ZZ)
  582. >>> R.dmp_discriminant(x**2*y + x*z + t)
  583. -4*y*t + z**2
  584. """
  585. if not u:
  586. return dup_discriminant(f, K)
  587. d, v = dmp_degree(f, u), u - 1
  588. if d <= 0:
  589. return dmp_zero(v)
  590. else:
  591. s = (-1)**((d*(d - 1)) // 2)
  592. c = dmp_LC(f, K)
  593. r = dmp_resultant(f, dmp_diff(f, 1, u, K), u, K)
  594. c = dmp_mul_ground(c, K(s), v, K)
  595. return dmp_quo(r, c, v, K)
  596. def _dup_rr_trivial_gcd(f, g, K):
  597. """Handle trivial cases in GCD algorithm over a ring. """
  598. if not (f or g):
  599. return [], [], []
  600. elif not f:
  601. if K.is_nonnegative(dup_LC(g, K)):
  602. return g, [], [K.one]
  603. else:
  604. return dup_neg(g, K), [], [-K.one]
  605. elif not g:
  606. if K.is_nonnegative(dup_LC(f, K)):
  607. return f, [K.one], []
  608. else:
  609. return dup_neg(f, K), [-K.one], []
  610. return None
  611. def _dup_ff_trivial_gcd(f, g, K):
  612. """Handle trivial cases in GCD algorithm over a field. """
  613. if not (f or g):
  614. return [], [], []
  615. elif not f:
  616. return dup_monic(g, K), [], [dup_LC(g, K)]
  617. elif not g:
  618. return dup_monic(f, K), [dup_LC(f, K)], []
  619. else:
  620. return None
  621. def _dmp_rr_trivial_gcd(f, g, u, K):
  622. """Handle trivial cases in GCD algorithm over a ring. """
  623. zero_f = dmp_zero_p(f, u)
  624. zero_g = dmp_zero_p(g, u)
  625. if_contain_one = dmp_one_p(f, u, K) or dmp_one_p(g, u, K)
  626. if zero_f and zero_g:
  627. return tuple(dmp_zeros(3, u, K))
  628. elif zero_f:
  629. if K.is_nonnegative(dmp_ground_LC(g, u, K)):
  630. return g, dmp_zero(u), dmp_one(u, K)
  631. else:
  632. return dmp_neg(g, u, K), dmp_zero(u), dmp_ground(-K.one, u)
  633. elif zero_g:
  634. if K.is_nonnegative(dmp_ground_LC(f, u, K)):
  635. return f, dmp_one(u, K), dmp_zero(u)
  636. else:
  637. return dmp_neg(f, u, K), dmp_ground(-K.one, u), dmp_zero(u)
  638. elif if_contain_one:
  639. return dmp_one(u, K), f, g
  640. elif query('USE_SIMPLIFY_GCD'):
  641. return _dmp_simplify_gcd(f, g, u, K)
  642. else:
  643. return None
  644. def _dmp_ff_trivial_gcd(f, g, u, K):
  645. """Handle trivial cases in GCD algorithm over a field. """
  646. zero_f = dmp_zero_p(f, u)
  647. zero_g = dmp_zero_p(g, u)
  648. if zero_f and zero_g:
  649. return tuple(dmp_zeros(3, u, K))
  650. elif zero_f:
  651. return (dmp_ground_monic(g, u, K),
  652. dmp_zero(u),
  653. dmp_ground(dmp_ground_LC(g, u, K), u))
  654. elif zero_g:
  655. return (dmp_ground_monic(f, u, K),
  656. dmp_ground(dmp_ground_LC(f, u, K), u),
  657. dmp_zero(u))
  658. elif query('USE_SIMPLIFY_GCD'):
  659. return _dmp_simplify_gcd(f, g, u, K)
  660. else:
  661. return None
  662. def _dmp_simplify_gcd(f, g, u, K):
  663. """Try to eliminate `x_0` from GCD computation in `K[X]`. """
  664. df = dmp_degree(f, u)
  665. dg = dmp_degree(g, u)
  666. if df > 0 and dg > 0:
  667. return None
  668. if not (df or dg):
  669. F = dmp_LC(f, K)
  670. G = dmp_LC(g, K)
  671. else:
  672. if not df:
  673. F = dmp_LC(f, K)
  674. G = dmp_content(g, u, K)
  675. else:
  676. F = dmp_content(f, u, K)
  677. G = dmp_LC(g, K)
  678. v = u - 1
  679. h = dmp_gcd(F, G, v, K)
  680. cff = [ dmp_quo(cf, h, v, K) for cf in f ]
  681. cfg = [ dmp_quo(cg, h, v, K) for cg in g ]
  682. return [h], cff, cfg
  683. def dup_rr_prs_gcd(f, g, K):
  684. """
  685. Computes polynomial GCD using subresultants over a ring.
  686. Returns ``(h, cff, cfg)`` such that ``a = gcd(f, g)``, ``cff = quo(f, h)``,
  687. and ``cfg = quo(g, h)``.
  688. Examples
  689. ========
  690. >>> from sympy.polys import ring, ZZ
  691. >>> R, x = ring("x", ZZ)
  692. >>> R.dup_rr_prs_gcd(x**2 - 1, x**2 - 3*x + 2)
  693. (x - 1, x + 1, x - 2)
  694. """
  695. result = _dup_rr_trivial_gcd(f, g, K)
  696. if result is not None:
  697. return result
  698. fc, F = dup_primitive(f, K)
  699. gc, G = dup_primitive(g, K)
  700. c = K.gcd(fc, gc)
  701. h = dup_subresultants(F, G, K)[-1]
  702. _, h = dup_primitive(h, K)
  703. c *= K.canonical_unit(dup_LC(h, K))
  704. h = dup_mul_ground(h, c, K)
  705. cff = dup_quo(f, h, K)
  706. cfg = dup_quo(g, h, K)
  707. return h, cff, cfg
  708. def dup_ff_prs_gcd(f, g, K):
  709. """
  710. Computes polynomial GCD using subresultants over a field.
  711. Returns ``(h, cff, cfg)`` such that ``a = gcd(f, g)``, ``cff = quo(f, h)``,
  712. and ``cfg = quo(g, h)``.
  713. Examples
  714. ========
  715. >>> from sympy.polys import ring, QQ
  716. >>> R, x = ring("x", QQ)
  717. >>> R.dup_ff_prs_gcd(x**2 - 1, x**2 - 3*x + 2)
  718. (x - 1, x + 1, x - 2)
  719. """
  720. result = _dup_ff_trivial_gcd(f, g, K)
  721. if result is not None:
  722. return result
  723. h = dup_subresultants(f, g, K)[-1]
  724. h = dup_monic(h, K)
  725. cff = dup_quo(f, h, K)
  726. cfg = dup_quo(g, h, K)
  727. return h, cff, cfg
  728. def dmp_rr_prs_gcd(f, g, u, K):
  729. """
  730. Computes polynomial GCD using subresultants over a ring.
  731. Returns ``(h, cff, cfg)`` such that ``a = gcd(f, g)``, ``cff = quo(f, h)``,
  732. and ``cfg = quo(g, h)``.
  733. Examples
  734. ========
  735. >>> from sympy.polys import ring, ZZ
  736. >>> R, x,y, = ring("x,y", ZZ)
  737. >>> f = x**2 + 2*x*y + y**2
  738. >>> g = x**2 + x*y
  739. >>> R.dmp_rr_prs_gcd(f, g)
  740. (x + y, x + y, x)
  741. """
  742. if not u:
  743. return dup_rr_prs_gcd(f, g, K)
  744. result = _dmp_rr_trivial_gcd(f, g, u, K)
  745. if result is not None:
  746. return result
  747. fc, F = dmp_primitive(f, u, K)
  748. gc, G = dmp_primitive(g, u, K)
  749. h = dmp_subresultants(F, G, u, K)[-1]
  750. c, _, _ = dmp_rr_prs_gcd(fc, gc, u - 1, K)
  751. _, h = dmp_primitive(h, u, K)
  752. h = dmp_mul_term(h, c, 0, u, K)
  753. unit = K.canonical_unit(dmp_ground_LC(h, u, K))
  754. if unit != K.one:
  755. h = dmp_mul_ground(h, unit, u, K)
  756. cff = dmp_quo(f, h, u, K)
  757. cfg = dmp_quo(g, h, u, K)
  758. return h, cff, cfg
  759. def dmp_ff_prs_gcd(f, g, u, K):
  760. """
  761. Computes polynomial GCD using subresultants over a field.
  762. Returns ``(h, cff, cfg)`` such that ``a = gcd(f, g)``, ``cff = quo(f, h)``,
  763. and ``cfg = quo(g, h)``.
  764. Examples
  765. ========
  766. >>> from sympy.polys import ring, QQ
  767. >>> R, x,y, = ring("x,y", QQ)
  768. >>> f = QQ(1,2)*x**2 + x*y + QQ(1,2)*y**2
  769. >>> g = x**2 + x*y
  770. >>> R.dmp_ff_prs_gcd(f, g)
  771. (x + y, 1/2*x + 1/2*y, x)
  772. """
  773. if not u:
  774. return dup_ff_prs_gcd(f, g, K)
  775. result = _dmp_ff_trivial_gcd(f, g, u, K)
  776. if result is not None:
  777. return result
  778. fc, F = dmp_primitive(f, u, K)
  779. gc, G = dmp_primitive(g, u, K)
  780. h = dmp_subresultants(F, G, u, K)[-1]
  781. c, _, _ = dmp_ff_prs_gcd(fc, gc, u - 1, K)
  782. _, h = dmp_primitive(h, u, K)
  783. h = dmp_mul_term(h, c, 0, u, K)
  784. h = dmp_ground_monic(h, u, K)
  785. cff = dmp_quo(f, h, u, K)
  786. cfg = dmp_quo(g, h, u, K)
  787. return h, cff, cfg
  788. HEU_GCD_MAX = 6
  789. def _dup_zz_gcd_interpolate(h, x, K):
  790. """Interpolate polynomial GCD from integer GCD. """
  791. f = []
  792. while h:
  793. g = h % x
  794. if g > x // 2:
  795. g -= x
  796. f.insert(0, g)
  797. h = (h - g) // x
  798. return f
  799. def dup_zz_heu_gcd(f, g, K):
  800. """
  801. Heuristic polynomial GCD in `Z[x]`.
  802. Given univariate polynomials `f` and `g` in `Z[x]`, returns
  803. their GCD and cofactors, i.e. polynomials ``h``, ``cff`` and ``cfg``
  804. such that::
  805. h = gcd(f, g), cff = quo(f, h) and cfg = quo(g, h)
  806. The algorithm is purely heuristic which means it may fail to compute
  807. the GCD. This will be signaled by raising an exception. In this case
  808. you will need to switch to another GCD method.
  809. The algorithm computes the polynomial GCD by evaluating polynomials
  810. f and g at certain points and computing (fast) integer GCD of those
  811. evaluations. The polynomial GCD is recovered from the integer image
  812. by interpolation. The final step is to verify if the result is the
  813. correct GCD. This gives cofactors as a side effect.
  814. Examples
  815. ========
  816. >>> from sympy.polys import ring, ZZ
  817. >>> R, x = ring("x", ZZ)
  818. >>> R.dup_zz_heu_gcd(x**2 - 1, x**2 - 3*x + 2)
  819. (x - 1, x + 1, x - 2)
  820. References
  821. ==========
  822. .. [1] [Liao95]_
  823. """
  824. result = _dup_rr_trivial_gcd(f, g, K)
  825. if result is not None:
  826. return result
  827. df = dup_degree(f)
  828. dg = dup_degree(g)
  829. gcd, f, g = dup_extract(f, g, K)
  830. if df == 0 or dg == 0:
  831. return [gcd], f, g
  832. f_norm = dup_max_norm(f, K)
  833. g_norm = dup_max_norm(g, K)
  834. B = K(2*min(f_norm, g_norm) + 29)
  835. x = max(min(B, 99*K.sqrt(B)),
  836. 2*min(f_norm // abs(dup_LC(f, K)),
  837. g_norm // abs(dup_LC(g, K))) + 4)
  838. for i in range(0, HEU_GCD_MAX):
  839. ff = dup_eval(f, x, K)
  840. gg = dup_eval(g, x, K)
  841. if ff and gg:
  842. h = K.gcd(ff, gg)
  843. cff = ff // h
  844. cfg = gg // h
  845. h = _dup_zz_gcd_interpolate(h, x, K)
  846. h = dup_primitive(h, K)[1]
  847. cff_, r = dup_div(f, h, K)
  848. if not r:
  849. cfg_, r = dup_div(g, h, K)
  850. if not r:
  851. h = dup_mul_ground(h, gcd, K)
  852. return h, cff_, cfg_
  853. cff = _dup_zz_gcd_interpolate(cff, x, K)
  854. h, r = dup_div(f, cff, K)
  855. if not r:
  856. cfg_, r = dup_div(g, h, K)
  857. if not r:
  858. h = dup_mul_ground(h, gcd, K)
  859. return h, cff, cfg_
  860. cfg = _dup_zz_gcd_interpolate(cfg, x, K)
  861. h, r = dup_div(g, cfg, K)
  862. if not r:
  863. cff_, r = dup_div(f, h, K)
  864. if not r:
  865. h = dup_mul_ground(h, gcd, K)
  866. return h, cff_, cfg
  867. x = 73794*x * K.sqrt(K.sqrt(x)) // 27011
  868. raise HeuristicGCDFailed('no luck')
  869. def _dmp_zz_gcd_interpolate(h, x, v, K):
  870. """Interpolate polynomial GCD from integer GCD. """
  871. f = []
  872. while not dmp_zero_p(h, v):
  873. g = dmp_ground_trunc(h, x, v, K)
  874. f.insert(0, g)
  875. h = dmp_sub(h, g, v, K)
  876. h = dmp_quo_ground(h, x, v, K)
  877. if K.is_negative(dmp_ground_LC(f, v + 1, K)):
  878. return dmp_neg(f, v + 1, K)
  879. else:
  880. return f
  881. def dmp_zz_heu_gcd(f, g, u, K):
  882. """
  883. Heuristic polynomial GCD in `Z[X]`.
  884. Given univariate polynomials `f` and `g` in `Z[X]`, returns
  885. their GCD and cofactors, i.e. polynomials ``h``, ``cff`` and ``cfg``
  886. such that::
  887. h = gcd(f, g), cff = quo(f, h) and cfg = quo(g, h)
  888. The algorithm is purely heuristic which means it may fail to compute
  889. the GCD. This will be signaled by raising an exception. In this case
  890. you will need to switch to another GCD method.
  891. The algorithm computes the polynomial GCD by evaluating polynomials
  892. f and g at certain points and computing (fast) integer GCD of those
  893. evaluations. The polynomial GCD is recovered from the integer image
  894. by interpolation. The evaluation process reduces f and g variable by
  895. variable into a large integer. The final step is to verify if the
  896. interpolated polynomial is the correct GCD. This gives cofactors of
  897. the input polynomials as a side effect.
  898. Examples
  899. ========
  900. >>> from sympy.polys import ring, ZZ
  901. >>> R, x,y, = ring("x,y", ZZ)
  902. >>> f = x**2 + 2*x*y + y**2
  903. >>> g = x**2 + x*y
  904. >>> R.dmp_zz_heu_gcd(f, g)
  905. (x + y, x + y, x)
  906. References
  907. ==========
  908. .. [1] [Liao95]_
  909. """
  910. if not u:
  911. return dup_zz_heu_gcd(f, g, K)
  912. result = _dmp_rr_trivial_gcd(f, g, u, K)
  913. if result is not None:
  914. return result
  915. gcd, f, g = dmp_ground_extract(f, g, u, K)
  916. f_norm = dmp_max_norm(f, u, K)
  917. g_norm = dmp_max_norm(g, u, K)
  918. B = K(2*min(f_norm, g_norm) + 29)
  919. x = max(min(B, 99*K.sqrt(B)),
  920. 2*min(f_norm // abs(dmp_ground_LC(f, u, K)),
  921. g_norm // abs(dmp_ground_LC(g, u, K))) + 4)
  922. for i in range(0, HEU_GCD_MAX):
  923. ff = dmp_eval(f, x, u, K)
  924. gg = dmp_eval(g, x, u, K)
  925. v = u - 1
  926. if not (dmp_zero_p(ff, v) or dmp_zero_p(gg, v)):
  927. h, cff, cfg = dmp_zz_heu_gcd(ff, gg, v, K)
  928. h = _dmp_zz_gcd_interpolate(h, x, v, K)
  929. h = dmp_ground_primitive(h, u, K)[1]
  930. cff_, r = dmp_div(f, h, u, K)
  931. if dmp_zero_p(r, u):
  932. cfg_, r = dmp_div(g, h, u, K)
  933. if dmp_zero_p(r, u):
  934. h = dmp_mul_ground(h, gcd, u, K)
  935. return h, cff_, cfg_
  936. cff = _dmp_zz_gcd_interpolate(cff, x, v, K)
  937. h, r = dmp_div(f, cff, u, K)
  938. if dmp_zero_p(r, u):
  939. cfg_, r = dmp_div(g, h, u, K)
  940. if dmp_zero_p(r, u):
  941. h = dmp_mul_ground(h, gcd, u, K)
  942. return h, cff, cfg_
  943. cfg = _dmp_zz_gcd_interpolate(cfg, x, v, K)
  944. h, r = dmp_div(g, cfg, u, K)
  945. if dmp_zero_p(r, u):
  946. cff_, r = dmp_div(f, h, u, K)
  947. if dmp_zero_p(r, u):
  948. h = dmp_mul_ground(h, gcd, u, K)
  949. return h, cff_, cfg
  950. x = 73794*x * K.sqrt(K.sqrt(x)) // 27011
  951. raise HeuristicGCDFailed('no luck')
  952. def dup_qq_heu_gcd(f, g, K0):
  953. """
  954. Heuristic polynomial GCD in `Q[x]`.
  955. Returns ``(h, cff, cfg)`` such that ``a = gcd(f, g)``,
  956. ``cff = quo(f, h)``, and ``cfg = quo(g, h)``.
  957. Examples
  958. ========
  959. >>> from sympy.polys import ring, QQ
  960. >>> R, x = ring("x", QQ)
  961. >>> f = QQ(1,2)*x**2 + QQ(7,4)*x + QQ(3,2)
  962. >>> g = QQ(1,2)*x**2 + x
  963. >>> R.dup_qq_heu_gcd(f, g)
  964. (x + 2, 1/2*x + 3/4, 1/2*x)
  965. """
  966. result = _dup_ff_trivial_gcd(f, g, K0)
  967. if result is not None:
  968. return result
  969. K1 = K0.get_ring()
  970. cf, f = dup_clear_denoms(f, K0, K1)
  971. cg, g = dup_clear_denoms(g, K0, K1)
  972. f = dup_convert(f, K0, K1)
  973. g = dup_convert(g, K0, K1)
  974. h, cff, cfg = dup_zz_heu_gcd(f, g, K1)
  975. h = dup_convert(h, K1, K0)
  976. c = dup_LC(h, K0)
  977. h = dup_monic(h, K0)
  978. cff = dup_convert(cff, K1, K0)
  979. cfg = dup_convert(cfg, K1, K0)
  980. cff = dup_mul_ground(cff, K0.quo(c, cf), K0)
  981. cfg = dup_mul_ground(cfg, K0.quo(c, cg), K0)
  982. return h, cff, cfg
  983. def dmp_qq_heu_gcd(f, g, u, K0):
  984. """
  985. Heuristic polynomial GCD in `Q[X]`.
  986. Returns ``(h, cff, cfg)`` such that ``a = gcd(f, g)``,
  987. ``cff = quo(f, h)``, and ``cfg = quo(g, h)``.
  988. Examples
  989. ========
  990. >>> from sympy.polys import ring, QQ
  991. >>> R, x,y, = ring("x,y", QQ)
  992. >>> f = QQ(1,4)*x**2 + x*y + y**2
  993. >>> g = QQ(1,2)*x**2 + x*y
  994. >>> R.dmp_qq_heu_gcd(f, g)
  995. (x + 2*y, 1/4*x + 1/2*y, 1/2*x)
  996. """
  997. result = _dmp_ff_trivial_gcd(f, g, u, K0)
  998. if result is not None:
  999. return result
  1000. K1 = K0.get_ring()
  1001. cf, f = dmp_clear_denoms(f, u, K0, K1)
  1002. cg, g = dmp_clear_denoms(g, u, K0, K1)
  1003. f = dmp_convert(f, u, K0, K1)
  1004. g = dmp_convert(g, u, K0, K1)
  1005. h, cff, cfg = dmp_zz_heu_gcd(f, g, u, K1)
  1006. h = dmp_convert(h, u, K1, K0)
  1007. c = dmp_ground_LC(h, u, K0)
  1008. h = dmp_ground_monic(h, u, K0)
  1009. cff = dmp_convert(cff, u, K1, K0)
  1010. cfg = dmp_convert(cfg, u, K1, K0)
  1011. cff = dmp_mul_ground(cff, K0.quo(c, cf), u, K0)
  1012. cfg = dmp_mul_ground(cfg, K0.quo(c, cg), u, K0)
  1013. return h, cff, cfg
  1014. def dup_inner_gcd(f, g, K):
  1015. """
  1016. Computes polynomial GCD and cofactors of `f` and `g` in `K[x]`.
  1017. Returns ``(h, cff, cfg)`` such that ``a = gcd(f, g)``,
  1018. ``cff = quo(f, h)``, and ``cfg = quo(g, h)``.
  1019. Examples
  1020. ========
  1021. >>> from sympy.polys import ring, ZZ
  1022. >>> R, x = ring("x", ZZ)
  1023. >>> R.dup_inner_gcd(x**2 - 1, x**2 - 3*x + 2)
  1024. (x - 1, x + 1, x - 2)
  1025. """
  1026. # XXX: This used to check for K.is_Exact but leads to awkward results when
  1027. # the domain is something like RR[z] e.g.:
  1028. #
  1029. # >>> g, p, q = Poly(1, x).cancel(Poly(51.05*x*y - 1.0, x))
  1030. # >>> g
  1031. # 1.0
  1032. # >>> p
  1033. # Poly(17592186044421.0, x, domain='RR[y]')
  1034. # >>> q
  1035. # Poly(898081097567692.0*y*x - 17592186044421.0, x, domain='RR[y]'))
  1036. #
  1037. # Maybe it would be better to flatten into multivariate polynomials first.
  1038. if K.is_RR or K.is_CC:
  1039. try:
  1040. exact = K.get_exact()
  1041. except DomainError:
  1042. return [K.one], f, g
  1043. f = dup_convert(f, K, exact)
  1044. g = dup_convert(g, K, exact)
  1045. h, cff, cfg = dup_inner_gcd(f, g, exact)
  1046. h = dup_convert(h, exact, K)
  1047. cff = dup_convert(cff, exact, K)
  1048. cfg = dup_convert(cfg, exact, K)
  1049. return h, cff, cfg
  1050. elif K.is_Field:
  1051. if K.is_QQ and query('USE_HEU_GCD'):
  1052. try:
  1053. return dup_qq_heu_gcd(f, g, K)
  1054. except HeuristicGCDFailed:
  1055. pass
  1056. return dup_ff_prs_gcd(f, g, K)
  1057. else:
  1058. if K.is_ZZ and query('USE_HEU_GCD'):
  1059. try:
  1060. return dup_zz_heu_gcd(f, g, K)
  1061. except HeuristicGCDFailed:
  1062. pass
  1063. return dup_rr_prs_gcd(f, g, K)
  1064. def _dmp_inner_gcd(f, g, u, K):
  1065. """Helper function for `dmp_inner_gcd()`. """
  1066. if not K.is_Exact:
  1067. try:
  1068. exact = K.get_exact()
  1069. except DomainError:
  1070. return dmp_one(u, K), f, g
  1071. f = dmp_convert(f, u, K, exact)
  1072. g = dmp_convert(g, u, K, exact)
  1073. h, cff, cfg = _dmp_inner_gcd(f, g, u, exact)
  1074. h = dmp_convert(h, u, exact, K)
  1075. cff = dmp_convert(cff, u, exact, K)
  1076. cfg = dmp_convert(cfg, u, exact, K)
  1077. return h, cff, cfg
  1078. elif K.is_Field:
  1079. if K.is_QQ and query('USE_HEU_GCD'):
  1080. try:
  1081. return dmp_qq_heu_gcd(f, g, u, K)
  1082. except HeuristicGCDFailed:
  1083. pass
  1084. return dmp_ff_prs_gcd(f, g, u, K)
  1085. else:
  1086. if K.is_ZZ and query('USE_HEU_GCD'):
  1087. try:
  1088. return dmp_zz_heu_gcd(f, g, u, K)
  1089. except HeuristicGCDFailed:
  1090. pass
  1091. return dmp_rr_prs_gcd(f, g, u, K)
  1092. def dmp_inner_gcd(f, g, u, K):
  1093. """
  1094. Computes polynomial GCD and cofactors of `f` and `g` in `K[X]`.
  1095. Returns ``(h, cff, cfg)`` such that ``a = gcd(f, g)``,
  1096. ``cff = quo(f, h)``, and ``cfg = quo(g, h)``.
  1097. Examples
  1098. ========
  1099. >>> from sympy.polys import ring, ZZ
  1100. >>> R, x,y, = ring("x,y", ZZ)
  1101. >>> f = x**2 + 2*x*y + y**2
  1102. >>> g = x**2 + x*y
  1103. >>> R.dmp_inner_gcd(f, g)
  1104. (x + y, x + y, x)
  1105. """
  1106. if not u:
  1107. return dup_inner_gcd(f, g, K)
  1108. J, (f, g) = dmp_multi_deflate((f, g), u, K)
  1109. h, cff, cfg = _dmp_inner_gcd(f, g, u, K)
  1110. return (dmp_inflate(h, J, u, K),
  1111. dmp_inflate(cff, J, u, K),
  1112. dmp_inflate(cfg, J, u, K))
  1113. def dup_gcd(f, g, K):
  1114. """
  1115. Computes polynomial GCD of `f` and `g` in `K[x]`.
  1116. Examples
  1117. ========
  1118. >>> from sympy.polys import ring, ZZ
  1119. >>> R, x = ring("x", ZZ)
  1120. >>> R.dup_gcd(x**2 - 1, x**2 - 3*x + 2)
  1121. x - 1
  1122. """
  1123. return dup_inner_gcd(f, g, K)[0]
  1124. def dmp_gcd(f, g, u, K):
  1125. """
  1126. Computes polynomial GCD of `f` and `g` in `K[X]`.
  1127. Examples
  1128. ========
  1129. >>> from sympy.polys import ring, ZZ
  1130. >>> R, x,y, = ring("x,y", ZZ)
  1131. >>> f = x**2 + 2*x*y + y**2
  1132. >>> g = x**2 + x*y
  1133. >>> R.dmp_gcd(f, g)
  1134. x + y
  1135. """
  1136. return dmp_inner_gcd(f, g, u, K)[0]
  1137. def dup_rr_lcm(f, g, K):
  1138. """
  1139. Computes polynomial LCM over a ring in `K[x]`.
  1140. Examples
  1141. ========
  1142. >>> from sympy.polys import ring, ZZ
  1143. >>> R, x = ring("x", ZZ)
  1144. >>> R.dup_rr_lcm(x**2 - 1, x**2 - 3*x + 2)
  1145. x**3 - 2*x**2 - x + 2
  1146. """
  1147. if not f or not g:
  1148. return dmp_zero(0)
  1149. fc, f = dup_primitive(f, K)
  1150. gc, g = dup_primitive(g, K)
  1151. c = K.lcm(fc, gc)
  1152. h = dup_quo(dup_mul(f, g, K),
  1153. dup_gcd(f, g, K), K)
  1154. u = K.canonical_unit(dup_LC(h, K))
  1155. return dup_mul_ground(h, c*u, K)
  1156. def dup_ff_lcm(f, g, K):
  1157. """
  1158. Computes polynomial LCM over a field in `K[x]`.
  1159. Examples
  1160. ========
  1161. >>> from sympy.polys import ring, QQ
  1162. >>> R, x = ring("x", QQ)
  1163. >>> f = QQ(1,2)*x**2 + QQ(7,4)*x + QQ(3,2)
  1164. >>> g = QQ(1,2)*x**2 + x
  1165. >>> R.dup_ff_lcm(f, g)
  1166. x**3 + 7/2*x**2 + 3*x
  1167. """
  1168. h = dup_quo(dup_mul(f, g, K),
  1169. dup_gcd(f, g, K), K)
  1170. return dup_monic(h, K)
  1171. def dup_lcm(f, g, K):
  1172. """
  1173. Computes polynomial LCM of `f` and `g` in `K[x]`.
  1174. Examples
  1175. ========
  1176. >>> from sympy.polys import ring, ZZ
  1177. >>> R, x = ring("x", ZZ)
  1178. >>> R.dup_lcm(x**2 - 1, x**2 - 3*x + 2)
  1179. x**3 - 2*x**2 - x + 2
  1180. """
  1181. if K.is_Field:
  1182. return dup_ff_lcm(f, g, K)
  1183. else:
  1184. return dup_rr_lcm(f, g, K)
  1185. def dmp_rr_lcm(f, g, u, K):
  1186. """
  1187. Computes polynomial LCM over a ring in `K[X]`.
  1188. Examples
  1189. ========
  1190. >>> from sympy.polys import ring, ZZ
  1191. >>> R, x,y, = ring("x,y", ZZ)
  1192. >>> f = x**2 + 2*x*y + y**2
  1193. >>> g = x**2 + x*y
  1194. >>> R.dmp_rr_lcm(f, g)
  1195. x**3 + 2*x**2*y + x*y**2
  1196. """
  1197. fc, f = dmp_ground_primitive(f, u, K)
  1198. gc, g = dmp_ground_primitive(g, u, K)
  1199. c = K.lcm(fc, gc)
  1200. h = dmp_quo(dmp_mul(f, g, u, K),
  1201. dmp_gcd(f, g, u, K), u, K)
  1202. return dmp_mul_ground(h, c, u, K)
  1203. def dmp_ff_lcm(f, g, u, K):
  1204. """
  1205. Computes polynomial LCM over a field in `K[X]`.
  1206. Examples
  1207. ========
  1208. >>> from sympy.polys import ring, QQ
  1209. >>> R, x,y, = ring("x,y", QQ)
  1210. >>> f = QQ(1,4)*x**2 + x*y + y**2
  1211. >>> g = QQ(1,2)*x**2 + x*y
  1212. >>> R.dmp_ff_lcm(f, g)
  1213. x**3 + 4*x**2*y + 4*x*y**2
  1214. """
  1215. h = dmp_quo(dmp_mul(f, g, u, K),
  1216. dmp_gcd(f, g, u, K), u, K)
  1217. return dmp_ground_monic(h, u, K)
  1218. def dmp_lcm(f, g, u, K):
  1219. """
  1220. Computes polynomial LCM of `f` and `g` in `K[X]`.
  1221. Examples
  1222. ========
  1223. >>> from sympy.polys import ring, ZZ
  1224. >>> R, x,y, = ring("x,y", ZZ)
  1225. >>> f = x**2 + 2*x*y + y**2
  1226. >>> g = x**2 + x*y
  1227. >>> R.dmp_lcm(f, g)
  1228. x**3 + 2*x**2*y + x*y**2
  1229. """
  1230. if not u:
  1231. return dup_lcm(f, g, K)
  1232. if K.is_Field:
  1233. return dmp_ff_lcm(f, g, u, K)
  1234. else:
  1235. return dmp_rr_lcm(f, g, u, K)
  1236. def dmp_content(f, u, K):
  1237. """
  1238. Returns GCD of multivariate coefficients.
  1239. Examples
  1240. ========
  1241. >>> from sympy.polys import ring, ZZ
  1242. >>> R, x,y, = ring("x,y", ZZ)
  1243. >>> R.dmp_content(2*x*y + 6*x + 4*y + 12)
  1244. 2*y + 6
  1245. """
  1246. cont, v = dmp_LC(f, K), u - 1
  1247. if dmp_zero_p(f, u):
  1248. return cont
  1249. for c in f[1:]:
  1250. cont = dmp_gcd(cont, c, v, K)
  1251. if dmp_one_p(cont, v, K):
  1252. break
  1253. if K.is_negative(dmp_ground_LC(cont, v, K)):
  1254. return dmp_neg(cont, v, K)
  1255. else:
  1256. return cont
  1257. def dmp_primitive(f, u, K):
  1258. """
  1259. Returns multivariate content and a primitive polynomial.
  1260. Examples
  1261. ========
  1262. >>> from sympy.polys import ring, ZZ
  1263. >>> R, x,y, = ring("x,y", ZZ)
  1264. >>> R.dmp_primitive(2*x*y + 6*x + 4*y + 12)
  1265. (2*y + 6, x + 2)
  1266. """
  1267. cont, v = dmp_content(f, u, K), u - 1
  1268. if dmp_zero_p(f, u) or dmp_one_p(cont, v, K):
  1269. return cont, f
  1270. else:
  1271. return cont, [ dmp_quo(c, cont, v, K) for c in f ]
  1272. def dup_cancel(f, g, K, include=True):
  1273. """
  1274. Cancel common factors in a rational function `f/g`.
  1275. Examples
  1276. ========
  1277. >>> from sympy.polys import ring, ZZ
  1278. >>> R, x = ring("x", ZZ)
  1279. >>> R.dup_cancel(2*x**2 - 2, x**2 - 2*x + 1)
  1280. (2*x + 2, x - 1)
  1281. """
  1282. return dmp_cancel(f, g, 0, K, include=include)
  1283. def dmp_cancel(f, g, u, K, include=True):
  1284. """
  1285. Cancel common factors in a rational function `f/g`.
  1286. Examples
  1287. ========
  1288. >>> from sympy.polys import ring, ZZ
  1289. >>> R, x,y = ring("x,y", ZZ)
  1290. >>> R.dmp_cancel(2*x**2 - 2, x**2 - 2*x + 1)
  1291. (2*x + 2, x - 1)
  1292. """
  1293. K0 = None
  1294. if K.is_Field and K.has_assoc_Ring:
  1295. K0, K = K, K.get_ring()
  1296. cq, f = dmp_clear_denoms(f, u, K0, K, convert=True)
  1297. cp, g = dmp_clear_denoms(g, u, K0, K, convert=True)
  1298. else:
  1299. cp, cq = K.one, K.one
  1300. _, p, q = dmp_inner_gcd(f, g, u, K)
  1301. if K0 is not None:
  1302. _, cp, cq = K.cofactors(cp, cq)
  1303. p = dmp_convert(p, u, K, K0)
  1304. q = dmp_convert(q, u, K, K0)
  1305. K = K0
  1306. p_neg = K.is_negative(dmp_ground_LC(p, u, K))
  1307. q_neg = K.is_negative(dmp_ground_LC(q, u, K))
  1308. if p_neg and q_neg:
  1309. p, q = dmp_neg(p, u, K), dmp_neg(q, u, K)
  1310. elif p_neg:
  1311. cp, p = -cp, dmp_neg(p, u, K)
  1312. elif q_neg:
  1313. cp, q = -cp, dmp_neg(q, u, K)
  1314. if not include:
  1315. return cp, cq, p, q
  1316. p = dmp_mul_ground(p, cp, u, K)
  1317. q = dmp_mul_ground(q, cq, u, K)
  1318. return p, q