def test_transform_with_sun_center(): sun_center = SkyCoord(0 * u.deg, 0 * u.deg, 0 * u.AU, frame=HeliographicStonyhurst(obstime="2001-01-01")) with transform_with_sun_center(): result1 = sun_center.transform_to( HeliographicStonyhurst(obstime="2001-02-01")) # The coordinate should stay pointing at Sun center assert_quantity_allclose(result1.lon, sun_center.lon) assert_quantity_allclose(result1.lat, sun_center.lat) assert_quantity_allclose(result1.radius, sun_center.radius) other = SkyCoord(10 * u.deg, 20 * u.deg, 1 * u.AU, frame=HeliographicStonyhurst(obstime="2001-01-01")) with transform_with_sun_center(): result2 = other.transform_to( HeliographicCarrington(observer='earth', obstime="2001-02-01")) # The coordinate should stay at the same latitude and the same distance from Sun center assert_quantity_allclose(result2.lat, other.lat) assert_quantity_allclose(result2.radius, other.radius)
def test_transform_with_sun_center_reset(): # This test sequence ensures that the context manager resets propoerly sun_center = SkyCoord(0*u.deg, 0*u.deg, 0*u.AU, frame=HeliographicStonyhurst(obstime="2001-01-01")) end_frame = HeliocentricInertial(obstime="2001-02-01") # Without the context manager, the coordinate should not point at Sun center result1 = sun_center.transform_to(end_frame) assert result1.lon != sun_center.lon assert result1.lat != sun_center.lat assert result1.distance != sun_center.radius # Using the context manager, the coordinate should point at Sun center with transform_with_sun_center(): result2 = sun_center.transform_to(end_frame) assert_quantity_allclose(result2.lon, sun_center.lon) assert_quantity_allclose(result2.lat, sun_center.lat) assert_quantity_allclose(result2.distance, sun_center.radius) # Exiting a nested context manager should not affect the outer context manager with transform_with_sun_center(): with transform_with_sun_center(): pass result2a = sun_center.transform_to(end_frame) assert_quantity_allclose(result2a.lon, result2.lon) assert_quantity_allclose(result2a.lat, result2.lat) assert_quantity_allclose(result2a.distance, result2.distance) # After the context manager, the coordinate should have the same result as the first transform result3 = sun_center.transform_to(end_frame) assert_quantity_allclose(result3.lon, result1.lon) assert_quantity_allclose(result3.lat, result1.lat) assert_quantity_allclose(result3.distance, result1.distance)
def test_consistency_with_rotatedsunframe(): old_observer = frames.HeliographicStonyhurst(10 * u.deg, 20 * u.deg, 1 * u.AU, obstime='2001-01-01') new_observer = frames.HeliographicStonyhurst(30 * u.deg, 40 * u.deg, 2 * u.AU, obstime='2001-01-08') hpc_coord = SkyCoord(100 * u.arcsec, 200 * u.arcsec, frame='helioprojective', observer=old_observer, obstime=old_observer.obstime) # Perform the differential rotation using solar_rotate_coordinate() result1 = solar_rotate_coordinate(hpc_coord, observer=new_observer) # Perform the differential rotation using RotatedSunFrame, with translational motion of the Sun # ignored using transform_with_sun_center() rsf_coord = RotatedSunFrame(base=hpc_coord, rotated_time=new_observer.obstime) with transform_with_sun_center(): result2 = rsf_coord.transform_to(result1.replicate_without_data()) assert_quantity_allclose(result1.Tx, result2.Tx) assert_quantity_allclose(result1.Ty, result2.Ty) assert_quantity_allclose(result1.distance, result2.distance)
def _warp_sun_coordinates(xy, smap, new_observer, **diff_rot_kwargs): """ This function takes pixel coordinates in the warped image (`xy`) and calculates the pixel locations of those pixels in the map. To do this it converts the input pixel coordinates to helioprojective coordinates as seen by new_observer, then transforms them to heliographic Stonyhurst, adds the differential rotation correction and then transforms them back to helioprojective coordinates as seen by the map observer and then calculates their corresponding pixel coordinates in the input map. This is an inverse function needed by `skimage.transform.warp`. Parameters ---------- xy : `numpy.ndarray` Pixel coordinates in the warped image. smap : `~sunpy.map.GenericMap` Original map that we want to transform. Returns ------- xy2 : `numpy.ndarray` Pixel coordinates in the map corresponding to the input pixels in the warped image. Notes ----- The translational motion of the Sun over the time interval will be ignored. See :func:`~sunpy.coordinates.transform_with_sun_center`. """ # Suppress NaN warnings in coordinate transforms with warnings.catch_warnings(): warnings.simplefilter('ignore') # The time interval between the new observer time and the map observation time. interval = (parse_time(new_observer.obstime) - parse_time(smap.date)).to(u.s) # We need to get the input pixel coordinates into the OUTPUT HPC frame. # To save us having to construct a WCS etc, we do the transformation # using the output map, and then replace the observer in place before # transforming to HGS. This is acceptable because the pixel -> world # transformation is independent of the observer. input_pixels = xy.T * u.pix map_coord = smap.pixel_to_world(*input_pixels) output_hpc_coords = SkyCoord(map_coord.Tx, map_coord.Ty, map_coord.distance, obstime=new_observer.obstime, observer=new_observer, frame=Helioprojective) heliographic_coordinate = output_hpc_coords.transform_to(HeliographicStonyhurst) # Compute the differential rotation. drot = diff_rot(interval, heliographic_coordinate.lat.to(u.degree), **diff_rot_kwargs) # The change in longitude is negative because we are mapping from the # new coordinates to the old. rotated_coord = SkyCoord(heliographic_coordinate.lon - drot, heliographic_coordinate.lat, heliographic_coordinate.radius, obstime=heliographic_coordinate.obstime, frame=HeliographicStonyhurst) with transform_with_sun_center(): # As seen from the map observer, which coordinates are behind the Sun. where_off_disk_from_map_observer = rotated_coord.transform_to( Heliocentric(observer=smap.observer_coordinate)).z.value < 0 # Re-project the pixels which are on disk back to location of the original observer coordinates_at_map_observer = rotated_coord.transform_to(smap.coordinate_frame) # Go back to pixel co-ordinates x2, y2 = smap.world_to_pixel(coordinates_at_map_observer) # Re-stack the data to make it correct output form xy2 = np.dstack([x2.T.value.flat, y2.T.value.flat])[0] # Set the off disk coordinates to NaN so they are not included in the output image. xy2[where_off_disk_from_map_observer.flat] = np.nan return xy2
def solar_rotate_coordinate(coordinate, observer=None, time=None, **diff_rot_kwargs): """ Given a coordinate on the Sun, calculate where that coordinate maps to as seen by a new observer at some later or earlier time, given that the input coordinate rotates according to the solar rotation profile. The amount of solar rotation is based on the amount of time between the observation time of the input coordinate and the observation time of the new observer. The new observer is specified in one of two ways, either using the "observer" or "time" keywords. If the "observer" keyword is set, it is used to specify the location of the new observer in space and time. The difference between the coordinate time and the new observer time is used to calculate the amount of solar rotation applied, and the location of the new observer in space is used to calculate where the rotated coordinate is as seen from the new observer. If the "time" keyword is set, it is used to specify the number of seconds to rotate the coordinate by. Note that using the "time" keyword assumes that the new observer is on the Earth. This may be a reasonable assumption depending on the application. Either the "observer" or "time" keyword must be specified, but both cannot be specified at the same time. Parameters ---------- coordinate : `~astropy.coordinates.SkyCoord` Any valid coordinate which is transformable to Heliographic Stonyhurst. observer : `~astropy.coordinates.BaseCoordinateFrame`, `~astropy.coordinates.SkyCoord`, None The location of the new observer in space and time (the observer must have an interpretable obstime property). time : `~astropy.time.Time`, `~astropy.time.TimeDelta`, `~astropy.units.Quantity`, None diff_rot_kwargs : `dict` Keyword arguments are passed on as keyword arguments to `~sunpy.physics.differential_rotation.diff_rot`. Note that the keyword "frame_time" is automatically set to the value "sidereal". Returns ------- coordinate : `~astropy.coordinates.SkyCoord` The locations of the input coordinates after the application of solar rotation as seen from the point-of-view of the new observer. Notes ----- The translational motion of the Sun over the time interval will be ignored. See :func:`~sunpy.coordinates.transform_with_sun_center`. Example ------- >>> import astropy.units as u >>> from astropy.coordinates import SkyCoord >>> from sunpy.coordinates import Helioprojective, get_body_heliographic_stonyhurst >>> from sunpy.physics.differential_rotation import solar_rotate_coordinate >>> from sunpy.time import parse_time >>> start_time = parse_time('2010-09-10 12:34:56') >>> c = SkyCoord(-570*u.arcsec, 120*u.arcsec, obstime=start_time, ... observer="earth", frame=Helioprojective) >>> solar_rotate_coordinate(c, time=start_time + 25*u.hr) # doctest: +SKIP <SkyCoord (Helioprojective: obstime=2010-09-11T13:34:56.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate (obstime=2010-09-11T13:34:56.000): (lon, lat, radius) in (deg, deg, AU) (-5.68434189e-14, 7.24318962, 1.00669016)>): (Tx, Ty, distance) in (arcsec, arcsec, AU) (-378.27830452, 105.70767875, 1.00245134)> >>> new_observer = get_body_heliographic_stonyhurst("earth", start_time + 6*u.day) >>> solar_rotate_coordinate(c, observer=new_observer) <SkyCoord (Helioprojective: obstime=2010-09-16T12:34:56.000, rsun=695700.0 km, observer=<HeliographicStonyhurst Coordinate (obstime=2010-09-16T12:34:56.000): (lon, lat, radius) in (deg, deg, AU) (2.65061438e-14, 7.18706547, 1.00534174)>): (Tx, Ty, distance) in (arcsec, arcsec, AU) (620.42567049, 126.13662663, 1.00185786)> """ # Check the input and create the new observer new_observer = _get_new_observer(coordinate.obstime, observer, time) # The keyword "frame_time" must be explicitly set to "sidereal" # when using this function. diff_rot_kwargs.update({"frame_time": "sidereal"}) # Calculate the interval between the start and end time interval = (new_observer.obstime - coordinate.obstime).to(u.s) # Ignore some invalid NaN comparisions within astropy # (fixed in astropy 4.0.1 https://github.com/astropy/astropy/pull/9843) with np.errstate(invalid='ignore'): # Compute Stonyhurst Heliographic co-ordinates - returns (longitude, # latitude). Points off the limb are returned as nan. heliographic_coordinate = coordinate.transform_to(HeliographicStonyhurst) # Compute the differential rotation drot = diff_rot(interval, heliographic_coordinate.lat.to(u.degree), **diff_rot_kwargs) # Rotate the input co-ordinate as seen by the original observer heliographic_rotated = SkyCoord(heliographic_coordinate.lon + drot, heliographic_coordinate.lat, heliographic_coordinate.radius, obstime=coordinate.obstime, frame=HeliographicStonyhurst) # Calculate where the rotated co-ordinate appears as seen by new observer # for the co-ordinate system of the input co-ordinate. The translational # motion of the Sun will be ignored for the transformation. frame_newobs = coordinate.frame.replicate_without_data(observer=new_observer, obstime=new_observer.obstime) with transform_with_sun_center(): return heliographic_rotated.transform_to(frame_newobs)