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

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

r""" 

Dense complex double vectors using a NumPy backend. 

  

EXAMPLES:: 

  

sage: v = vector(CDF,[(1,-1), (2,pi), (3,5)]) 

sage: v 

(1.0 - 1.0*I, 2.0 + 3.141592653589793*I, 3.0 + 5.0*I) 

sage: type(v) 

<type 'sage.modules.vector_complex_double_dense.Vector_complex_double_dense'> 

sage: parent(v) 

Vector space of dimension 3 over Complex Double Field 

sage: v[0] = 5 

sage: v 

(5.0, 2.0 + 3.141592653589793*I, 3.0 + 5.0*I) 

sage: loads(dumps(v)) == v 

True 

  

TESTS:: 

  

sage: v = vector(CDF, [2, 2]) 

sage: v - v 

(0.0, 0.0) 

sage: (v - v).norm() 

0.0 

  

AUTHORS: 

  

-- Jason Grout, Oct 2008: switch to NumPy backend, factored out 

Vector_double_dense class 

""" 

  

#***************************************************************************** 

# Copyright (C) 2008 Jason Grout <jason-sage@creativetrax.com> 

# 

# This program is free software: you can redistribute it and/or modify 

# it under the terms of the GNU General Public License as published by 

# the Free Software Foundation, either version 2 of the License, or 

# (at your option) any later version. 

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

#***************************************************************************** 

from __future__ import absolute_import 

  

from sage.rings.complex_double import CDF 

  

cimport numpy 

  

  

cdef class Vector_complex_double_dense(Vector_double_dense): 

""" 

Vectors over the Complex Double Field. These are supposed to be 

fast vector operations using C doubles. Most operations are 

implemented using numpy which will call the underlying BLAS, if 

needed, on the system. 

  

EXAMPLES: 

sage: v = vector(CDF,[(1,-1), (2,pi), (3,5)]) 

sage: v 

(1.0 - 1.0*I, 2.0 + 3.141592653589793*I, 3.0 + 5.0*I) 

sage: v*v # rel tol 1e-15 

-21.86960440108936 + 40.56637061435917*I 

""" 

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

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

self._python_dtype = complex 

self._numpy_dtypeint = numpy.NPY_CDOUBLE 

# TODO: Make ComplexDoubleElement instead of CDF for speed 

self._sage_dtype = CDF 

self.__create_vector__() 

return 

  

def __reduce__(self): 

""" 

Pickling 

  

EXAMPLES: 

sage: a = vector(CDF, range(9)) 

sage: loads(dumps(a)) == a 

True 

""" 

return (unpickle_v1, (self._parent, self.list(), self._degree, self._is_mutable)) 

  

  

# For backwards compatibility, we must keep the function unpickle_v0 

def unpickle_v0(parent, entries, degree): 

""" 

Create a complex double vector containing the entries. 

  

EXAMPLES: 

sage: v = vector(CDF, [1,2,3]) 

sage: w = sage.modules.vector_complex_double_dense.unpickle_v0(v.parent(), list(v), v.degree()) 

sage: v == w 

True 

""" 

return unpickle_v1(parent, entries, degree) 

  

def unpickle_v1(parent, entries, degree, is_mutable=None): 

""" 

Create a complex double vector with the given parent, entries, 

degree, and mutability. 

  

EXAMPLES: 

sage: v = vector(CDF, [1,2,3]) 

sage: w = sage.modules.vector_complex_double_dense.unpickle_v1(v.parent(), list(v), v.degree(), v.is_mutable()) 

sage: v == w 

True 

""" 

cdef Vector_complex_double_dense v = Vector_complex_double_dense(parent, entries) 

if is_mutable is not None: 

v._is_mutable = is_mutable 

return v