.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "examples/visualize/flowvel_streamlines.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_visualize_flowvel_streamlines.py: Flow velocities and streamlines =============================== In this section we will plot flow velocities and streamlines for some model results. .. GENERATED FROM PYTHON SOURCE LINES 8-13 .. code-block:: default import matplotlib.pyplot as plt import numpy as np import xarray as xr .. GENERATED FROM PYTHON SOURCE LINES 14-15 We'll start with the usual imports .. GENERATED FROM PYTHON SOURCE LINES 15-18 .. code-block:: default import imod .. GENERATED FROM PYTHON SOURCE LINES 20-21 Load and unpack the data .. GENERATED FROM PYTHON SOURCE LINES 21-31 .. code-block:: default ds_fluxes = imod.data.fluxes() ds_fluxes = ds_fluxes.isel(time=-1) ds_fluxes lower = ds_fluxes["bdgflf"] right = ds_fluxes["bdgfrf"] front = ds_fluxes["bdgfff"] heads = ds_fluxes["head"] .. GENERATED FROM PYTHON SOURCE LINES 32-39 Calculating flow velocity ------------------------- The imod-python function imod.evaluate.flow_velocity() computes flow velocities in m/d based on the budget results (bdgflf - flow lower face, bdgfrf - flow right face and bdgfff - flow front face). To apply this function, we first need to define a top_bot array. .. GENERATED FROM PYTHON SOURCE LINES 39-45 .. code-block:: default top_bot = xr.full_like(lower, 1.0) top_bot["top"] = top_bot["z"] - 0.5 * top_bot["dz"] top_bot["bot"] = top_bot["z"] + 0.5 * top_bot["dz"] top_bot .. raw:: html
<xarray.DataArray 'bdgflf' (layer: 49, y: 218, x: 248)>
    array([[[1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            ...,
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.]],

           [[1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            ...,
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.]],

           [[1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            ...,
    ...
            ...,
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.]],

           [[1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            ...,
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.]],

           [[1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            ...,
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.],
            [1., 1., 1., ..., 1., 1., 1.]]], dtype=float32)
    Coordinates:
      * x        (x) float64 9.095e+04 9.105e+04 9.115e+04 ... 1.156e+05 1.156e+05
      * y        (y) float64 4.676e+05 4.674e+05 4.674e+05 ... 4.46e+05 4.458e+05
        dx       float64 100.0
        dy       float64 -100.0
        time     datetime64[ns] 1979-12-31
      * layer    (layer) int32 1 2 3 4 5 6 7 8 9 10 ... 41 42 43 44 45 46 47 48 49
        dz       (layer) float64 -2.0 -2.0 -2.0 -2.0 ... -10.0 -10.0 -10.0 -10.0
        z        (layer) float64 9.0 7.0 5.0 3.0 1.0 ... -215.0 -225.0 -235.0 -245.0
        top      (layer) float64 10.0 8.0 6.0 4.0 ... -210.0 -220.0 -230.0 -240.0
        bot      (layer) float64 8.0 6.0 4.0 2.0 0.0 ... -220.0 -230.0 -240.0 -250.0


.. GENERATED FROM PYTHON SOURCE LINES 46-48 Next we'll calculate the velocities and plot a cross-section of the vertical velocity at for y = 450050 .. GENERATED FROM PYTHON SOURCE LINES 48-57 .. code-block:: default fig, ax = plt.subplots() vx, vy, vz = imod.evaluate.flow_velocity( front=front, lower=lower, right=right, top_bot=top_bot, porosity=0.3 ) vz.sel(y=450050.0, method="nearest").plot(cmap="RdYlBu", yincrease=False) plt.title("Vz") .. image-sg:: /examples/visualize/images/sphx_glr_flowvel_streamlines_001.png :alt: Vz :srcset: /examples/visualize/images/sphx_glr_flowvel_streamlines_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none /opt/conda/envs/imod/lib/python3.10/site-packages/xarray/core/dataarray.py:885: FutureWarning: elementwise comparison failed; returning scalar instead, but in the future will perform elementwise comparison return key in self.data Text(0.5, 1.0, 'Vz') .. GENERATED FROM PYTHON SOURCE LINES 58-64 Quiver plot ----------- It is also possible to make quiver plots for a cross section defined by two pairs of coordinates. We will first define arrays indicating the start and end location of the cross section to be evaluated. .. GENERATED FROM PYTHON SOURCE LINES 64-68 .. code-block:: default start = np.array([97132.710, 457177.928]) end = np.array([103736.517, 457215.557]) .. GENERATED FROM PYTHON SOURCE LINES 69-70 Using the function ``imod.evaluate.quiver_line()`` considering the starting and ending points .. GENERATED FROM PYTHON SOURCE LINES 70-73 .. code-block:: default u, v = imod.evaluate.quiver_line(right, front, lower, start, end) .. GENERATED FROM PYTHON SOURCE LINES 74-75 Adding top and bottom information to the previously obtained arrays .. GENERATED FROM PYTHON SOURCE LINES 75-81 .. code-block:: default u["top"] = u["z"] + 0.5 * u["dz"] u["bottom"] = u["z"] - 0.5 * u["dz"] v["top"] = v["z"] + 0.5 * v["dz"] v["bottom"] = v["z"] - 0.5 * v["dz"] .. GENERATED FROM PYTHON SOURCE LINES 82-84 Defining a cross section that shows the heads for the same location where the start and end points were defined, to use as a background image for the plot .. GENERATED FROM PYTHON SOURCE LINES 84-89 .. code-block:: default cross_section = imod.select.cross_section_line(heads, start, end) cross_section["top"] = cross_section["z"] - 0.5 * cross_section["dz"] cross_section["bottom"] = cross_section["z"] + 0.5 * cross_section["dz"] .. GENERATED FROM PYTHON SOURCE LINES 90-91 Ploting the cross section .. GENERATED FROM PYTHON SOURCE LINES 91-100 .. code-block:: default colors = "magma" levels = np.arange(-8, 0.0, 0.5) skip = (slice(None, None, 2), slice(None, None, 2)) fig, ax = plt.subplots() fig, ax = imod.visualize.cross_section( cross_section, colors=colors, levels=levels, ax=ax, fig=fig ) imod.visualize.quiver(u[skip], v[skip], ax, kwargs_quiver={"color": "k"}) .. image-sg:: /examples/visualize/images/sphx_glr_flowvel_streamlines_002.png :alt: flowvel streamlines :srcset: /examples/visualize/images/sphx_glr_flowvel_streamlines_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 101-106 Streamline function ------------------- This function shows the streamlines for a line cross section through a 3D flow field. We will use the previously created arrays that indicate the start and end location of the cross section and apply them in the function imod.evaluate.streamfunction_line() .. GENERATED FROM PYTHON SOURCE LINES 106-109 .. code-block:: default streamfunction = imod.evaluate.streamfunction_line(right, front, start, end) .. GENERATED FROM PYTHON SOURCE LINES 110-116 The previous array contains the streamfunction projected on the cross-section defined by provided linestring, with new dimension “s” along the cross-section. The cellsizes along “s” are given in the “ds” coordinate. The streamline function can be plotted using imod.visualize.streamfunction(), but first we need to define the 'top' and 'bottom' of the layers in the array Note: By default the steamlines are plotted in white. .. GENERATED FROM PYTHON SOURCE LINES 116-126 .. code-block:: default streamfunction["bottom"] = streamfunction["z"] + streamfunction["dz"] streamfunction["top"] = streamfunction["z"] fig, ax = plt.subplots() fig, ax = imod.visualize.cross_section( cross_section, colors=colors, levels=levels, ax=ax, fig=fig ) imod.visualize.streamfunction(streamfunction, ax=ax, n_streamlines=10) .. image-sg:: /examples/visualize/images/sphx_glr_flowvel_streamlines_003.png :alt: flowvel streamlines :srcset: /examples/visualize/images/sphx_glr_flowvel_streamlines_003.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 1.392 seconds) .. _sphx_glr_download_examples_visualize_flowvel_streamlines.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: flowvel_streamlines.py ` .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: flowvel_streamlines.ipynb ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_