File: test_representation_methods.py

package info (click to toggle)
python-astropy 1.3-8~bpo8%2B2
  • links: PTS, VCS
  • area: main
  • in suites: jessie-backports
  • size: 44,292 kB
  • sloc: ansic: 160,360; python: 137,322; sh: 11,493; lex: 7,638; yacc: 4,956; xml: 1,796; makefile: 474; cpp: 364
file content (204 lines) | stat: -rw-r--r-- 9,461 bytes parent folder | download | duplicates (2)
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
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
# Licensed under a 3-clause BSD style license - see LICENSE.rst

# TEST_UNICODE_LITERALS
import numpy as np

from ... import units as u
from .. import SphericalRepresentation, Longitude, Latitude
from ...tests.helper import pytest
from ...utils.compat.numpycompat import NUMPY_LT_1_9
from ...utils.compat.numpy import broadcast_to as np_broadcast_to


class TestManipulation():
    """Manipulation of Representation shapes.

    Checking that attributes are manipulated correctly.

    Even more exhaustive tests are done in time.tests.test_methods
    """

    def setup(self):
        lon = Longitude(np.arange(0, 24, 4), u.hourangle)
        lat = Latitude(np.arange(-90, 91, 30), u.deg)
        # With same-sized arrays
        self.s0 = SphericalRepresentation(
            lon[:, np.newaxis] * np.ones(lat.shape),
            lat * np.ones(lon.shape)[:, np.newaxis],
            np.ones(lon.shape + lat.shape) * u.kpc)
        # With unequal arrays -> these will be broadcast.
        self.s1 = SphericalRepresentation(lon[:, np.newaxis], lat, 1. * u.kpc)
        # For completeness on some tests, also a cartesian one
        self.c0 = self.s0.to_cartesian()

    def test_ravel(self):
        s0_ravel = self.s0.ravel()
        assert type(s0_ravel) is type(self.s0)
        assert s0_ravel.shape == (self.s0.size,)
        assert np.all(s0_ravel.lon == self.s0.lon.ravel())
        assert np.may_share_memory(s0_ravel.lon, self.s0.lon)
        assert np.may_share_memory(s0_ravel.lat, self.s0.lat)
        assert np.may_share_memory(s0_ravel.distance, self.s0.distance)
        # Since s1 was broadcast, ravel needs to make a copy.
        s1_ravel = self.s1.ravel()
        assert type(s1_ravel) is type(self.s1)
        assert s1_ravel.shape == (self.s1.size,)
        assert np.all(s1_ravel.lon == self.s1.lon.ravel())
        assert not np.may_share_memory(s1_ravel.lat, self.s1.lat)

    def test_copy(self):
        s0_copy = self.s0.copy()
        assert s0_copy.shape == self.s0.shape
        assert np.all(s0_copy.lon == self.s0.lon)
        assert np.all(s0_copy.lat == self.s0.lat)
        # Check copy was made of internal data.
        assert not np.may_share_memory(s0_copy.distance, self.s0.distance)

    def test_flatten(self):
        s0_flatten = self.s0.flatten()
        assert s0_flatten.shape == (self.s0.size,)
        assert np.all(s0_flatten.lon == self.s0.lon.flatten())
        # Flatten always copies.
        assert not np.may_share_memory(s0_flatten.distance, self.s0.distance)
        s1_flatten = self.s1.flatten()
        assert s1_flatten.shape == (self.s1.size,)
        assert np.all(s1_flatten.lon == self.s1.lon.flatten())
        assert not np.may_share_memory(s1_flatten.lat, self.s1.lat)

    def test_transpose(self):
        s0_transpose = self.s0.transpose()
        assert s0_transpose.shape == (7, 6)
        assert np.all(s0_transpose.lon == self.s0.lon.transpose())
        assert np.may_share_memory(s0_transpose.distance, self.s0.distance)
        s1_transpose = self.s1.transpose()
        assert s1_transpose.shape == (7, 6)
        assert np.all(s1_transpose.lat == self.s1.lat.transpose())
        assert np.may_share_memory(s1_transpose.lat, self.s1.lat)
        # Only one check on T, since it just calls transpose anyway.
        # Doing it on the CartesianRepresentation just for variety's sake.
        c0_T = self.c0.T
        assert c0_T.shape == (7, 6)
        assert np.all(c0_T.x == self.c0.x.T)
        assert np.may_share_memory(c0_T.y, self.c0.y)

    def test_diagonal(self):
        s0_diagonal = self.s0.diagonal()
        assert s0_diagonal.shape == (6,)
        assert np.all(s0_diagonal.lat == self.s0.lat.diagonal())
        if not NUMPY_LT_1_9:
            assert np.may_share_memory(s0_diagonal.lat, self.s0.lat)

    def test_swapaxes(self):
        s1_swapaxes = self.s1.swapaxes(0, 1)
        assert s1_swapaxes.shape == (7, 6)
        assert np.all(s1_swapaxes.lat == self.s1.lat.swapaxes(0, 1))
        assert np.may_share_memory(s1_swapaxes.lat, self.s1.lat)

    def test_reshape(self):
        s0_reshape = self.s0.reshape(2, 3, 7)
        assert s0_reshape.shape == (2, 3, 7)
        assert np.all(s0_reshape.lon == self.s0.lon.reshape(2, 3, 7))
        assert np.all(s0_reshape.lat == self.s0.lat.reshape(2, 3, 7))
        assert np.all(s0_reshape.distance == self.s0.distance.reshape(2, 3, 7))
        assert np.may_share_memory(s0_reshape.lon, self.s0.lon)
        assert np.may_share_memory(s0_reshape.lat, self.s0.lat)
        assert np.may_share_memory(s0_reshape.distance, self.s0.distance)
        s1_reshape = self.s1.reshape(3, 2, 7)
        assert s1_reshape.shape == (3, 2, 7)
        assert np.all(s1_reshape.lat == self.s1.lat.reshape(3, 2, 7))
        assert np.may_share_memory(s1_reshape.lat, self.s1.lat)
        # For reshape(3, 14), copying is necessary for lon, lat, but not for d
        s1_reshape2 = self.s1.reshape(3, 14)
        assert s1_reshape2.shape == (3, 14)
        assert np.all(s1_reshape2.lon == self.s1.lon.reshape(3, 14))
        assert not np.may_share_memory(s1_reshape2.lon, self.s1.lon)
        assert s1_reshape2.distance.shape == (3, 14)
        assert np.may_share_memory(s1_reshape2.distance, self.s1.distance)

    def test_shape_setting(self):
        # Shape-setting should be on the object itself, since copying removes
        # zero-strides due to broadcasting.  We reset the objects at the end.
        self.s0.shape = (2, 3, 7)
        assert self.s0.shape == (2, 3, 7)
        assert self.s0.lon.shape == (2, 3, 7)
        assert self.s0.lat.shape == (2, 3, 7)
        assert self.s0.distance.shape == (2, 3, 7)
        # this works with the broadcasting.
        self.s1.shape = (2, 3, 7)
        assert self.s1.shape == (2, 3, 7)
        assert self.s1.lon.shape == (2, 3, 7)
        assert self.s1.lat.shape == (2, 3, 7)
        assert self.s1.distance.shape == (2, 3, 7)
        assert self.s1.distance.strides == (0, 0, 0)
        # but this one does not.
        oldshape = self.s1.shape
        with pytest.raises(AttributeError):
            self.s1.shape = (42,)
        assert self.s1.shape == oldshape
        assert self.s1.lon.shape == oldshape
        assert self.s1.lat.shape == oldshape
        assert self.s1.distance.shape == oldshape
        # Finally, a more complicated one that checks that things get reset
        # properly if it is not the first component that fails.
        s2 = SphericalRepresentation(self.s1.lon.copy(), self.s1.lat,
                                     self.s1.distance, copy=False)
        assert 0 not in s2.lon.strides
        assert 0 in s2.lat.strides
        with pytest.raises(AttributeError):
            s2.shape = (42,)
        assert s2.shape == oldshape
        assert s2.lon.shape == oldshape
        assert s2.lat.shape == oldshape
        assert s2.distance.shape == oldshape
        assert 0 not in s2.lon.strides
        assert 0 in s2.lat.strides
        self.setup()

    def test_squeeze(self):
        s0_squeeze = self.s0.reshape(3, 1, 2, 1, 7).squeeze()
        assert s0_squeeze.shape == (3, 2, 7)
        assert np.all(s0_squeeze.lat == self.s0.lat.reshape(3, 2, 7))
        assert np.may_share_memory(s0_squeeze.lat, self.s0.lat)

    def test_add_dimension(self):
        s0_adddim = self.s0[:, np.newaxis, :]
        assert s0_adddim.shape == (6, 1, 7)
        assert np.all(s0_adddim.lon == self.s0.lon[:, np.newaxis, :])
        assert np.may_share_memory(s0_adddim.lat, self.s0.lat)

    def test_take(self):
        s0_take = self.s0.take((5, 2))
        assert s0_take.shape == (2,)
        assert np.all(s0_take.lon == self.s0.lon.take((5, 2)))

    def test_broadcast_to(self):
        s0_broadcast = self.s0._apply(np_broadcast_to, (3, 6, 7), subok=True)
        assert type(s0_broadcast) is type(self.s0)
        assert s0_broadcast.shape == (3, 6, 7)
        assert np.all(s0_broadcast.lon == self.s0.lon)
        assert np.all(s0_broadcast.lat == self.s0.lat)
        assert np.all(s0_broadcast.distance == self.s0.distance)
        assert np.may_share_memory(s0_broadcast.lon, self.s0.lon)
        assert np.may_share_memory(s0_broadcast.lat, self.s0.lat)
        assert np.may_share_memory(s0_broadcast.distance, self.s0.distance)
        s1_broadcast = self.s1._apply(np_broadcast_to, shape=(3, 6, 7),
                                      subok=True)
        assert s1_broadcast.shape == (3, 6, 7)
        assert np.all(s1_broadcast.lat == self.s1.lat)
        assert np.all(s1_broadcast.lon == self.s1.lon)
        assert np.all(s1_broadcast.distance == self.s1.distance)
        assert s1_broadcast.distance.shape == (3, 6, 7)
        assert np.may_share_memory(s1_broadcast.lat, self.s1.lat)
        assert np.may_share_memory(s1_broadcast.lon, self.s1.lon)
        assert np.may_share_memory(s1_broadcast.distance, self.s1.distance)

        # A final test that "may_share_memory" equals "does_share_memory"
        # Do this on a copy, to keep self.s0 unchanged.
        sc = self.s0.copy()
        assert not np.may_share_memory(sc.lon, self.s0.lon)
        assert not np.may_share_memory(sc.lat, self.s0.lat)
        sc_broadcast = sc._apply(np_broadcast_to, (3, 6, 7), subok=True)
        assert np.may_share_memory(sc_broadcast.lon, sc.lon)
        # Can only write to copy, not to broadcast version.
        sc.lon[0, 0] = 22. * u.hourangle
        assert np.all(sc_broadcast.lon[:, 0, 0] == 22. * u.hourangle)