Package sage :: Package gsl :: Module dft :: Class IndexedSequence
[hide private]
[frames] | no frames]

Class IndexedSequence

source code

                      object --+    
                               |    
structure.sage_object.SageObject --+
                                   |
                                  IndexedSequence

Instance Methods [hide private]
 
__init__(self, L, index_object)
\code{index_object} must be a SAGE object with an _iter_ method containing the same number of elements as self, which is a list of elements taken from a field.
source code
 
dict(self) source code
 
list(self) source code
 
base_ring(self)
This just returns the common parent R of the N list elements.
source code
 
index_object(self) source code
 
_repr_(self)
Implements print method.
source code
 
plot_histogram(self)
Plots the histogram plot of the sequence, which is assumed to be real or from a finite field, with a real indexing set I coercible into RR.
source code
 
plot(self)
Plots the points of the sequence, whose elements are assumed to be real or from a finite field, with a real indexing set I = range(len(self)).
source code
 
dft(self, chi=<function <lambda> at 0x4f22848>)
Implements a discrete Fourier transform "over QQ" using exact N-th roots of unity.
source code
 
idft(self)
Implements a discrete inverse Fourier transform.
source code
 
dct(self)
Implements a discrete Cosine transform...
source code
 
dst(self)
Implements a discrete Sine transform...
source code
 
convolution(self, other)
Convolves two sequences of the same length (automatically expands the shortest one by extending it by 0 if they have different lengths).
source code
 
convolution_periodic(self, other)
Convolves two collections indexed by a range(...) of the same length (automatically expands the shortest one by extending it by 0 if they have different lengths).
source code
 
__mul__(self, other)
Implements scalar multiplication (on the right).
source code
 
__eq__(self, other)
Implements boolean ==.
source code
 
fft(self)
Wraps the gsl FastFourierTransform.forward in fft.pyx.
source code
 
ifft(self)
Implements the gsl FastFourierTransform.inverse in fft.pyx.
source code
 
dwt(self, other='haar', wavelet_k=2)
Wraps the gsl WaveletTransform.forward in dwt.pyx (written by Johua Kantor).
source code
 
idwt(self, other='haar', wavelet_k=2)
Implements the gsl WaveletTransform.backward in dwt.pyx.
source code

Inherited from structure.sage_object.SageObject: __hash__, __new__, __repr__, _axiom_, _axiom_init_, _gap_, _gap_init_, _gp_, _gp_init_, _interface_, _interface_init_, _interface_is_cached_, _kash_, _kash_init_, _macaulay2_, _macaulay2_init_, _magma_, _magma_init_, _maple_, _maple_init_, _mathematica_, _mathematica_init_, _maxima_, _maxima_init_, _octave_, _octave_init_, _pari_, _pari_init_, _r_init_, _sage_, _singular_, _singular_init_, category, db, dump, dumps, rename, reset_name, save, version

Inherited from object: __delattr__, __getattribute__, __reduce__, __reduce_ex__, __setattr__, __str__

Properties [hide private]

Inherited from object: __class__

Method Details [hide private]

__init__(self, L, index_object)
(Constructor)

source code 

\code{index_object} must be a SAGE object with an _iter_ method
containing the same number of elements as self, which is a
list of elements taken from a field.

EXAMPLES:
    sage: J = range(10)
    sage: A = [1/10 for j in J]
    sage: s = IndexedSequence(A,J)
    sage: s
    Indexed sequence: [1/10, 1/10, 1/10, 1/10, 1/10, 1/10, 1/10, 1/10, 1/10, 1/10]
        indexed by [0, 1, 2, 3, 4, 5, 6, 7, 8, 9]
    sage: s.dict()
    {0: 1/10,
     1: 1/10,
     2: 1/10,
     3: 1/10,
     4: 1/10,
     5: 1/10,
     6: 1/10,
     7: 1/10,
     8: 1/10,
     9: 1/10}
    sage: s.list()
    [1/10, 1/10, 1/10, 1/10, 1/10, 1/10, 1/10, 1/10, 1/10, 1/10]
    sage: s.index_object()
    [0, 1, 2, 3, 4, 5, 6, 7, 8, 9]
    sage: s.base_ring()
    Rational Field

Overrides: object.__init__

base_ring(self)

source code 

This just returns the common parent R of the N list
elements. In some applications (say, when computing the
discrete Fourier transform, dft), it is more accurate to think
of the base_ring as the group ring QQ(zeta_N)[R].

_repr_(self)

source code 

Implements print method.

EXAMPLES:
    sage: A = [ZZ(i) for i in range(3)]
    sage: I = range(3)
    sage: s = IndexedSequence(A,I)
    sage: s
    Indexed sequence: [0, 1, 2]
     indexed by [0, 1, 2]
    sage: print s
    Indexed sequence: [0, 1, 2]
     indexed by [0, 1, 2]
    sage: I = GF(3)
    sage: A = [i^2 for i in I]
    sage: s = IndexedSequence(A,I)
    sage: s
    Indexed sequence: [0, 1, 1]
     indexed by Finite Field of size 3

plot_histogram(self)

source code 

Plots the histogram plot of the sequence, which is assumed to be real
or from a finite field, with a real indexing set I coercible into RR.

EXAMPLES:
    sage: J = range(3)
    sage: A = [ZZ(i^2)+1 for i in J]
    sage: s = IndexedSequence(A,J)
    sage: P = s.plot_histogram()

Now type show(P) to view this in a browser.

plot(self)

source code 

Plots the points of the sequence, whose elements are assumed
to be real or from a finite field, with a real indexing set I
= range(len(self)).

EXAMPLES:
    sage: I = range(3)
    sage: A = [ZZ(i^2)+1 for i in I]
    sage: s = IndexedSequence(A,I)
    sage: P = s.plot()

Now type show(P) to view this in a browser.

Overrides: structure.sage_object.SageObject.plot

dft(self, chi=<function <lambda> at 0x4f22848>)

source code 

Implements a discrete Fourier transform "over QQ" using exact
N-th roots of unity.

EXAMPLES:
    sage: J = range(6)
    sage: A = [ZZ(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: s.dft(lambda x:x^2)
    Indexed sequence: [6, 0, 0, 6, 0, 0]
     indexed by [0, 1, 2, 3, 4, 5]
    sage: s.dft()
    Indexed sequence: [6, 0, 0, 0, 0, 0]
     indexed by [0, 1, 2, 3, 4, 5]
    sage: G = SymmetricGroup(3)
    sage: J = G.conjugacy_classes_representatives()
    sage: s = IndexedSequence([1,2,3],J) # 1,2,3 are the values of a class fcn on G
    sage: s.dft()   # the "scalar-valued Fourier transform" of this class fcn
    Indexed sequence: [8, 2, 2]
     indexed by [(), (1,2), (1,2,3)]
    sage: J = AbelianGroup(2,[2,3],names='ab')
    sage: s = IndexedSequence([1,2,3,4,5,6],J)
    sage: s.dft()   # the precision of output is somewhat random and architecture dependent.
    Indexed sequence: [21.0000000000000, -2.99999999999997 - 1.73205080756885*I, -2.99999999999999 + 1.73205080756888*I, -9.00000000000000 + 0.0000000000000485744257349999*I, -0.00000000000000976996261670137 - 0.0000000000000159872115546022*I, -0.00000000000000621724893790087 - 0.0000000000000106581410364015*I]                        
        indexed by Multiplicative Abelian Group isomorphic to C2 x C3
    sage: J = CyclicPermutationGroup(6)
    sage: s = IndexedSequence([1,2,3,4,5,6],J)
    sage: s.dft()   # the precision of output is somewhat random and architecture dependent.
    Indexed sequence: [21.0000000000000, -2.99999999999997 - 1.73205080756885*I, -2.99999999999999 + 1.73205080756888*I, -9.00000000000000 + 0.0000000000000485744257349999*I, -0.00000000000000976996261670137 - 0.0000000000000159872115546022*I, -0.00000000000000621724893790087 - 0.0000000000000106581410364015*I]
        indexed by Cyclic group of order 6 as a permutation group
    sage: p = 7; J = range(p); A = [kronecker_symbol(j,p) for j in J]
    sage: s = IndexedSequence(A,J)
    sage: Fs = s.dft()
    sage: c = Fs.list()[1]; [x/c for x in Fs.list()]; s.list()
    [0, 1, 1, -1, 1, -1, -1]
    [0, 1, 1, -1, 1, -1, -1]            
    
The DFT of the values of the quadratic residue symbol is itself, up to
a constant factor (denoted c on the last line above).

TODO: Read the parent of the elements of S; if QQ or CC leave as
is; if AbelianGroup, use abelian_group_dual; if some other
implemented Group (permutation, matrix), call .characters()
and test if the index list is the set of conjugacy classes.

idft(self)

source code 

Implements a discrete inverse Fourier transform. Only works over QQ.

EXAMPLES:
    sage: J = range(5)
    sage: A = [ZZ(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: fs = s.dft(); fs
    Indexed sequence: [5, 0, 0, 0, 0]
        indexed by [0, 1, 2, 3, 4]
    sage: it = fs.idft(); it
    Indexed sequence: [1, 1, 1, 1, 1]
        indexed by [0, 1, 2, 3, 4]
    sage: it == s
    True

dct(self)

source code 

Implements a discrete Cosine transform

EXAMPLES:
    sage: J = range(5)
    sage: A = [exp(-2*pi*i*I/5) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: s.dct()    # discrete cosine   (random low bits)
    Indexed sequence: [2.50000000000000 - 0.000000000000000111022302462516*I, 2.50000000000000 - 0.000000000000000111022302462516*I, 2.50000000000000 - 0.000000000000000111022302462516*I, 2.50000000000000 - 0.000000000000000111022302462516*I, 2.50000000000000 - 0.000000000000000111022302462516*I]
       indexed by [0, 1, 2, 3, 4]

dst(self)

source code 

Implements a discrete Sine transform

EXAMPLES:
    sage: J = range(5)
    sage: I = CC.0; pi = CC(pi)
    sage: A = [exp(-2*pi*i*I/5) for i in J]
    sage: s = IndexedSequence(A,J)

    sage: s.dst()        # discrete sine
    Indexed sequence: [1.11022302462516e-16 - 2.50000000000000*I, 1.11022302462516e-16 - 2.50000000000000*I, 1.11022302462516e-16 - 2.50000000000000*I, 1.11022302462516e-16 - 2.50000000000000*I, 1.11022302462516e-16 - 2.50000000000000*I]
        indexed by [0, 1, 2, 3, 4]

convolution(self, other)

source code 

Convolves two sequences of the same length (automatically expands
the shortest one by extending it by 0 if they have different lengths).
If {a_n} and {b_n} are sequences of length N (n=0,1,...,N-1), extended 
by zero for all n in ZZ, then the convolution is

         c_j = \sum_{i=0}^{N-1} a_ib_{j-i}.

INPUT: 
    self, other   --  a collection of elements of a ring with
                      index set a finite abelian group (under +)

OUTPUT:
    self*other -- the Dirichlet convolution

EXAMPLES:
    sage: J = range(5)
    sage: A = [ZZ(1) for i in J]
    sage: B = [ZZ(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: t = IndexedSequence(B,J)
    sage: s.convolution(t)
    [1, 2, 3, 4, 5, 4, 3, 2, 1]

AUTHOR: David Joyner (9-2006)

convolution_periodic(self, other)

source code 

Convolves two collections indexed by a range(...) of the same length (automatically 
expands the shortest one by extending it by 0 if they have different lengths).
If {a_n} and {b_n} are sequences of length N (n=0,1,...,N-1), extended 
periodically for all n in ZZ, then the convolution is

         c_j = \sum_{i=0}^{N-1} a_ib_{j-i}.

INPUT: 
    self, other   --  a sequence of elements of CC, RR or GF(q)

OUTPUT:
    self*other -- the Dirichlet convolution

EXAMPLES:
    sage: I = range(5)
    sage: A = [ZZ(1) for i in I]
    sage: B = [ZZ(1) for i in I]
    sage: s = IndexedSequence(A,I)
    sage: t = IndexedSequence(B,I)       
    sage: s.convolution_periodic(t)
    [5, 5, 5, 5, 5, 5, 5, 5, 5]

AUTHOR: David Joyner (9-2006)

__mul__(self, other)

source code 

Implements scalar multiplication (on the right).

EXAMPLES:
    sage: J = range(5)
    sage: A = [ZZ(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: s.base_ring()
    Integer Ring
    sage: t = s*(1/3); t; t.base_ring()
    Indexed sequence: [1/3, 1/3, 1/3, 1/3, 1/3]
        indexed by [0, 1, 2, 3, 4]
    Rational Field            

__eq__(self, other)
(Equality operator)

source code 

Implements boolean ==.

EXAMPLES:
    sage: J = range(5)
    sage: A = [ZZ(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: t = s*(1/3)
    sage: t*3==s
    1

WARNING: ** elements are considered different if they differ
by 10^(-8), which is pretty arbitrary -- use with CAUTION!! **

fft(self)

source code 

Wraps the gsl FastFourierTransform.forward in fft.pyx.  If the
length is a power of 2 then this automatically uses the radix2
method. If the number of sample points in the input is a power
of 2 then the wrapper for the GSL function
gsl_fft_complex_radix2_forward is automatically called.
Otherwise, gsl_fft_complex_forward is used.

EXAMPLES:
    sage: J = range(5)
    sage: A = [RR(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: t = s.fft(); t
    Indexed sequence: [5.00000000000000, 0, 0, 0, 0]
     indexed by [0, 1, 2, 3, 4]

ifft(self)

source code 

Implements the gsl FastFourierTransform.inverse in fft.pyx.
If the number of sample points in the input is a power of 2 
then the wrapper for the GSL function 
gsl_fft_complex_radix2_inverse is automatically called.
Otherwise, gsl_fft_complex_inverse is used.

EXAMPLES:
    sage: J = range(5)
    sage: A = [RR(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: t = s.fft(); t
    Indexed sequence: [5.00000000000000, 0, 0, 0, 0]
       indexed by [0, 1, 2, 3, 4]
    sage: t.ifft()
    Indexed sequence: [1.00000000000000, 1.00000000000000, 1.00000000000000, 1.00000000000000, 1.00000000000000]
        indexed by [0, 1, 2, 3, 4]            
    sage: t.ifft() == s
    1

dwt(self, other='haar', wavelet_k=2)

source code 

Wraps the gsl WaveletTransform.forward in dwt.pyx (written
by Johua Kantor). Assumes the length of the sample is a power of 2. 
Uses the GSL function gsl_wavelet_transform_forward.

other -- the wavelet_type:   the name of the type of wavelet, 
                             valid choices are: 
                             'daubechies','daubechies_centered',
                             'haar' (default),'haar_centered',
                             'bspline', and 'bspline_centered'. 

wavelet_k -- For daubechies wavelets, wavelet_k specifies a 
             daubechie wavelet with k/2 vanishing moments. 
             k = 4,6,...,20 for k even are the only ones implemented. 
             For Haar wavelets, wavelet_k must be 2. 
             For bspline wavelets, 
             wavelet_k = 103,105,202,204,206,208,301,305, 307,309 
             will give biorthogonal B-spline wavelets of order (i,j) where
             wavelet_k=100*i+j. 

The wavelet transform uses J=log_2(n) levels.

EXAMPLES:
    sage: J = range(8)
    sage: A = [RR(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: t = s.dwt()
    sage: t        # slightly random output
    Indexed sequence: [2.82842712474999, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000]
        indexed by [0, 1, 2, 3, 4, 5, 6, 7]

idwt(self, other='haar', wavelet_k=2)

source code 

Implements the gsl WaveletTransform.backward in dwt.pyx.
other must be an element of 
 {"haar", "daubechies", "daubechies_centered",
  "haar_centered", "bspline", "bspline_centered"}. 
Assumes the length of the sample is a power of 2. Uses the          
GSL function gsl_wavelet_transform_backward.

INPUT:
    other -- the wavelet_type:   the name of the type of wavelet, 
                         valid choices are: 
                         'daubechies','daubechies_centered',
                         'haar' (default),'haar_centered',
                         'bspline', and 'bspline_centered'. 

    wavelet_k -- For daubechies wavelets, wavelet_k specifies a 
                 daubechie wavelet with k/2 vanishing moments. 
                 k = 4,6,...,20 for k even are the only ones implemented. 
                 For Haar wavelets, wavelet_k must be 2. 
                 For bspline wavelets, 
                 wavelet_k = 103,105,202,204,206,208,301,305, 307,309 
                 will give biorthogonal B-spline wavelets of order (i,j) where
                 wavelet_k=100*i+j. 


EXAMPLES:
    sage: J = range(8)
    sage: A = [RR(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: t = s.dwt()
    sage: t            # random arch dependent output
    Indexed sequence: [2.82842712474999, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000]
        indexed by [0, 1, 2, 3, 4, 5, 6, 7]            
    sage: t.idwt()                  # random arch dependent output
    Indexed sequence: [1.00000000000000, 1.00000000000000, 1.00000000000000, 1.00000000000000, 1.00000000000000, 1.00000000000000, 1.00000000000000, 1.00000000000000]
        indexed by [0, 1, 2, 3, 4, 5, 6, 7]            
    sage: t.idwt() == s
    True
    sage: J = range(16)
    sage: A = [RR(1) for i in J]
    sage: s = IndexedSequence(A,J)
    sage: t = s.dwt("bspline", 103)
    sage: t   # random arch dependent output
    Indexed sequence: [4.00000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000]
        indexed by [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15]
    sage: t.idwt("bspline", 103) == s
    True