.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "examples/selection.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_examples_selection.py: Select unstructured data ======================== Xarray has flexible tools for label based selection, in the form of ``.sel`` and ``.isel`` for index selection. This works well for structured data since the orthogonality of the x and y axes is reflected in the axes of the underlying arrays. This orthogonality does not exist for unstructured grids, as the data for all faces cannot be stored in a two-dimensional array and is stored in a one-dimensional array instead. Xugrid provides tools for convenient spatial selection, primarily via the ``.ugrid.sel`` method; its behavior is comparable to xarray's ``.sel`` method. The ``.ugrid.sel`` method should only be used for selection in the x or y dimension. Selections along other dimension (such as time) should be performed by xarray's ``.sel`` instead (without the ``ugrid`` accessor). The examples below demonstrate the various ways to select data. Imports ------- The following imports suffice for the examples. .. GENERATED FROM PYTHON SOURCE LINES 27-33 .. code-block:: Python import matplotlib.pyplot as plt import numpy as np import shapely import xugrid as xu .. GENERATED FROM PYTHON SOURCE LINES 34-36 We will take a look at a sample dataset: a triangular grid with the surface elevation of the Netherlands. .. GENERATED FROM PYTHON SOURCE LINES 36-40 .. code-block:: Python uda = xu.data.elevation_nl() uda.ugrid.plot(vmin=-20, vmax=90, cmap="terrain") .. image-sg:: /examples/images/sphx_glr_selection_001.png :alt: selection :srcset: /examples/images/sphx_glr_selection_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 41-57 We will start by demonstrating the behavior of ``.ugrid.sel``. This method takes several types of arguments, like its xarray equivalent. The return type and shape of the selection operation depends on the argument given. ========== =========== Selection Result type ========== =========== Subset xugrid Point xarray Line xarray ========== =========== Grid subset selection --------------------- A subset of the unstructured grid is returned by using slices without a step: .. GENERATED FROM PYTHON SOURCE LINES 57-61 .. code-block:: Python subset = uda.ugrid.sel(x=slice(100_000.0, 200_000.0), y=slice(450_000.0, 550_000.0)) subset.ugrid.plot(vmin=-20, vmax=90, cmap="terrain") .. image-sg:: /examples/images/sphx_glr_selection_002.png :alt: selection :srcset: /examples/images/sphx_glr_selection_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 62-64 The default arguments of ``x`` and ``y`` are: ``slice(None, None)``. In such a case the entire grid is returned. .. GENERATED FROM PYTHON SOURCE LINES 64-68 .. code-block:: Python subset = uda.ugrid.sel() subset.ugrid.plot(vmin=-20, vmax=90, cmap="terrain") .. image-sg:: /examples/images/sphx_glr_selection_003.png :alt: selection :srcset: /examples/images/sphx_glr_selection_003.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 69-75 .. note:: ``None`` in a Python slice can be interpreted as "from the start" or "up to and including the end". This means we can easily select along a single dimension: .. GENERATED FROM PYTHON SOURCE LINES 75-79 .. code-block:: Python subset = uda.ugrid.sel(x=slice(100_000.0, 200_000.0)) subset.ugrid.plot(vmin=-20, vmax=90, cmap="terrain", aspect=1, size=5) .. image-sg:: /examples/images/sphx_glr_selection_004.png :alt: selection :srcset: /examples/images/sphx_glr_selection_004.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 80-81 Or, using ``None`` if we only care about the start: .. GENERATED FROM PYTHON SOURCE LINES 81-85 .. code-block:: Python subset = uda.ugrid.sel(x=slice(100_000.0, None)) subset.ugrid.plot(vmin=-20, vmax=90, cmap="terrain", aspect=1, size=5) .. image-sg:: /examples/images/sphx_glr_selection_005.png :alt: selection :srcset: /examples/images/sphx_glr_selection_005.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 86-94 Point selection --------------- Since point data can be represented as an ordinary xarray DataArray with x and y coordinates, all point selection result in xarray DataArrays rather than UgridDataArrays with an associated unstructured grid topology. We will use a utility function to show what is selected on the map: .. GENERATED FROM PYTHON SOURCE LINES 94-103 .. code-block:: Python def show_point_selection(uda, da): _, ax = plt.subplots() uda.ugrid.plot(ax=ax, vmin=-20, vmax=90, cmap="terrain") ax.scatter(da["mesh2d_x"], da["mesh2d_y"], color="red") ax.set_aspect(1.0) .. GENERATED FROM PYTHON SOURCE LINES 104-105 Two values will select a point: .. GENERATED FROM PYTHON SOURCE LINES 105-110 .. code-block:: Python da = uda.ugrid.sel(x=150_000.0, y=463_000.0) show_point_selection(uda, da) da .. image-sg:: /examples/images/sphx_glr_selection_006.png :alt: selection :srcset: /examples/images/sphx_glr_selection_006.png :class: sphx-glr-single-img .. raw:: html
<xarray.DataArray 'elevation' (mesh2d_nFaces: 1)> Size: 4B
    [1 values with dtype=float32]
    Coordinates:
      * mesh2d_nFaces  (mesh2d_nFaces) int64 8B 4221
        mesh2d_face_x  (mesh2d_nFaces) float64 8B ...
        mesh2d_face_y  (mesh2d_nFaces) float64 8B ...
        mesh2d_index   (mesh2d_nFaces) int64 8B 0
        mesh2d_x       (mesh2d_nFaces) float64 8B 1.5e+05
        mesh2d_y       (mesh2d_nFaces) float64 8B 4.63e+05
    Attributes:
        unit:     m NAP


.. GENERATED FROM PYTHON SOURCE LINES 111-114 Multiple values are broadcasted against each other ("outer indexing"). If we select by three x values and two y values, the result is a collection of six points: .. GENERATED FROM PYTHON SOURCE LINES 114-119 .. code-block:: Python da = uda.ugrid.sel(x=[125_000.0, 150_000.0, 175_000.0], y=[400_000.0, 465_000.0]) show_point_selection(uda, da) da .. image-sg:: /examples/images/sphx_glr_selection_007.png :alt: selection :srcset: /examples/images/sphx_glr_selection_007.png :class: sphx-glr-single-img .. raw:: html
<xarray.DataArray 'elevation' (mesh2d_nFaces: 6)> Size: 24B
    [6 values with dtype=float32]
    Coordinates:
      * mesh2d_nFaces  (mesh2d_nFaces) int64 48B 1905 1356 372 3057 4198 4113
        mesh2d_face_x  (mesh2d_nFaces) float64 48B ...
        mesh2d_face_y  (mesh2d_nFaces) float64 48B ...
        mesh2d_index   (mesh2d_nFaces) int64 48B 0 1 2 3 4 5
        mesh2d_x       (mesh2d_nFaces) float64 48B 1.25e+05 1.5e+05 ... 1.75e+05
        mesh2d_y       (mesh2d_nFaces) float64 48B 4e+05 4e+05 ... 4.65e+05 4.65e+05
    Attributes:
        unit:     m NAP


.. GENERATED FROM PYTHON SOURCE LINES 120-121 To select points without broadcasting, use ``.ugrid.sel_points`` instead: .. GENERATED FROM PYTHON SOURCE LINES 121-128 .. code-block:: Python da = uda.ugrid.sel_points( x=[125_000.0, 150_000.0, 175_000.0], y=[400_000.0, 430_000.0, 465_000.0] ) show_point_selection(uda, da) da .. image-sg:: /examples/images/sphx_glr_selection_008.png :alt: selection :srcset: /examples/images/sphx_glr_selection_008.png :class: sphx-glr-single-img .. raw:: html
<xarray.DataArray 'elevation' (mesh2d_nFaces: 3)> Size: 12B
    [3 values with dtype=float32]
    Coordinates:
      * mesh2d_nFaces  (mesh2d_nFaces) int64 24B 1905 2675 4113
        mesh2d_face_x  (mesh2d_nFaces) float64 24B ...
        mesh2d_face_y  (mesh2d_nFaces) float64 24B ...
        mesh2d_index   (mesh2d_nFaces) int64 24B 0 1 2
        mesh2d_x       (mesh2d_nFaces) float64 24B 1.25e+05 1.5e+05 1.75e+05
        mesh2d_y       (mesh2d_nFaces) float64 24B 4e+05 4.3e+05 4.65e+05
    Attributes:
        unit:     m NAP


.. GENERATED FROM PYTHON SOURCE LINES 129-130 We can sample points along a line as well by providing slices **with** a step: .. GENERATED FROM PYTHON SOURCE LINES 130-135 .. code-block:: Python da = uda.ugrid.sel(x=slice(100_000.0, 200_000.0, 10_000.0), y=465_000.0) show_point_selection(uda, da) da .. image-sg:: /examples/images/sphx_glr_selection_009.png :alt: selection :srcset: /examples/images/sphx_glr_selection_009.png :class: sphx-glr-single-img .. raw:: html
<xarray.DataArray 'elevation' (mesh2d_nFaces: 10)> Size: 40B
    [10 values with dtype=float32]
    Coordinates:
      * mesh2d_nFaces  (mesh2d_nFaces) int64 80B 3263 3225 3144 ... 2941 3135 4256
        mesh2d_face_x  (mesh2d_nFaces) float64 80B ...
        mesh2d_face_y  (mesh2d_nFaces) float64 80B ...
        mesh2d_index   (mesh2d_nFaces) int64 80B 0 1 2 3 4 5 6 7 8 9
        mesh2d_x       (mesh2d_nFaces) float64 80B 1e+05 1.1e+05 ... 1.8e+05 1.9e+05
        mesh2d_y       (mesh2d_nFaces) float64 80B 4.65e+05 4.65e+05 ... 4.65e+05
    Attributes:
        unit:     m NAP


.. GENERATED FROM PYTHON SOURCE LINES 136-137 Two slices with a step results in broadcasting: .. GENERATED FROM PYTHON SOURCE LINES 137-144 .. code-block:: Python da = uda.ugrid.sel( x=slice(100_000.0, 200_000.0, 10_000.0), y=slice(400_000.0, 500_000.0, 10_000.0) ) show_point_selection(uda, da) da .. image-sg:: /examples/images/sphx_glr_selection_010.png :alt: selection :srcset: /examples/images/sphx_glr_selection_010.png :class: sphx-glr-single-img .. raw:: html
<xarray.DataArray 'elevation' (mesh2d_nFaces: 100)> Size: 400B
    [100 values with dtype=float32]
    Coordinates:
      * mesh2d_nFaces  (mesh2d_nFaces) int64 800B 1867 1882 1890 ... 4121 4077 1581
        mesh2d_face_x  (mesh2d_nFaces) float64 800B ...
        mesh2d_face_y  (mesh2d_nFaces) float64 800B ...
        mesh2d_index   (mesh2d_nFaces) int64 800B 0 1 2 3 4 5 ... 94 95 96 97 98 99
        mesh2d_x       (mesh2d_nFaces) float64 800B 1e+05 1.1e+05 ... 1.9e+05
        mesh2d_y       (mesh2d_nFaces) float64 800B 4e+05 4e+05 ... 4.9e+05 4.9e+05
    Attributes:
        unit:     m NAP


.. GENERATED FROM PYTHON SOURCE LINES 145-146 As well as a slice with a step and multiple values: .. GENERATED FROM PYTHON SOURCE LINES 146-151 .. code-block:: Python da = uda.ugrid.sel(x=slice(100_000.0, 200_000.0, 10_000.0), y=[400_000.0, 430_000.0]) show_point_selection(uda, da) da .. image-sg:: /examples/images/sphx_glr_selection_011.png :alt: selection :srcset: /examples/images/sphx_glr_selection_011.png :class: sphx-glr-single-img .. raw:: html
<xarray.DataArray 'elevation' (mesh2d_nFaces: 20)> Size: 80B
    [20 values with dtype=float32]
    Coordinates:
      * mesh2d_nFaces  (mesh2d_nFaces) int64 160B 1867 1882 1890 ... 2889 226 3035
        mesh2d_face_x  (mesh2d_nFaces) float64 160B ...
        mesh2d_face_y  (mesh2d_nFaces) float64 160B ...
        mesh2d_index   (mesh2d_nFaces) int64 160B 0 1 2 3 4 5 ... 14 15 16 17 18 19
        mesh2d_x       (mesh2d_nFaces) float64 160B 1e+05 1.1e+05 ... 1.9e+05
        mesh2d_y       (mesh2d_nFaces) float64 160B 4e+05 4e+05 ... 4.3e+05 4.3e+05
    Attributes:
        unit:     m NAP


.. GENERATED FROM PYTHON SOURCE LINES 152-163 Line selection -------------- Since line data can be represented as an ordinary xarray DataArray with x and y coordinates, all line selection result in xarray DataArrays rather than UgridDataArrays with an associated unstructured grid topology. Line selection is performed by finding all faces that are intersected by the line. We start by defining a utility to show the selection again: .. GENERATED FROM PYTHON SOURCE LINES 163-178 .. code-block:: Python def show_line_selection(uda, da, line_x=None, line_y=None): _, (ax0, ax1) = plt.subplots(ncols=2, figsize=(10, 5)) uda.ugrid.plot(ax=ax0, vmin=-20, vmax=90, cmap="terrain") da.plot(ax=ax1, x="mesh2d_s") if line_x is None: ax0.axhline(line_y, color="red") elif line_y is None: ax0.axvline(line_x, color="red") else: ax0.plot(line_x, line_y, color="red") ax0.set_aspect(1.0) .. GENERATED FROM PYTHON SOURCE LINES 179-181 A single value for either x or y in ``.ugrid.sel`` will select values along a line: .. GENERATED FROM PYTHON SOURCE LINES 181-185 .. code-block:: Python da = uda.ugrid.sel(y=465_000.0) show_line_selection(uda, da, line_y=465_000.0) .. image-sg:: /examples/images/sphx_glr_selection_012.png :alt: selection :srcset: /examples/images/sphx_glr_selection_012.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 186-188 Line segments that are not axis aligned can be selected with ``.ugrid.intersect_line``: .. GENERATED FROM PYTHON SOURCE LINES 188-192 .. code-block:: Python da = uda.ugrid.intersect_line(start=(60_000.0, 400_000.0), end=(190_000.0, 475_000.0)) show_line_selection(uda, da, (60_000.0, 190_000.0), (400_000.0, 475_000.0)) .. image-sg:: /examples/images/sphx_glr_selection_013.png :alt: selection :srcset: /examples/images/sphx_glr_selection_013.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 193-194 Linestrings can be selected with ``.ugrid.intersect_linestring``: .. GENERATED FROM PYTHON SOURCE LINES 194-207 .. code-block:: Python linestring = shapely.geometry.LineString( [ (60_000.0, 400_000.0), (190_000.0, 400_000.0), (120_000.0, 575_000.0), (250_000.0, 575_000.0), ] ) da = uda.ugrid.intersect_linestring(linestring) show_line_selection(uda, da, *shapely.get_coordinates(linestring).T) .. image-sg:: /examples/images/sphx_glr_selection_014.png :alt: selection :srcset: /examples/images/sphx_glr_selection_014.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 208-209 This will work for any type of shapely line: .. GENERATED FROM PYTHON SOURCE LINES 209-214 .. code-block:: Python ring = shapely.geometry.Point(155_000.0, 463_000).buffer(50_000.0).exterior da = uda.ugrid.intersect_linestring(ring) show_line_selection(uda, da, *shapely.get_coordinates(ring).T) .. image-sg:: /examples/images/sphx_glr_selection_015.png :alt: selection :srcset: /examples/images/sphx_glr_selection_015.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 215-221 Index selection --------------- We may also use ordinary index selection to create a subset. This does not require the ``.ugrid`` accessor. For example, to take only the first thousands faces: .. GENERATED FROM PYTHON SOURCE LINES 221-225 .. code-block:: Python subset = uda.isel(mesh2d_nFaces=np.arange(1000)) subset.ugrid.plot(vmin=-20, vmax=90, cmap="terrain", aspect=1, size=5) .. image-sg:: /examples/images/sphx_glr_selection_016.png :alt: selection :srcset: /examples/images/sphx_glr_selection_016.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 226-238 For a 2D topology, selecting faces by an index always results in a valid topology. However, selecting by node or edge does not give a guarantee that the result forms a valid 2D topology: e.g. if we only select two nodes, or only two edges from a face, the result cannot form a valid 2D face. To avoid generating invalid topologies, xugrid always checks whether the result of a selection results in a valid 2D topology and raises an error if the result is invalid. In general, index selection should only be performed on the "core" dimension of the UGRID topology. This is the edge dimension for 1D topologies, and the face dimension for 2D topologies. .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 3.324 seconds) .. _sphx_glr_download_examples_selection.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: selection.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: selection.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: selection.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_