def TelenVanBarel(initial_poly_list, accuracy = 1.e-10): """Uses Telen and VanBarels matrix reduction method to find a vector basis for the system of polynomials. Parameters -------- initial_poly_list: list The polynomials in the system we are solving. run_checks : bool If True, checks will be run to make sure TVB works and if it doesn't an S-polynomial will be searched for to fix it. accuracy: float How small we want a number to be before assuming it is zero. Returns ----------- basisDict : dict A dictionary of terms not in the vector basis a matrixes of things in the vector basis that the term can be reduced to. VB : numpy array The terms in the vector basis, each row being a term. degree : int The degree of the Macaualy matrix that was constructed. """ power = is_power(initial_poly_list) dim = initial_poly_list[0].dim poly_coeff_list = [] degree = find_degree(initial_poly_list) #This sorting is required for fast matrix construction. Ascending should be False. initial_poly_list = sort_polys_by_degree(initial_poly_list, ascending = False) """This is the first construction option, simple monomial multiplication.""" for poly in initial_poly_list: poly_coeff_list = add_polys(degree, poly, poly_coeff_list) """This is the second construction option, it uses the fancy triangle method that is faster but less stable.""" #for deg in reversed(range(min([poly.degree for poly in initial_poly_list]), degree+1)): # poly_coeff_list += deg_d_polys(initial_poly_list, deg, dim) #Creates the matrix for either of the above two methods. Comment out if using the third method. matrix, matrix_terms, cuts = create_matrix(poly_coeff_list, degree, dim) """This is the thrid matrix construction option, it uses the permutation arrays.""" #if power: # matrix, matrix_terms, cuts = createMatrixFast(initial_poly_list, degree, dim) #else: # matrix, matrix_terms, cuts = construction(initial_poly_list, degree, dim) matrix, matrix_terms = rrqr_reduceTelenVanBarel2(matrix, matrix_terms, cuts, accuracy = accuracy) #matrix, matrix_terms = rrqr_reduceTelenVanBarelFullRank(matrix, matrix_terms, cuts, accuracy = accuracy) height = matrix.shape[0] matrix[:,height:] = solve_triangular(matrix[:,:height],matrix[:,height:]) matrix[:,:height] = np.eye(height) VB = matrix_terms[height:] basisDict = makeBasisDict(matrix, matrix_terms, VB, power) return basisDict, VB, degree
def solve(polys, verbose=False, sim_diag="Telen"): ''' Finds the roots of a list of multivariate polynomials by simultaneously diagonalizing a multiplication matrices created with the TVB method. Parameters ---------- polys : list of polynomial objects Polynomials to find the common roots of. verbose : bool Print information about how the roots are computed. returns ------- roots : numpy array The common roots of the polynomials. Each row is a root. ''' polys = match_poly_dimensions(polys) poly_type = is_power(polys, return_string = True) dim = polys[0].dim #Reduce the Macaulay Matrix TVB-style to generate a basis for C[]/I and # a dictionary that expresses other monomials in terms of that basis basisDict, VB = TVB_MacaulayReduction(polys, accuracy = 1.e-10, verbose=verbose) if verbose: print('Basis for C[]/I\n', VB) print('Dictionary which represents non-basis monomials in terms of the basis\n', basisDict) #See what happened to dimension of C[]/I. Did we loose roots? len_VB = len(VB) degrees = [poly.degree for poly in polys] max_number_of_roots = np.prod(degrees) if len_VB < max_number_of_roots: raise MacaulayError('Roots were lost during the Macaulay Reduction') if len_VB > max_number_of_roots: raise MacaulayError('Dimension of C[]/I is too large') #make Mx1, ..., Mxn mult_matrices = np.zeros((dim, len_VB, len_VB)) for i in range(1, dim+1): mult_matrices[i-1] = Mxi_Matrix(i, basisDict, VB, dim, poly_type) if verbose: print('The Multiplication Matrices:\n', mult_matrices) #simultaneously diagonalize and return roots if sim_diag == "Telen": #Determine if the system contains only homogenous polynomials homogenous = all([is_homogenous(poly) for poly in polys]) #Find the roots roots = Telen_sim_diag(mult_matrices, verbose=verbose, homogenous=homogenous).T else: roots = sim_diag(mult_matrices, verbose=verbose).T if verbose: print("Roots:\n", roots) return roots
def newton_polish(polys, root, niter=100, tol=1e-5): """ Perform Newton's method on a system of N polynomials in M variables. Parameters ---------- polys : list A list of polynomial objects of the same type (MultiPower or MultiCheb). root : ndarray An initial guess for Newton's method, intended to be a candidate root from root_finder. niter : int A maximum number of iterations of Newton's method. tol : float Tolerance for convergence of Newton's method. Returns ------- x1 : ndarray The terminal point of Newton's method, an estimation for a root of the system """ poly_type = is_power(polys, return_string=True) m = len(polys) dim = max(poly.dim for poly in polys) f_x = np.empty(m, dtype="complex_") jac = np.empty((m, dim), dtype="complex_") def f(x): #f_x = np.empty(m,dtype="complex_") for i, poly in enumerate(polys): f_x[i] = poly(x) return f_x def Df(x): #jac = np.empty((m,dim),dtype="complex_") for i, poly in enumerate(polys): jac[i] = poly.grad(x) return jac i = 0 x0, x1 = root, root while True: if i == niter: break delta = np.linalg.solve(Df(x0), -f(x0)) norm = np.linalg.norm(delta) x1 = delta + x0 if norm < tol or norm > .1: break x0 = x1 i += 1 return x1
def new_macaulay(initial_poly_list, global_accuracy=1.e-10): """ Accepts a list of polynomials and use them to construct a Macaulay matrix. The matrix is constucted one degree at a time. parameters -------- initial_poly_list: list Polynomials for Macaulay construction. global_accuracy : float Round-off parameter: values within global_accuracy of zero are rounded to zero. Defaults to 1e-10. Returns ------- final_polys : list Reduced Macaulay matrix that can be passed into the root finder. """ Power = is_power(initial_poly_list) poly_coeff_list = list() dim = initial_poly_list[0].dim degree = max([i.degree for i in initial_poly_list]) total_degree = find_degree(initial_poly_list) for i in initial_poly_list: poly_coeff_list = add_polys(degree, i, poly_coeff_list) initial_poly_list = list() polys_old = create_reduce(poly_coeff_list, Power, False, 'matrix') poly_coeff_list = list() for i in polys_old: poly_coeff_list.append(i) polys_new = list() mons = mon_combos(np.zeros(dim, dtype=int), 1) mons = mons[1:] while degree < total_degree: poly_coeff_list, polys_old = add_polys_to_deg_x( degree, poly_coeff_list, polys_old, mons, Power) print("Reducing large matrix") poly_coeff_list = create_reduce(poly_coeff_list, Power, False, "matrix") degree += 1 print("Final reduction") final_polys = create_reduce(poly_coeff_list, Power, True) return final_polys
def Macaulay(initial_poly_list, global_accuracy=1.e-10): """ Accepts a list of polynomials and use them to construct a Macaulay matrix. parameters -------- initial_poly_list: list Polynomials for Macaulay construction. global_accuracy : float Round-off parameter: values within global_accuracy of zero are rounded to zero. Defaults to 1e-10. Returns ------- final_polys : list Reduced Macaulay matrix that can be passed into the root finder. """ Power = is_power(initial_poly_list) poly_coeff_list = [] degree = find_degree(initial_poly_list) for i in initial_poly_list: poly_coeff_list = add_polys(degree, i, poly_coeff_list) matrix, matrix_terms = create_matrix(poly_coeff_list) plt.matshow(matrix) plt.show() #rrqr_reduce2 and rrqr_reduce same pretty matched on stability, though I feel like 2 should be better. matrix = utils.rrqr_reduce2(matrix, global_accuracy=global_accuracy) # here matrix = clean_zeros_from_matrix(matrix) non_zero_rows = np.sum(np.abs(matrix), axis=1) != 0 matrix = matrix[non_zero_rows, :] #Only keeps the non_zero_polymonials matrix = triangular_solve(matrix) matrix = clean_zeros_from_matrix(matrix) #The other reduction option. I thought it would be really stable but seems to be the worst of the three. #matrix = matrixReduce(matrix, triangular_solve = True, global_accuracy = global_accuracy) rows = get_good_rows(matrix, matrix_terms) final_polys = get_polys_from_matrix(matrix, matrix_terms, rows, Power) return final_polys
def F4(polys, reducedGroebner=True, accuracy=1.e-10, phi=True): '''Uses the F4 algorithm to find a Groebner Basis. Parameters ---------- polys : list The polynomails for which a Groebner basis is computed. reducedGroebner : bool Defaults to True. If True then a reduced Groebner Basis is found. If False, just a Groebner Basis is found. accuracy : float Defaults to 1.e-10. What is determined to be zero in the code. Returns ------- groebner_basis : list The polynomials in the Groebner Basis. ''' power = is_power(polys) old_polys = list() new_polys = polys polys_were_added = True while len(new_polys) > 0: old_polys, new_polys, matrix_polys = sort_reducible_polys( old_polys, new_polys) matrix_polys = add_phi_to_matrix(old_polys, new_polys, matrix_polys, phi=phi) matrix_polys, matrix_terms = add_r_to_matrix(matrix_polys, old_polys + new_polys) matrix, matrix_terms = create_matrix(matrix_polys, matrix_terms) old_polys += new_polys new_polys = get_new_polys(matrix, matrix_terms, accuracy=accuracy, power=power) groebner_basis = old_polys if reducedGroebner: groebner_basis = reduce_groebner_basis(groebner_basis, power) return groebner_basis
def roots(polys, method='Groebner'): ''' Finds the roots of the given list of polynomials. Parameters ---------- polys : list of polynomial objects Polynomials to find the common roots of. method : string The root finding method to be used. Can be either 'Groebner', 'Macaulay', or 'TVB'. returns ------- list of numpy arrays the common roots of the polynomials ''' polys = match_poly_dimensions(polys) # Determine polynomial type poly_type = is_power(polys, return_string=True) if method == 'TVB' or method == 'new_TVB': try: m_f, var_dict = TVBMultMatrix(polys, poly_type, method) except TVBError as e: if str( e ) == "Doesn't have all x^n's on diagonal. Do linear transformation": raise e ''' #Optionally have it do F4 instead in thise case. warnings.warn("TVB method failed. Trying F4 instead. \ Error message from TVB is - {}".format(e), InstabilityWarning) method = 'Groebner' GB, m_f, var_dict = groebnerMultMatrix(polys, poly_type, method) ''' elif str(e) == 'Polys are non-zero dimensional': return -1 else: raise e else: GB, m_f, var_dict = groebnerMultMatrix(polys, poly_type, method) # both TVBMultMatrix and groebnerMultMatrix will return m_f as # -1 if the ideal is not zero dimensional or if there are no roots if type(m_f) == int: return -1 # Get list of indexes of single variables and store vars that were not # in the vector space basis. dim = max(f.dim for f in polys) var_list = get_var_list(dim) var_indexes = [-1] * dim vars_not_in_basis = {} for i in range(len(var_list)): var = var_list[i] # x_i if var in var_dict: var_indexes[i] = var_dict[var] else: # maps the position in the root to its variable vars_not_in_basis[i] = var vnib = False if len(vars_not_in_basis) != 0: if method == 'TVB': print("This isn't working yet...") return -1 vnib = True # Get left eigenvectors e = np.linalg.eig(m_f.T) eig = e[1] num_vectors = eig.shape[1] eig_vectors = [eig[:, i] for i in range(num_vectors)] # columns of eig roots = [] for v in eig_vectors: if v[var_dict[tuple(0 for i in range(dim))]] == 0: continue root = np.zeros(dim, dtype=complex) # This will always work because var_indexes and root have the # same length - dim - and var_indexes has the variables in the # order they should be in the root for i in range(dim): x_i_pos = var_indexes[i] if x_i_pos != -1: root[i] = v[x_i_pos] / v[var_dict[tuple(0 for i in range(dim))]] if vnib: # Go through the indexes of variables not in the basis in # decreasing order. It must be done in decreasing order for the # roots to be calculated correctly, since the vars with lower # indexes depend on the ones with higher indexes for pos in list(vars_not_in_basis.keys())[::-1]: GB_poly = _get_poly_with_LT(vars_not_in_basis[pos], GB) var_value = GB_poly(root) * -1 root[pos] = var_value roots.append(root) #roots.append(newton_polish(polys,root)) return roots
def solve(polys, MSmatrix=0, eigvals=True, verbose=False): ''' Finds the roots of the given list of polynomials. Parameters ---------- polys : list of polynomial objects Polynomials to find the common roots of. MSmatrix : int Controls which Moller-Stetter matrix is constructed For a univariate polynomial, the options are: 0 (default) -- The companion or colleague matrix, rotated 180 degrees 1 -- The unrotated companion or colleague matrix -1 -- The inverse of the companion or colleague matrix For a multivariate polynomial, the options are: 0 (default) -- The Moller-Stetter matrix of a random polynomial Some positive integer i <= dimension -- The Moller-Stetter matrix of x_i, where variables are index from x1, ..., xn Some negative integer i >= -dimension -- The Moller-Stetter matrix of x_i-inverse eigvals : bool Whether to compute roots of univariate polynomials from eigenvalues (True) or eigenvectors (False). Roots of multivariate polynomials are always comptued from eigenvectors verbose : bool Prints information about how the roots are computed. returns ------- roots : numpy array The common roots of the polynomials. Each row is a root. ''' polys = match_poly_dimensions(polys) # Determine polynomial type and dimension of the system poly_type = is_power(polys, return_string=True) dim = polys[0].dim if dim == 1: if len(polys) == 1: return oneD.solve(polys[0], MSmatrix=MSmatrix, eigvals=eigvals, verbose=verbose) else: zeros = np.unique( oneD.solve(polys[0], MSmatrix=MSmatrix, eigvals=eigvals, verbose=verbose)) #Finds the roots of each succesive polynomial and checks which roots are common. for poly in polys[1:]: if len(zeros) == 0: break zeros2 = np.unique( oneD.solve(poly, MSmatrix=MSmatrix, eigvals=eigvals, verbose=verbose)) common = list() tol = 1.e-10 for zero in zeros2: spot = np.where(np.abs(zeros - zero) < tol) if len(spot[0]) > 0: common.append(zero) zeros = common return zeros else: if MSmatrix < 0: return division(polys, verbose=verbose, divisor_var=-MSmatrix - 1) else: return multiplication(polys, verbose=verbose, MSmatrix=MSmatrix)
def TVB_MacaulayReduction(initial_poly_list, accuracy = 1.e-10, verbose=False): ''' Reduces the Macaulay matrix to find a vector basis for the system of polynomials. Parameters -------- polys: list The polynomials in the system we are solving. accuracy: float How small we want a number to be before assuming it is zero. Returns ----------- basisDict : dict A dictionary of terms not in the vector basis a matrixes of things in the vector basis that the term can be reduced to. VB : nparray The terms in the vector basis, each list being a term. ''' power = is_power(initial_poly_list) dim = initial_poly_list[0].dim poly_coeff_list = [] degree = find_degree(initial_poly_list) initial_poly_list = sort_polys_by_degree(initial_poly_list, ascending = False) for poly in initial_poly_list: poly_coeff_list = add_polys(degree, poly, poly_coeff_list) #Creates the matrix for either of the above two methods. Comment out if using the third method. matrix, matrix_terms, cuts = create_matrix(poly_coeff_list, degree, dim) if verbose: np.set_printoptions(suppress=False, linewidth=200) print('\nStarting Macaulay Matrix\n', matrix) print('\nColumns in Macaulay Matrix\nFirst element in tuple is degree of x monomial, Second element is degree of y monomial\n', matrix_terms) print('\nLocation of Cuts in the Macaulay Matrix into [ Mb | M1* | M2* ]\n', cuts) #First QR reduction #If bottom left is zero only does the first QR reduction on top part of matrix (for speed). Otherwise does it on the whole thing if np.allclose(matrix[cuts[0]:,:cuts[0]], 0): #RRQR reduces A and D without pivoting, sticking the result in it's place and multiplying the rest of the matrix by Q.T C1,matrix[:cuts[0],:cuts[0]] = qr_multiply(matrix[:,:cuts[0]], matrix[:,cuts[0]:].T, mode = 'right') matrix[:cuts[0],cuts[0]:] = C1.T C1 = 0 #check if there are zeros along the diagonal of R1 if any(np.isclose(np.diag(matrix[:,:cuts[0]]),0, rtol=accuracy)): raise MacaulayError("R1 IS NOT FULL RANK") #set small values to zero before backsolving matrix[np.isclose(matrix, 0, rtol=accuracy)] = 0 matrix[:cuts[0],cuts[0]:] = solve_triangular(matrix[:cuts[0],:cuts[0]],matrix[:cuts[0],cuts[0]:]) matrix[:cuts[0],:cuts[0]] = np.eye(cuts[0]) matrix[cuts[0]:,cuts[0]:] -= (matrix[cuts[0]:,:cuts[0]])@matrix[:cuts[0],cuts[0]:] #? else: #RRQR reduces A and D without pivoting, sticking the result in it's place. Q1,matrix[:,:cuts[0]] = qr(matrix[:,:cuts[0]]) #check if there are zeros along the diagonal of R1 if any(np.isclose(np.diag(matrix[:,:cuts[0]]),0, rtol=accuracy)): raise MacaulayError("R1 IS NOT FULL RANK") #Multiplying the rest of the matrix by Q.T matrix[:,cuts[0]:] = Q1.T@matrix[:,cuts[0]:] Q1 = 0 #Get rid of Q1 for memory purposes. #Second QR reduction on all the rest of the matrix matrix[cuts[0]:,cuts[0]:],P = qr(matrix[cuts[0]:,cuts[0]:], mode = 'r', pivoting = True) #Shift around the top right columns matrix[:cuts[0],cuts[0]:] = matrix[:cuts[0],cuts[0]:][:,P] #Shift around the columns labels matrix_terms[cuts[0]:] = matrix_terms[cuts[0]:][P] P = 0 #set small values to zero matrix[np.isclose(matrix, 0, rtol=accuracy)] = 0 #eliminate zero rows from the bottom of the matrix. Zero rows above #nonzero elements are not eliminated. This saves time since Macaulay matrices #we deal with are only zero at the very bottom #matrix = row_swap_matrix(matrix) #not for TVB's method needed bc QRP is on the whole matrix for row in matrix[::-1]: if np.allclose(row, 0): matrix = matrix[:-1] else: break height = matrix.shape[0] matrix[:,height:] = solve_triangular(matrix[:,:height],matrix[:,height:]) matrix[:,:height] = np.eye(height) #return np.vstack((matrix[:,height:].T,np.eye(height))), matrix_terms if verbose: np.set_printoptions(suppress=True, linewidth=200) print("\nFinal Macaulay Matrix\n", matrix) print("\nColumns in Macaulay Matrix\n", matrix_terms) VB = matrix_terms[height:] basisDict = makeBasisDict(matrix, matrix_terms, VB, power) return basisDict, VB
def MacaulayReduction(initial_poly_list, max_number_of_roots, accuracy=1.e-10, verbose=False): """Reduces the Macaulay matrix to find a vector basis for the system of polynomials. Parameters -------- initial_poly_list: list The polynomials in the system we are solving. accuracy: float How small we want a number to be before assuming it is zero. Returns ----------- basisDict : dict A dictionary of terms not in the vector basis a matrixes of things in the vector basis that the term can be reduced to. VB : numpy array The terms in the vector basis, each row being a term. """ power = is_power(initial_poly_list) dim = initial_poly_list[0].dim poly_coeff_list = [] degree = find_degree(initial_poly_list) #This sorting is required for fast matrix construction. Ascending should be False. initial_poly_list = sort_polys_by_degree(initial_poly_list, ascending=False) """This is the first construction option, simple monomial multiplication.""" for poly in initial_poly_list: poly_coeff_list = add_polys(degree, poly, poly_coeff_list) """This is the second construction option, it uses the fancy triangle method that is faster but less stable.""" #for deg in reversed(range(min([poly.degree for poly in initial_poly_list]), degree+1)): # poly_coeff_list += deg_d_polys(initial_poly_list, deg, dim) #Creates the matrix for either of the above two methods. Comment out if using the third method. #try: matrix, matrix_terms, cuts = create_matrix(poly_coeff_list, degree, dim) if verbose: np.set_printoptions(suppress=False, linewidth=200) print('\nStarting Macaulay Matrix\n', matrix) print( '\nColumns in Macaulay Matrix\nFirst element in tuple is degree of x, Second element is degree of y\n', matrix_terms) print( '\nLocation of Cuts in the Macaulay Matrix into [ Mb | M1* | M2* ]\n', cuts) """This is the thrid matrix construction option, it uses the permutation arrays.""" #if power: # matrix, matrix_terms, cuts = createMatrixFast(initial_poly_list, degree, dim) #else: # matrix, matrix_terms, cuts = construction(initial_poly_list, degree, dim) #If bottom left is zero only does the first QR reduction on top part of matrix (for speed). Otherwise does it on the whole thing if np.allclose(matrix[cuts[0]:, :cuts[0]], 0): matrix, matrix_terms = rrqr_reduceMacaulay2(matrix, matrix_terms, cuts, max_number_of_roots, accuracy=accuracy) else: matrix, matrix_terms = rrqr_reduceMacaulay(matrix, matrix_terms, cuts, max_number_of_roots, accuracy=accuracy) #Make there are enough rows in the reduced Macaulay matrix, i.e. didn't loose a row assert matrix.shape[0] >= matrix.shape[1] - max_number_of_roots #matrix, matrix_terms = rrqr_reduceMacaulayFullRank(matrix, matrix_terms, cuts, number_of_roots, accuracy = accuracy) height = matrix.shape[0] matrix[:, height:] = solve_triangular(matrix[:, :height], matrix[:, height:]) # except Exception as e: # if str(e)[:46] == 'singular matrix: resolution failed at diagonal': # matrix, matrix_terms, cuts = create_matrix(poly_coeff_list, degree+1, dim) # if verbose: # np.set_printoptions(suppress=False, linewidth=200) # print('\nNew, Bigger Macaulay Matrix\n', matrix) # print('\nColumns in Macaulay Matrix\nFirst element in tuple is degree of x, Second element is degree of y\n', matrix_terms) # print('\nLocation of Cuts in the Macaulay Matrix into [ Mb | M1* | M2* ]\n', cuts) # # """This is the thrid matrix construction option, it uses the permutation arrays.""" # #if power: # # matrix, matrix_terms, cuts = createMatrixFast(initial_poly_list, degree, dim) # #else: # # matrix, matrix_terms, cuts = construction(initial_poly_list, degree, dim) # # #If bottom left is zero only does the first QR reduction on top part of matrix (for speed). Otherwise does it on the whole thing # if np.allclose(matrix[cuts[0]:,:cuts[0]], 0): # matrix, matrix_terms = rrqr_reduceMacaulay2(matrix, matrix_terms, cuts, max_number_of_roots, accuracy = accuracy) # else: # matrix, matrix_terms = rrqr_reduceMacaulay(matrix, matrix_terms, cuts, max_number_of_roots, accuracy = accuracy) # # #Make there are enough rows in the reduced Macaulay matrix, i.e. didn't loose a row # assert matrix.shape[0] >= matrix.shape[1] - max_number_of_roots # # #matrix, matrix_terms = rrqr_reduceMacaulayFullRank(matrix, matrix_terms, cuts, number_of_roots, accuracy = accuracy) # height = matrix.shape[0] # matrix[:,height:] = solve_triangular(matrix[:,:height],matrix[:,height:]) # else: # raise e #Make there are enough rows in the reduced Macaulay matrix, i.e. didn't loose a row assert matrix.shape[0] >= matrix.shape[1] - max_number_of_roots matrix[:, :height] = np.eye(height) #return np.vstack((matrix[:,height:].T,np.eye(height))), matrix_terms if verbose: np.set_printoptions(suppress=True, linewidth=200) print("\nFinal Macaulay Matrix\n", matrix) print("\nColumns in Macaulay Matrix\n", matrix_terms) VB = matrix_terms[height:] #plt.plot(matrix_terms[:,0],matrix_terms[:,1],'kx') #plt.plot(VB[:,0],VB[:,1],'r.') basisDict = makeBasisDict(matrix, matrix_terms, VB, power) return basisDict, VB
def multiplication(polys, verbose=False, MSmatrix=0, rotate=False): ''' Finds the roots of the given list of multidimensional polynomials using a multiplication matrix. Parameters ---------- polys : list of polynomial objects Polynomials to find the common roots of. MSmatrix : int Controls which Moller-Stetter matrix is constructed. The options are: 0 (default) -- The Moller-Stetter matrix of a random polynomial Some positive integer i < dimension -- The Moller-Stetter matrix of x_i verbose : bool Prints information about how the roots are computed. returns ------- roots : numpy array The common roots of the polynomials. Each row is a root. ''' poly_type = is_power(polys, return_string=True) dim = polys[0].dim if MSmatrix not in list(range(dim + 1)): raise ValueError( 'MSmatrix must be 0 (random polynomial), or the index of a variable' ) #By Bezout's Theorem. Useful for making sure that the reduced Macaulay Matrix is as we expect degrees = [poly.degree for poly in polys] max_number_of_roots = np.prod(degrees) m_f, var_dict = MSMultMatrix(polys, poly_type, max_number_of_roots, verbose=verbose, MSmatrix=MSmatrix) if rotate: #rotate multiplication matrix 180 degrees m_f = np.rot90(m_f, 2) if verbose: print("\nM_f:\n", m_f[::-1, ::-1]) # both MSMultMatrix and will return m_f as # -1 if the ideal is not zero dimensional or if there are no roots if type(m_f) == int: return -1 # Get list of indexes of single variables and store vars that were not # in the vector space basis. var_spots = list() spot = np.zeros(dim) for i in range(dim): spot[i] = 1 if not rotate: var_spots.append(var_dict[tuple(spot)]) else: #if m_f is rotated 180, the eigenvectors are backwards var_spots.append(m_f.shape[0] - 1 - var_dict[tuple(spot)]) spot[i] = 0 # Get left eigenvectors vals, vecs = np.linalg.eig(m_f.T) if verbose: print('\nLeft Eigenvectors (as rows)\n', vecs.T) print('\nEigenvals\n', vals) zeros_spot = var_dict[tuple(0 for i in range(dim))] if rotate: #if m_f is rotate 180, the eigenvectors are backwards zeros_spot = m_f.shape[0] - 1 - zeros_spot #vecs = vecs[:,np.abs(vecs[zeros_spot]) > 1.e-10] if verbose: print('\nVariable Spots in the Vector\n', var_spots) print('\nEigeinvecs at the Variable Spots:\n', vecs[var_spots]) print('\nConstant Term Spot in the Vector\n', zeros_spot) print('\nEigeinvecs at the Constant Term\n', vecs[zeros_spot]) roots = vecs[var_spots] / vecs[zeros_spot] #Checks that the algorithm finds the correct number of roots with Bezout's Theorem assert roots.shape[ 1] <= max_number_of_roots, "Found too many roots" #Check if too many roots #if roots.shape[1] < max_number_of_roots: # warnings.warn('Expected ' + str(max_number_of_roots) # + " roots, Found " + str(roots.shape[1]) , Warning) # print("Number of Roots Lost:", max_number_of_roots - roots.shape[1]) return roots.T
def division(polys, get_divvar_coord_from_eigval=False, divisor_var=0, tol=1.e-12, verbose=False, polish=False): '''Calculates the common zeros of polynomials using a division matrix. Parameters -------- polys: MultiCheb Polynomials The polynomials for which the common roots are found. divisor_var : int What variable is being divided by. 0 is x, 1 is y, etc. Defaults to x. get_divvar_coord_from_eigval: bool Whether the divisor_var-coordinate of the roots is calculated from the eigenvalue or eigenvector. More stable to use the eigenvector. Defaults to false Returns ----------- zeros : numpy array The common roots of the polynomials. Each row is a root. ''' #This first section creates the Macaulay Matrix with the monomials that don't have #the divisor variable in the first columns. power = is_power(polys) dim = polys[0].dim #By Bezout's Theorem. Useful for making sure that the reduced Macaulay Matrix is as we expect degrees = [poly.degree for poly in polys] max_number_of_roots = np.prod(degrees) matrix_degree = np.sum(poly.degree for poly in polys) - len(polys) + 1 poly_coeff_list = [] for poly in polys: poly_coeff_list = add_polys(matrix_degree, poly, poly_coeff_list) matrix, matrix_terms, cuts = create_matrix(poly_coeff_list, matrix_degree, dim, divisor_var, get_divvar_coord_from_eigval) if verbose: np.set_printoptions(suppress=False, linewidth=200) print('\nStarting Macaulay Matrix\n', matrix) print( '\nColumns in Macaulay Matrix\nFirst element in tuple is degree of x monomial, Second element is degree of y monomial \n', matrix_terms) print( '\nLocation of Cuts in the Macaulay Matrix into [ Mb | M1* | M2* ]\n', cuts) #If bottom left is zero only does the first QR reduction on top part of matrix (for speed). Otherwise does it on the whole thing if np.allclose(matrix[cuts[0]:, :cuts[0]], 0): matrix, matrix_terms = rrqr_reduceMacaulay2(matrix, matrix_terms, cuts, max_number_of_roots, tol) else: matrix, matrix_terms = rrqr_reduceMacaulay(matrix, matrix_terms, cuts, max_number_of_roots, tol) rows, columns = matrix.shape #Make there are enough rows in the reduced Macaulay matrix, i.e. didn't loose a row #This should be a valid assert statement but the max_number_of_roots won't be the product of dimensions for non homogenous polynomials #assert rows >= columns - max_number_of_roots VB = matrix_terms[matrix.shape[0]:] matrix = np.hstack( (np.eye(rows), solve_triangular(matrix[:, :rows], matrix[:, rows:]))) if verbose: np.set_printoptions(suppress=True, linewidth=200) print("\nFinal Macaulay Matrix\n", matrix) print("\nColumns in Macaulay Matrix\n", matrix_terms) #------------> chebyshev if not power: #Builds the inverse matrix. The terms are the vector basis as well as y^k/x terms for all k. Reducing #this matrix allows the y^k/x terms to be reduced back into the vector basis. inverses = matrix_terms[np.where(matrix_terms[:, divisor_var] == 0)[0]] inverses[:, divisor_var] = -np.ones(inverses.shape[0], dtype='int') inv_matrix_terms = np.vstack((inverses, VB)) inv_matrix = np.zeros([len(inverses), len(inv_matrix_terms)]) #A bunch of different dictionaries are used below for speed purposes and to prevent repeat calculations. #A dictionary of term in inv_matrix_terms to their spot in inv_matrix_terms. inv_spot_dict = dict() spot = 0 for term in inv_matrix_terms: inv_spot_dict[tuple(term)] = spot spot += 1 #A dictionary of terms on the diagonal to their reduction in the vector basis. diag_reduction_dict = dict() for i in range(matrix.shape[0]): term = matrix_terms[i] diag_reduction_dict[tuple(term)] = matrix[i][-len(VB):] #A dictionary of terms to the terms in their quotient when divided by x. divisor_terms_dict = dict() for term in matrix_terms: divisor_terms_dict[tuple(term)] = get_divisor_terms( term, divisor_var) #A dictionary of terms to their quotient when divided by x. term_divide_dict = dict() for term in matrix_terms[-len(VB):]: term_divide_dict[tuple(term)] = divide_term( term, inv_matrix_terms, inv_spot_dict, diag_reduction_dict, len(VB), divisor_terms_dict) #Builds the inv_matrix by dividing the rows of matrix by x. for i in range(cuts[0]): inv_matrix[i] = divide_row(matrix[i][-len(VB):], matrix_terms[-len(VB):], term_divide_dict, len(inv_matrix_terms)) spot = matrix_terms[i] spot[divisor_var] -= 1 inv_matrix[i][inv_spot_dict[tuple(spot)]] += 1 #Reduces the inv_matrix to solve for the y^k/x terms in the vector basis. Q, R = qr(inv_matrix) inv_solutions = np.hstack((np.eye(R.shape[0]), solve_triangular(R[:, :R.shape[0]], R[:, R.shape[0]:]))) #A dictionary of term in the vector basis to their spot in the vector basis. VB_spot_dict = dict() spot = 0 for row in VB: VB_spot_dict[tuple(row)] = spot spot += 1 #A dictionary of terms of type y^k/x to their reduction in the vector basis. inv_reduction_dict = dict() for i in range(len(inv_solutions)): inv_reduction_dict[tuple( inv_matrix_terms[i])] = inv_solutions[i][len(inv_solutions):] #Builds the division matrix and finds the eigenvalues and eigenvectors. division_matrix = build_division_matrix(VB, VB_spot_dict, diag_reduction_dict, inv_reduction_dict, divisor_terms_dict) #<---------end Chebyshev else: #--------->Power basisDict = makeBasisDict(matrix, matrix_terms, VB) #Dictionary of terms in the vector basis their spots in the matrix. VBdict = {} spot = 0 for row in VB: VBdict[tuple(row)] = spot spot += 1 # Build division matrix division_matrix = np.zeros((len(VB), len(VB))) for i in range(VB.shape[0]): var = np.zeros(dim) var[divisor_var] = 1 term = tuple(VB[i] - var) if term in VBdict: division_matrix[VBdict[term]][i] += 1 else: division_matrix[:, i] -= basisDict[term] #<----------end Power vals, vecs = eig(division_matrix.T) if verbose: print("\nDivision Matrix\n", np.round(division_matrix[::-1, ::-1], 2)) print("\nLeft Eigenvectors (as rows)\n", vecs.T) #Calculates the zeros, the x values from the eigenvalues and the y values from the eigenvectors. zeros = list() for i in range(len(vals)): if abs(vecs[-1][i]) < 1.e-3: #This root has magnitude greater than 1, will possibly generate a false root due to instability continue root = np.zeros(dim, dtype=complex) if get_divvar_coord_from_eigval: for spot in range(0, divisor_var): root[spot] = vecs[-(2 + spot)][i] / vecs[-1][i] for spot in range(divisor_var + 1, dim): root[spot] = vecs[-(1 + spot)][i] / vecs[-1][i] else: for spot in range(0, divisor_var): root[spot] = vecs[-(3 + spot)][i] / vecs[-1][i] for spot in range(divisor_var + 1, dim): root[spot] = vecs[-(2 + spot)][i] / vecs[-1][i] if get_divvar_coord_from_eigval: #If get_divvar_coord_from_eigval, use the eigenval to calulate that coordinate root[divisor_var] = 1 / vals[i] elif power: #If using the eigenvector and it's in power form, normal calculation of coordinate root[divisor_var] = vecs[-(2)][i] / vecs[-1][i] else: #If using the eigenvector and it's in cheb form, have to use quadratic formula to calculate coordinate. vecval = vecs[-(2)][i] / vecs[-1][ i] #vecval = cT2(x)/cT1(x) = c(2x^2-1)/cx = (2x^2-1)/x root1 = root.copy() root1[divisor_var] = (vecval + np.sqrt(vecval**2 + 8)) / 4 root2 = root.copy() root2[divisor_var] = (vecval - np.sqrt(vecval**2 + 8)) / 4 if sum(np.abs(poly(root1)) for poly in polys) < sum( np.abs(poly(root2)) for poly in polys): root = root1 else: root = root2 if polish: root = newton_polish(polys, root, tol=tol) zeros.append(root) zeros = np.array(zeros) #Checks that the algorithm finds the correct number of roots with Bezout's Theorem assert zeros.shape[ 0] <= max_number_of_roots, "Found too many roots" #Check if too many roots #if zeros.shape[0] < max_number_of_roots: # warnings.warn('Expected ' + str(max_number_of_roots) # + " roots, Found " + str(zeros.shape[0]) , Warning) # print("Number of Roots Lost:", max_number_of_roots - zeros.shape[0]) return zeros