Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

""" 

Dense matrices over the Complex Double Field using NumPy 

  

EXAMPLES:: 

  

sage: b = Mat(CDF,2,3).basis() 

sage: b[0,0] 

[1.0 0.0 0.0] 

[0.0 0.0 0.0] 

  

We deal with the case of zero rows or zero columns:: 

  

sage: m = MatrixSpace(CDF,0,3) 

sage: m.zero_matrix() 

[] 

  

TESTS:: 

  

sage: a = matrix(CDF,2,[i+(4-i)*I for i in range(4)], sparse=False) 

sage: TestSuite(a).run() 

sage: Mat(CDF,0,0).zero_matrix().inverse() 

[] 

  

AUTHORS: 

  

- Jason Grout (2008-09): switch to NumPy backend 

  

- Josh Kantor 

  

- William Stein: many bug fixes and touch ups. 

""" 

  

############################################################################## 

# Copyright (C) 2004,2005,2006 Joshua Kantor <kantor.jm@gmail.com> 

# Distributed under the terms of the GNU General Public License (GPL) 

# The full text of the GPL is available at: 

# http://www.gnu.org/licenses/ 

############################################################################## 

from __future__ import absolute_import 

  

from sage.rings.complex_double import CDF 

  

cimport numpy as cnumpy 

  

numpy=None 

  

cdef class Matrix_complex_double_dense(Matrix_double_dense): 

""" 

Class that implements matrices over the real double field. These 

are supposed to be fast matrix operations using C doubles. Most 

operations are implemented using numpy which will call the 

underlying BLAS on the system. 

  

EXAMPLES:: 

  

sage: m = Matrix(CDF, [[1,2*I],[3+I,4]]) 

sage: m**2 

[-1.0 + 6.0*I 10.0*I] 

[15.0 + 5.0*I 14.0 + 6.0*I] 

sage: n= m^(-1); n # abs tol 1e-15 

[ 0.3333333333333333 + 0.3333333333333333*I 0.16666666666666669 - 0.16666666666666666*I] 

[-0.16666666666666666 - 0.3333333333333333*I 0.08333333333333331 + 0.08333333333333333*I] 

  

To compute eigenvalues the use the functions ``left_eigenvectors`` or 

``right_eigenvectors``:: 

  

sage: p,e = m.right_eigenvectors() 

  

the result of eigen is a pair (p,e), where p is a list of 

eigenvalues and the e is a matrix whose columns are the 

eigenvectors. 

  

To solve a linear system Ax = b where A = [[1,2] and b = [5,6] 

[3,4]] 

  

:: 

  

sage: b = vector(CDF,[5,6]) 

sage: m.solve_right(b) # abs tol 1e-14 

(2.6666666666666665 + 0.6666666666666669*I, -0.3333333333333333 - 1.1666666666666667*I) 

  

See the commands qr, lu, and svd for QR, LU, and singular value 

decomposition. 

""" 

def __cinit__(self, parent, entries, copy, coerce): 

global numpy 

if numpy is None: 

import numpy 

self._numpy_dtype = numpy.dtype('complex128') 

self._python_dtype = complex 

self._numpy_dtypeint = cnumpy.NPY_CDOUBLE 

# TODO: Make ComplexDoubleElement instead of CDF for speed 

self._sage_dtype = CDF 

self.__create_matrix__() 

return