.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "user-guide\07-time-discretization.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_user-guide_07-time-discretization.py: Model time discretization ========================= iMOD Python provides nice functionality to discretize your models into stress periods, depending on the timesteps you assigned your boundary conditions. This functionality is activated with the ``create_time_discretization()`` method. .. GENERATED FROM PYTHON SOURCE LINES 12-27 Basics ------ To demonstrate the ``create_time_discretization()`` method, we first have to create a Model object. In this case we'll use a Modflow 6 simulation, but the ``imod.wq.SeawatModel`` and imod.flow.ImodflowModel also support this. .. note:: The simulation created in this example misses some mandatory packes, so is not able to write a functional simulation. It is purely intended to describe the workings of the ``time_discretization`` method. For a fully functional Modflow 6 model, see the examples in :ref:`mf6-introduction`. Wel'll start off with the usual imports: .. GENERATED FROM PYTHON SOURCE LINES 27-34 .. code-block:: Python import numpy as np import pandas as pd import xarray as xr import imod .. GENERATED FROM PYTHON SOURCE LINES 35-37 We can discretize a simulation as follows, this creates a TimeDiscretization object under the key ``"time_discretization"``. .. GENERATED FROM PYTHON SOURCE LINES 37-45 .. code-block:: Python simulation = imod.mf6.Modflow6Simulation("example") simulation.create_time_discretization( additional_times=["2000-01-01", "2000-01-02", "2000-01-04"] ) simulation["time_discretization"] .. raw:: html
TimeDiscretization
<xarray.Dataset> Size: 48B
    Dimensions:              (time: 2)
    Coordinates:
      * time                 (time) datetime64[ns] 16B 2000-01-01 2000-01-02
    Data variables:
        timestep_duration    (time) float64 16B 1.0 2.0
        n_timesteps          int64 8B 1
        timestep_multiplier  float64 8B 1.0


.. GENERATED FROM PYTHON SOURCE LINES 46-47 To view the data inside TimeDiscretization object: .. GENERATED FROM PYTHON SOURCE LINES 47-50 .. code-block:: Python simulation["time_discretization"].dataset .. raw:: html
<xarray.Dataset> Size: 48B
    Dimensions:              (time: 2)
    Coordinates:
      * time                 (time) datetime64[ns] 16B 2000-01-01 2000-01-02
    Data variables:
        timestep_duration    (time) float64 16B 1.0 2.0
        n_timesteps          int64 8B 1
        timestep_multiplier  float64 8B 1.0


.. GENERATED FROM PYTHON SOURCE LINES 51-56 Notice that even though we specified three points in time, only two timesteps are included in the ``time`` coordinate, this is because Modflow requires a start time and a duration of each stress period. iMOD Python therefore uses three points in time to compute two stress periods with a duration. .. GENERATED FROM PYTHON SOURCE LINES 56-59 .. code-block:: Python simulation["time_discretization"].dataset["timestep_duration"] .. raw:: html
<xarray.DataArray 'timestep_duration' (time: 2)> Size: 16B
    array([1., 2.])
    Coordinates:
      * time     (time) datetime64[ns] 16B 2000-01-01 2000-01-02


.. GENERATED FROM PYTHON SOURCE LINES 60-69 These two stress periods use their respective start time in their ``time`` coordinate. Boundary Conditions ------------------- The ``create_time_discretization`` method becomes especially useful if we add boundary conditions to our groundwater model. We'll first still have to initialize a groundwater flow model though: .. GENERATED FROM PYTHON SOURCE LINES 69-72 .. code-block:: Python gwf_model = imod.mf6.GroundwaterFlowModel() .. GENERATED FROM PYTHON SOURCE LINES 73-75 We'll create a dummy grid, which is used to provide an appropriate grid to the boundary condition packages. .. GENERATED FROM PYTHON SOURCE LINES 75-84 .. code-block:: Python template_data = np.ones((2, 1, 2, 2)) layer = [1] y = [1.5, 0.5] x = [0.5, 1.5] dims = ("time", "layer", "y", "x") .. GENERATED FROM PYTHON SOURCE LINES 85-86 Next, we can assign a Constant Head boundary with two timesteps: .. GENERATED FROM PYTHON SOURCE LINES 86-98 .. code-block:: Python chd_times = pd.to_datetime(["2000-01-01", "2000-01-04"]) chd_data = xr.DataArray( data=template_data, dims=dims, coords={"time": chd_times, "layer": layer, "y": y, "x": x}, ) gwf_model["chd"] = imod.mf6.ConstantHead(head=chd_data) gwf_model["chd"].dataset .. raw:: html
<xarray.Dataset> Size: 139B
    Dimensions:        (time: 2, layer: 1, y: 2, x: 2)
    Coordinates:
      * time           (time) datetime64[ns] 16B 2000-01-01 2000-01-04
      * layer          (layer) int64 8B 1
      * y              (y) float64 16B 1.5 0.5
      * x              (x) float64 16B 0.5 1.5
    Data variables:
        head           (time, layer, y, x) float64 64B 1.0 1.0 1.0 ... 1.0 1.0 1.0
        print_input    bool 1B False
        print_flows    bool 1B False
        save_flows     bool 1B False
        observations   object 8B None
        repeat_stress  object 8B None


.. GENERATED FROM PYTHON SOURCE LINES 99-101 We'll also assign a Recharge boundary with two timesteps, which differ from the ConstantHead boundary: .. GENERATED FROM PYTHON SOURCE LINES 101-113 .. code-block:: Python rch_times = pd.to_datetime(["2000-01-01", "2000-01-02"]) rch_data = xr.DataArray( data=template_data, dims=dims, coords={"time": rch_times, "layer": layer, "y": y, "x": x}, ) gwf_model["rch"] = imod.mf6.Recharge(rate=rch_data) gwf_model["rch"].dataset .. raw:: html
<xarray.Dataset> Size: 140B
    Dimensions:        (time: 2, layer: 1, y: 2, x: 2)
    Coordinates:
      * time           (time) datetime64[ns] 16B 2000-01-01 2000-01-02
      * layer          (layer) int64 8B 1
      * y              (y) float64 16B 1.5 0.5
      * x              (x) float64 16B 0.5 1.5
    Data variables:
        rate           (time, layer, y, x) float64 64B 1.0 1.0 1.0 ... 1.0 1.0 1.0
        print_input    bool 1B False
        print_flows    bool 1B False
        save_flows     bool 1B False
        observations   object 8B None
        repeat_stress  object 8B None
        fixed_cell     bool 1B False


.. GENERATED FROM PYTHON SOURCE LINES 114-117 We can now let iMOD Python figure out how the simulation's time should be discretized. It is important that we provide an endtime, otherwise the duration of the last stress period cannot be determined: .. GENERATED FROM PYTHON SOURCE LINES 117-127 .. code-block:: Python endtime = pd.to_datetime(["2000-01-06"]) simulation_bc = imod.mf6.Modflow6Simulation("example_bc") simulation_bc["gwf_1"] = gwf_model simulation_bc.create_time_discretization(additional_times=endtime) simulation_bc["time_discretization"].dataset .. raw:: html
<xarray.Dataset> Size: 64B
    Dimensions:              (time: 3)
    Coordinates:
      * time                 (time) datetime64[ns] 24B 2000-01-01 ... 2000-01-04
    Data variables:
        timestep_duration    (time) float64 24B 1.0 2.0 2.0
        n_timesteps          int64 8B 1
        timestep_multiplier  float64 8B 1.0


.. GENERATED FROM PYTHON SOURCE LINES 128-140 Notice that iMOD Python figured out that the two boundary conditions, both with two timesteps, should lead to three stress periods! Specifying extra settings ------------------------- The ``TimeDiscretization`` package also supports other settings, like the amount of timesteps which Modflow should use within a stress period, as well as a timestep multiplier, to gradually increase the timesteps modflow uses within a stress period. This can be useful when boundary conditions change very abruptly between stress periods. These settings are set by accessing their respective variables in the ``dataset``. .. GENERATED FROM PYTHON SOURCE LINES 140-150 .. code-block:: Python times = simulation_bc["time_discretization"].dataset.coords["time"] simulation_bc["time_discretization"].dataset["timestep_multiplier"] = 1.5 simulation_bc["time_discretization"].dataset["n_timesteps"] = xr.DataArray( data=[2, 4, 4], dims=("time",), coords={"time": times} ) simulation_bc["time_discretization"].dataset .. raw:: html
<xarray.Dataset> Size: 80B
    Dimensions:              (time: 3)
    Coordinates:
      * time                 (time) datetime64[ns] 24B 2000-01-01 ... 2000-01-04
    Data variables:
        timestep_duration    (time) float64 24B 1.0 2.0 2.0
        n_timesteps          (time) int64 24B 2 4 4
        timestep_multiplier  float64 8B 1.5


.. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.141 seconds) .. _sphx_glr_download_user-guide_07-time-discretization.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: 07-time-discretization.ipynb <07-time-discretization.ipynb>` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: 07-time-discretization.py <07-time-discretization.py>` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: 07-time-discretization.zip <07-time-discretization.zip>` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_