
.. DO NOT EDIT.
.. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY.
.. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE:
.. "examples/plot_diags_prior_eddies.py"
.. LINE NUMBERS ARE GIVEN BELOW.

.. only:: html

    .. note::
        :class: sphx-glr-download-link-note

        :ref:`Go to the end <sphx_glr_download_examples_plot_diags_prior_eddies.py>`
        to download the full example code.

.. rst-class:: sphx-glr-example-title

.. _sphx_glr_examples_plot_diags_prior_eddies.py:


Diagnostics prior to eddy detection
===================================

.. GENERATED FROM PYTHON SOURCE LINES 8-12

Initialisations
-----------------

Import needed stuff.

.. GENERATED FROM PYTHON SOURCE LINES 12-23

.. code-block:: Python

    import cmocean as cm
    import matplotlib.pyplot as plt
    import xarray as xr

    from shoot.dyn import get_geos, get_lnam, get_okuboweiss
    from shoot.grid import get_dx_dy
    from shoot.plot import create_map, pcarr
    from shoot.samples import get_sample_file

    xr.set_options(display_style="text")





.. rst-class:: sphx-glr-script-out

 .. code-block:: none


    <xarray.core.options.set_options object at 0x70269976bf90>



.. GENERATED FROM PYTHON SOURCE LINES 24-25

Read data

.. GENERATED FROM PYTHON SOURCE LINES 25-30

.. code-block:: Python

    root_path = "OBS/SATELLITE/jan2024_ionian_sea_duacs.nc"
    path = get_sample_file(root_path)
    ds = xr.open_dataset(path).isel(time=0)
    # ds = ds.sel(longitude=slice(24, 34), latitude=slice(31, 35))








.. GENERATED FROM PYTHON SOURCE LINES 31-32

Compute grid metrics only once

.. GENERATED FROM PYTHON SOURCE LINES 32-34

.. code-block:: Python

    dx, dy = get_dx_dy(ds)








.. GENERATED FROM PYTHON SOURCE LINES 35-39

Geostrophic currents
--------------------

ADT

.. GENERATED FROM PYTHON SOURCE LINES 39-43

.. code-block:: Python

    fig, ax = create_map(ds.longitude, ds.latitude, figsize=(8, 5))
    ds.adt.plot(ax=ax, cmap="cmo.balance", add_colorbar=False, transform=pcarr)
    ds.adt.plot.contour(ax=ax, transform=pcarr, colors="k", levels=20)




.. image-sg:: /examples/images/sphx_glr_plot_diags_prior_eddies_001.png
   :alt: time = 2024-01-01
   :srcset: /examples/images/sphx_glr_plot_diags_prior_eddies_001.png
   :class: sphx-glr-single-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none


    <cartopy.mpl.contour.GeoContourSet object at 0x7026992a4e90>



.. GENERATED FROM PYTHON SOURCE LINES 44-45

Compute gestrophic current

.. GENERATED FROM PYTHON SOURCE LINES 45-47

.. code-block:: Python

    ugeos, vgeos = get_geos(ds.adt, dx=dx, dy=dy)








.. GENERATED FROM PYTHON SOURCE LINES 48-49

Compare them to dataset currents

.. GENERATED FROM PYTHON SOURCE LINES 49-65

.. code-block:: Python

    fig0, ax0 = create_map(ds.longitude, ds.latitude, figsize=(8, 5))
    kwqv = dict(units="dots", width=1, scale_units="dots", scale=1 / 20, transform=pcarr)
    ds.plot.quiver(x="longitude", y="latitude", u="ugos", v="vgos", color="k", ax=ax0, label="dataset", **kwqv)
    ax0.quiver(
        ds.longitude.values,
        ds.latitude.values,
        ugeos.values,
        vgeos.values,
        color="tab:orange",
        label="computed",
        **kwqv,
    )
    plt.title("Geostrophic currents")
    plt.legend()





.. image-sg:: /examples/images/sphx_glr_plot_diags_prior_eddies_002.png
   :alt: Geostrophic currents
   :srcset: /examples/images/sphx_glr_plot_diags_prior_eddies_002.png
   :class: sphx-glr-single-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none


    <matplotlib.legend.Legend object at 0x7026992b41d0>



.. GENERATED FROM PYTHON SOURCE LINES 66-72

Local normalized angular momentum
----------------------------------

The normalized angular momentum is computed at the center of 2D scanning window.

Window in km

.. GENERATED FROM PYTHON SOURCE LINES 72-74

.. code-block:: Python

    window_lnam = 50








.. GENERATED FROM PYTHON SOURCE LINES 75-76

Local normalized angular momentum

.. GENERATED FROM PYTHON SOURCE LINES 76-78

.. code-block:: Python

    lnam = get_lnam(ds.ugos, ds.vgos, window_lnam, dx=dx, dy=dy)








.. GENERATED FROM PYTHON SOURCE LINES 79-80

Plot

.. GENERATED FROM PYTHON SOURCE LINES 80-85

.. code-block:: Python

    fig1, ax1 = create_map(ds.longitude, ds.latitude, figsize=(8, 5))
    lnam.plot(cmap=cm.cm.diff, ax=ax1, add_colorbar=False, transform=pcarr)
    plt.title(f"Local angular momentum [{window_lnam}km]")





.. image-sg:: /examples/images/sphx_glr_plot_diags_prior_eddies_003.png
   :alt: Local angular momentum [50km]
   :srcset: /examples/images/sphx_glr_plot_diags_prior_eddies_003.png
   :class: sphx-glr-single-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none


    Text(0.5, 1.0, 'Local angular momentum [50km]')



.. GENERATED FROM PYTHON SOURCE LINES 86-90

Okubo-Weiss
-----------

Diagnostic

.. GENERATED FROM PYTHON SOURCE LINES 90-93

.. code-block:: Python

    ow = get_okuboweiss(ds.ugos, ds.vgos, dx=dx, dy=dy)









.. GENERATED FROM PYTHON SOURCE LINES 94-95

Plot

.. GENERATED FROM PYTHON SOURCE LINES 95-99

.. code-block:: Python

    fig2, ax2 = create_map(ds.longitude, ds.latitude, figsize=(8, 5))
    ow.plot(cmap=cm.cm.balance, ax=ax2, add_colorbar=False, transform=pcarr)
    ow.plot.contour(levels=[0], colors="k", transform=pcarr, ax=ax2)
    plt.title("Okubo-Weiss")



.. image-sg:: /examples/images/sphx_glr_plot_diags_prior_eddies_004.png
   :alt: Okubo-Weiss
   :srcset: /examples/images/sphx_glr_plot_diags_prior_eddies_004.png
   :class: sphx-glr-single-img


.. rst-class:: sphx-glr-script-out

 .. code-block:: none


    Text(0.5, 1.0, 'Okubo-Weiss')




.. rst-class:: sphx-glr-timing

   **Total running time of the script:** (0 minutes 1.217 seconds)


.. _sphx_glr_download_examples_plot_diags_prior_eddies.py:

.. only:: html

  .. container:: sphx-glr-footer sphx-glr-footer-example

    .. container:: sphx-glr-download sphx-glr-download-jupyter

      :download:`Download Jupyter notebook: plot_diags_prior_eddies.ipynb <plot_diags_prior_eddies.ipynb>`

    .. container:: sphx-glr-download sphx-glr-download-python

      :download:`Download Python source code: plot_diags_prior_eddies.py <plot_diags_prior_eddies.py>`

    .. container:: sphx-glr-download sphx-glr-download-zip

      :download:`Download zipped: plot_diags_prior_eddies.zip <plot_diags_prior_eddies.zip>`


.. only:: html

 .. rst-class:: sphx-glr-signature

    `Gallery generated by Sphinx-Gallery <https://sphinx-gallery.github.io>`_
