def numerical_points(self, F=None, **kwds): """ Return some or all numerical approximations of rational points of a projective scheme. This is for dimension 0 subschemes only and the points are determined through a groebner calculation over the base ring and then numerically approximating the roots of the resulting polynomials. If the base ring is a number field, the embedding into ``F`` must be known. INPUT: ``F`` - numerical ring kwds: - ``point_tolerance`` - positive real number (optional, default=10^(-10)). For numerically inexact fields, two points are considered the same if their coordinates are within tolerance. - ``zero_tolerance`` - positive real number (optional, default=10^(-10)). For numerically inexact fields, points are on the subscheme if they satisfy the equations to within tolerance. OUTPUT: A list of points in the ambient space. .. WARNING:: For numerically inexact fields the list of points returned may contain repeated or be missing points due to tolerance. EXAMPLES:: sage: P.<x,y,z> = ProjectiveSpace(QQ, 2) sage: E = P.subscheme([y^3 - x^3 - x*z^2, x*y*z]) sage: L = E(QQ).numerical_points(F=RR); L [(0.000000000000000 : 0.000000000000000 : 1.00000000000000), (1.00000000000000 : 1.00000000000000 : 0.000000000000000)] sage: L[0].codomain() Projective Space of dimension 2 over Real Field with 53 bits of precision :: sage: S.<a> = QQ[] sage: K.<v> = NumberField(a^5 - 7, embedding=CC((7)**(1/5))) sage: P.<x,y,z> = ProjectiveSpace(K,2) sage: X = P.subscheme([x^2 - v^2*z^2, y-v*z]) sage: len(X(K).numerical_points(F=CDF)) 2 :: sage: P.<x1, x2, x3> = ProjectiveSpace(QQ, 2) sage: E = P.subscheme([3000*x1^50 + 9875643*x2^2*x3^48 + 12334545*x2^50, x1 + x2]) sage: len(E(P.base_ring()).numerical_points(F=CDF, zero_tolerance=1e-6)) 49 TESTS:: sage: P.<x,y,z> = ProjectiveSpace(QQ, 2) sage: E = P.subscheme([y^3 - x^3 - x*z^2, x*y*z]) sage: E(QQ).numerical_points(F=CDF, point_tolerance=-1) Traceback (most recent call last): ... ValueError: tolerance must be positive :: sage: P.<x,y,z> = ProjectiveSpace(QQ, 2) sage: E = P.subscheme([y^3 - x^3 - x*z^2, x*y*z]) sage: E(QQ).numerical_points(F=CC, zero_tolerance=-1) Traceback (most recent call last): ... ValueError: tolerance must be positive :: sage: P.<x,y,z> = ProjectiveSpace(QQ, 2) sage: E = P.subscheme([y^3 - x^3 - x*z^2, x*y*z]) sage: E(QQ).numerical_points(F=QQbar) Traceback (most recent call last): ... TypeError: F must be a numerical field """ from sage.schemes.projective.projective_space import is_ProjectiveSpace if F is None: F = CC if F not in Fields() or not hasattr(F, 'precision'): raise TypeError('F must be a numerical field') X = self.codomain() if X.base_ring() not in NumberFields(): raise TypeError('base ring must be a number field') PP = X.ambient_space().change_ring(F) if not is_ProjectiveSpace(X) and X.base_ring() in Fields(): #Then it must be a subscheme dim_ideal = X.defining_ideal().dimension() if dim_ideal < 1: # no points return [] if dim_ideal == 1: # if X zero-dimensional pt_tol = RR(kwds.pop('point_tolerance', 10**(-10))) zero_tol = RR(kwds.pop('zero_tolerance', 10**(-10))) if pt_tol <= 0 or zero_tol <= 0: raise ValueError("tolerance must be positive") rat_points = set() PS = X.ambient_space() N = PS.dimension_relative() BR = X.base_ring() #need a lexicographic ordering for elimination R = PolynomialRing(BR, N + 1, PS.variable_names(), order='lex') RF = R.change_ring(F) I = R.ideal(X.defining_polynomials()) #Determine the points through elimination #This is much faster than using the I.variety() function on each affine chart. for k in range(N + 1): #create the elimination ideal for the kth affine patch G = I.substitute({R.gen(k): 1}).groebner_basis() G = [RF(g) for g in G] if G != [1]: P = {} #keep track that we know the kth coordinate is 1 P.update({RF.gen(k): 1}) points = [P] #work backwards from solving each equation for the possible #values of the next coordinate for i in range(len(G) - 1, -1, -1): new_points = [] good = 0 for P in points: #substitute in our dictionary entry that has the values #of coordinates known so far. This results in a single #variable polynomial (by elimination) L = G[i].substitute(P) if len(RF(L).variables()) == 1: for pol in L.univariate_polynomial().roots( ring=F, multiplicities=False): r = L.variables()[0] varindex = RF.gens().index(r) P.update({RF.gen(varindex): pol}) new_points.append(copy(P)) good = 1 else: new_points.append(P) good = 1 if good: points = new_points #the dictionary entries now have values for all coordinates #they are approximate solutions to the equations #make them into projective points polys = [ g.change_ring(F) for g in X.defining_polynomials() ] for i in range(len(points)): if len(points[i]) == N + 1: S = PP([ points[i][RF.gen(j)] for j in range(N + 1) ]) S.normalize_coordinates() if all(g(list(S)) < zero_tol for g in polys): rat_points.add(S) # remove duplicate element using tolerance #since they are normalized we can just compare coefficients dupl_points = list(rat_points) for i in range(len(dupl_points)): u = dupl_points[i] for j in range(i + 1, len(dupl_points)): v = dupl_points[j] if all((u[k] - v[k]).abs() < pt_tol for k in range(len(u))): rat_points.remove(u) break rat_points = sorted(rat_points) return rat_points raise NotImplementedError( 'numerical approximation of points only for dimension 0 subschemes' )
def numerical_points(self, F=None, **kwds): """ Return some or all numerical approximations of rational points of an affine scheme. This is for dimension 0 subschemes only and the points are determined through a groebner calculation over the base ring and then numerically approximating the roots of the resulting polynomials. If the base ring is a number field, the embedding into ``F`` must be known. INPUT: ``F`` - numerical ring kwds: - ``zero_tolerance`` - positive real number (optional, default=10^(-10)). For numerically inexact fields, points are on the subscheme if they satisfy the equations to within tolerance. OUTPUT: A list of points in the ambient space. .. WARNING:: For numerically inexact fields the list of points returned may contain repeated or be missing points due to tolerance. EXAMPLES:: sage: K.<v> = QuadraticField(3) sage: A.<x,y> = AffineSpace(K, 2) sage: X = A.subscheme([x^3 - v^2*y, y - v*x^2 + 3]) sage: L = X(K).numerical_points(F=RR); L # abs tol 1e-14 [(-1.18738247880014, -0.558021142104134), (1.57693558184861, 1.30713548084184), (4.80659931965815, 37.0162574656220)] sage: L[0].codomain() Affine Space of dimension 2 over Real Field with 53 bits of precision :: sage: A.<x,y> = AffineSpace(QQ, 2) sage: X = A.subscheme([y^2 - x^2 - 3*x, x^2 - 10*y]) sage: len(X(QQ).numerical_points(F=ComplexField(100))) 4 :: sage: A.<x1, x2> = AffineSpace(QQ, 2) sage: E = A.subscheme([30*x1^100 + 1000*x2^2 + 2000*x1*x2 + 1, x1 + x2]) sage: len(E(A.base_ring()).numerical_points(F=CDF, zero_tolerance=1e-9)) 100 TESTS:: sage: A.<x,y> = AffineSpace(QQ, 2) sage: X = A.subscheme([y^2 - x^2 - 3*x, x^2 - 10*y]) sage: X(QQ).numerical_points(F=QQ) Traceback (most recent call last): ... TypeError: F must be a numerical field :: sage: A.<x,y> = AffineSpace(QQ, 2) sage: X = A.subscheme([y^2 - x^2 - 3*x, x^2 - 10*y]) sage: X(QQ).numerical_points(F=CC, zero_tolerance=-1) Traceback (most recent call last): ... ValueError: tolerance must be positive """ from sage.schemes.affine.affine_space import is_AffineSpace if F is None: F = CC if not F in Fields() or not hasattr(F, 'precision'): raise TypeError('F must be a numerical field') X = self.codomain() if X.base_ring() not in NumberFields(): raise TypeError('base ring must be a number field') AA = X.ambient_space().change_ring(F) if not is_AffineSpace(X) and X.base_ring() in Fields(): # Then X must be a subscheme dim_ideal = X.defining_ideal().dimension() if dim_ideal != 0: # no points return [] else: return [] # if X zero-dimensional zero_tol = RR(kwds.pop('zero_tolerance', 10**(-10))) if zero_tol <= 0: raise ValueError("tolerance must be positive") rat_points = [] PS = X.ambient_space() N = PS.dimension_relative() BR = X.base_ring() # need a lexicographic ordering for elimination R = PolynomialRing(BR, N, PS.gens(), order='lex') RF = R.change_ring(F) I = R.ideal(X.defining_polynomials()) I0 = R.ideal(0) # Determine the points through elimination This is much faster # than using the I.variety() function on each affine chart. G = I.groebner_basis() G = [RF(g) for g in G] if G != [1]: P = {} points = [P] # work backwards from solving each equation for the possible # values of the next coordinate for g in reversed(G): new_points = [] good = False for P in points: # substitute in our dictionary entry that has the # values of coordinates known so far. This results # in a single variable polynomial (by elimination) L = g.substitute(P) if len(RF(L).variables()) == 1: r = L.variables()[0] var = RF.gen(RF.gens().index(r)) for pol in L.univariate_polynomial().roots(ring=F, multiplicities=False): P[var] = pol new_points.append(copy(P)) good = True else: new_points.append(P) good = True if good: points = new_points # the dictionary entries now have values for all # coordinates they are the rational solutions to the # equations make them into affine points polys = [g.change_ring(F) for g in X.defining_polynomials()] for P in points: if len(P) == N: S = AA([P[R.gen(j)] for j in range(N)]) if all([g(list(S)) < zero_tol for g in polys]): rat_points.append(S) rat_points = sorted(rat_points) return rat_points