.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/ex_2d.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_2d.py: **************** 2D interpolation **************** Interpolation of a two-dimensional regular grid. Bivariate ######### Perform a :py:func:`bivariate ` interpolation of gridded data points. .. GENERATED FROM PYTHON SOURCE LINES 14-23 .. 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 24-34 The first step is to load the data into memory and create the interpolator object: .. note :: An exception will be thrown if the constructor is not able to determine which axes are the longitudes and latitudes. You can force the data to be read by specifying on the longitude and latitude axes the respective ``degrees_east`` and ``degrees_north`` attribute ``units``. If your grid does not contain geodetic coordinates, set the ``geodetic`` option of the constructor to ``False``. .. GENERATED FROM PYTHON SOURCE LINES 34-37 .. code-block:: Python ds = pyinterp.tests.load_grid2d() interpolator = pyinterp.backends.xarray.Grid2D(ds.mss) .. GENERATED FROM PYTHON SOURCE LINES 38-43 We will then build the coordinates on which we want to interpolate our grid: .. note:: The coordinates used for interpolation are shifted to avoid using the points of the bivariate function. .. GENERATED FROM PYTHON SOURCE LINES 43-47 .. code-block:: Python mx, my = numpy.meshgrid(numpy.arange(-180, 180, 1) + 1 / 3.0, numpy.arange(-89, 89, 1) + 1 / 3.0, indexing='ij') .. GENERATED FROM PYTHON SOURCE LINES 48-50 The grid is :py:meth:`interpolated ` to the desired coordinates: .. GENERATED FROM PYTHON SOURCE LINES 50-55 .. code-block:: Python mss = interpolator.bivariate(coords={ 'lon': mx.ravel(), 'lat': my.ravel() }).reshape(mx.shape) .. GENERATED FROM PYTHON SOURCE LINES 56-57 Let's visualize the original grid and the result of the interpolation. .. GENERATED FROM PYTHON SOURCE LINES 57-85 .. code-block:: Python fig = matplotlib.pyplot.figure(figsize=(10, 8)) ax1 = fig.add_subplot( 211, projection=cartopy.crs.PlateCarree(central_longitude=180)) lons, lats = numpy.meshgrid(ds.lon, ds.lat, indexing='ij') pcm = ax1.pcolormesh(lons, lats, ds.mss.T, cmap='jet', shading='auto', transform=cartopy.crs.PlateCarree(), vmin=-0.1, vmax=0.1) ax1.coastlines() ax1.set_title('Original MSS') ax2 = fig.add_subplot(212, projection=cartopy.crs.PlateCarree()) pcm = ax2.pcolormesh(mx, my, mss, cmap='jet', shading='auto', transform=cartopy.crs.PlateCarree(), vmin=-0.1, vmax=0.1) ax2.coastlines() ax2.set_title('Bilinear Interpolated MSS') fig.colorbar(pcm, ax=[ax1, ax2], shrink=0.8) fig.show() .. image-sg:: /auto_examples/images/sphx_glr_ex_2d_001.png :alt: Original MSS, Bilinear Interpolated MSS :srcset: /auto_examples/images/sphx_glr_ex_2d_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 86-110 Values can be interpolated with several methods: *bilinear*, *nearest*, and *inverse distance weighting*. Distance calculations, if necessary, are calculated using the `Haversine formula `_. Bicubic ####### To interpolate data points on a regular two-dimensional grid. The interpolated surface is smoother than the corresponding surfaces obtained by bilinear interpolation. .. warning:: When using this interpolator, pay attention to the undefined values. Because as long as the calculation window uses an indefinite point, the interpolator will compute indeterminate values. In other words, this interpolator increases the area covered by the masked values. To avoid this behavior, it is necessary to :doc:`pre-process ` the grid to delete undefined values. The interpolation :py:meth:`bicubic ` function has more parameters to define the data frame used by the spline functions and how to process the edges of the regional grids: .. GENERATED FROM PYTHON SOURCE LINES 110-119 .. code-block:: Python mss = interpolator.bicubic( coords={ 'lon': mx.ravel(), 'lat': my.ravel() }, nx=3, ny=3, ).reshape(mx.shape) .. GENERATED FROM PYTHON SOURCE LINES 120-131 .. warning:: The grid provided must have strictly increasing axes to meet the specifications of the interpolation. When building the grid, specify the ``increasing_axes`` option to flip the decreasing axes and the grid automatically. For example: .. code:: python interpolator = pyinterp.backends.xarray.Grid2D( ds.mss, increasing_axes=True) .. GENERATED FROM PYTHON SOURCE LINES 131-157 .. code-block:: Python fig = matplotlib.pyplot.figure(figsize=(10, 8)) ax1 = fig.add_subplot( 211, projection=cartopy.crs.PlateCarree(central_longitude=180)) pcm = ax1.pcolormesh(lons, lats, ds.mss.T, cmap='jet', shading='auto', transform=cartopy.crs.PlateCarree(), vmin=-0.1, vmax=0.1) ax1.coastlines() ax1.set_title('Original MSS') ax2 = fig.add_subplot(212, projection=cartopy.crs.PlateCarree()) pcm = ax2.pcolormesh(mx, my, mss, cmap='jet', shading='auto', transform=cartopy.crs.PlateCarree(), vmin=-0.1, vmax=0.1) ax2.coastlines() ax2.set_title('Bicubic Interpolated MSS') fig.colorbar(pcm, ax=[ax1, ax2], shrink=0.8) fig.show() .. image-sg:: /auto_examples/images/sphx_glr_ex_2d_002.png :alt: Original MSS, Bicubic Interpolated MSS :srcset: /auto_examples/images/sphx_glr_ex_2d_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 3.335 seconds) .. _sphx_glr_download_auto_examples_ex_2d.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_2d.ipynb :alt: Launch binder :width: 150 px .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: ex_2d.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: ex_2d.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: ex_2d.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_