Coverage for C:\src\imod-python\imod\mf6\wel.py: 93%
227 statements
« prev ^ index » next coverage.py v7.5.1, created at 2024-05-08 14:15 +0200
« prev ^ index » next coverage.py v7.5.1, created at 2024-05-08 14:15 +0200
1from __future__ import annotations
3import warnings
4from typing import Any, Optional, Tuple, Union
6import cftime
7import numpy as np
8import numpy.typing as npt
9import pandas as pd
10import xarray as xr
11import xugrid as xu
13import imod
14from imod.logging import init_log_decorator
15from imod.mf6.boundary_condition import (
16 BoundaryCondition,
17 DisStructuredBoundaryCondition,
18 DisVerticesBoundaryCondition,
19)
20from imod.mf6.interfaces.ipointdatapackage import IPointDataPackage
21from imod.mf6.mf6_wel_adapter import Mf6Wel
22from imod.mf6.package import Package
23from imod.mf6.utilities.dataset import remove_inactive
24from imod.mf6.utilities.regrid import RegridderType
25from imod.mf6.validation import validation_pkg_error_message
26from imod.mf6.write_context import WriteContext
27from imod.prepare import assign_wells
28from imod.prepare.layer import create_layered_top
29from imod.schemata import (
30 AnyNoDataSchema,
31 DTypeSchema,
32 EmptyIndexesSchema,
33 ValidationError,
34)
35from imod.select.points import points_indices, points_values
36from imod.typing import GridDataArray
37from imod.typing.grid import is_spatial_2D, ones_like
38from imod.util.structured import values_within_range
41def _assign_dims(arg: Any) -> Tuple | xr.DataArray:
42 is_da = isinstance(arg, xr.DataArray)
43 if is_da and "time" in arg.coords:
44 if arg.ndim != 2:
45 raise ValueError("time varying variable: must be 2d")
46 if arg.dims[0] != "time":
47 arg = arg.transpose()
48 da = xr.DataArray(
49 data=arg.values, coords={"time": arg["time"]}, dims=["time", "index"]
50 )
51 return da
52 elif is_da:
53 return "index", arg.values
54 else:
55 return "index", arg
58def mask_2D(package: Well, domain_2d: GridDataArray) -> Well:
59 point_active = points_values(domain_2d, x=package.x, y=package.y)
61 is_inside_exterior = point_active == 1
62 selection = package.dataset.loc[{"index": is_inside_exterior}]
64 cls = type(package)
65 new = cls.__new__(cls)
66 new.dataset = selection
67 return new
70class Well(BoundaryCondition, IPointDataPackage):
71 """
72 Agnostic WEL package, which accepts x, y and a top and bottom of the well screens.
74 This package can be written to any provided model grid.
75 Any number of WEL Packages can be specified for a single groundwater flow model.
76 https://water.usgs.gov/water-resources/software/MODFLOW-6/mf6io_6.0.4.pdf#page=63
78 Parameters
79 ----------
81 y: float or list of floats
82 is the y location of the well.
83 x: float or list of floats
84 is the x location of the well.
85 screen_top: float or list of floats
86 is the top of the well screen.
87 screen_bottom: float or list of floats
88 is the bottom of the well screen.
89 rate: float, list of floats or xr.DataArray
90 is the volumetric well rate. A positive value indicates well
91 (injection) and a negative value indicates discharge (extraction) (q).
92 If provided as DataArray, an ``"index"`` dimension is required and an
93 optional ``"time"`` dimension and coordinate specify transient input.
94 In the latter case, it is important that dimensions are in the order:
95 ``("time", "index")``
96 concentration: array of floats (xr.DataArray, optional)
97 if this flow package is used in simulations also involving transport, then this array is used
98 as the concentration for inflow over this boundary.
99 concentration_boundary_type: ({"AUX", "AUXMIXED"}, optional)
100 if this flow package is used in simulations also involving transport, then this keyword specifies
101 how outflow over this boundary is computed.
102 id: list of Any, optional
103 assign an identifier code to each well. if not provided, one will be generated
104 Must be convertible to string, and unique entries.
105 minimum_k: float, optional
106 on creating point wells, no point wells will be placed in cells with a lower horizontal conductivity than this
107 minimum_thickness: float, optional
108 on creating point wells, no point wells will be placed in cells with a lower thickness than this
109 print_input: ({True, False}, optional)
110 keyword to indicate that the list of well information will be written to
111 the listing file immediately after it is read.
112 Default is False.
113 print_flows: ({True, False}, optional)
114 Indicates that the list of well flow rates will be printed to the
115 listing file for every stress period time step in which "BUDGET PRINT"
116 is specified in Output Control. If there is no Output Control option
117 and PRINT FLOWS is specified, then flow rates are printed for the last
118 time step of each stress period.
119 Default is False.
120 save_flows: ({True, False}, optional)
121 Indicates that well flow terms will be written to the file specified
122 with "BUDGET FILEOUT" in Output Control.
123 Default is False.
124 observations: [Not yet supported.]
125 Default is None.
126 validate: {True, False}
127 Flag to indicate whether the package should be validated upon
128 initialization. This raises a ValidationError if package input is
129 provided in the wrong manner. Defaults to True.
130 repeat_stress: Optional[xr.DataArray] of datetimes
131 Used to repeat data for e.g. repeating stress periods such as
132 seasonality without duplicating the values. The DataArray should have
133 dimensions ``("repeat", "repeat_items")``. The ``repeat_items``
134 dimension should have size 2: the first value is the "key", the second
135 value is the "value". For the "key" datetime, the data of the "value"
136 datetime will be used. Can also be set with a dictionary using the
137 ``set_repeat_stress`` method.
139 Examples
140 ---------
142 >>> screen_top = [0.0, 0.0]
143 >>> screen_bottom = [-2.0, -2.0]
144 >>> y = [83.0, 77.0]
145 >>> x = [81.0, 82.0]
146 >>> rate = [1.0, 1.0]
148 >>> imod.mf6.Well(x, y, screen_top, screen_bottom, rate)
150 For a transient well:
152 >>> weltimes = pd.date_range("2000-01-01", "2000-01-03")
154 >>> rate_factor_time = xr.DataArray([0.5, 1.0], coords={"time": weltimes}, dims=("time",))
155 >>> rate_transient = rate_factor_time * xr.DataArray(rate, dims=("index",))
157 >>> imod.mf6.Well(x, y, screen_top, screen_bottom, rate_transient)
158 """
160 @property
161 def x(self) -> npt.NDArray[np.float64]:
162 return self.dataset["x"].values
164 @property
165 def y(self) -> npt.NDArray[np.float64]:
166 return self.dataset["y"].values
168 _pkg_id = "wel"
170 _auxiliary_data = {"concentration": "species"}
171 _init_schemata = {
172 "screen_top": [DTypeSchema(np.floating)],
173 "screen_bottom": [DTypeSchema(np.floating)],
174 "y": [DTypeSchema(np.floating)],
175 "x": [DTypeSchema(np.floating)],
176 "rate": [DTypeSchema(np.floating)],
177 "concentration": [DTypeSchema(np.floating)],
178 }
179 _write_schemata = {
180 "screen_top": [AnyNoDataSchema(), EmptyIndexesSchema()],
181 "screen_bottom": [AnyNoDataSchema(), EmptyIndexesSchema()],
182 "y": [AnyNoDataSchema(), EmptyIndexesSchema()],
183 "x": [AnyNoDataSchema(), EmptyIndexesSchema()],
184 "rate": [AnyNoDataSchema(), EmptyIndexesSchema()],
185 "concentration": [AnyNoDataSchema(), EmptyIndexesSchema()],
186 }
188 _regrid_method: dict[str, Tuple[RegridderType, str]] = {}
190 @init_log_decorator()
191 def __init__(
192 self,
193 x: list[float],
194 y: list[float],
195 screen_top: list[float],
196 screen_bottom: list[float],
197 rate: list[float] | xr.DataArray,
198 concentration: Optional[list[float] | xr.DataArray] = None,
199 concentration_boundary_type="aux",
200 id: Optional[list[Any]] = None,
201 minimum_k: float = 0.1,
202 minimum_thickness: float = 1.0,
203 print_input: bool = False,
204 print_flows: bool = False,
205 save_flows: bool = False,
206 observations=None,
207 validate: bool = True,
208 repeat_stress: Optional[xr.DataArray] = None,
209 ):
210 if id is None:
211 id = [str(i) for i in range(len(x))]
212 else:
213 set_id = set(id)
214 if len(id) != len(set_id):
215 raise ValueError("id's must be unique")
216 id = [str(i) for i in id]
217 dict_dataset = {
218 "screen_top": _assign_dims(screen_top),
219 "screen_bottom": _assign_dims(screen_bottom),
220 "y": _assign_dims(y),
221 "x": _assign_dims(x),
222 "rate": _assign_dims(rate),
223 "id": _assign_dims(id),
224 "minimum_k": minimum_k,
225 "minimum_thickness": minimum_thickness,
226 "print_input": print_input,
227 "print_flows": print_flows,
228 "save_flows": save_flows,
229 "observations": observations,
230 "repeat_stress": repeat_stress,
231 "concentration": concentration,
232 "concentration_boundary_type": concentration_boundary_type,
233 }
234 super().__init__(dict_dataset)
235 # Set index as coordinate
236 index_coord = np.arange(self.dataset.dims["index"])
237 self.dataset = self.dataset.assign_coords(index=index_coord)
238 self._validate_init_schemata(validate)
240 @classmethod
241 def is_grid_agnostic_package(cls) -> bool:
242 return True
244 def clip_box(
245 self,
246 time_min: Optional[cftime.datetime | np.datetime64 | str] = None,
247 time_max: Optional[cftime.datetime | np.datetime64 | str] = None,
248 layer_min: Optional[int] = None,
249 layer_max: Optional[int] = None,
250 x_min: Optional[float] = None,
251 x_max: Optional[float] = None,
252 y_min: Optional[float] = None,
253 y_max: Optional[float] = None,
254 top: Optional[GridDataArray] = None,
255 bottom: Optional[GridDataArray] = None,
256 ) -> Package:
257 """
258 Clip a package by a bounding box (time, layer, y, x).
260 The well package doesn't use the layer attribute to describe its depth and length.
261 Instead, it uses the screen_top and screen_bottom parameters which corresponds with
262 the z-coordinates of the top and bottom of the well. To go from a layer_min and
263 layer_max to z-values used for clipping the well a top and bottom array have to be
264 provided as well.
266 Slicing intervals may be half-bounded, by providing None:
268 * To select 500.0 <= x <= 1000.0:
269 ``clip_box(x_min=500.0, x_max=1000.0)``.
270 * To select x <= 1000.0: ``clip_box(x_min=None, x_max=1000.0)``
271 or ``clip_box(x_max=1000.0)``.
272 * To select x >= 500.0: ``clip_box(x_min = 500.0, x_max=None.0)``
273 or ``clip_box(x_min=1000.0)``.
275 Parameters
276 ----------
277 time_min: optional
278 time_max: optional
279 layer_min: optional, int
280 layer_max: optional, int
281 x_min: optional, float
282 x_max: optional, float
283 y_min: optional, float
284 y_max: optional, float
285 top: optional, GridDataArray
286 bottom: optional, GridDataArray
287 state_for_boundary: optional, GridDataArray
289 Returns
290 -------
291 sliced : Package
292 """
293 if (layer_max or layer_min) and (top is None or bottom is None):
294 raise ValueError(
295 "When clipping by layer both the top and bottom should be defined"
296 )
298 if top is not None:
299 # Bug in mypy when using unions in isInstance
300 if not isinstance(top, GridDataArray) or "layer" not in top.coords: # type: ignore
301 top = create_layered_top(bottom, top)
303 # The super method will select in the time dimension without issues.
304 new = super().clip_box(time_min=time_min, time_max=time_max)
306 ds = new.dataset
308 z_max = self._find_well_value_at_layer(ds, top, layer_max)
309 z_min = self._find_well_value_at_layer(ds, bottom, layer_min)
311 if z_max is not None:
312 ds["screen_top"] = ds["screen_top"].clip(None, z_max)
313 if z_min is not None:
314 ds["screen_bottom"] = ds["screen_bottom"].clip(z_min, None)
316 # Initiate array of True with right shape to deal with case no spatial
317 # selection needs to be done.
318 in_bounds = np.full(ds.dims["index"], True)
319 # Select all variables along "index" dimension
320 in_bounds &= values_within_range(ds["x"], x_min, x_max)
321 in_bounds &= values_within_range(ds["y"], y_min, y_max)
322 in_bounds &= values_within_range(ds["screen_top"], z_min, z_max)
323 in_bounds &= values_within_range(ds["screen_bottom"], z_min, z_max)
324 # remove wells where the screen bottom and top are the same
325 in_bounds &= abs(ds["screen_bottom"] - ds["screen_top"]) > 1e-5
326 # Replace dataset with reduced dataset based on booleans
327 new.dataset = ds.loc[{"index": in_bounds}]
329 return new
331 @staticmethod
332 def _find_well_value_at_layer(
333 well_dataset: xr.Dataset, grid: GridDataArray, layer: Optional[int]
334 ):
335 value = None if layer is None else grid.isel(layer=layer)
337 # if value is a grid select the values at the well locations and drop the dimensions
338 if (value is not None) and is_spatial_2D(value):
339 value = imod.select.points_values(
340 value,
341 x=well_dataset["x"].values,
342 y=well_dataset["y"].values,
343 out_of_bounds="ignore",
344 ).drop_vars(lambda x: x.coords)
346 return value
348 def write(
349 self,
350 pkgname: str,
351 globaltimes: Union[list[np.datetime64], np.ndarray],
352 write_context: WriteContext,
353 ):
354 raise NotImplementedError(
355 "To write a wel package first convert it to a MF6 well using to_mf6_pkg."
356 )
358 def __create_wells_df(self) -> pd.DataFrame:
359 wells_df = self.dataset.to_dataframe()
360 wells_df = wells_df.rename(
361 columns={
362 "screen_top": "top",
363 "screen_bottom": "bottom",
364 }
365 )
367 return wells_df
369 def __create_assigned_wells(
370 self,
371 wells_df: pd.DataFrame,
372 active: GridDataArray,
373 top: GridDataArray,
374 bottom: GridDataArray,
375 k: GridDataArray,
376 minimum_k: float,
377 minimum_thickness: float,
378 ):
379 # Ensure top, bottom & k
380 # are broadcasted to 3d grid
381 like = ones_like(active)
382 bottom = like * bottom
383 top_2d = (like * top).sel(layer=1)
384 top_3d = bottom.shift(layer=1).fillna(top_2d)
386 k = like * k
388 index_names = wells_df.index.names
390 # Unset multi-index, because assign_wells cannot deal with
391 # multi-indices which is returned by self.dataset.to_dataframe() in
392 # case of a "time" and "species" coordinate.
393 wells_df = wells_df.reset_index()
395 wells_assigned = assign_wells(
396 wells_df, top_3d, bottom, k, minimum_thickness, minimum_k, True
397 )
398 # Set multi-index again
399 wells_assigned = wells_assigned.set_index(index_names).sort_index()
401 return wells_assigned
403 def __create_dataset_vars(
404 self, wells_assigned: pd.DataFrame, wells_df: pd.DataFrame, cellid: xr.DataArray
405 ) -> xr.Dataset:
406 """
407 Create dataset with all variables (rate, concentration), with a similar shape as the cellids.
408 """
409 data_vars = ["rate"]
410 if "concentration" in wells_assigned.columns:
411 data_vars.append("concentration")
413 ds_vars = wells_assigned[data_vars].to_xarray()
414 # "rate" variable in conversion from multi-indexed DataFrame to xarray
415 # DataArray results in duplicated values for "rate" along dimension
416 # "species". Select first species to reduce this again.
417 index_names = wells_df.index.names
418 if "species" in index_names:
419 ds_vars["rate"] = ds_vars["rate"].isel(species=0)
421 # Carefully rename the dimension and set coordinates
422 d_rename = {"index": "ncellid"}
423 ds_vars = ds_vars.rename_dims(**d_rename).rename_vars(**d_rename)
424 ds_vars = ds_vars.assign_coords(**{"ncellid": cellid.coords["ncellid"].values})
426 return ds_vars
428 def __create_cellid(self, wells_assigned: pd.DataFrame, active: xr.DataArray):
429 like = ones_like(active)
431 # Groupby index and select first, to unset any duplicate records
432 # introduced by the multi-indexed "time" dimension.
433 df_for_cellid = wells_assigned.groupby("index").first()
434 d_for_cellid = df_for_cellid[["x", "y", "layer"]].to_dict("list")
436 return self.__derive_cellid_from_points(like, **d_for_cellid)
438 @staticmethod
439 def __derive_cellid_from_points(
440 dst_grid: GridDataArray,
441 x: list,
442 y: list,
443 layer: list,
444 ) -> GridDataArray:
445 """
446 Create DataArray with Modflow6 cell identifiers based on x, y coordinates
447 in a dataframe. For structured grid this DataArray contains 3 columns:
448 ``layer, row, column``. For unstructured grids, this contains 2 columns:
449 ``layer, cell2d``.
450 See also: https://water.usgs.gov/water-resources/software/MODFLOW-6/mf6io_6.4.0.pdf#page=35
452 Note
453 ----
454 The "layer" coordinate should already be provided in the dataframe.
455 To determine the layer coordinate based on screen depts, look at
456 :func:`imod.prepare.wells.assign_wells`.
458 Parameters
459 ----------
460 dst_grid: {xr.DataArray, xu.UgridDataArray}
461 Destination grid to map the points to based on their x and y coordinates.
462 x: {list, np.array}
463 array-like with x-coordinates
464 y: {list, np.array}
465 array-like with y-coordinates
466 layer: {list, np.array}
467 array-like with layer-coordinates
469 Returns
470 -------
471 cellid : xr.DataArray
472 2D DataArray with a ``ncellid`` rows and 3 to 2 columns, depending
473 on whether on a structured or unstructured grid."""
475 # Find indices belonging to x, y coordinates
476 indices_cell2d = points_indices(dst_grid, out_of_bounds="ignore", x=x, y=y)
477 # Convert cell2d indices from 0-based to 1-based.
478 indices_cell2d = {dim: index + 1 for dim, index in indices_cell2d.items()}
479 # Prepare layer indices, for later concatenation
481 if isinstance(dst_grid, xu.UgridDataArray):
482 indices_layer = xr.DataArray(
483 layer, coords=indices_cell2d["mesh2d_nFaces"].coords
484 )
485 face_dim = dst_grid.ugrid.grid.face_dimension
486 indices_cell2d_dims = [face_dim]
487 cell2d_coords = ["cell2d"]
488 else:
489 indices_layer = xr.DataArray(layer, coords=indices_cell2d["x"].coords)
490 indices_cell2d_dims = ["y", "x"]
491 cell2d_coords = ["row", "column"]
493 # Prepare cellid array of the right shape.
494 cellid_ls = [indices_layer] + [
495 indices_cell2d[dim] for dim in indices_cell2d_dims
496 ]
497 cellid = xr.concat(cellid_ls, dim="nmax_cellid")
498 # Rename generic dimension name "index" to ncellid.
499 cellid = cellid.rename(index="ncellid")
500 # Put dimensions in right order after concatenation.
501 cellid = cellid.transpose("ncellid", "nmax_cellid")
502 # Assign extra coordinate names.
503 coords = {
504 "nmax_cellid": ["layer"] + cell2d_coords,
505 "x": ("ncellid", x),
506 "y": ("ncellid", y),
507 }
508 cellid = cellid.assign_coords(**coords)
510 return cellid
512 def render(self, directory, pkgname, globaltimes, binary):
513 raise NotImplementedError(
514 f"{self.__class__.__name__} is a grid-agnostic package and does not have a render method. To render the package, first convert to a Modflow6 package by calling pkg.to_mf6_pkg()"
515 )
517 def to_mf6_pkg(
518 self,
519 active: GridDataArray,
520 top: GridDataArray,
521 bottom: GridDataArray,
522 k: GridDataArray,
523 validate: bool = False,
524 is_partitioned: bool = False,
525 ) -> Mf6Wel:
526 """
527 Write package to Modflow 6 package.
529 Based on the model grid and top and bottoms, cellids are determined.
530 When well screens hit multiple layers, groundwater extractions are
531 distributed based on layer transmissivities. Wells located in inactive
532 cells are removed.
534 Note
535 ----
536 The well distribution based on transmissivities assumes confined
537 aquifers. If wells fall dry (and the rate distribution has to be
538 recomputed at runtime), it is better to use the Multi-Aquifer Well
539 package.
541 Parameters
542 ----------
543 is_partitioned: bool
544 validate: bool
545 Run validation before converting
546 active: {xarry.DataArray, xugrid.UgridDataArray}
547 Grid with active cells.
548 top: {xarry.DataArray, xugrid.UgridDataArray}
549 Grid with top of model layers.
550 bottom: {xarry.DataArray, xugrid.UgridDataArray}
551 Grid with bottom of model layers.
552 k: {xarry.DataArray, xugrid.UgridDataArray}
553 Grid with hydraulic conductivities.
554 Returns
555 -------
556 Mf6Wel
557 Object with wells as list based input.
558 """
559 if validate:
560 errors = self._validate(self._write_schemata)
561 if len(errors) > 0:
562 message = validation_pkg_error_message(errors)
563 raise ValidationError(message)
565 minimum_k = self.dataset["minimum_k"].item()
566 minimum_thickness = self.dataset["minimum_thickness"].item()
568 wells_df = self.__create_wells_df()
569 wells_assigned = self.__create_assigned_wells(
570 wells_df, active, top, bottom, k, minimum_k, minimum_thickness
571 )
573 nwells_df = len(wells_df["id"].unique())
574 nwells_assigned = (
575 0 if wells_assigned.empty else len(wells_assigned["id"].unique())
576 )
578 if nwells_df == 0:
579 raise ValueError("No wells were assigned in package. None were present.")
581 if not is_partitioned and nwells_df != nwells_assigned:
582 raise ValueError(
583 "One or more well(s) are completely invalid due to minimum conductivity and thickness constraints."
584 )
586 ds = xr.Dataset()
587 ds["cellid"] = self.__create_cellid(wells_assigned, active)
589 ds_vars = self.__create_dataset_vars(wells_assigned, wells_df, ds["cellid"])
590 ds = ds.assign(**ds_vars.data_vars)
592 ds = remove_inactive(ds, active)
593 ds["save_flows"] = self["save_flows"].values[()]
594 ds["print_flows"] = self["print_flows"].values[()]
595 ds["print_input"] = self["print_input"].values[()]
597 return Mf6Wel(**ds.data_vars)
599 def mask(self, domain: GridDataArray) -> Well:
600 """
601 Mask wells based on two-dimensional domain. For three-dimensional
602 masking: Wells falling in inactive cells are automatically removed in
603 the call to write to Modflow 6 package. You can verify this by calling
604 the ``to_mf6_pkg`` method.
605 """
607 # Drop layer coordinate if present, otherwise a layer coordinate is assigned
608 # which causes conflicts downstream when assigning wells and deriving
609 # cellids.
610 domain_2d = domain.isel(layer=0, drop=True, missing_dims="ignore").drop_vars(
611 "layer", errors="ignore"
612 )
613 return mask_2D(self, domain_2d)
615 def get_regrid_methods(self) -> Optional[dict[str, Tuple[RegridderType, str]]]:
616 return self._regrid_method
619class WellDisStructured(DisStructuredBoundaryCondition):
620 """
621 WEL package for structured discretization (DIS) models .
622 Any number of WEL Packages can be specified for a single groundwater flow model.
623 https://water.usgs.gov/water-resources/software/MODFLOW-6/mf6io_6.0.4.pdf#page=63
625 .. warning::
626 This class is deprecated and will be deleted in a future release.
627 Consider changing your code to use the ``imod.mf6.Well`` package.
629 Parameters
630 ----------
631 layer: list of int
632 Model layer in which the well is located.
633 row: list of int
634 Row in which the well is located.
635 column: list of int
636 Column in which the well is located.
637 rate: float or list of floats
638 is the volumetric well rate. A positive value indicates well
639 (injection) and a negative value indicates discharge (extraction) (q).
640 concentration: array of floats (xr.DataArray, optional)
641 if this flow package is used in simulations also involving transport, then this array is used
642 as the concentration for inflow over this boundary.
643 concentration_boundary_type: ({"AUX", "AUXMIXED"}, optional)
644 if this flow package is used in simulations also involving transport, then this keyword specifies
645 how outflow over this boundary is computed.
646 print_input: ({True, False}, optional)
647 keyword to indicate that the list of well information will be written to
648 the listing file immediately after it is read.
649 Default is False.
650 print_flows: ({True, False}, optional)
651 Indicates that the list of well flow rates will be printed to the
652 listing file for every stress period time step in which "BUDGET PRINT"
653 is specified in Output Control. If there is no Output Control option
654 and PRINT FLOWS is specified, then flow rates are printed for the last
655 time step of each stress period.
656 Default is False.
657 save_flows: ({True, False}, optional)
658 Indicates that well flow terms will be written to the file specified
659 with "BUDGET FILEOUT" in Output Control.
660 Default is False.
661 observations: [Not yet supported.]
662 Default is None.
663 validate: {True, False}
664 Flag to indicate whether the package should be validated upon
665 initialization. This raises a ValidationError if package input is
666 provided in the wrong manner. Defaults to True.
667 repeat_stress: Optional[xr.DataArray] of datetimes
668 Used to repeat data for e.g. repeating stress periods such as
669 seasonality without duplicating the values. The DataArray should have
670 dimensions ``("repeat", "repeat_items")``. The ``repeat_items``
671 dimension should have size 2: the first value is the "key", the second
672 value is the "value". For the "key" datetime, the data of the "value"
673 datetime will be used. Can also be set with a dictionary using the
674 ``set_repeat_stress`` method.
675 """
677 _pkg_id = "wel"
678 _period_data = ("layer", "row", "column", "rate")
679 _keyword_map = {}
680 _template = DisStructuredBoundaryCondition._initialize_template(_pkg_id)
681 _auxiliary_data = {"concentration": "species"}
683 _init_schemata = {
684 "layer": [DTypeSchema(np.integer)],
685 "row": [DTypeSchema(np.integer)],
686 "column": [DTypeSchema(np.integer)],
687 "rate": [DTypeSchema(np.floating)],
688 "concentration": [DTypeSchema(np.floating)],
689 }
691 _write_schemata = {}
693 @init_log_decorator()
694 def __init__(
695 self,
696 layer,
697 row,
698 column,
699 rate,
700 concentration=None,
701 concentration_boundary_type="aux",
702 print_input=False,
703 print_flows=False,
704 save_flows=False,
705 observations=None,
706 validate: bool = True,
707 repeat_stress=None,
708 ):
709 dict_dataset = {
710 "layer": _assign_dims(layer),
711 "row": _assign_dims(row),
712 "column": _assign_dims(column),
713 "rate": _assign_dims(rate),
714 "print_input": print_input,
715 "print_flows": print_flows,
716 "save_flows": save_flows,
717 "observations": observations,
718 "repeat_stress": repeat_stress,
719 "concentration": concentration,
720 "concentration_boundary_type": concentration_boundary_type,
721 }
722 super().__init__(dict_dataset)
723 self._validate_init_schemata(validate)
725 warnings.warn(
726 f"{self.__class__.__name__} is deprecated and will be removed in the v1.0 release."
727 "Please adapt your code to use the imod.mf6.Well package",
728 DeprecationWarning,
729 )
731 def clip_box(
732 self,
733 time_min: Optional[cftime.datetime | np.datetime64 | str] = None,
734 time_max: Optional[cftime.datetime | np.datetime64 | str] = None,
735 layer_min: Optional[int] = None,
736 layer_max: Optional[int] = None,
737 x_min: Optional[float] = None,
738 x_max: Optional[float] = None,
739 y_min: Optional[float] = None,
740 y_max: Optional[float] = None,
741 top: Optional[GridDataArray] = None,
742 bottom: Optional[GridDataArray] = None,
743 ) -> Package:
744 """
745 Clip a package by a bounding box (time, layer, y, x).
747 Slicing intervals may be half-bounded, by providing None:
749 * To select 500.0 <= x <= 1000.0:
750 ``clip_box(x_min=500.0, x_max=1000.0)``.
751 * To select x <= 1000.0: ``clip_box(x_min=None, x_max=1000.0)``
752 or ``clip_box(x_max=1000.0)``.
753 * To select x >= 500.0: ``clip_box(x_min = 500.0, x_max=None.0)``
754 or ``clip_box(x_min=1000.0)``.
756 Parameters
757 ----------
758 time_min: optional
759 time_max: optional
760 layer_min: optional, int
761 layer_max: optional, int
762 x_min: optional, float
763 x_max: optional, float
764 y_min: optional, float
765 y_max: optional, float
766 top: optional, GridDataArray
767 bottom: optional, GridDataArray
768 state_for_boundary: optional, GridDataArray
770 Returns
771 -------
772 sliced : Package
773 """
774 # TODO: include x and y values.
775 for arg in (
776 layer_min,
777 layer_max,
778 x_min,
779 x_max,
780 y_min,
781 y_max,
782 ):
783 if arg is not None:
784 raise NotImplementedError("Can only clip_box in time for Well packages")
786 # The super method will select in the time dimension without issues.
787 new = super().clip_box(time_min=time_min, time_max=time_max)
788 return new
791class WellDisVertices(DisVerticesBoundaryCondition):
792 """
793 WEL package for discretization by vertices (DISV) models. Any number of WEL
794 Packages can be specified for a single groundwater flow model.
795 https://water.usgs.gov/water-resources/software/MODFLOW-6/mf6io_6.0.4.pdf#page=63
797 .. warning::
798 This class is deprecated and will be deleted in a future release.
799 Consider changing your code to use the ``imod.mf6.Well`` package.
801 Parameters
802 ----------
803 layer: list of int
804 Modellayer in which the well is located.
805 cell2d: list of int
806 Cell in which the well is located.
807 rate: float or list of floats
808 is the volumetric well rate. A positive value indicates well (injection)
809 and a negative value indicates discharge (extraction) (q).
810 concentration: array of floats (xr.DataArray, optional)
811 if this flow package is used in simulations also involving transport,
812 then this array is used as the concentration for inflow over this
813 boundary.
814 concentration_boundary_type: ({"AUX", "AUXMIXED"}, optional)
815 if this flow package is used in simulations also involving transport,
816 then this keyword specifies how outflow over this boundary is computed.
817 print_input: ({True, False}, optional)
818 keyword to indicate that the list of well information will be written to
819 the listing file immediately after it is read. Default is False.
820 print_flows: ({True, False}, optional)
821 Indicates that the list of well flow rates will be printed to the
822 listing file for every stress period time step in which "BUDGET PRINT"
823 is specified in Output Control. If there is no Output Control option and
824 PRINT FLOWS is specified, then flow rates are printed for the last time
825 step of each stress period. Default is False.
826 save_flows: ({True, False}, optional)
827 Indicates that well flow terms will be written to the file specified
828 with "BUDGET FILEOUT" in Output Control. Default is False.
829 observations: [Not yet supported.]
830 Default is None.
831 validate: {True, False}
832 Flag to indicate whether the package should be validated upon
833 initialization. This raises a ValidationError if package input is
834 provided in the wrong manner. Defaults to True.
835 """
837 _pkg_id = "wel"
838 _period_data = ("layer", "cell2d", "rate")
839 _keyword_map = {}
840 _template = DisVerticesBoundaryCondition._initialize_template(_pkg_id)
841 _auxiliary_data = {"concentration": "species"}
843 _init_schemata = {
844 "layer": [DTypeSchema(np.integer)],
845 "cell2d": [DTypeSchema(np.integer)],
846 "rate": [DTypeSchema(np.floating)],
847 "concentration": [DTypeSchema(np.floating)],
848 }
850 _write_schemata = {}
852 @init_log_decorator()
853 def __init__(
854 self,
855 layer,
856 cell2d,
857 rate,
858 concentration=None,
859 concentration_boundary_type="aux",
860 print_input=False,
861 print_flows=False,
862 save_flows=False,
863 observations=None,
864 validate: bool = True,
865 ):
866 dict_dataset = {
867 "layer": _assign_dims(layer),
868 "cell2d": _assign_dims(cell2d),
869 "rate": _assign_dims(rate),
870 "print_input": print_input,
871 "print_flows": print_flows,
872 "save_flows": save_flows,
873 "observations": observations,
874 "concentration": concentration,
875 "concentration_boundary_type": concentration_boundary_type,
876 }
877 super().__init__(dict_dataset)
878 self._validate_init_schemata(validate)
880 warnings.warn(
881 f"{self.__class__.__name__} is deprecated and will be removed in the v1.0 release."
882 "Please adapt your code to use the imod.mf6.Well package",
883 DeprecationWarning,
884 )
886 def clip_box(
887 self,
888 time_min: Optional[cftime.datetime | np.datetime64 | str] = None,
889 time_max: Optional[cftime.datetime | np.datetime64 | str] = None,
890 layer_min: Optional[int] = None,
891 layer_max: Optional[int] = None,
892 x_min: Optional[float] = None,
893 x_max: Optional[float] = None,
894 y_min: Optional[float] = None,
895 y_max: Optional[float] = None,
896 top: Optional[GridDataArray] = None,
897 bottom: Optional[GridDataArray] = None,
898 ) -> Package:
899 """
900 Clip a package by a bounding box (time, layer, y, x).
902 Slicing intervals may be half-bounded, by providing None:
904 * To select 500.0 <= x <= 1000.0:
905 ``clip_box(x_min=500.0, x_max=1000.0)``.
906 * To select x <= 1000.0: ``clip_box(x_min=None, x_max=1000.0)``
907 or ``clip_box(x_max=1000.0)``.
908 * To select x >= 500.0: ``clip_box(x_min = 500.0, x_max=None.0)``
909 or ``clip_box(x_min=1000.0)``.
911 Parameters
912 ----------
913 time_min: optional
914 time_max: optional
915 layer_min: optional, int
916 layer_max: optional, int
917 x_min: optional, float
918 x_max: optional, float
919 y_min: optional, float
920 y_max: optional, float
921 top: optional, GridDataArray
922 bottom: optional, GridDataArray
923 state_for_boundary: optional, GridDataArray
925 Returns
926 -------
927 clipped: Package
928 """
929 # TODO: include x and y values.
930 for arg in (
931 layer_min,
932 layer_max,
933 x_min,
934 x_max,
935 y_min,
936 y_max,
937 ):
938 if arg is not None:
939 raise NotImplementedError("Can only clip_box in time for Well packages")
941 # The super method will select in the time dimension without issues.
942 new = super().clip_box(time_min=time_min, time_max=time_max)
943 return new