1# -*- coding: utf-8 -*-
2# Licensed under a 3-clause BSD style license - see LICENSE.rst
3
4
5import numpy as np
6import pytest
7
8from astropy import units as u
9from astropy.coordinates import transformations as t
10from astropy.coordinates.builtin_frames import ICRS, FK5, FK4, FK4NoETerms, Galactic, AltAz, HCRS
11from astropy.coordinates import representation as r
12from astropy.coordinates.baseframe import frame_transform_graph
13from astropy.coordinates.matrix_utilities import rotation_matrix
14from astropy.tests.helper import assert_quantity_allclose as assert_allclose
15from astropy.time import Time
16from astropy.units import allclose as quantity_allclose
17from astropy.utils.exceptions import AstropyWarning
18
19
20# Coordinates just for these tests.
21class TCoo1(ICRS):
22    pass
23
24
25class TCoo2(ICRS):
26    pass
27
28
29class TCoo3(ICRS):
30    pass
31
32
33def test_transform_classes():
34    """
35    Tests the class-based/OO syntax for creating transforms
36    """
37    def tfun(c, f):
38        return f.__class__(ra=c.ra, dec=c.dec)
39
40    _ = t.FunctionTransform(tfun, TCoo1, TCoo2,
41                            register_graph=frame_transform_graph)
42
43    c1 = TCoo1(ra=1*u.radian, dec=0.5*u.radian)
44    c2 = c1.transform_to(TCoo2())
45    assert_allclose(c2.ra.radian, 1)
46    assert_allclose(c2.dec.radian, 0.5)
47
48    def matfunc(coo, fr):
49        return [[1, 0, 0],
50                [0, coo.ra.degree, 0],
51                [0, 0, 1]]
52    trans2 = t.DynamicMatrixTransform(matfunc, TCoo1, TCoo2)
53    trans2.register(frame_transform_graph)
54
55    c3 = TCoo1(ra=1*u.deg, dec=2*u.deg)
56    c4 = c3.transform_to(TCoo2())
57
58    assert_allclose(c4.ra.degree, 1)
59    assert_allclose(c4.ra.degree, 1)
60
61    # be sure to unregister the second one - no need for trans1 because it
62    # already got unregistered when trans2 was created.
63    trans2.unregister(frame_transform_graph)
64
65
66def test_transform_decos():
67    """
68    Tests the decorator syntax for creating transforms
69    """
70    c1 = TCoo1(ra=1*u.deg, dec=2*u.deg)
71
72    @frame_transform_graph.transform(t.FunctionTransform, TCoo1, TCoo2)
73    def trans(coo1, f):
74        return TCoo2(ra=coo1.ra, dec=coo1.dec * 2)
75
76    c2 = c1.transform_to(TCoo2())
77    assert_allclose(c2.ra.degree, 1)
78    assert_allclose(c2.dec.degree, 4)
79
80    c3 = TCoo1(r.CartesianRepresentation(x=1*u.pc, y=1*u.pc, z=2*u.pc))
81
82    @frame_transform_graph.transform(t.StaticMatrixTransform, TCoo1, TCoo2)
83    def matrix():
84        return [[2, 0, 0],
85                [0, 1, 0],
86                [0, 0, 1]]
87
88    c4 = c3.transform_to(TCoo2())
89
90    assert_allclose(c4.cartesian.x, 2*u.pc)
91    assert_allclose(c4.cartesian.y, 1*u.pc)
92    assert_allclose(c4.cartesian.z, 2*u.pc)
93
94
95def test_shortest_path():
96    class FakeTransform:
97        def __init__(self, pri):
98            self.priority = pri
99
100    g = t.TransformGraph()
101
102    # cheating by adding graph elements directly that are not classes - the
103    # graphing algorithm still works fine with integers - it just isn't a valid
104    # TransformGraph
105
106    # the graph looks is a down-going diamond graph with the lower-right slightly
107    # heavier and a cycle from the bottom to the top
108    # also, a pair of nodes isolated from 1
109
110    g._graph[1][2] = FakeTransform(1)
111    g._graph[1][3] = FakeTransform(1)
112    g._graph[2][4] = FakeTransform(1)
113    g._graph[3][4] = FakeTransform(2)
114    g._graph[4][1] = FakeTransform(5)
115
116    g._graph[5][6] = FakeTransform(1)
117
118    path, d = g.find_shortest_path(1, 2)
119    assert path == [1, 2]
120    assert d == 1
121    path, d = g.find_shortest_path(1, 3)
122    assert path == [1, 3]
123    assert d == 1
124    path, d = g.find_shortest_path(1, 4)
125    print('Cached paths:', g._shortestpaths)
126    assert path == [1, 2, 4]
127    assert d == 2
128
129    # unreachable
130    path, d = g.find_shortest_path(1, 5)
131    assert path is None
132    assert d == float('inf')
133
134    path, d = g.find_shortest_path(5, 6)
135    assert path == [5, 6]
136    assert d == 1
137
138
139def test_sphere_cart():
140    """
141    Tests the spherical <-> cartesian transform functions
142    """
143    from astropy.utils import NumpyRNGContext
144    from astropy.coordinates import spherical_to_cartesian, cartesian_to_spherical
145
146    x, y, z = spherical_to_cartesian(1, 0, 0)
147    assert_allclose(x, 1)
148    assert_allclose(y, 0)
149    assert_allclose(z, 0)
150
151    x, y, z = spherical_to_cartesian(0, 1, 1)
152    assert_allclose(x, 0)
153    assert_allclose(y, 0)
154    assert_allclose(z, 0)
155
156    x, y, z = spherical_to_cartesian(5, 0, np.arcsin(4. / 5.))
157    assert_allclose(x, 3)
158    assert_allclose(y, 4)
159    assert_allclose(z, 0)
160
161    r, lat, lon = cartesian_to_spherical(0, 1, 0)
162    assert_allclose(r, 1)
163    assert_allclose(lat, 0 * u.deg)
164    assert_allclose(lon, np.pi / 2 * u.rad)
165
166    # test round-tripping
167    with NumpyRNGContext(13579):
168        x, y, z = np.random.randn(3, 5)
169
170    r, lat, lon = cartesian_to_spherical(x, y, z)
171    x2, y2, z2 = spherical_to_cartesian(r, lat, lon)
172
173    assert_allclose(x, x2)
174    assert_allclose(y, y2)
175    assert_allclose(z, z2)
176
177
178def test_transform_path_pri():
179    """
180    This checks that the transformation path prioritization works by
181    making sure the ICRS -> Gal transformation always goes through FK5
182    and not FK4.
183    """
184    frame_transform_graph.invalidate_cache()
185    tpath, td = frame_transform_graph.find_shortest_path(ICRS, Galactic)
186    assert tpath == [ICRS, FK5, Galactic]
187    assert td == 2
188
189    # but direct from FK4 to Galactic should still be possible
190    tpath, td = frame_transform_graph.find_shortest_path(FK4, Galactic)
191    assert tpath == [FK4, FK4NoETerms, Galactic]
192    assert td == 2
193
194
195def test_obstime():
196    """
197    Checks to make sure observation time is
198    accounted for at least in FK4 <-> ICRS transformations
199    """
200    b1950 = Time('B1950')
201    j1975 = Time('J1975')
202
203    fk4_50 = FK4(ra=1*u.deg, dec=2*u.deg, obstime=b1950)
204    fk4_75 = FK4(ra=1*u.deg, dec=2*u.deg, obstime=j1975)
205
206    icrs_50 = fk4_50.transform_to(ICRS())
207    icrs_75 = fk4_75.transform_to(ICRS())
208
209    # now check that the resulting coordinates are *different* - they should be,
210    # because the obstime is different
211    assert icrs_50.ra.degree != icrs_75.ra.degree
212    assert icrs_50.dec.degree != icrs_75.dec.degree
213
214# ------------------------------------------------------------------------------
215# Affine transform tests and helpers:
216
217# just acting as a namespace
218
219
220class transfunc:
221    rep = r.CartesianRepresentation(np.arange(3)*u.pc)
222    dif = r.CartesianDifferential(*np.arange(3, 6)*u.pc/u.Myr)
223    rep0 = r.CartesianRepresentation(np.zeros(3)*u.pc)
224
225    @classmethod
226    def both(cls, coo, fr):
227        # exchange x <-> z and offset
228        M = np.array([[0., 0., 1.],
229                      [0., 1., 0.],
230                      [1., 0., 0.]])
231        return M, cls.rep.with_differentials(cls.dif)
232
233    @classmethod
234    def just_matrix(cls, coo, fr):
235        # exchange x <-> z and offset
236        M = np.array([[0., 0., 1.],
237                      [0., 1., 0.],
238                      [1., 0., 0.]])
239        return M, None
240
241    @classmethod
242    def no_matrix(cls, coo, fr):
243        return None, cls.rep.with_differentials(cls.dif)
244
245    @classmethod
246    def no_pos(cls, coo, fr):
247        return None, cls.rep0.with_differentials(cls.dif)
248
249    @classmethod
250    def no_vel(cls, coo, fr):
251        return None, cls.rep
252
253
254@pytest.mark.parametrize('transfunc', [transfunc.both, transfunc.no_matrix,
255                                       transfunc.no_pos, transfunc.no_vel,
256                                       transfunc.just_matrix])
257@pytest.mark.parametrize('rep', [
258    r.CartesianRepresentation(5, 6, 7, unit=u.pc),
259    r.CartesianRepresentation(5, 6, 7, unit=u.pc,
260                              differentials=r.CartesianDifferential(8, 9, 10,
261                                                                    unit=u.pc/u.Myr)),
262    r.CartesianRepresentation(5, 6, 7, unit=u.pc,
263                              differentials=r.CartesianDifferential(8, 9, 10,
264                                                                    unit=u.pc/u.Myr))
265     .represent_as(r.CylindricalRepresentation, r.CylindricalDifferential)
266])
267def test_affine_transform_succeed(transfunc, rep):
268    c = TCoo1(rep)
269
270    # compute expected output
271    M, offset = transfunc(c, TCoo2)
272
273    _rep = rep.to_cartesian()
274    diffs = dict([(k, diff.represent_as(r.CartesianDifferential, rep))
275                  for k, diff in rep.differentials.items()])
276    expected_rep = _rep.with_differentials(diffs)
277
278    if M is not None:
279        expected_rep = expected_rep.transform(M)
280
281    expected_pos = expected_rep.without_differentials()
282    if offset is not None:
283        expected_pos = expected_pos + offset.without_differentials()
284
285    expected_vel = None
286    if c.data.differentials:
287        expected_vel = expected_rep.differentials['s']
288
289        if offset and offset.differentials:
290            expected_vel = (expected_vel + offset.differentials['s'])
291
292    # register and do the transformation and check against expected
293    trans = t.AffineTransform(transfunc, TCoo1, TCoo2)
294    trans.register(frame_transform_graph)
295
296    c2 = c.transform_to(TCoo2())
297
298    assert quantity_allclose(c2.data.to_cartesian().xyz,
299                             expected_pos.to_cartesian().xyz)
300
301    if expected_vel is not None:
302        diff = c2.data.differentials['s'].to_cartesian(base=c2.data)
303        assert quantity_allclose(diff.xyz, expected_vel.d_xyz)
304
305    trans.unregister(frame_transform_graph)
306
307
308# these should fail
309def transfunc_invalid_matrix(coo, fr):
310    return np.eye(4), None
311
312# Leaving this open in case we want to add more functions to check for failures
313
314
315@pytest.mark.parametrize('transfunc', [transfunc_invalid_matrix])
316def test_affine_transform_fail(transfunc):
317    diff = r.CartesianDifferential(8, 9, 10, unit=u.pc/u.Myr)
318    rep = r.CartesianRepresentation(5, 6, 7, unit=u.pc, differentials=diff)
319    c = TCoo1(rep)
320
321    # register and do the transformation and check against expected
322    trans = t.AffineTransform(transfunc, TCoo1, TCoo2)
323    trans.register(frame_transform_graph)
324
325    with pytest.raises(ValueError):
326        c.transform_to(TCoo2())
327
328    trans.unregister(frame_transform_graph)
329
330
331def test_too_many_differentials():
332    dif1 = r.CartesianDifferential(*np.arange(3, 6)*u.pc/u.Myr)
333    dif2 = r.CartesianDifferential(*np.arange(3, 6)*u.pc/u.Myr**2)
334    rep = r.CartesianRepresentation(np.arange(3)*u.pc,
335                                    differentials={'s': dif1, 's2': dif2})
336
337    with pytest.raises(ValueError):
338        c = TCoo1(rep)
339
340    # register and do the transformation and check against expected
341    trans = t.AffineTransform(transfunc.both, TCoo1, TCoo2)
342    trans.register(frame_transform_graph)
343
344    # Check that if frame somehow gets through to transformation, multiple
345    # differentials are caught
346    c = TCoo1(rep.without_differentials())
347    c._data = c._data.with_differentials({'s': dif1, 's2': dif2})
348    with pytest.raises(ValueError):
349        c.transform_to(TCoo2())
350
351    trans.unregister(frame_transform_graph)
352
353# A matrix transform of a unit spherical with differentials should work
354
355
356@pytest.mark.parametrize('rep', [
357    r.UnitSphericalRepresentation(lon=15*u.degree, lat=-11*u.degree,
358        differentials=r.SphericalDifferential(d_lon=15*u.mas/u.yr,
359                                              d_lat=11*u.mas/u.yr,
360                                              d_distance=-110*u.km/u.s)),
361    r.UnitSphericalRepresentation(lon=15*u.degree, lat=-11*u.degree,
362        differentials={'s': r.RadialDifferential(d_distance=-110*u.km/u.s)}),
363    r.SphericalRepresentation(lon=15*u.degree, lat=-11*u.degree,
364                              distance=150*u.pc,
365        differentials={'s': r.RadialDifferential(d_distance=-110*u.km/u.s)})
366])
367def test_unit_spherical_with_differentials(rep):
368
369    c = TCoo1(rep)
370
371    # register and do the transformation and check against expected
372    trans = t.AffineTransform(transfunc.just_matrix, TCoo1, TCoo2)
373    trans.register(frame_transform_graph)
374    c2 = c.transform_to(TCoo2())
375
376    assert 's' in rep.differentials
377    assert isinstance(c2.data.differentials['s'],
378                      rep.differentials['s'].__class__)
379
380    if isinstance(rep.differentials['s'], r.RadialDifferential):
381        assert c2.data.differentials['s'] is rep.differentials['s']
382
383    trans.unregister(frame_transform_graph)
384
385    # should fail if we have to do offsets
386    trans = t.AffineTransform(transfunc.both, TCoo1, TCoo2)
387    trans.register(frame_transform_graph)
388
389    with pytest.raises(TypeError):
390        c.transform_to(TCoo2())
391
392    trans.unregister(frame_transform_graph)
393
394
395def test_vel_transformation_obstime_err():
396    # TODO: replace after a final decision on PR #6280
397    from astropy.coordinates.sites import get_builtin_sites
398
399    diff = r.CartesianDifferential([.1, .2, .3]*u.km/u.s)
400    rep = r.CartesianRepresentation([1, 2, 3]*u.au, differentials=diff)
401
402    loc = get_builtin_sites()['example_site']
403
404    aaf = AltAz(obstime='J2010', location=loc)
405    aaf2 = AltAz(obstime=aaf.obstime + 3*u.day, location=loc)
406    aaf3 = AltAz(obstime=aaf.obstime + np.arange(3)*u.day, location=loc)
407    aaf4 = AltAz(obstime=aaf.obstime, location=loc)
408
409    aa = aaf.realize_frame(rep)
410
411    with pytest.raises(NotImplementedError) as exc:
412        aa.transform_to(aaf2)
413    assert 'cannot transform' in exc.value.args[0]
414
415    with pytest.raises(NotImplementedError) as exc:
416        aa.transform_to(aaf3)
417    assert 'cannot transform' in exc.value.args[0]
418
419    aa.transform_to(aaf4)
420
421    aa.transform_to(ICRS())
422
423
424def test_function_transform_with_differentials():
425    def tfun(c, f):
426        return f.__class__(ra=c.ra, dec=c.dec)
427
428    _ = t.FunctionTransform(tfun, TCoo3, TCoo2,
429                            register_graph=frame_transform_graph)
430
431    t3 = TCoo3(ra=1*u.deg, dec=2*u.deg, pm_ra_cosdec=1*u.marcsec/u.yr,
432               pm_dec=1*u.marcsec/u.yr,)
433
434    with pytest.warns(AstropyWarning, match=r'.*they have been dropped.*') as w:
435        t3.transform_to(TCoo2())
436    assert len(w) == 1
437
438
439def test_frame_override_component_with_attribute():
440    """
441    It was previously possible to define a frame with an attribute with the
442    same name as a component. We don't want to allow this!
443    """
444    from astropy.coordinates.baseframe import BaseCoordinateFrame
445    from astropy.coordinates.attributes import Attribute
446
447    class BorkedFrame(BaseCoordinateFrame):
448        ra = Attribute(default=150)
449        dec = Attribute(default=150)
450
451    def trans_func(coo1, f):
452        pass
453
454    trans = t.FunctionTransform(trans_func, BorkedFrame, ICRS)
455    with pytest.raises(ValueError) as exc:
456        trans.register(frame_transform_graph)
457
458    assert ('BorkedFrame' in exc.value.args[0] and
459            "'ra'" in exc.value.args[0] and
460            "'dec'" in exc.value.args[0])
461
462
463def test_static_matrix_combine_paths():
464    """
465    Check that combined staticmatrixtransform matrices provide the same
466    transformation as using an intermediate transformation.
467
468    This is somewhat of a regression test for #7706
469    """
470    from astropy.coordinates.baseframe import BaseCoordinateFrame
471    from astropy.coordinates.matrix_utilities import rotation_matrix
472
473    class AFrame(BaseCoordinateFrame):
474        default_representation = r.SphericalRepresentation
475        default_differential = r.SphericalCosLatDifferential
476
477    t1 = t.StaticMatrixTransform(rotation_matrix(30.*u.deg, 'z'),
478                                 ICRS, AFrame)
479    t1.register(frame_transform_graph)
480    t2 = t.StaticMatrixTransform(rotation_matrix(30.*u.deg, 'z').T,
481                                 AFrame, ICRS)
482    t2.register(frame_transform_graph)
483
484    class BFrame(BaseCoordinateFrame):
485        default_representation = r.SphericalRepresentation
486        default_differential = r.SphericalCosLatDifferential
487
488    t3 = t.StaticMatrixTransform(rotation_matrix(30.*u.deg, 'x'),
489                                 ICRS, BFrame)
490    t3.register(frame_transform_graph)
491    t4 = t.StaticMatrixTransform(rotation_matrix(30.*u.deg, 'x').T,
492                                 BFrame, ICRS)
493    t4.register(frame_transform_graph)
494
495    c = Galactic(123*u.deg, 45*u.deg)
496    c1 = c.transform_to(BFrame())  # direct
497    c2 = c.transform_to(AFrame()).transform_to(BFrame())  # thru A
498    c3 = c.transform_to(ICRS()).transform_to(BFrame())  # thru ICRS
499
500    assert quantity_allclose(c1.lon, c2.lon)
501    assert quantity_allclose(c1.lat, c2.lat)
502
503    assert quantity_allclose(c1.lon, c3.lon)
504    assert quantity_allclose(c1.lat, c3.lat)
505
506    for t_ in [t1, t2, t3, t4]:
507        t_.unregister(frame_transform_graph)
508
509
510def test_multiple_aliases():
511    from astropy.coordinates.baseframe import BaseCoordinateFrame
512
513    # Define a frame with multiple aliases
514    class MultipleAliasesFrame(BaseCoordinateFrame):
515        name = ['alias_1', 'alias_2']
516        default_representation = r.SphericalRepresentation
517
518    def tfun(c, f):
519        return f.__class__(lon=c.lon, lat=c.lat)
520
521    # Register a transform
522    graph = t.TransformGraph()
523    _ = t.FunctionTransform(tfun, MultipleAliasesFrame, MultipleAliasesFrame,
524                            register_graph=graph)
525
526    # Test that both aliases have been added to the transform graph
527    assert graph.lookup_name('alias_1') == MultipleAliasesFrame
528    assert graph.lookup_name('alias_2') == MultipleAliasesFrame
529
530    # Test that both aliases appear in the graphviz DOT format output
531    dotstr = graph.to_dot_graph()
532    assert '`alias_1`\\n`alias_2`' in dotstr
533
534
535def test_remove_transform_and_unregister():
536    def tfun(c, f):
537        f.__class__(ra=c.ra, dec=c.dec)
538
539    # Register transforms
540    graph = t.TransformGraph()
541    ftrans1 = t.FunctionTransform(tfun, TCoo1, TCoo1, register_graph=graph)
542    ftrans2 = t.FunctionTransform(tfun, TCoo2, TCoo2, register_graph=graph)
543    _ = t.FunctionTransform(tfun, TCoo1, TCoo2, register_graph=graph)
544
545    # Confirm that the frames are part of the graph
546    assert TCoo1 in graph.frame_set
547    assert TCoo2 in graph.frame_set
548
549    # Use all three ways to remove a transform
550
551    # Remove the only transform with TCoo2 as the "from" frame
552    ftrans2.unregister(graph)
553    # TCoo2 should still be part of the graph because it is the "to" frame of a transform
554    assert TCoo2 in graph.frame_set
555
556    # Remove the remaining transform that involves TCoo2
557    graph.remove_transform(TCoo1, TCoo2, None)
558    # Now TCoo2 should not be part of the graph
559    assert TCoo2 not in graph.frame_set
560
561    # Remove the remaining  transform that involves TCoo1
562    graph.remove_transform(None, None, ftrans1)
563    # Now TCoo1 should not be part of the graph
564    assert TCoo1 not in graph.frame_set
565
566
567def test_remove_transform_errors():
568    def tfun(c, f):
569        return f.__class__(ra=c.ra, dec=c.dec)
570
571    graph = t.TransformGraph()
572    _ = t.FunctionTransform(tfun, TCoo1, TCoo1, register_graph=graph)
573
574    # Test bad calls to remove_transform
575
576    with pytest.raises(ValueError):
577        graph.remove_transform(None, TCoo1, None)
578
579    with pytest.raises(ValueError):
580        graph.remove_transform(TCoo1, None, None)
581
582    with pytest.raises(ValueError):
583        graph.remove_transform(None, None, None)
584
585    with pytest.raises(ValueError):
586        graph.remove_transform(None, None, 1)
587
588    with pytest.raises(ValueError):
589        graph.remove_transform(TCoo1, TCoo1, 1)
590
591
592def test_impose_finite_difference_dt():
593    class H1(HCRS):
594        pass
595
596    class H2(HCRS):
597        pass
598
599    class H3(HCRS):
600        pass
601
602    graph = t.TransformGraph()
603    tfun = lambda c, f: f.__class__(ra=c.ra, dec=c.dec)
604
605    # Set up a number of transforms with different time steps
606    old_dt = 1*u.min
607    transform1 = t.FunctionTransformWithFiniteDifference(tfun, H1, H1, register_graph=graph,
608                                                         finite_difference_dt=old_dt)
609    transform2 = t.FunctionTransformWithFiniteDifference(tfun, H2, H2, register_graph=graph,
610                                                         finite_difference_dt=old_dt * 2)
611    transform3 = t.FunctionTransformWithFiniteDifference(tfun, H2, H3, register_graph=graph,
612                                                         finite_difference_dt=old_dt * 3)
613
614    # Check that all of the transforms have the same new time step
615    new_dt = 1*u.yr
616    with graph.impose_finite_difference_dt(new_dt):
617        assert transform1.finite_difference_dt == new_dt
618        assert transform2.finite_difference_dt == new_dt
619        assert transform3.finite_difference_dt == new_dt
620
621    # Check that all of the original time steps have been restored
622    assert transform1.finite_difference_dt == old_dt
623    assert transform2.finite_difference_dt == old_dt * 2
624    assert transform3.finite_difference_dt == old_dt * 3
625
626
627@pytest.mark.parametrize("first, second, check",
628                         [((rotation_matrix(30*u.deg), None),
629                           (rotation_matrix(45*u.deg), None),
630                           (rotation_matrix(75*u.deg), None)),
631                          ((rotation_matrix(30*u.deg), r.CartesianRepresentation([1, 0, 0])),
632                           (rotation_matrix(45*u.deg), None),
633                           (rotation_matrix(75*u.deg), r.CartesianRepresentation([1/np.sqrt(2), -1/np.sqrt(2), 0]))),
634                          ((rotation_matrix(30*u.deg), None),
635                           (rotation_matrix(45*u.deg), r.CartesianRepresentation([0, 0, 1])),
636                           (rotation_matrix(75*u.deg), r.CartesianRepresentation([0, 0, 1]))),
637                          ((rotation_matrix(30*u.deg), r.CartesianRepresentation([1, 0, 0])),
638                           (rotation_matrix(45*u.deg), r.CartesianRepresentation([0, 0, 1])),
639                           (rotation_matrix(75*u.deg), r.CartesianRepresentation([1/np.sqrt(2), -1/np.sqrt(2), 1]))),
640                          ((rotation_matrix(30*u.deg), r.CartesianRepresentation([1, 2 ,3])),
641                           (None, r.CartesianRepresentation([4, 5, 6])),
642                           (rotation_matrix(30*u.deg), r.CartesianRepresentation([5, 7, 9]))),
643                          ((None, r.CartesianRepresentation([1, 2, 3])),
644                           (rotation_matrix(45*u.deg), r.CartesianRepresentation([4, 5, 6])),
645                           (rotation_matrix(45*u.deg), r.CartesianRepresentation([3/np.sqrt(2)+4, 1/np.sqrt(2)+5, 9]))),
646                          ((None, r.CartesianRepresentation([1, 2, 3])),
647                           (None, r.CartesianRepresentation([4, 5, 6])),
648                           (None, r.CartesianRepresentation([5, 7, 9]))),
649                          ((rotation_matrix(30*u.deg), r.CartesianRepresentation([1, 0, 0])),
650                           (None, None),
651                           (rotation_matrix(30*u.deg), r.CartesianRepresentation([1, 0, 0]))),
652                          ((None, None),
653                           (rotation_matrix(45*u.deg), r.CartesianRepresentation([0, 0, 1])),
654                           (rotation_matrix(45*u.deg), r.CartesianRepresentation([0, 0, 1]))),
655                          ((None, None),
656                           (None, None),
657                           (None, None))])
658def test_combine_affine_params(first, second, check):
659    result = t._combine_affine_params(first, second)
660    if check[0] is None:
661        assert result[0] is None
662    else:
663        assert_allclose(result[0], check[0])
664    if check[1] is None:
665        assert result[1] is None
666    else:
667        assert_allclose(result[1].xyz, check[1].xyz)
668