Ticket #11595: trac_11595-exact-eigenspaces.patch
File trac_11595-exact-eigenspaces.patch, 32.9 KB (added by , 11 years ago) |
---|
-
sage/calculus/wester.py
# HG changeset patch # User Rob Beezer <beezer@ups.edu> # Date 1311566250 25200 # Node ID d5a48f59037fce6106dcdd8bac96df5da8ac5281 # Parent 8532a2ad1e558cbc91ddaaa6b7cc79956dd1e8ba 11595: update exact eigenspace routines diff --git a/sage/calculus/wester.py b/sage/calculus/wester.py
a b 468 468 469 469 sage: # (YES) Find the eigenvalues of a 3x3 integer matrix. 470 470 sage: m = matrix(QQ, 3, [5,-3,-7, -2,1,2, 2,-3,-4]) 471 sage: m.eigenspaces ()471 sage: m.eigenspaces_left() 472 472 [ 473 473 (3, Vector space of degree 3 and dimension 1 over Rational Field 474 474 User basis matrix: -
sage/graphs/generic_graph.py
diff --git a/sage/graphs/generic_graph.py b/sage/graphs/generic_graph.py
a b 14242 14242 M = self.kirchhoff_matrix() 14243 14243 else: 14244 14244 M = self.adjacency_matrix() 14245 return M.right_eigenspaces(algebraic_multiplicity=False) 14245 # could pass format='all' to get QQbar eigenvalues and eigenspaces 14246 # which would be a change in default behavior 14247 return M.right_eigenspaces(format='galois', algebraic_multiplicity=False) 14246 14248 14247 14249 ### Automorphism and isomorphism 14248 14250 -
sage/matrix/matrix2.pyx
diff --git a/sage/matrix/matrix2.pyx b/sage/matrix/matrix2.pyx
a b 4361 4361 if f.degree() == n: 4362 4362 break 4363 4363 return f 4364 4365 def eigenspaces(self, var='a', even_if_inexact=None): 4364 4365 # Deprecated Aug 2008 Trac #3794 4366 # Removed July 2011 4367 # def eigenspaces(self, var='a', even_if_inexact=None): 4368 4369 def eigenspaces_left(self, format='all', var='a', algebraic_multiplicity=False): 4366 4370 r""" 4367 Deprecated: Instead of eigenspaces, use eigenspaces_left 4368 """ 4369 # sage.misc.misc.deprecation("Use eigenspaces_left") 4370 if even_if_inexact is not None: 4371 sage.misc.misc.deprecation("The 'even_if_inexact' parameter is deprecated; a warning will be issued if called over an inexact ring.") 4372 return self.eigenspaces_left(var=var) 4373 4374 4375 4376 def eigenspaces_left(self, var='a', algebraic_multiplicity=False): 4377 r""" 4378 Compute left eigenspaces of a matrix. 4379 4371 Compute the left eigenspaces of a matrix. 4372 4373 Note that ``eigenspaces_left()`` and ``left_eigenspaces()`` 4374 are identical methods. Here "left" refers to the eigenvectors 4375 being placed to the left of the matrix. 4376 4377 INPUT: 4378 4379 - ``self`` - a square matrix over an exact field. For inexact 4380 matrices consult the numerical or symbolic matrix classes. 4381 4382 - ``format`` - default: 'all' 4383 - 'all' - attempts to create every eigenspace. This will 4384 always be possible for matrices with rational entries. 4385 - 'galois' - for each irreducible factor of the characteristic 4386 polynomial, a single eigenspace will be output for a 4387 single root/eigenvalue for the irreducible factor. 4388 4389 - ``var`` - default: 'a' - variable name used to 4390 represent elements of the root field of each 4391 irreducible factor of the characteristic polynomial. 4392 If var='a', then the root fields will be in terms of 4393 a0, a1, a2, ...., where the numbering runs across all 4394 the irreducible factors of the characteristic polynomial, 4395 even for linear factors. 4396 4397 - ``algebraic_multiplicity`` - default: False - whether or 4398 not to include the algebraic multiplicity of each eigenvalue 4399 in the output. See the discussion below. 4400 4401 OUTPUT: 4402 4380 4403 If algebraic_multiplicity=False, return a list of pairs (e, V) 4381 where e runs through all eigenvalues (up to Galois conjugation) of 4382 this matrix, and V is the corresponding left eigenspace. 4383 4384 If algebraic_multiplicity=True, return a list of pairs (e, V, n) 4404 where e is an eigenvalue of the matrix, and V is the corresponding 4405 left eigenspace. For Galois conjugates of eigenvalues, there 4406 may be just one representative eigenspace, depending on the 4407 ``format`` keyword. 4408 4409 If algebraic_multiplicity=True, return a list of triples (e, V, n) 4385 4410 where e and V are as above and n is the algebraic multiplicity of 4386 the eigenvalue. If the eigenvalues are given symbolically, as roots 4387 of an irreducible factor of the characteristic polynomial, then the 4388 algebraic multiplicity returned is the multiplicity of each 4389 conjugate eigenvalue. 4390 4391 The eigenspaces are returned sorted by the corresponding 4392 characteristic polynomials, where polynomials are sorted in 4393 dictionary order starting with constant terms. 4394 4395 INPUT: 4396 4397 4398 - ``var`` - variable name used to represent elements 4399 of the root field of each irreducible factor of the characteristic 4400 polynomial I.e., if var='a', then the root fields will be in terms 4401 of a0, a1, a2, ..., ak. 4402 4403 4411 the eigenvalue. 4412 4404 4413 .. warning:: 4405 4414 4406 4415 Uses a somewhat naive algorithm (simply factors the 4407 4416 characteristic polynomial and computes kernels directly 4408 4417 over the extension field). 4409 4418 4410 TODO: 4411 4412 Maybe implement the better algorithm that is in 4413 dual_eigenvector in sage/modular/hecke/module.py. 4414 4415 EXAMPLES: We compute the left eigenspaces of a `3\times 3` 4416 rational matrix. 4417 4418 :: 4419 4419 EXAMPLES: 4420 4421 We compute the left eigenspaces of a `3\times 3` 4422 rational matrix. First, we request `all` of the eigenvalues, 4423 so the results are in the field of algebraic numbers, `QQbar`. 4424 Then we request just one eigenspace per irreducible factor of 4425 the characteristic polynomial with the `galois` keyword. :: 4426 4420 4427 sage: A = matrix(QQ,3,3,range(9)); A 4421 4428 [0 1 2] 4422 4429 [3 4 5] 4423 4430 [6 7 8] 4424 sage: es = A.eigenspaces_left(); es 4431 sage: es = A.eigenspaces_left(format='all'); es 4432 [ 4433 (0, Vector space of degree 3 and dimension 1 over Rational Field 4434 User basis matrix: 4435 [ 1 -2 1]), 4436 (-1.348469228349535?, Vector space of degree 3 and dimension 1 over Algebraic Field 4437 User basis matrix: 4438 [ 1 0.3101020514433644? -0.3797958971132713?]), 4439 (13.34846922834954?, Vector space of degree 3 and dimension 1 over Algebraic Field 4440 User basis matrix: 4441 [ 1 1.289897948556636? 1.579795897113272?]) 4442 ] 4443 4444 sage: es = A.eigenspaces_left(format='galois'); es 4425 4445 [ 4426 4446 (0, Vector space of degree 3 and dimension 1 over Rational Field 4427 4447 User basis matrix: … … 4430 4450 User basis matrix: 4431 4451 [ 1 1/15*a1 + 2/5 2/15*a1 - 1/5]) 4432 4452 ] 4433 sage: es = A.eigenspaces_left( algebraic_multiplicity=True); es4453 sage: es = A.eigenspaces_left(format='galois', algebraic_multiplicity=True); es 4434 4454 [ 4435 4455 (0, Vector space of degree 3 and dimension 1 over Rational Field 4436 4456 User basis matrix: … … 4443 4463 sage: delta = e*v - v*A 4444 4464 sage: abs(abs(delta)) < 1e-10 4445 4465 True 4446 4447 The same computation, but with implicit base change to a field ::4448 4466 4467 The same computation, but with implicit base change to a field. :: 4468 4449 4469 sage: A = matrix(ZZ,3,range(9)); A 4450 4470 [0 1 2] 4451 4471 [3 4 5] 4452 4472 [6 7 8] 4453 sage: A.eigenspaces_left( )4473 sage: A.eigenspaces_left(format='galois') 4454 4474 [ 4455 4475 (0, Vector space of degree 3 and dimension 1 over Rational Field 4456 4476 User basis matrix: … … 4459 4479 User basis matrix: 4460 4480 [ 1 1/15*a1 + 2/5 2/15*a1 - 1/5]) 4461 4481 ] 4462 4482 4463 4483 We compute the left eigenspaces of the matrix of the Hecke operator 4464 `T_2` on level 43 modular symbols. 4465 4466 :: 4467 4484 `T_2` on level 43 modular symbols, both with all eigenvalues (the default) 4485 and with one subspace per factor. :: 4486 4468 4487 sage: A = ModularSymbols(43).T(2).matrix(); A 4469 4488 [ 3 0 0 0 0 0 -1] 4470 4489 [ 0 -2 1 0 0 0 0] … … 4488 4507 User basis matrix: 4489 4508 [ 0 1 0 1 -1 1 -1] 4490 4509 [ 0 0 1 0 -1 2 -1], 2), 4510 (-1.414213562373095?, Vector space of degree 7 and dimension 2 over Algebraic Field 4511 User basis matrix: 4512 [ 0 1 0 -1 0.4142135623730951? 1 -1] 4513 [ 0 0 1 0 -1 0 2.414213562373095?], 2), 4514 (1.414213562373095?, Vector space of degree 7 and dimension 2 over Algebraic Field 4515 User basis matrix: 4516 [ 0 1 0 -1 -2.414213562373095? 1 -1] 4517 [ 0 0 1 0 -1 0 -0.4142135623730951?], 2) 4518 ] 4519 sage: A.eigenspaces_left(format='galois', algebraic_multiplicity=True) 4520 [ 4521 (3, Vector space of degree 7 and dimension 1 over Rational Field 4522 User basis matrix: 4523 [ 1 0 1/7 0 -1/7 0 -2/7], 1), 4524 (-2, Vector space of degree 7 and dimension 2 over Rational Field 4525 User basis matrix: 4526 [ 0 1 0 1 -1 1 -1] 4527 [ 0 0 1 0 -1 2 -1], 2), 4491 4528 (a2, Vector space of degree 7 and dimension 2 over Number Field in a2 with defining polynomial x^2 - 2 4492 4529 User basis matrix: 4493 4530 [ 0 1 0 -1 -a2 - 1 1 -1] 4494 4531 [ 0 0 1 0 -1 0 -a2 + 1], 2) 4495 4532 ] 4496 4497 Next we compute the left eigenspaces over the finite field of order 4498 11:: 4499 4533 4534 Next we compute the left eigenspaces over the finite field of order 11. :: 4535 4500 4536 sage: A = ModularSymbols(43, base_ring=GF(11), sign=1).T(2).matrix(); A 4501 4537 [ 3 9 0 0] 4502 4538 [ 0 9 0 1] … … 4506 4542 Finite Field of size 11 4507 4543 sage: A.charpoly() 4508 4544 x^4 + 10*x^3 + 3*x^2 + 2*x + 1 4509 sage: A.eigenspaces_left( var = 'beta')4545 sage: A.eigenspaces_left(format='galois', var = 'beta') 4510 4546 [ 4511 4547 (9, Vector space of degree 4 and dimension 1 over Finite Field of size 11 4512 4548 User basis matrix: … … 4518 4554 User basis matrix: 4519 4555 [ 0 1 0 5*beta2 + 10]) 4520 4556 ] 4521 4522 TESTS: 4523 4524 Warnings are issued if the generic algorithm is used over 4525 inexact fields. Garbage may result in these cases because of 4526 numerical precision issues. 4527 4528 :: 4529 4530 sage: R=RealField(30) 4531 sage: M=matrix(R,2,[2,1,1,1]) 4532 sage: M.eigenspaces_left() # random output from numerical issues 4557 4558 This method is only applicable to exact matrices. 4559 The "eigenmatrix" routines for matrices with double-precision 4560 floating-point entries (``RDF``, ``CDF``) are the best 4561 alternative. There are also "eigenmatrix" routines for 4562 matrices with symbolic entries. :: 4563 4564 sage: A = matrix(QQ, 3, 3, range(9)) 4565 sage: A.change_ring(RR).eigenspaces_left() 4566 Traceback (most recent call last): 4567 ... 4568 NotImplementedError: eigenspaces cannot be computed reliably for inexact rings such as Real Field with 53 bits of precision, 4569 consult numerical or symbolic matrix classes for other options 4570 4571 sage: em = A.change_ring(RDF).eigenmatrix_left() 4572 sage: eigenvalues = em[0]; eigenvalues 4573 [ 13.3484692... 0 0] 4574 [ 0 -1.34846922... 0] 4575 [ 0 0 -6.2265089...e-16] 4576 sage: eigenvectors = em[1]; eigenvectors 4577 [ 0.440242867... 0.567868371... 0.695493875...] 4578 [ 0.897878732... 0.278434036... -0.341010658...] 4579 [ 0.408248290... -0.816496580... 0.408248290...] 4580 4581 sage: x, y = var('x y') 4582 sage: S = matrix([[x, y], [y, 3*x^2]]) 4583 sage: em = S.eigenmatrix_left() 4584 sage: eigenvalues = em[0]; eigenvalues 4585 [-1/2*sqrt(9*x^4 - 6*x^3 + x^2 + 4*y^2) + 3/2*x^2 + 1/2*x 0] 4586 [ 0 3/2*x^2 + 1/2*x + 1/2*sqrt(9*x^4 - 6*x^3 + x^2 + 4*y^2)] 4587 sage: eigenvectors = em[1]; eigenvectors 4588 [ 1 1/2*(3*x^2 - x - sqrt(9*x^4 - 6*x^3 + x^2 + 4*y^2))/y] 4589 [ 1 -1/2*(x - sqrt(9*x^4 - 6*x^3 + x^2 + 4*y^2) - 3*x^2)/y] 4590 4591 A request for ``'all'`` the eigenvalues, when it is not 4592 possible, will raise an error. Using the ``'galois'`` 4593 format option is more likely to be successful. :: 4594 4595 sage: F.<b> = FiniteField(11^2) 4596 sage: A = matrix(F, [[b + 1, b + 1], [10*b + 4, 5*b + 4]]) 4597 sage: A.eigenspaces_left(format='all') 4598 Traceback (most recent call last): 4599 ... 4600 NotImplementedError: unable to construct eigenspaces for eigenvalues outside the base field, 4601 try the keyword option: format='galois' 4602 4603 sage: A.eigenspaces_left(format='galois') 4533 4604 [ 4534 ( 2.6180340, Vector space of degree 2 and dimension 0 over Real Field with 30 bits of precision4605 (a0, Vector space of degree 2 and dimension 1 over Univariate Quotient Polynomial Ring in a0 over Finite Field in b of size 11^2 with modulus x^2 + (5*b + 6)*x + 8*b + 10 4535 4606 User basis matrix: 4536 []), 4537 (0.38196601, Vector space of degree 2 and dimension 0 over Real Field with 30 bits of precision 4538 User basis matrix: 4539 []) 4607 [ 1 6*b*a0 + 3*b + 1]) 4540 4608 ] 4541 sage: R=ComplexField(30) 4542 sage: N=matrix(R,2,[2,1,1,1]) 4543 sage: N.eigenspaces_left() # random output from numerical issues 4544 [ 4545 (2.6180340, Vector space of degree 2 and dimension 0 over Complex Field with 30 bits of precision 4546 User basis matrix: 4547 []), 4548 (0.38196601, Vector space of degree 2 and dimension 0 over Complex Field with 30 bits of precision 4549 User basis matrix: 4550 []) 4551 ] 4552 """ 4553 x = self.fetch('eigenspaces_left') 4609 4610 TESTS:: 4611 4612 sage: B = matrix(QQ, 2, 3, range(6)) 4613 sage: B.eigenspaces_left() 4614 Traceback (most recent call last): 4615 ... 4616 TypeError: matrix must be square, not 2 x 3 4617 4618 sage: B = matrix(QQ, 4, 4, range(16)) 4619 sage: B.eigenspaces_left(format='junk') 4620 Traceback (most recent call last): 4621 ... 4622 ValueError: format keyword must be 'all' or 'galois', not junk 4623 4624 sage: B.eigenspaces_left(algebraic_multiplicity='garbage') 4625 Traceback (most recent call last): 4626 ... 4627 ValueError: algebraic_multiplicity keyword must be True or False 4628 4629 """ 4630 if not format in ['all', 'galois']: 4631 raise ValueError("format keyword must be 'all' or 'galois', not {0}".format(format)) 4632 if not algebraic_multiplicity in [True, False]: 4633 raise ValueError('algebraic_multiplicity keyword must be True or False'.format(algebraic_multiplicity)) 4634 if not self.is_square(): 4635 raise TypeError('matrix must be square, not {0} x {1}'.format(self.nrows(), self.ncols())) 4636 if not self.base_ring().is_exact(): 4637 msg = ("eigenspaces cannot be computed reliably for inexact rings such as {0},\n", 4638 "consult numerical or symbolic matrix classes for other options") 4639 raise NotImplementedError(''.join(msg).format(self.base_ring())) 4640 4641 key = 'eigenspaces_left_' + format + '{0}'.format(var) 4642 x = self.fetch(key) 4554 4643 if not x is None: 4555 4644 if algebraic_multiplicity: 4556 4645 return x 4557 4646 else: 4558 return Sequence([(e[0],e[1]) for e in x], cr=True) 4559 4560 if not self.base_ring().is_exact(): 4561 from warnings import warn 4562 warn("Using generic algorithm for an inexact ring, which will probably give incorrect results due to numerical precision issues.") 4563 4564 4565 4566 # minpoly is rarely implemented and is unreliable (leading to hangs) via linbox when implemented 4567 # as of 2007-03-25. 4568 #try: 4569 # G = self.minpoly().factor() # can be computed faster when available. 4570 #except NotImplementedError: 4571 G = self.fcp() # factored charpoly of self. 4647 return Sequence([(e[0],e[1]) for e in x], cr=True, check=False) 4648 4649 # Possible improvements: 4650 # algorithm for dual_eigenvector in sage/modular/hecke/module.py 4651 # use minpoly when algebraic_multiplicity=False 4652 # as of 2007-03-25 minpoly is unreliable via linbox 4653 4654 import sage.rings.qqbar 4655 import sage.categories.homset 4656 G = self.fcp() # factored characteristic polynomial 4572 4657 V = [] 4573 i = 04658 i = -1 # variable name index, increments for each eigenvalue 4574 4659 for h, e in G: 4660 i = i + 1 4575 4661 if h.degree() == 1: 4576 4662 alpha = -h[0]/h[1] 4577 4663 F = alpha.parent() 4578 4664 if F != self.base_ring(): 4579 4665 self = self.change_ring(F) 4580 4666 A = self - alpha 4667 W = A.kernel() 4668 V.append((alpha, W.ambient_module().span_of_basis(W.basis()), e)) 4581 4669 else: 4582 F = h.root_field(' %s%s'%(var,i))4670 F = h.root_field('{0}{1}'.format(var,i)) 4583 4671 alpha = F.gen(0) 4584 4672 A = self.change_ring(F) - alpha 4585 W = A.kernel() 4586 i = i + 1 4587 V.append((alpha, W.ambient_module().span_of_basis(W.basis()), e)) 4588 V = Sequence(V, cr=True) 4589 self.cache('eigenspaces_left', V) 4673 W = A.kernel() 4674 WB = W.basis() 4675 if format == 'galois': 4676 V.append((alpha, W.ambient_module().span_of_basis(WB), e)) 4677 elif format == 'all': 4678 try: 4679 alpha_conj = alpha.galois_conjugates(sage.rings.qqbar.QQbar) 4680 except (AttributeError, TypeError): 4681 msg = ("unable to construct eigenspaces for eigenvalues outside the base field,\n" 4682 "try the keyword option: format='galois'") 4683 raise NotImplementedError(''.join(msg)) 4684 for ev in alpha_conj: 4685 m = sage.categories.homset.hom(alpha.parent(), ev.parent(), ev) 4686 space = (ev.parent())**self.nrows() 4687 evec_list = [(space)([m(i) for i in v]) for v in WB] 4688 V.append((ev, space.span_of_basis(evec_list, already_echelonized=True), e)) 4689 V = Sequence(V, cr=True, check=False) 4690 self.cache(key, V) 4590 4691 if algebraic_multiplicity: 4591 4692 return V 4592 4693 else: 4593 return Sequence([(e[0],e[1]) for e in V], cr=True) 4594 4595 def eigenspaces_right(self, var='a', algebraic_multiplicity=False): 4694 return Sequence([(e[0],e[1]) for e in V], cr=True, check=False) 4695 4696 left_eigenspaces = eigenspaces_left 4697 4698 def eigenspaces_right(self, format='all', var='a', algebraic_multiplicity=False): 4596 4699 r""" 4597 Compute right eigenspaces of a matrix. 4598 4700 Compute the right eigenspaces of a matrix. 4701 4702 Note that ``eigenspaces_right()`` and ``right_eigenspaces()`` 4703 are identical methods. Here "right" refers to the eigenvectors 4704 being placed to the right of the matrix. 4705 4706 INPUT: 4707 4708 - ``self`` - a square matrix over an exact field. For inexact 4709 matrices consult the numerical or symbolic matrix classes. 4710 4711 - ``format`` - default: 'all' 4712 - 'all' - attempts to create every eigenspace. This will 4713 always be possible for matrices with rational entries. 4714 - 'galois' - for each irreducible factor of the characteristic 4715 polynomial, a single eigenspace will be output for a 4716 single root/eigenvalue for the irreducible factor. 4717 4718 - ``var`` - default: 'a' - variable name used to 4719 represent elements of the root field of each 4720 irreducible factor of the characteristic polynomial. 4721 If var='a', then the root fields will be in terms of 4722 a0, a1, a2, ...., where the numbering runs across all 4723 the irreducible factors of the characteristic polynomial, 4724 even for linear factors. 4725 4726 - ``algebraic_multiplicity`` - default: False - whether or 4727 not to include the algebraic multiplicity of each eigenvalue 4728 in the output. See the discussion below. 4729 4730 OUTPUT: 4731 4599 4732 If algebraic_multiplicity=False, return a list of pairs (e, V) 4600 where e runs through all eigenvalues (up to Galois conjugation) of 4601 this matrix, and V is the corresponding right eigenspace. 4602 4603 If algebraic_multiplicity=True, return a list of pairs (e, V, n) 4733 where e is an eigenvalue of the matrix, and V is the corresponding 4734 left eigenspace. For Galois conjugates of eigenvalues, there 4735 may be just one representative eigenspace, depending on the 4736 ``format`` keyword. 4737 4738 If algebraic_multiplicity=True, return a list of triples (e, V, n) 4604 4739 where e and V are as above and n is the algebraic multiplicity of 4605 the eigenvalue. If the eigenvalues are given symbolically, as roots 4606 of an irreducible factor of the characteristic polynomial, then the 4607 algebraic multiplicity returned is the multiplicity of each 4608 conjugate eigenvalue. 4609 4610 The eigenspaces are returned sorted by the corresponding 4611 characteristic polynomials, where polynomials are sorted in 4612 dictionary order starting with constant terms. 4613 4614 INPUT: 4615 4616 4617 - ``var`` - variable name used to represent elements 4618 of the root field of each irreducible factor of the characteristic 4619 polynomial I.e., if var='a', then the root fields will be in terms 4620 of a0, a1, a2, ..., ak. 4621 4622 4740 the eigenvalue. 4741 4623 4742 .. warning:: 4624 4743 4625 4744 Uses a somewhat naive algorithm (simply factors the 4626 4745 characteristic polynomial and computes kernels directly 4627 4746 over the extension field). 4628 4747 4629 TODO: Maybe implement the better algorithm that4630 is in dual_eigenvector in sage/modular/hecke/module.py. 4631 4632 EXAMPLES: We compute the right eigenspaces of a `3\times 3`4633 rational matrix.4634 4635 ::4636 4637 sage: A = matrix(QQ, 3,3,range(9)); A4748 EXAMPLES: 4749 4750 Right eigenspaces are computed from the left eigenspaces of the 4751 transpose of the matrix. As such, there is a greater collection 4752 of illustrative examples at the :meth:`eigenspaces_left`. 4753 4754 We compute the right eigenspaces of a `3\times 3` rational matrix. :: 4755 4756 sage: A = matrix(QQ, 3 ,3, range(9)); A 4638 4757 [0 1 2] 4639 4758 [3 4 5] 4640 4759 [6 7 8] 4641 sage: es = A.eigenspaces_right(); es 4760 sage: A.eigenspaces_right() 4761 [ 4762 (0, Vector space of degree 3 and dimension 1 over Rational Field 4763 User basis matrix: 4764 [ 1 -2 1]), 4765 (-1.348469228349535?, Vector space of degree 3 and dimension 1 over Algebraic Field 4766 User basis matrix: 4767 [ 1 0.1303061543300932? -0.7393876913398137?]), 4768 (13.34846922834954?, Vector space of degree 3 and dimension 1 over Algebraic Field 4769 User basis matrix: 4770 [ 1 3.069693845669907? 5.139387691339814?]) 4771 ] 4772 sage: es = A.eigenspaces_right(format='galois'); es 4642 4773 [ 4643 4774 (0, Vector space of degree 3 and dimension 1 over Rational Field 4644 4775 User basis matrix: … … 4647 4778 User basis matrix: 4648 4779 [ 1 1/5*a1 + 2/5 2/5*a1 - 1/5]) 4649 4780 ] 4650 sage: es = A.eigenspaces_right( algebraic_multiplicity=True); es4781 sage: es = A.eigenspaces_right(format='galois', algebraic_multiplicity=True); es 4651 4782 [ 4652 4783 (0, Vector space of degree 3 and dimension 1 over Rational Field 4653 4784 User basis matrix: … … 4660 4791 sage: delta = v*e - A*v 4661 4792 sage: abs(abs(delta)) < 1e-10 4662 4793 True 4663 4794 4664 4795 The same computation, but with implicit base change to a field:: 4665 4666 sage: A = matrix(ZZ, 3,range(9)); A4796 4797 sage: A = matrix(ZZ, 3, range(9)); A 4667 4798 [0 1 2] 4668 4799 [3 4 5] 4669 4800 [6 7 8] 4670 sage: A.eigenspaces_right( )4801 sage: A.eigenspaces_right(format='galois') 4671 4802 [ 4672 4803 (0, Vector space of degree 3 and dimension 1 over Rational Field 4673 4804 User basis matrix: … … 4676 4807 User basis matrix: 4677 4808 [ 1 1/5*a1 + 2/5 2/5*a1 - 1/5]) 4678 4809 ] 4679 4680 TESTS: Warnings are issued if the generic algorithm is used over 4681 inexact fields. Garbage may result in these cases because of 4682 numerical precision issues. 4683 4684 :: 4685 4686 sage: R=RealField(30) 4687 sage: M=matrix(R,2,[2,1,1,1]) 4688 sage: M.eigenspaces_right() # random output from numerical issues 4689 [(2.6180340, 4690 Vector space of degree 2 and dimension 0 over Real Field with 30 bits of precision 4691 User basis matrix: 4692 []), 4693 (0.38196601, 4694 Vector space of degree 2 and dimension 0 over Real Field with 30 bits of precision 4695 User basis matrix: 4696 [])] 4697 sage: R=ComplexField(30) 4698 sage: N=matrix(R,2,[2,1,1,1]) 4699 sage: N.eigenspaces_right() # random output from numerical issues 4700 [(2.6180340, 4701 Vector space of degree 2 and dimension 0 over Complex Field with 30 bits of precision 4702 User basis matrix: 4703 []), 4704 (0.38196601, 4705 Vector space of degree 2 and dimension 0 over Complex Field with 30 bits of precision 4706 User basis matrix: 4707 [])] 4708 """ 4709 return self.transpose().eigenspaces_left(var=var, algebraic_multiplicity=algebraic_multiplicity) 4810 4811 This method is only applicable to exact matrices. 4812 The "eigenmatrix" routines for matrices with double-precision 4813 floating-point entries (``RDF``, ``CDF``) are the best 4814 alternative. There are also "eigenmatrix" routines for 4815 matrices with symbolic entries. :: 4816 4817 sage: B = matrix(RR, 3, 3, range(9)) 4818 sage: B.eigenspaces_right() 4819 Traceback (most recent call last): 4820 ... 4821 NotImplementedError: eigenspaces cannot be computed reliably for inexact rings such as Real Field with 53 bits of precision, 4822 consult numerical or symbolic matrix classes for other options 4823 4824 sage: em = B.change_ring(RDF).eigenmatrix_right() 4825 sage: eigenvalues = em[0]; eigenvalues 4826 [ 13.3484692... 0 0] 4827 [ 0 -1.34846922... 0] 4828 [ 0 0 -8.86256604...e-16] 4829 sage: eigenvectors = em[1]; eigenvectors 4830 [ 0.164763817... 0.799699663... 0.408248290...] 4831 [ 0.505774475... 0.104205787... -0.816496580...] 4832 [ 0.846785134... -0.591288087... 0.408248290...] 4833 4834 sage: x, y = var('x y') 4835 sage: S = matrix([[x, y], [y, 3*x^2]]) 4836 sage: em = S.eigenmatrix_right() 4837 sage: eigenvalues = em[0]; eigenvalues 4838 [-1/2*sqrt(9*x^4 - 6*x^3 + x^2 + 4*y^2) + 3/2*x^2 + 1/2*x 0] 4839 [ 0 3/2*x^2 + 1/2*x + 1/2*sqrt(9*x^4 - 6*x^3 + x^2 + 4*y^2)] 4840 sage: eigenvectors = em[1]; eigenvectors 4841 [ 1 1] 4842 [ 1/2*(3*x^2 - x - sqrt(9*x^4 - 6*x^3 + x^2 + 4*y^2))/y -1/2*(x - sqrt(9*x^4 - 6*x^3 + x^2 + 4*y^2) - 3*x^2)/y] 4843 4844 TESTS:: 4845 4846 sage: B = matrix(QQ, 2, 3, range(6)) 4847 sage: B.eigenspaces_right() 4848 Traceback (most recent call last): 4849 ... 4850 TypeError: matrix must be square, not 2 x 3 4851 4852 sage: B = matrix(QQ, 4, 4, range(16)) 4853 sage: B.eigenspaces_right(format='junk') 4854 Traceback (most recent call last): 4855 ... 4856 ValueError: format keyword must be 'all' or 'galois', not junk 4857 4858 sage: B.eigenspaces_right(algebraic_multiplicity='garbage') 4859 Traceback (most recent call last): 4860 ... 4861 ValueError: algebraic_multiplicity keyword must be True or False 4862 """ 4863 if not format in ['all', 'galois']: 4864 raise ValueError("format keyword must be 'all' or 'galois', not {0}".format(format)) 4865 if not algebraic_multiplicity in [True, False]: 4866 raise ValueError('algebraic_multiplicity keyword must be True or False'.format(algebraic_multiplicity)) 4867 if not self.is_square(): 4868 raise TypeError('matrix must be square, not {0} x {1}'.format(self.nrows(), self.ncols())) 4869 if not self.base_ring().is_exact(): 4870 msg = ("eigenspaces cannot be computed reliably for inexact rings such as {0},\n", 4871 "consult numerical or symbolic matrix classes for other options") 4872 raise NotImplementedError(''.join(msg).format(self.base_ring())) 4873 4874 key = 'eigenspaces_right_' + format + '{0}'.format(var) 4875 x = self.fetch(key) 4876 if not x is None: 4877 if algebraic_multiplicity: 4878 return x 4879 else: 4880 return Sequence([(e[0],e[1]) for e in x], cr=True, check=False) 4881 4882 V = self.transpose().eigenspaces_left(format=format, var=var, algebraic_multiplicity=True) 4883 4884 self.cache(key, V) 4885 if algebraic_multiplicity: 4886 return V 4887 else: 4888 return Sequence([(e[0],e[1]) for e in V], cr=True, check=False) 4710 4889 4711 4890 right_eigenspaces = eigenspaces_right 4712 4891 … … 4883 5062 V = [] 4884 5063 from sage.rings.qqbar import QQbar 4885 5064 from sage.categories.homset import hom 4886 eigenspaces = self.eigenspaces_left( algebraic_multiplicity=True)5065 eigenspaces = self.eigenspaces_left(format='galois', algebraic_multiplicity=True) 4887 5066 evec_list=[] 4888 5067 n = self._nrows 4889 5068 evec_eval_list = [] -
sage/modules/free_module.py
diff --git a/sage/modules/free_module.py b/sage/modules/free_module.py
a b 816 816 817 817 :: 818 818 819 sage: v = matrix(RDF, 3, range(9)).eigenspaces ()[0][1].basis()[0]819 sage: v = matrix(RDF, 3, range(9)).eigenspaces_left()[0][1].basis()[0] 820 820 sage: v.complex_vector() 821 821 (...0.440242867..., ...0.567868371..., ...0.695493875...) 822 822 sage: v.complex_vector().parent().echelonized_basis_matrix()[0] * v[0] -
sage/quadratic_forms/quadratic_form__local_representation_conditions.py
diff --git a/sage/quadratic_forms/quadratic_form__local_representation_conditions.py b/sage/quadratic_forms/quadratic_form__local_representation_conditions.py
a b 143 143 ## Compute the local conditions at the real numbers (i.e. "p = infinity") 144 144 ## ---------------------------------------------------------------------- 145 145 M = Q.matrix() 146 E = M.eigenspaces ()146 E = M.eigenspaces_left() 147 147 M_eigenvalues = [E[i][0] for i in range(len(E))] 148 148 149 149 pos_flag = infinity