.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "generated/gallery/oceanography/plot_atlantic_profiles.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_gallery_oceanography_plot_atlantic_profiles.py: Oceanographic Profiles and T-S Diagrams ======================================= This example demonstrates how to plot vertical profiles of different variables in the same axes, and how to make a scatter plot of two variables. There is an oceanographic theme but the same techniques are equally applicable to atmospheric or other kinds of data. The data used are profiles of potential temperature and salinity in the Equatorial and South Atlantic, output from an ocean model. The y-axis of the first plot produced will be automatically inverted due to the presence of the attribute positive=down on the depth coordinate. This means depth values intuitively increase downward on the y-axis. .. GENERATED FROM PYTHON SOURCE LINES 18-93 .. rst-class:: sphx-glr-horizontal * .. image-sg:: /generated/gallery/oceanography/images/sphx_glr_plot_atlantic_profiles_001.png :alt: plot atlantic profiles :srcset: /generated/gallery/oceanography/images/sphx_glr_plot_atlantic_profiles_001.png :class: sphx-glr-multi-img * .. image-sg:: /generated/gallery/oceanography/images/sphx_glr_plot_atlantic_profiles_002.png :alt: plot atlantic profiles :srcset: /generated/gallery/oceanography/images/sphx_glr_plot_atlantic_profiles_002.png :class: sphx-glr-multi-img .. code-block:: Python import matplotlib.pyplot as plt import iris import iris.iterate import iris.plot as iplt def main(): # Load the gridded temperature and salinity data. fname = iris.sample_data_path("atlantic_profiles.nc") cubes = iris.load(fname) (theta,) = cubes.extract("sea_water_potential_temperature") (salinity,) = cubes.extract("sea_water_practical_salinity") # Extract profiles of temperature and salinity from a particular point in # the southern portion of the domain, and limit the depth of the profile # to 1000m. lon_cons = iris.Constraint(longitude=330.5) lat_cons = iris.Constraint(latitude=lambda lat: -10 < lat < -9) depth_cons = iris.Constraint(depth=lambda d: d <= 1000) theta_1000m = theta.extract(depth_cons & lon_cons & lat_cons) salinity_1000m = salinity.extract(depth_cons & lon_cons & lat_cons) # Plot these profiles on the same set of axes. Depth is automatically # recognised as a vertical coordinate and placed on the y-axis. # The first plot is in the default axes. We'll use the same color for the # curve and its axes/tick labels. plt.figure(figsize=(5, 6)) temperature_color = (0.3, 0.4, 0.5) ax1 = plt.gca() iplt.plot( theta_1000m, linewidth=2, color=temperature_color, alpha=0.75, ) ax1.set_xlabel("Potential Temperature / K", color=temperature_color) ax1.set_ylabel("Depth / m") for ticklabel in ax1.get_xticklabels(): ticklabel.set_color(temperature_color) # To plot salinity in the same axes we use twiny(). We'll use a different # color to identify salinity. salinity_color = (0.6, 0.1, 0.15) ax2 = plt.gca().twiny() iplt.plot( salinity_1000m, linewidth=2, color=salinity_color, alpha=0.75, ) ax2.set_xlabel("Salinity / PSU", color=salinity_color) for ticklabel in ax2.get_xticklabels(): ticklabel.set_color(salinity_color) plt.tight_layout() iplt.show() # Now plot a T-S diagram using scatter. We'll use all the profiles here, # and each point will be coloured according to its depth. plt.figure(figsize=(6, 6)) depth_values = theta.coord("depth").points for s, t in iris.iterate.izip(salinity, theta, coords="depth"): iplt.scatter(s, t, c=depth_values, marker="+", cmap="RdYlBu_r") ax = plt.gca() ax.set_xlabel("Salinity / PSU") ax.set_ylabel("Potential Temperature / K") cb = plt.colorbar(orientation="horizontal") cb.set_label("Depth / m") plt.tight_layout() iplt.show() if __name__ == "__main__": main() .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.752 seconds) .. _sphx_glr_download_generated_gallery_oceanography_plot_atlantic_profiles.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_atlantic_profiles.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_atlantic_profiles.py ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_