From 3a649e4b39172dc29ad79a54fa8f1f9284a0cd21 Mon Sep 17 00:00:00 2001 From: =?UTF-8?q?Fr=C3=A9d=C3=A9ric=20Chapoton?= Date: Sat, 13 Sep 2025 20:16:40 +0200 Subject: [PATCH] add Newton polytopes for Laurent polynomials --- .../rings/polynomial/laurent_polynomial.pyx | 24 +++++++++++++++++++ .../polynomial/laurent_polynomial_mpair.pyx | 23 ++++++++++++++++++ 2 files changed, 47 insertions(+) diff --git a/src/sage/rings/polynomial/laurent_polynomial.pyx b/src/sage/rings/polynomial/laurent_polynomial.pyx index 81c110ebd94..d17cb3cc630 100644 --- a/src/sage/rings/polynomial/laurent_polynomial.pyx +++ b/src/sage/rings/polynomial/laurent_polynomial.pyx @@ -913,6 +913,30 @@ cdef class LaurentPolynomial_univariate(LaurentPolynomial): """ return [i + self.__n for i in self.__u.exponents()] + def newton_polytope(self): + r""" + Return the Newton polytope of this Laurent polynomial. + + EXAMPLES:: + + sage: R. = LaurentPolynomialRing(QQ) + sage: f = 1 + x + 33 * x^-3 + sage: P = f.newton_polytope(); P # needs sage.geometry.polyhedron + A 1-dimensional polyhedron in ZZ^1 defined as the convex hull of 2 vertices + + TESTS:: + + sage: R. = LaurentPolynomialRing(QQ) + sage: R(0).newton_polytope() # needs sage.geometry.polyhedron + The empty polyhedron in ZZ^0 + sage: R(1).newton_polytope() # needs sage.geometry.polyhedron + A 0-dimensional polyhedron in ZZ^1 defined as the convex hull of 1 vertex + """ + from sage.geometry.polyhedron.constructor import Polyhedron + from sage.rings.integer_ring import ZZ + return Polyhedron(vertices=[(e,) for e in self.exponents()], + base_ring=ZZ) + def __setitem__(self, n, value): """ EXAMPLES:: diff --git a/src/sage/rings/polynomial/laurent_polynomial_mpair.pyx b/src/sage/rings/polynomial/laurent_polynomial_mpair.pyx index 67446b49e25..6aa37df3611 100644 --- a/src/sage/rings/polynomial/laurent_polynomial_mpair.pyx +++ b/src/sage/rings/polynomial/laurent_polynomial_mpair.pyx @@ -1325,6 +1325,29 @@ cdef class LaurentPolynomial_mpair(LaurentPolynomial): # Find the minimal valuation of x by checking each term return Integer(min(e[i] for e in self.exponents())) + def newton_polytope(self): + r""" + Return the Newton polytope of this Laurent polynomial. + + EXAMPLES:: + + sage: R. = LaurentPolynomialRing(QQ) + sage: f = 1 + x*y + y**2 + 33 * x^-3 + sage: P = f.newton_polytope(); P # needs sage.geometry.polyhedron + A 2-dimensional polyhedron in ZZ^2 defined as the convex hull of 4 vertices + + TESTS:: + + sage: R. = LaurentPolynomialRing(QQ) + sage: R(0).newton_polytope() # needs sage.geometry.polyhedron + The empty polyhedron in ZZ^0 + sage: R(1).newton_polytope() # needs sage.geometry.polyhedron + A 0-dimensional polyhedron in ZZ^2 defined as the convex hull of 1 vertex + """ + from sage.geometry.polyhedron.constructor import Polyhedron + from sage.rings.integer_ring import ZZ + return Polyhedron(vertices=self.exponents(), base_ring=ZZ) + def has_inverse_of(self, i): """ INPUT: