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

.. only:: html

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

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

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

.. _sphx_glr_advanced_examples_plot_slope_methods.py:


Slope and aspect methods
========================

Calculating terrain attributes—not only slope and aspect but also curvatures—requires estimating the
elevation derivatives of the surface. xDEM offers three different ways to
calculate elevation derivatives, which can result in slightly different results.

Here is an example of how to generate the two with each method, and understand their differences.

See also the :ref:`terrain-attributes` feature page.

**References:** `Horn (1981) <https://ieeexplore.ieee.org/document/1456186>`_, `Zevenbergen and Thorne (1987) <http://dx.doi.org/10.1002/esp.3290120107>`_, `Florinsky (2009) <https://doi.org/10.1080/13658810802527499>`_.

.. GENERATED FROM PYTHON SOURCE LINES 15-21

.. code-block:: Python


    import matplotlib.pyplot as plt
    import numpy as np

    import xdem








.. GENERATED FROM PYTHON SOURCE LINES 22-23

We open example data.

.. GENERATED FROM PYTHON SOURCE LINES 23-45

.. code-block:: Python

    dem = xdem.DEM(xdem.examples.get_path("longyearbyen_ref_dem"))


    def plot_attribute(attribute, cmap, label=None, vlim=None):

        if vlim is not None:
            if isinstance(vlim, (int, np.integer, float, np.floating)):
                vlims = {"vmin": -vlim, "vmax": vlim}
            elif len(vlim) == 2:
                vlims = {"vmin": vlim[0], "vmax": vlim[1]}
        else:
            vlims = {}

        attribute.plot(cmap=cmap, cbar_title=label, **vlims)

        plt.xticks([])
        plt.yticks([])
        plt.tight_layout()

        plt.show()









.. GENERATED FROM PYTHON SOURCE LINES 46-47

Slope with method of Horn (1981).

.. GENERATED FROM PYTHON SOURCE LINES 49-50

.. note::  (GDAL default), based on a refined approximation of the gradient (page 18, bottom left, and pages 20-21).

.. GENERATED FROM PYTHON SOURCE LINES 50-55

.. code-block:: Python


    slope_horn = xdem.terrain.slope(dem, surface_fit="Horn")

    plot_attribute(slope_horn, "Reds", "Slope of Horn (1981) (°)")




.. image-sg:: /advanced_examples/images/sphx_glr_plot_slope_methods_001.png
   :alt: plot slope methods
   :srcset: /advanced_examples/images/sphx_glr_plot_slope_methods_001.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 56-57

Slope with method of Zevenbergen and Thorne (1987), Equation 13.

.. GENERATED FROM PYTHON SOURCE LINES 57-62

.. code-block:: Python


    slope_zevenberg = xdem.terrain.slope(dem, surface_fit="ZevenbergThorne")

    plot_attribute(slope_zevenberg, "Reds", "Slope of Zevenberg and Thorne (1987) (°)")




.. image-sg:: /advanced_examples/images/sphx_glr_plot_slope_methods_002.png
   :alt: plot slope methods
   :srcset: /advanced_examples/images/sphx_glr_plot_slope_methods_002.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 63-64

Slope with method of Florinsky (2009).

.. GENERATED FROM PYTHON SOURCE LINES 64-70

.. code-block:: Python


    slope_florinsky = xdem.terrain.slope(dem, surface_fit="Florinsky")

    plot_attribute(slope_florinsky, "Reds", "Slope of Florinsky (2009) (°)")





.. image-sg:: /advanced_examples/images/sphx_glr_plot_slope_methods_003.png
   :alt: plot slope methods
   :srcset: /advanced_examples/images/sphx_glr_plot_slope_methods_003.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 71-72

We can compute the difference between the different slope computations - for instance, here between with Horn and Zevenberg methods.

.. GENERATED FROM PYTHON SOURCE LINES 72-77

.. code-block:: Python


    diff_slope = slope_horn - slope_zevenberg

    plot_attribute(diff_slope, "RdYlBu", "Slope of Horn (1981) minus\n slope of Zevenberg and Thorne (1987) (°)", vlim=3)




.. image-sg:: /advanced_examples/images/sphx_glr_plot_slope_methods_004.png
   :alt: plot slope methods
   :srcset: /advanced_examples/images/sphx_glr_plot_slope_methods_004.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 78-80

The differences are negative, implying that the method of Horn always provides flatter slopes.
Additionally, they seem to occur in places of high curvatures. We verify this by plotting the maximal curvature.

.. GENERATED FROM PYTHON SOURCE LINES 80-85

.. code-block:: Python


    maxc = xdem.terrain.max_curvature(dem)

    plot_attribute(maxc, "RdYlBu", "Maximal curvature (100 m $^{-1}$)", vlim=2)




.. image-sg:: /advanced_examples/images/sphx_glr_plot_slope_methods_005.png
   :alt: plot slope methods
   :srcset: /advanced_examples/images/sphx_glr_plot_slope_methods_005.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 86-88

We quantify the relationship by computing the median of slope differences in bins of curvatures, and plot the
result. We define custom bins for curvature, due to its skewed distribution.

.. GENERATED FROM PYTHON SOURCE LINES 88-106

.. code-block:: Python


    df_bin = xdem.spatialstats.nd_binning(
        values=diff_slope[:],
        list_var=[maxc[:]],
        list_var_names=["maxc"],
        list_var_bins=30,
        statistics=[np.nanmedian, "count"],
    )

    xdem.spatialstats.plot_1d_binning(
        df_bin,
        var_name="maxc",
        statistic_name="nanmedian",
        label_var="Maximal absolute curvature (100 m$^{-1}$)",
        label_statistic="Slope of Horn (1981) minus\n " "slope of Zevenberg and Thorne (1987) (°)",
    )





.. image-sg:: /advanced_examples/images/sphx_glr_plot_slope_methods_006.png
   :alt: plot slope methods
   :srcset: /advanced_examples/images/sphx_glr_plot_slope_methods_006.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 107-109

We perform the same exercise to analyze the differences in terrain aspect. We compute the difference modulo 360°,
to account for the circularity of aspect.

.. GENERATED FROM PYTHON SOURCE LINES 109-120

.. code-block:: Python


    aspect_horn = xdem.terrain.aspect(dem)
    aspect_zevenberg = xdem.terrain.aspect(dem, method="ZevenbergThorne")

    diff_aspect = aspect_horn - aspect_zevenberg
    diff_aspect_mod = np.minimum(diff_aspect % 360, 360 - diff_aspect % 360)

    plot_attribute(
        diff_aspect_mod, "Spectral", "Aspect of Horn (1981) minus\n aspect of Zevenberg and Thorne (1987) (°)", vlim=[0, 90]
    )




.. image-sg:: /advanced_examples/images/sphx_glr_plot_slope_methods_007.png
   :alt: plot slope methods
   :srcset: /advanced_examples/images/sphx_glr_plot_slope_methods_007.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 121-124

Same as for slope, differences in aspect seem to coincide with high curvature areas. We observe also observe large
differences for areas with nearly flat slopes, owing to the high sensitivity of orientation estimation
for flat terrain.

.. GENERATED FROM PYTHON SOURCE LINES 126-127

.. note:: The default aspect for a 0° slope is 180°, as in GDAL.


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

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


.. _sphx_glr_download_advanced_examples_plot_slope_methods.py:

.. only:: html

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

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

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

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

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

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

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


.. only:: html

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

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