.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/grid/ex_3d.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_grid_ex_3d.py: .. _example_3d_interpolation: 3D Interpolation ================ This example demonstrates how to perform 3D interpolation on a regular grid. The pyinterp library supports both trivariate and bicubic interpolation for 3D data. .. GENERATED FROM PYTHON SOURCE LINES 10-19 .. code-block:: Python import cartopy.crs import matplotlib import matplotlib.pyplot import numpy import pyinterp import pyinterp.backends.xarray import pyinterp.tests .. GENERATED FROM PYTHON SOURCE LINES 20-29 Trivariate Interpolation ------------------------ Trivariate interpolation extends bivariate interpolation to three dimensions. It performs a bilinear interpolation on the 2D spatial plane (longitude and latitude) and then a linear interpolation on the third dimension (in this case, time). First, we load the 3D dataset and create the interpolator object. .. GENERATED FROM PYTHON SOURCE LINES 29-32 .. code-block:: Python ds = pyinterp.tests.load_grid3d() interpolator = pyinterp.backends.xarray.Grid3D(ds.tcw) .. GENERATED FROM PYTHON SOURCE LINES 33-42 Next, we define the coordinates for interpolation. To avoid interpolating at the exact grid points, we introduce a slight shift. .. note:: When working with a time axis, ensure that the date units are consistent between the grid and the interpolation coordinates. The :py:meth:`pyinterp.TemporalAxis.safe_cast` method can help manage date conversions and prevent inconsistencies. .. GENERATED FROM PYTHON SOURCE LINES 42-48 .. code-block:: Python mx, my, mz = numpy.meshgrid(numpy.arange(-180, 180, 0.25) + 1 / 3.0, numpy.arange(-80, 80, 0.25) + 1 / 3.0, numpy.array(['2002-07-02T15:00:00'], dtype='datetime64'), indexing='ij') .. GENERATED FROM PYTHON SOURCE LINES 49-50 Now, we perform the trivariate interpolation. .. GENERATED FROM PYTHON SOURCE LINES 50-56 .. code-block:: Python trivariate = interpolator.trivariate({ 'longitude': mx.ravel(), 'latitude': my.ravel(), 'time': mz.ravel() }) .. GENERATED FROM PYTHON SOURCE LINES 57-68 Bicubic Interpolation on a 3D Grid ---------------------------------- For smoother results, you can use bicubic interpolation for the spatial dimensions, followed by a linear interpolation on the third dimension. .. note:: Bicubic interpolation requires that the grid axes are strictly increasing. If your latitudes are in descending order, you can set the `increasing_axes` parameter to ``True`` to automatically flip them. .. GENERATED FROM PYTHON SOURCE LINES 68-71 .. code-block:: Python interpolator = pyinterp.backends.xarray.Grid3D(ds.data_vars['tcw'], increasing_axes=True) .. GENERATED FROM PYTHON SOURCE LINES 72-73 We then perform the bicubic interpolation. .. GENERATED FROM PYTHON SOURCE LINES 73-79 .. code-block:: Python bicubic = interpolator.bicubic({ 'longitude': mx.ravel(), 'latitude': my.ravel(), 'time': mz.ravel() }) .. GENERATED FROM PYTHON SOURCE LINES 80-82 To visualize the results, we reshape the output arrays and extract the longitude and latitude coordinates. .. GENERATED FROM PYTHON SOURCE LINES 82-87 .. code-block:: Python trivariate = trivariate.reshape(mx.shape).squeeze(axis=2) bicubic = bicubic.reshape(mx.shape).squeeze(axis=2) lons = mx[:, 0].squeeze() lats = my[0, :].squeeze() .. GENERATED FROM PYTHON SOURCE LINES 88-89 Finally, let's plot the results of both trivariate and bicubic interpolation. .. GENERATED FROM PYTHON SOURCE LINES 89-115 .. code-block:: Python fig = matplotlib.pyplot.figure(figsize=(10, 8)) fig.subplots_adjust(left=0.05, right=0.95, top=0.95, bottom=0.05, hspace=0.25) ax1 = fig.add_subplot( 211, projection=cartopy.crs.PlateCarree(central_longitude=180)) ax1.set_extent([-180, 180, -90, 90], crs=cartopy.crs.PlateCarree()) pcm = ax1.pcolormesh(lons, lats, trivariate.T, cmap='jet', shading='auto', transform=cartopy.crs.PlateCarree()) ax1.coastlines() ax1.set_title('Trivariate Interpolation') ax2 = fig.add_subplot( 212, projection=cartopy.crs.PlateCarree(central_longitude=180)) ax2.set_extent([-180, 180, -90, 90], crs=cartopy.crs.PlateCarree()) pcm = ax2.pcolormesh(lons, lats, bicubic.T, cmap='jet', shading='auto', transform=cartopy.crs.PlateCarree()) ax2.coastlines() ax2.set_title('Bicubic Interpolation') fig.colorbar(pcm, ax=[ax1, ax2], shrink=0.8) .. image-sg:: /auto_examples/grid/images/sphx_glr_ex_3d_001.png :alt: Trivariate Interpolation, Bicubic Interpolation :srcset: /auto_examples/grid/images/sphx_glr_ex_3d_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 116-118 Same plot as above, but zoomed into a specific region to better highlight the differences between the two interpolation methods. .. GENERATED FROM PYTHON SOURCE LINES 118-148 .. code-block:: Python fig = matplotlib.pyplot.figure(figsize=(5, 8)) fig.subplots_adjust(left=0.05, right=0.95, top=0.95, bottom=0.05, hspace=0.25) ax1 = fig.add_subplot( 211, projection=cartopy.crs.PlateCarree(central_longitude=180)) pcm = ax1.pcolormesh(lons, lats, trivariate.T, cmap='jet', shading='auto', transform=cartopy.crs.PlateCarree(), vmin=0, vmax=80) ax1.coastlines() ax1.set_extent([80, 170, -45, 30], crs=cartopy.crs.PlateCarree()) ax1.set_title('Trilinear') ax2 = fig.add_subplot( 212, projection=cartopy.crs.PlateCarree(central_longitude=180)) pcm = ax2.pcolormesh(lons, lats, bicubic.T, cmap='jet', shading='auto', transform=cartopy.crs.PlateCarree(), vmin=0, vmax=80) ax2.coastlines() ax2.set_extent([80, 170, -45, 30], crs=cartopy.crs.PlateCarree()) ax2.set_title('Spline & Linear in time') fig.colorbar(pcm, ax=[ax1, ax2], shrink=0.8) .. image-sg:: /auto_examples/grid/images/sphx_glr_ex_3d_002.png :alt: Trilinear, Spline & Linear in time :srcset: /auto_examples/grid/images/sphx_glr_ex_3d_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 10.236 seconds) .. _sphx_glr_download_auto_examples_grid_ex_3d.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/grid/ex_3d.ipynb :alt: Launch binder :width: 150 px .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: ex_3d.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: ex_3d.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: ex_3d.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_