.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/orbit/ex_orbit.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. or to run this example in your browser via Binder .. rst-class:: sphx-glr-example-title .. _sphx_glr_auto_examples_orbit_ex_orbit.py: ******************* Orbit Interpolation ******************* This example demonstrates how to use ``pyinterp`` to interpolate satellite orbits. The library can propagate orbit ephemerides from a template file that contains satellite positions for a single orbit cycle. This is useful for simulations and other applications that require accurate orbit propagation over time. To begin, we will load the orbit ephemeris from a test file. .. GENERATED FROM PYTHON SOURCE LINES 15-75 .. code-block:: Python import pathlib import cartopy.crs as ccrs import matplotlib.pyplot as plt import numpy as np import pandas as pd import pyinterp import pyinterp.tests def load_test_ephemeris( filename: pathlib.Path, ) -> tuple[float, np.ndarray, np.ndarray, np.ndarray, np.timedelta64]: """Loads the ephemeris from a text file. Args: filename: Name of the file to be loaded. Returns: A tuple containing the height of the orbit, the ephemeris, and the duration of the cycle. """ with open(filename, encoding='utf-8') as stream: lines = stream.readlines() def to_dict(comments) -> dict[str, float]: """Returns a dictionary describing the parameters of the orbit.""" result = {} for item in comments: if not item.startswith('#'): raise ValueError('Comments must start with #') key, value = item[1:].split('=') result[key.strip()] = float(value) return result # The first two lines are the header and contain the height and the # duration of the cycle in fractional days. settings = to_dict(lines[:2]) del lines[:2] # The rest of the lines are the ephemeris. ephemeris = np.loadtxt( lines, delimiter=' ', dtype={ 'names': ('time', 'longitude', 'latitude', 'height'), 'formats': ('f8', 'f8', 'f8', 'f8'), }, ) return ( settings['height'], ephemeris['longitude'], ephemeris['latitude'], ephemeris['time'].astype('timedelta64[s]'), np.timedelta64(int(settings['cycle_duration'] * 86400.0 * 1e9), 'ns'), ) .. GENERATED FROM PYTHON SOURCE LINES 76-79 Loading the Ephemeris ===================== First, we set the path to the test file and load the ephemeris data. .. GENERATED FROM PYTHON SOURCE LINES 79-82 .. code-block:: Python swot_calval_ephemeris_path = pyinterp.tests.swot_calval_ephemeris_path() ephemeris = load_test_ephemeris(swot_calval_ephemeris_path) .. GENERATED FROM PYTHON SOURCE LINES 83-87 Calculating Orbit Properties ============================ With the ephemeris loaded, we can compute the orbit properties using the :py:func:`pyinterp.calculate_orbit` function. .. GENERATED FROM PYTHON SOURCE LINES 87-89 .. code-block:: Python orbit = pyinterp.calculate_orbit(*ephemeris) .. GENERATED FROM PYTHON SOURCE LINES 90-93 The resulting :py:class:`pyinterp.Orbit` object provides methods to access various orbit properties. For example, you can get the number of passes per cycle: .. GENERATED FROM PYTHON SOURCE LINES 93-95 .. code-block:: Python print(f'Passes per cycle: {orbit.passes_per_cycle()}') .. rst-class:: sphx-glr-script-out .. code-block:: none Passes per cycle: 28 .. GENERATED FROM PYTHON SOURCE LINES 96-97 You can also get the cycle duration and the orbit duration: .. GENERATED FROM PYTHON SOURCE LINES 97-100 .. code-block:: Python print(f'Cycle duration: {orbit.cycle_duration().astype("m8[ms]").item()}') print(f'Orbit duration: {orbit.orbit_duration().astype("m8[ms]").item()}') .. rst-class:: sphx-glr-script-out .. code-block:: none Cycle duration: 23:25:23.172000 Orbit duration: 1:40:23.083000 .. GENERATED FROM PYTHON SOURCE LINES 101-102 The duration of a specific pass can also be retrieved: .. GENERATED FROM PYTHON SOURCE LINES 102-104 .. code-block:: Python print(f'Pass 2 duration: {orbit.pass_duration(2).astype("m8[ms]").item()}') .. rst-class:: sphx-glr-script-out .. code-block:: none Pass 2 duration: 0:55:44.350000 .. GENERATED FROM PYTHON SOURCE LINES 105-110 Encoding and Decoding Pass Numbers ================================== A utility function is provided to compute an absolute pass number from a relative pass number and a cycle number. This is useful for storing pass information in a database or for indexing passes in a file. .. GENERATED FROM PYTHON SOURCE LINES 110-114 .. code-block:: Python absolute_pass_number = orbit.encode_absolute_pass_number(cycle_number=11, pass_number=2) print(f'Absolute pass number: {absolute_pass_number}') .. rst-class:: sphx-glr-script-out .. code-block:: none Absolute pass number: 282 .. GENERATED FROM PYTHON SOURCE LINES 115-116 You can decode the absolute pass number to get the cycle and pass numbers: .. GENERATED FROM PYTHON SOURCE LINES 116-120 .. code-block:: Python cycle_number, pass_number = orbit.decode_absolute_pass_number( absolute_pass_number) print(f'Cycle: {cycle_number}, Pass: {pass_number}') .. rst-class:: sphx-glr-script-out .. code-block:: none Cycle: 11, Pass: 2 .. GENERATED FROM PYTHON SOURCE LINES 121-137 Interpolating the Orbit ======================= The next step is to interpolate the orbit ephemerides over time to get the satellite positions for a given relative pass number. .. note:: You can iterate over the relative pass numbers for a given time period using the :py:meth:`pyinterp.Orbit.iterate` method. .. code-block:: python for cycle_number, pass_number, first_location_date in orbit.iterate( start_date, end_date): ... .. GENERATED FROM PYTHON SOURCE LINES 137-145 .. code-block:: Python pass_ = pyinterp.calculate_pass(2, orbit) assert pass_ is not None pd.DataFrame({ 'time': pass_.time, 'lon_nadir': pass_.lon_nadir, 'lat_nadir': pass_.lat_nadir, }) .. raw:: html
time lon_nadir lat_nadir
0 0 days 00:25:53.799409419 -29.393586 -77.497916
1 0 days 00:25:54.519992305 -29.310591 -77.497883
2 0 days 00:25:55.241146762 -29.227597 -77.497825
3 0 days 00:25:55.962772250 -29.144604 -77.497741
4 0 days 00:25:56.684960867 -29.061612 -77.497632
... ... ... ...
9953 0 days 01:21:25.850882766 146.412747 77.727877
9954 0 days 01:21:28.310619399 146.497271 77.727982
9955 0 days 01:21:30.770464201 146.581797 77.728059
9956 0 days 01:21:33.230367978 146.666324 77.728106
9957 0 days 01:21:35.690281545 146.750851 77.728124

9958 rows × 3 columns



.. GENERATED FROM PYTHON SOURCE LINES 146-155 The :py:class:`pyinterp.Pass` object contains the satellite's nadir positions for the given pass, including: * Nadir longitude (in degrees) * Nadir latitude (in degrees) * Time for each position (as a ``numpy.datetime64`` array) * Along-track distance (in meters) It also provides the coordinates of the satellite at the equator: .. GENERATED FROM PYTHON SOURCE LINES 155-157 .. code-block:: Python print(f'Equator coordinates: {pass_.equator_coordinates}') .. rst-class:: sphx-glr-script-out .. code-block:: none Equator coordinates: EquatorCoordinates(longitude=59.20360095529902, time=np.timedelta64(3087382968771,'ns')) .. GENERATED FROM PYTHON SOURCE LINES 158-168 .. note:: The ``pass_`` variable will be ``None`` if the pass number is outside the bounding box defined during the instantiation of the :py:class:`pyinterp.Orbit` object. Calculating the Swath ===================== Finally, we can calculate the satellite positions over a swath for a given pass. .. GENERATED FROM PYTHON SOURCE LINES 168-170 .. code-block:: Python swath = pyinterp.calculate_swath(pass_) .. GENERATED FROM PYTHON SOURCE LINES 171-175 The :py:class:`pyinterp.Swath` object contains the properties of the pass, similar to the :py:class:`pyinterp.Pass` object. Additionally, it includes the coordinates of the satellite over the swath for each location on the nadir track. .. GENERATED FROM PYTHON SOURCE LINES 175-181 .. code-block:: Python pd.DataFrame({ 'time': swath.time, 'lon_nadir': swath.lon_nadir, 'lat_nadir': swath.lat_nadir, }) .. raw:: html
time lon_nadir lat_nadir
0 0 days 00:25:53.799409419 -29.393586 -77.497916
1 0 days 00:25:54.519992305 -29.310591 -77.497883
2 0 days 00:25:55.241146762 -29.227597 -77.497825
3 0 days 00:25:55.962772250 -29.144604 -77.497741
4 0 days 00:25:56.684960867 -29.061612 -77.497632
... ... ... ...
9953 0 days 01:21:25.850882766 146.412747 77.727877
9954 0 days 01:21:28.310619399 146.497271 77.727982
9955 0 days 01:21:30.770464201 146.581797 77.728059
9956 0 days 01:21:33.230367978 146.666324 77.728106
9957 0 days 01:21:35.690281545 146.750851 77.728124

9958 rows × 3 columns



.. GENERATED FROM PYTHON SOURCE LINES 182-185 The following DataFrame shows the longitude and latitude coordinates of the satellite for the first two lines of the swath. The index ``x_ac`` represents the across-track distance in meters. .. GENERATED FROM PYTHON SOURCE LINES 185-197 .. code-block:: Python df = pd.DataFrame( { 'lon_0': swath.lon[0, :], 'lon_1': swath.lon[1, :], 'lat_0': swath.lat[0, :], 'lat_1': swath.lat[1, :], }, index=swath.x_ac[0, :], ) df.index.name = 'x_ac' df .. raw:: html
lon_0 lon_1 lat_0 lat_1
x_ac
-70000.0 -29.396648 -29.317550 -76.870888 -76.870857
-68000.0 -29.396565 -29.317360 -76.888803 -76.888772
-66000.0 -29.396481 -29.317170 -76.906719 -76.906688
-64000.0 -29.396397 -29.316979 -76.924634 -76.924603
-62000.0 -29.396312 -29.316788 -76.942550 -76.942519
... ... ... ... ...
62000.0 -29.390618 -29.303828 -78.053251 -78.053217
64000.0 -29.390518 -29.303599 -78.071164 -78.071130
66000.0 -29.390417 -29.303370 -78.089078 -78.089044
68000.0 -29.390317 -29.303140 -78.106991 -78.106957
70000.0 -29.390216 -29.302910 -78.124905 -78.124871

70 rows × 4 columns



.. GENERATED FROM PYTHON SOURCE LINES 198-201 Visualizing the Swath ===================== We can now plot the satellite positions over the swath. .. GENERATED FROM PYTHON SOURCE LINES 201-226 .. code-block:: Python fig, (ax1, ax2) = plt.subplots(1, 2, figsize=(20, 10), subplot_kw={'projection': ccrs.PlateCarree()}) fig.subplots_adjust(left=0.05, right=0.95, top=0.95, bottom=0.05, wspace=0.1) # Zoomed plot ax1.set_extent([58, 70, 12, 32], crs=ccrs.PlateCarree()) ax1.plot(swath.lon[::4, ::4].ravel(), swath.lat[::4, ::4].ravel(), 'b.', markersize=1) ax1.plot(swath.lon_nadir, swath.lat_nadir, 'r.', markersize=0.5) ax1.set_title('Satellite positions - Zoomed') ax1.coastlines() ax1.gridlines(draw_labels=True) # Full swath plot ax2.plot(swath.lon.ravel(), swath.lat.ravel(), 'b.', markersize=1) ax2.plot(swath.lon_nadir, swath.lat_nadir, 'r.', markersize=0.5) ax2.set_title('Satellite positions - Full Swath') ax2.coastlines() ax2.gridlines(draw_labels=True) ax1.set_aspect('auto') ax2.set_aspect('auto') .. image-sg:: /auto_examples/orbit/images/sphx_glr_ex_orbit_001.png :alt: Satellite positions - Zoomed, Satellite positions - Full Swath :srcset: /auto_examples/orbit/images/sphx_glr_ex_orbit_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 1.291 seconds) .. _sphx_glr_download_auto_examples_orbit_ex_orbit.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: binder-badge .. image:: images/binder_badge_logo.svg :target: https://mybinder.org/v2/gh/CNES/pangeo-pyinterp/master?urlpath=lab/tree/notebooks/auto_examples/orbit/ex_orbit.ipynb :alt: Launch binder :width: 150 px .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: ex_orbit.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: ex_orbit.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: ex_orbit.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_