.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/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_ex_3d.py: **************** 3D interpolation **************** Interpolation of a three-dimensional regular grid. Trivariate ========== The :py:func:`trivariate ` interpolation allows obtaining values at arbitrary points in a 3D space of a function defined on a grid. This method performs a bilinear interpolation in 2D space by considering the axes of longitude and latitude of the grid, then performs a linear interpolation in the third dimension. Its interface is similar to the :py:func:`bivariate ` class except for a third axis, which is handled by this object. .. note:: When using a time axis, care must be taken to use the same unit of dates, between the axis defined and the dates supplied during interpolation. The function :py:meth:`pyinterp.TemporalAxis.safe_cast` automates this task and will warn you if there is an inconsistency during the date conversion. .. GENERATED FROM PYTHON SOURCE LINES 27-36 .. 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 37-39 The first step is to load the data into memory and create the interpolator object: .. GENERATED FROM PYTHON SOURCE LINES 39-42 .. code-block:: Python ds = pyinterp.tests.load_grid3d() interpolator = pyinterp.backends.xarray.Grid3D(ds.tcw) .. GENERATED FROM PYTHON SOURCE LINES 43-56 We will build a new grid that will be used to build a new interpolated grid. .. note :: The coordinates used for interpolation are shifted to avoid using the points of the trivariate function. .. warning :: When using a time axis, care must be taken to use the same unit of dates, between the axis defined and the dates supplied during interpolation. The function :py:meth:`pyinterp.TemporalAxis.safe_cast` automates this task and will warn you if there is an inconsistency during the date conversion. .. GENERATED FROM PYTHON SOURCE LINES 56-62 .. 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 63-65 We interpolate our grid using a :py:meth:`classical `: .. GENERATED FROM PYTHON SOURCE LINES 65-71 .. code-block:: Python trivariate = interpolator.trivariate({ 'longitude': mx.ravel(), 'latitude': my.ravel(), 'time': mz.ravel() }) .. GENERATED FROM PYTHON SOURCE LINES 72-79 Bicubic on 3D grid ================== Used grid organizes the latitudes in descending order. We ask our constructor to flip this axis in order to correctly evaluate the bicubic interpolation from this 3D cube (only necessary to perform a bicubic interpolation). .. GENERATED FROM PYTHON SOURCE LINES 79-82 .. code-block:: Python interpolator = pyinterp.backends.xarray.Grid3D(ds.data_vars['tcw'], increasing_axes=True) .. GENERATED FROM PYTHON SOURCE LINES 83-86 We interpolate our grid using a :py:meth:`bicubic ` interpolation in space followed by a linear interpolation in the temporal axis: .. GENERATED FROM PYTHON SOURCE LINES 86-92 .. code-block:: Python bicubic = interpolator.bicubic({ 'longitude': mx.ravel(), 'latitude': my.ravel(), 'time': mz.ravel() }) .. GENERATED FROM PYTHON SOURCE LINES 93-94 We transform our result cubes into a matrix. .. GENERATED FROM PYTHON SOURCE LINES 94-99 .. 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 100-101 Let's visualize our results. .. GENERATED FROM PYTHON SOURCE LINES 101-131 .. code-block:: Python fig = matplotlib.pyplot.figure(figsize=(5, 8)) 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) fig.show() .. image-sg:: /auto_examples/images/sphx_glr_ex_3d_001.png :alt: Trilinear, Spline & Linear in time :srcset: /auto_examples/images/sphx_glr_ex_3d_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 5.071 seconds) .. _sphx_glr_download_auto_examples_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/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 `_