
.. DO NOT EDIT.
.. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY.
.. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE:
.. "examples/plot_croco_section.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_croco_section.py>`
        to download the full example code or to run this example in your browser via Binder.

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

.. _sphx_glr_examples_plot_croco_section.py:


Interpolate a meridional section of CROCO outputs to regular depths
===================================================================

In this notebook, we show:

* how to compute the depths from s-coordinates,
* how to easily find the name of variables and coordinates,
* how to interpolate a 3D field with varying depths to regular depths,
* how to compute the mixed layer depth from temperature.

.. GENERATED FROM PYTHON SOURCE LINES 16-20

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

Import needed modules.

.. GENERATED FROM PYTHON SOURCE LINES 20-32

.. code-block:: Python


    import numpy as np
    import xarray as xr
    import matplotlib.pyplot as plt
    import cmocean  # noqa
    import xoa
    from xoa.regrid import regrid1d
    from xoa.thermdyn import mixed_layer_depth
    import xoa.meta as xmeta

    xr.set_options(display_style="text")





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

 .. code-block:: none


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



.. GENERATED FROM PYTHON SOURCE LINES 33-34

Register the :meth:`xarray.Dataset.decode_sigma` callable accessor.

.. GENERATED FROM PYTHON SOURCE LINES 34-37

.. code-block:: Python


    xoa.register_accessors(decode_sigma=True)








.. GENERATED FROM PYTHON SOURCE LINES 38-40

The :ref:`xoa <accessors>` accessor is also registered by default, and give access
to most of the functionalities of the other accessors.

.. GENERATED FROM PYTHON SOURCE LINES 42-43

Register the internal CROCO naming specifications

.. GENERATED FROM PYTHON SOURCE LINES 43-48

.. code-block:: Python


    croco_cfg_file = xoa.get_meta_config_file("croco")
    print(croco_cfg_file)
    xmeta.register_meta_specs(croco_cfg_file)





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

 .. code-block:: none

    /home/docs/checkouts/readthedocs.org/user_builds/xoa/conda/develop/lib/python3.14/site-packages/xoa/meta/configs/croco.cfg




.. GENERATED FROM PYTHON SOURCE LINES 49-55

You can provide your own configuration file.

In this way, the :mod:`xoa.meta` module will recognise the
CROCO netcdf names.
It would be equivalent to force loading name specs with
`xoa.meta.set_meta_specs(xoa.get_meta_config_file("croco"))`.

.. GENERATED FROM PYTHON SOURCE LINES 57-60

Read the model
--------------
This sample is a meridional extraction of a full 3D CROCO output.

.. GENERATED FROM PYTHON SOURCE LINES 60-65

.. code-block:: Python

    sample_file = xoa.get_data_sample("MODELS/CROCO/SOUTH-AFRICA/croco.south-africa.meridional.nc")
    print(sample_file)
    ds = xr.open_dataset(sample_file)
    print(ds)





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

 .. code-block:: none

    /home/docs/.cache/xoa/MODELS/CROCO/SOUTH-AFRICA/croco.south-africa.meridional.nc
    <xarray.Dataset> Size: 48kB
    Dimensions:     (time: 1, s_w: 33, eta_rho: 56, xi_rho: 1, s_rho: 32, auxil: 4,
                     xi_u: 1, eta_v: 55)
    Coordinates: (12/13)
      * time        (time) float64 8B 2.592e+06
      * s_w         (s_w) float32 132B -1.0 -0.9688 -0.9375 ... -0.0625 -0.03125 0.0
      * eta_rho     (eta_rho) float32 224B 6.0 7.0 8.0 9.0 ... 58.0 59.0 60.0 61.0
      * xi_rho      (xi_rho) float32 4B 131.0
        lat_rho     (eta_rho, xi_rho) float32 224B ...
        lat_v       (eta_v, xi_rho) float32 220B ...
        ...          ...
        lon_v       (eta_v, xi_rho) float32 220B ...
      * s_rho       (s_rho) float32 128B -0.9844 -0.9531 ... -0.04688 -0.01562
      * xi_u        (xi_u) float32 4B 131.5
        lat_u       (eta_rho, xi_u) float32 224B ...
        lon_u       (eta_rho, xi_u) float32 224B ...
      * eta_v       (eta_v) float32 220B 6.5 7.5 8.5 9.5 ... 57.5 58.5 59.5 60.5
    Dimensions without coordinates: auxil
    Data variables: (12/23)
        AKt         (time, s_w, eta_rho, xi_rho) float32 7kB ...
        Cs_r        (s_rho) float32 128B ...
        Cs_w        (s_w) float32 132B ...
        Vtransform  float32 4B ...
        angle       (eta_rho, xi_rho) float32 224B ...
        el          float32 4B ...
        ...          ...
        time_step   (time, auxil) int32 16B ...
        u           (time, s_rho, eta_rho, xi_u) float32 7kB ...
        v           (time, s_rho, eta_v, xi_rho) float32 7kB ...
        w           (time, s_rho, eta_rho, xi_rho) float32 7kB ...
        xl          float32 4B ...
        zeta        (time, eta_rho, xi_rho) float32 224B ...
    Attributes: (12/57)
        type:           ROMS history file
        title:          BENGUELA TEST MODEL
        date:           
        rst_file:       CROCO_FILES/croco_rst.nc
        his_file:       CROCO_FILES/croco_his.nc
        avg_file:       CROCO_FILES/croco_avg.nc
        ...             ...
        v_sponge:       0.0
        sponge_expl:    Sponge parameters : extent (m) & viscosity (m2.s-1)
        SRCS:           main.F step.F read_inp.F timers_roms.F init_scalars.F ini...
        CPP-options:    REGIONAL BENGUELA_VHR MPI OBC_EAST OBC_WEST OBC_NORTH OBC...
        history:        Tue Mar 31 16:26:24 2020: ncks -O -d time,30 -d xi_rho,13...
        NCO:            4.4.2




.. GENERATED FROM PYTHON SOURCE LINES 66-77

Compute depths from s-coordinates
---------------------------------

Decode the dataset according to the CF conventions:

1. Find sigma terms
2. Compute depths
3. Assign depths as coordinates

Note that the :meth:`xarray.Dataset.decode_sigma` callable accessor
calls the :func:`xoa.sigma.decode_cf_sigma` function.

.. GENERATED FROM PYTHON SOURCE LINES 77-81

.. code-block:: Python


    ds = ds.decode_sigma()
    print(ds.depth)





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

 .. code-block:: none

    <xarray.DataArray 'depth' (time: 1, s_rho: 32, eta_rho: 56, xi_rho: 1)> Size: 14kB
    array([[[[-4.47877110e+03],
             [-4.45341971e+03],
             [-4.43332704e+03],
             ...,
             [-7.38942221e+01],
             [-7.34145067e+01],
             [-7.34091375e+01]],

            [[-4.10998820e+03],
             [-4.08677643e+03],
             [-4.06837911e+03],
             ...,
             [-7.04949884e+01],
             [-7.00423304e+01],
             [-7.00372638e+01]],

            [[-3.70399833e+03],
             [-3.68314549e+03],
             [-3.66661727e+03],
             ...,
    ...
             ...,
             [-4.29529272e+00],
             [-4.27484202e+00],
             [-4.27461272e+00]],

            [[-1.00794918e+01],
             [-1.00480588e+01],
             [-1.00177861e+01],
             ...,
             [-2.57388434e+00],
             [-2.56165088e+00],
             [-2.56151371e+00]],

            [[-3.17876268e+00],
             [-3.15406668e+00],
             [-3.12913607e+00],
             ...,
             [-8.56901999e-01],
             [-8.52836074e-01],
             [-8.52790485e-01]]]], shape=(1, 32, 56, 1))
    Coordinates:
      * time     (time) float64 8B 2.592e+06
      * s_rho    (s_rho) float32 128B -0.9844 -0.9531 -0.9219 ... -0.04688 -0.01562
      * eta_rho  (eta_rho) float32 224B 6.0 7.0 8.0 9.0 10.0 ... 58.0 59.0 60.0 61.0
      * xi_rho   (xi_rho) float32 4B 131.0
        lat_rho  (eta_rho, xi_rho) float32 224B -37.67 -37.61 ... -34.03 -33.96
        lon_rho  (eta_rho, xi_rho) float32 224B 18.83 18.83 18.83 ... 18.83 18.83
        depth    (time, s_rho, eta_rho, xi_rho) float64 14kB -4.479e+03 ... -0.8528
    Attributes:
        long_name:      Ocean s roms coordinate at rho point
        Vtransform:     2
        standard_name:  ocean_layer_depth
        units:          m




.. GENERATED FROM PYTHON SOURCE LINES 82-89

Find coordinate names from CF conventions
-----------------------------------------

The `depth` was assigned as coordinates at the previous stage.
We use the :ref:`xoa <accessors.dataset>` accessor to easily access the temperature, latitude and depth arrays.
The default configuration exposes shortcuts for some variables and coordinates
as shown in :metasec:`accessors`.

.. GENERATED FROM PYTHON SOURCE LINES 89-94

.. code-block:: Python


    temp = ds.xoa.temp.squeeze()
    temp = temp.where(temp != 0)  # convert zeros to nans
    lat_name = temp.xoa.lat.name








.. GENERATED FROM PYTHON SOURCE LINES 95-101

Interpolate at regular depths
-----------------------------

We interpolate the temperature array from irregular to regular depths.

Let's create the output depths.

.. GENERATED FROM PYTHON SOURCE LINES 101-106

.. code-block:: Python


    depth = xr.DataArray(
        np.linspace(ds.depth.values.min(), ds.depth.values.max(), 1000), name="depth", dims="depth"
    )








.. GENERATED FROM PYTHON SOURCE LINES 107-108

Let's interpolate the temperature.

.. GENERATED FROM PYTHON SOURCE LINES 108-111

.. code-block:: Python


    tempz = regrid1d(temp, depth, extrap="top")








.. GENERATED FROM PYTHON SOURCE LINES 112-118

Compute the mixed layer depths
-------------------------------

The mixed layer depths are computed here as the depth at which the temperature
is `deltatemp` lower than the surface temperature,
thanks to the :func:`xoa.thermdyn.mixed_layer_depth` function.

.. GENERATED FROM PYTHON SOURCE LINES 118-123

.. code-block:: Python


    deltatemp = 0.2
    mld = -mixed_layer_depth(temp, deltatemp=deltatemp)
    mldz = -mixed_layer_depth(temp, deltatemp=deltatemp)








.. GENERATED FROM PYTHON SOURCE LINES 124-126

Plots
-----

.. GENERATED FROM PYTHON SOURCE LINES 128-129

Plot the full section.

.. GENERATED FROM PYTHON SOURCE LINES 129-139

.. code-block:: Python


    fig, axs = plt.subplots(ncols=2, sharex=True, sharey=True, figsize=(10, 4))
    kw = dict(levels=np.arange(0, 23))
    temp.plot.contourf(lat_name, "depth", cmap="cmo.thermal", ax=axs[0], **kw)
    temp.plot.contour(lat_name, "depth", colors='k', linewidths=0.3, ax=axs[0], **kw)
    tempz.plot.contourf(lat_name, "depth", cmap="cmo.thermal", ax=axs[1], **kw)
    tempz.plot.contour(lat_name, "depth", colors='k', linewidths=0.3, ax=axs[1], **kw)
    axs[0].set_title("Original")
    axs[1].set_title("Interpolated")




.. image-sg:: /examples/images/sphx_glr_plot_croco_section_001.png
   :alt: Original, Interpolated
   :srcset: /examples/images/sphx_glr_plot_croco_section_001.png
   :class: sphx-glr-single-img


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

 .. code-block:: none

    /home/docs/checkouts/readthedocs.org/user_builds/xoa/conda/develop/lib/python3.14/site-packages/xarray/plot/accessor.py:663: FutureWarning: Using positional arguments is deprecated for plot methods, use keyword arguments instead.
      return dataarray_plot.contourf(self._da, *args, **kwargs)
    /home/docs/checkouts/readthedocs.org/user_builds/xoa/conda/develop/lib/python3.14/site-packages/xarray/plot/accessor.py:542: FutureWarning: Using positional arguments is deprecated for plot methods, use keyword arguments instead.
      return dataarray_plot.contour(self._da, *args, **kwargs)
    /home/docs/checkouts/readthedocs.org/user_builds/xoa/conda/develop/lib/python3.14/site-packages/xarray/plot/accessor.py:663: FutureWarning: Using positional arguments is deprecated for plot methods, use keyword arguments instead.
      return dataarray_plot.contourf(self._da, *args, **kwargs)
    /home/docs/checkouts/readthedocs.org/user_builds/xoa/conda/develop/lib/python3.14/site-packages/xarray/plot/accessor.py:542: FutureWarning: Using positional arguments is deprecated for plot methods, use keyword arguments instead.
      return dataarray_plot.contour(self._da, *args, **kwargs)

    Text(0.5, 1.0, 'Interpolated')



.. GENERATED FROM PYTHON SOURCE LINES 140-141

Plot a zoom near the surface and add the mixed layer depth isoline.

.. GENERATED FROM PYTHON SOURCE LINES 141-154

.. code-block:: Python


    fig, axs = plt.subplots(ncols=2, sharex=True, sharey=True, figsize=(10, 3))
    kw = dict(levels=np.arange(0, 23))
    temp.plot.contourf(lat_name, "depth", cmap="cmo.thermal", ax=axs[0], **kw)
    temp.plot.contour(lat_name, "depth", colors='k', linewidths=0.3, ax=axs[0], **kw)
    mld.plot.line(x=lat_name, color="k", linewidth=2, linestyle="--", ax=axs[0])
    tempz.plot.contourf(lat_name, "depth", cmap="cmo.thermal", ax=axs[1], **kw)
    tempz.plot.contour(lat_name, "depth", colors='k', linewidths=0.3, ax=axs[1], **kw)
    mldz.plot.line(x=lat_name, color="k", linewidth=2, linestyle="--", ax=axs[1])
    axs[0].set_ylim(-300, 0)
    axs[0].set_title("Original")
    axs[1].set_title("Interpolated")




.. image-sg:: /examples/images/sphx_glr_plot_croco_section_002.png
   :alt: Original, Interpolated
   :srcset: /examples/images/sphx_glr_plot_croco_section_002.png
   :class: sphx-glr-single-img


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

 .. code-block:: none


    Text(0.5, 1.0, 'Interpolated')



.. GENERATED FROM PYTHON SOURCE LINES 155-156

Et voilà!


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

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


.. _sphx_glr_download_examples_plot_croco_section.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/shom-fr/xoa/master?urlpath=lab/tree/notebooks/examples/plot_croco_section.ipynb
        :alt: Launch binder
        :width: 150 px

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

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

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

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

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

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


.. only:: html

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

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