.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "generated/examples/coordinates/plot_obs-planning.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code .. rst-class:: sphx-glr-example-title .. _sphx_glr_generated_examples_coordinates_plot_obs-planning.py: =================================================================== Determining and plotting the altitude/azimuth of a celestial object =================================================================== This example demonstrates coordinate transformations and the creation of visibility curves to assist with observing run planning. In this example, we make a `~astropy.coordinates.SkyCoord` instance for M33. The altitude-azimuth coordinates are then found using `astropy.coordinates.EarthLocation` and `astropy.time.Time` objects. This example is meant to demonstrate the capabilities of the `astropy.coordinates` package. For more convenient and/or complex observation planning, consider the `astroplan `_ package. *By: Erik Tollerud, Kelle Cruz* *License: BSD* .. GENERATED FROM PYTHON SOURCE LINES 27-35 Let's suppose you are planning to visit picturesque Bear Mountain State Park in New York, USA. You're bringing your telescope with you (of course), and someone told you M33 is a great target to observe there. You happen to know you're free at 11:00 pm local time, and you want to know if it will be up. Astropy can answer that. Import numpy and matplotlib. For the latter, use a nicer set of plot parameters and set up support for plotting/converting quantities. .. GENERATED FROM PYTHON SOURCE LINES 35-45 .. code-block:: Python import matplotlib.pyplot as plt import numpy as np from astropy.visualization import astropy_mpl_style, quantity_support plt.style.use(astropy_mpl_style) quantity_support() .. rst-class:: sphx-glr-script-out .. code-block:: none .MplQuantityConverter object at 0x7f46cb4735d0> .. GENERATED FROM PYTHON SOURCE LINES 46-48 Import the packages necessary for finding coordinates and making coordinate transformations .. GENERATED FROM PYTHON SOURCE LINES 48-53 .. code-block:: Python import astropy.units as u from astropy.coordinates import AltAz, EarthLocation, SkyCoord from astropy.time import Time .. GENERATED FROM PYTHON SOURCE LINES 54-58 `astropy.coordinates.SkyCoord.from_name` uses Simbad to resolve object names and retrieve coordinates. Get the coordinates of M33: .. GENERATED FROM PYTHON SOURCE LINES 58-61 .. code-block:: Python m33 = SkyCoord.from_name('M33') .. GENERATED FROM PYTHON SOURCE LINES 62-64 Use `astropy.coordinates.EarthLocation` to provide the location of Bear Mountain and set the time to 11pm EDT on 2012 July 12: .. GENERATED FROM PYTHON SOURCE LINES 64-69 .. code-block:: Python bear_mountain = EarthLocation(lat=41.3*u.deg, lon=-74*u.deg, height=390*u.m) utcoffset = -4*u.hour # Eastern Daylight Time time = Time('2012-7-12 23:00:00') - utcoffset .. GENERATED FROM PYTHON SOURCE LINES 70-76 `astropy.coordinates.EarthLocation.get_site_names` and `~astropy.coordinates.EarthLocation.get_site_names` can be used to get locations of major observatories. Use `astropy.coordinates` to find the Alt, Az coordinates of M33 at as observed from Bear Mountain at 11pm on 2012 July 12. .. GENERATED FROM PYTHON SOURCE LINES 76-80 .. code-block:: Python m33altaz = m33.transform_to(AltAz(obstime=time,location=bear_mountain)) print(f"M33's Altitude = {m33altaz.alt:.2}") .. rst-class:: sphx-glr-script-out .. code-block:: none M33's Altitude = 0.13 deg .. GENERATED FROM PYTHON SOURCE LINES 81-87 This is helpful since it turns out M33 is barely above the horizon at this time. It's more informative to find M33's airmass over the course of the night. Find the alt,az coordinates of M33 at 100 times evenly spaced between 10pm and 7am EDT: .. GENERATED FROM PYTHON SOURCE LINES 87-94 .. code-block:: Python midnight = Time('2012-7-13 00:00:00') - utcoffset delta_midnight = np.linspace(-2, 10, 100)*u.hour frame_July13night = AltAz(obstime=midnight+delta_midnight, location=bear_mountain) m33altazs_July13night = m33.transform_to(frame_July13night) .. GENERATED FROM PYTHON SOURCE LINES 95-96 convert alt, az to airmass with `~astropy.coordinates.AltAz.secz` attribute: .. GENERATED FROM PYTHON SOURCE LINES 96-99 .. code-block:: Python m33airmasss_July13night = m33altazs_July13night.secz .. GENERATED FROM PYTHON SOURCE LINES 100-101 Plot the airmass as a function of time: .. GENERATED FROM PYTHON SOURCE LINES 101-109 .. code-block:: Python plt.plot(delta_midnight, m33airmasss_July13night) plt.xlim(-2, 10) plt.ylim(1, 4) plt.xlabel('Hours from EDT Midnight') plt.ylabel('Airmass [Sec(z)]') plt.show() .. image-sg:: /generated/examples/coordinates/images/sphx_glr_plot_obs-planning_001.png :alt: plot obs planning :srcset: /generated/examples/coordinates/images/sphx_glr_plot_obs-planning_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 110-112 Use `~astropy.coordinates.get_sun` to find the location of the Sun at 1000 evenly spaced times between noon on July 12 and noon on July 13: .. GENERATED FROM PYTHON SOURCE LINES 112-121 .. code-block:: Python from astropy.coordinates import get_sun delta_midnight = np.linspace(-12, 12, 1000)*u.hour times_July12_to_13 = midnight + delta_midnight frame_July12_to_13 = AltAz(obstime=times_July12_to_13, location=bear_mountain) sunaltazs_July12_to_13 = get_sun(times_July12_to_13).transform_to(frame_July12_to_13) .. GENERATED FROM PYTHON SOURCE LINES 122-125 Do the same with `~astropy.coordinates.get_body` to find when the moon is up. Be aware that this will need to download a 10MB file from the internet to get a precise location of the moon. .. GENERATED FROM PYTHON SOURCE LINES 125-131 .. code-block:: Python from astropy.coordinates import get_body moon_July12_to_13 = get_body("moon", times_July12_to_13) moonaltazs_July12_to_13 = moon_July12_to_13.transform_to(frame_July12_to_13) .. GENERATED FROM PYTHON SOURCE LINES 132-133 Find the alt,az coordinates of M33 at those same times: .. GENERATED FROM PYTHON SOURCE LINES 133-136 .. code-block:: Python m33altazs_July12_to_13 = m33.transform_to(frame_July12_to_13) .. GENERATED FROM PYTHON SOURCE LINES 137-139 Make a beautiful figure illustrating nighttime and the altitudes of M33 and the Sun over that time: .. GENERATED FROM PYTHON SOURCE LINES 139-157 .. code-block:: Python plt.plot(delta_midnight, sunaltazs_July12_to_13.alt, color='r', label='Sun') plt.plot(delta_midnight, moonaltazs_July12_to_13.alt, color=[0.75]*3, ls='--', label='Moon') plt.scatter(delta_midnight, m33altazs_July12_to_13.alt, c=m33altazs_July12_to_13.az.value, label='M33', lw=0, s=8, cmap='viridis') plt.fill_between(delta_midnight, 0*u.deg, 90*u.deg, sunaltazs_July12_to_13.alt < -0*u.deg, color='0.5', zorder=0) plt.fill_between(delta_midnight, 0*u.deg, 90*u.deg, sunaltazs_July12_to_13.alt < -18*u.deg, color='k', zorder=0) plt.colorbar().set_label('Azimuth [deg]') plt.legend(loc='upper left') plt.xlim(-12*u.hour, 12*u.hour) plt.xticks((np.arange(13)*2-12)*u.hour) plt.ylim(0*u.deg, 90*u.deg) plt.xlabel('Hours from EDT Midnight') plt.ylabel('Altitude [deg]') plt.show() .. image-sg:: /generated/examples/coordinates/images/sphx_glr_plot_obs-planning_002.png :alt: plot obs planning :srcset: /generated/examples/coordinates/images/sphx_glr_plot_obs-planning_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 3.159 seconds) .. _sphx_glr_download_generated_examples_coordinates_plot_obs-planning.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_obs-planning.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_obs-planning.py ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_