pygetm.domain module
- class pygetm.domain.Domain(nx: int, ny: int, *, lon: ArrayLike | None = None, lat: ArrayLike | None = None, x: ArrayLike | None = None, y: ArrayLike | None = None, coordinate_type: CoordinateType | None = None, mask: ArrayLike | None = 1, H: ArrayLike | None = None, z0: ArrayLike | None = 0.0, f: ArrayLike | None = None, periodic_x: bool = False, periodic_y: bool = False, comm: Comm | None = None, logger: Logger | None = None)[source]
Bases:
objectCreate new domain.
- Parameters:
nx – number of tracer points in x-direction
ny – number of tracer points in y-direction
lon – longitude (°East)
lat – latitude (°North)
x – x coordinate (m)
y – y coordinate (m)
coordinate_type – preferred coordinate type for plots and output
mask – initial mask (0: land, 1: water)
H – initial bathymetric depth. This is the distance between the bottom and some arbitrary depth reference (m, positive if bottom lies below the depth reference). Typically the depth reference is mean sea level. Points with NaN or masked values will be marked as land in the mask.
z0 – minimum hydrodynamic bottom roughness (m)
f – Coriolis parameter (rad s-1). If not provided, it will be calculated from latitude. In that case argument lat must be provided.
periodic_x – use periodic boundary in x-direction (left == right)
periodic_y – use periodic boundary in y-direction (top == bottom)
comm – MPI communicator that comprises all processes that should get access to the domain
logger – logger for diagnostic messages
- property H: ndarray | None
bathymetric depth (m)
This is the distance between the bottom and some arbitrary depth reference (m, positive if bottom lies below the depth reference). Typically the depth reference is mean sea level.
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].It can be changed by assigning directly to the attribute or to slices of it. If assigning directly, values defined on the T grid
(ny x nx), the X grid(ny+1 x nx+1), or the supergrid(ny*2+1 x nx*2+1)are accepted, as are scalars. Assigned values are then interpolated to the supergrid. When assigning directly, masked or NaN values will be marked as land in the mask.This attribute is None on non-root MPI nodes or if H was not provided yet.
- property area: ndarray | None
grid cell area (m²)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes.
- cfl_check(z: float = 0.0, return_location: bool = False) float | tuple[float, int, int, float][source]
Determine maximum time step (s) for depth-integrated equations
- Parameters:
z – maximum surface elevation (m)
return_location – whether to also return the location and depth that determined the maximum step
Note: this returns global indices for the T grid, not the supergrid
- contains(x: float, y: float, coordinate_type: CoordinateType | None = None) bool[source]
Determine whether the domain contains the specified point.
- Parameters:
x – x coordinate
y – y coordinate
coordinate_type – coordinate system of the provided coordinates
- Returns:
True if the point falls within the domain, False otherwise
- create_grids(nz: int | None, halox: int, haloy: int, fields: Mapping[str, Array] | None = None, tiling: Tiling | None = None, input_manager=None, velocity_grids: int = 0, t_postfix: str = '') Grid[source]
- create_tiling(**kwargs) Tiling[source]
Create tiling object representing the domain decomposition.
- Parameters:
**kwargs – Additional arguments to pass to
parallel.Tiling.autodetect()- Returns:
tiling
- property dx: ndarray | None
grid cell length in x-direction (m)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes.
- property dy: ndarray | None
grid cell length in y-direction (m)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes.
- extract_rectangle(istart: int, istop: int, jstart: int, jstop: int) Domain[source]
Return a copy of the domain corresponding to the specified rectangle. Indices must be provided for the T grid. Negative indices can be used to indicate positions relative to the end of the domain. For example, istart=1, istop=-1, jstart=1, jstop=-1 removes the outermost strip of cells, and thus shrinks the domain by two cells in both x- and y-direction.
- Parameters:
istart – lower x index (first that is included)
istop – upper x index (first that is EXcluded)
jstart – lower y index (first that is included)
jstop – upper y index (first that is EXcluded)
- Returns:
Domain corresponding to the specified rectangle
- property f: ndarray | None
Coriolis parameter (rad s-1)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes or if
fwas not provided at domain creation. In the latter case, it will be calculated fromlat.
- get_final_mask() ndarray | None[source]
Infer masks for U, V, X points from T point mask. This closes interfaces and corners that neighbor one or more dry points. Additionally, it sets the mask (values 2,3,4) within and alongside open boundaries.
- get_rx0(zmin: float = 0.0, Dmin: float = 0.0) tuple[ndarray, ndarray][source]
Calculates the slope factor
rx0as defined in https://doi.org/10.1016/j.ocemod.2009.03.009At interfaces of shallow points (the depth of at least one wet neighbor being less than
Dmin), the surface elevation is clipped to the minimum value that still keeps both points wet.- Parameters:
zmin – minimum surface elevation (m)
Dmin – minimum depth (m)
- Returns:
a tuple with slope factors at U and V points, positioned at
[1::2, 2:-2:2]and[2:-2:2, 1::2], respectively
- property lat: ndarray | None
latitude (°North)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes or if
latwas not provided at domain creation (for instance, for Cartesian domains with prescribed Coriolis parameter).
- limit_velocity_depth(critical_depth: float = inf)[source]
Decrease bathymetric depth of velocity (U, V) points to the minimum of the bathymetric depth of both neighboring T points, wherever one of these two points is shallower than the specified critical depth.
- Parameters:
critical_depth – neighbor depth at which the limiting starts. If either neighbor (T grid) is shallower than this value, the depth of velocity point (U or V grid) is restricted.
- property lon: ndarray | None
longitude (°East)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes or if
lonwas not provided at domain creation (for instance, for Cartesian domains).
- property mask: ndarray | None
water)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].It can be changed by assigning directly to the attribute or to slices of it. If assigning directly, values defined on the T grid
(ny x nx), the X grid(ny+1 x nx+1), or the supergrid(ny*2+1 x nx*2+1)are accepted, as are scalars. Assigned values are then interpolated to the supergrid.This attribute is None on non-root MPI nodes.
- Type:
land-sea mask (0
- Type:
land, 1
- mask_indices(istart: int, istop: int, jstart: int, jstop: int, mask_value: int = 0)[source]
Mask all points that fall within the specified rectangle. Indices must be provided for the T grid, and thus range between 0 and nx for i, and between 0 and ny for j.
- Parameters:
istart – lower x index (first that is included)
istop – upper x index (first that is EXcluded)
jstart – lower y index (first that is included)
jstop – upper y index (first that is EXcluded)
- mask_rectangle(xmin: float | None = None, xmax: float | None = None, ymin: float | None = None, ymax: float | None = None, mask_value: int = 0, coordinate_type: CoordinateType | None = None)[source]
Mask all points that fall within the specified rectangle.
- Parameters:
xmin – lower native x coordinate of the rectangle to mask (default: left boundary of the domain)
xmax – upper native x coordinate of the rectangle to mask (default: right boundary of the domain)
ymin – lower native y coordinate of the rectangle to mask (default: bottom boundary of the domain)
ymax – upper native y coordinate of the rectangle to mask (default: top boundary of the domain)
Coordinates will be interpreted as longitude, latitude if the domain is configured as spherical; otherwise they will be interpreted as Cartesian x and y (m).
- mask_shallow(minimum_depth: float)[source]
Mask all points shallower less the specified value.
- Parameters:
minimum_depth – minimum bathmetric depth
H; points that are shallower will be masked
- mask_subbasins(nkeep: int = 1)[source]
Identify all separate basins (each a collection of unmasked cell centers connected via top/bottom/left/right interfaces), and mask all but the largest one(s).
- Parameters:
nkeep – number of basins to keep
- property max_rx0: float
//doi.org/10.1016/j.ocemod.2009.03.009
- Type:
Maximum slope factor rx0 as defined in https
- nearest_point(x: float, y: float, *, coordinate_type: CoordinateType | None = None, valid_cell_types: Iterable[CellType] = (CellType.ACTIVE,)) tuple[int, int][source]
- plot(field: ndarray | None = None, *, fig: matplotlib.figure.Figure | None = None, show_bathymetry: bool = True, show_mask: bool = False, show_mesh: bool = False, show_rivers: bool = True, show_open_boundaries: bool = True, show_subdomains: bool = False, editable: bool = False, coordinate_type: CoordinateType | None = None, tiling: Tiling | None = None, label: str | None = None, cmap: None | matplotlib.colors.Colormap | str = None) matplotlib.figure.Figure | None[source]
Plot the domain, optionally including bathymetric depth, mesh and river positions.
- Parameters:
field – 2D field to plot. it must be defined on the supergrid.
fig –
matplotlib.figure.Figureinstance to plot to. If not provided, a new figure is created.show_bathymetry – show bathymetry as color map
show_mask – show mask as color map (this disables
show_bathymetry)show_mesh – show model grid
show_rivers – show rivers with position and name
show_subdomains – show subdomain decompositon
editable – allow interactive selection of rectangular regions in the domain plot that are subsequently masked out
coordinate_type – coordinates to use for x and y axes (x/y, lon/lat, or i/j)
tiling – subdomain decomposition. This must be provided if show_subdomains is set
label – label for the colorbar of the plotted field
cmap – colormap to use
- Returns:
matplotlib.figure.Figureinstance for processes with rank 0, otherwiseNone
- property rotation: ndarray | None
grid rotation with respect to true North (rad)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes, or if the y-axis points to true North in every point (i.e, rotation is zero everywhere).
- smooth(rx0: float = 0.2) ndarray | None[source]
Smooth bathymetry by reducing slope factor to specified maximum.
The criterion to be satisfied at every interface between wet points is:
abs(H1 - H2) - rx0 * H1 - rx0 * H2 < 0
Handling abs as described at https://lpsolve.sourceforge.net/5.1/absolute.htm, we obtain the two criteria:
( 1-rx0) * H1 + (-1-rx0) * H2 < 0 (-1-rx0) * H1 + ( 1-rx0) * H2 < 0
Splitting into original bathymetry H and correction H’, and rearranging:
( 1-rx0)H1’ + (-1-rx0)H2’ < -( 1-rx0)H1 - (-1-rx0)H2 (-1-rx0)H1’ + ( 1-rx0)H2’ < -(-1-rx0)H1 - ( 1-rx0)H2
H1’ and H2’ are the corrections to the bathymetry, estimated by linear programming as described in https://doi.org/10.1016/j.ocemod.2009.03.009
- Parameters:
rx0 – maximum slope factor
- Returns:
bathymetry corrections (m). These are defined at T points, that is, at
[1::2, 1::2]
- to_xarray() Dataset[source]
Convert domain to
xarray.Dataset. The result can be loaded into a domain withfrom_xarray().- Returns:
xarray Dataset with domain data
- property x: ndarray | None
x coordinate (m)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes or if
xwas not provided at domain creation (for instance, for spherical domains).
- property y: ndarray | None
y coordinate (m)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].This attribute is None on non-root MPI nodes or if
ywas not provided at domain creation (for instance, for spherical domains).
- property z0: ndarray | None
minimum hydrodynamic bottom roughness (m)
It is defined on the supergrid and thus has shape
(ny*2+1, nx*2+1). Cell centers (T points) are at[1::2, 1::2], interfaces at[1::2, ::2](U points) and[::2, 1::2](V points), corners (X points) at[::2, ::2].It can be changed by assigning directly to the attribute or to slices of it. If assigning directly, values defined on the T grid
(ny x nx), the X grid(ny+1 x nx+1), or the supergrid(ny*2+1 x nx*2+1)are accepted, as are scalars. Assigned values are then interpolated to the supergrid.This attribute is None on non-root MPI nodes
- pygetm.domain.centers_to_supergrid_1d(source: ArrayLike, *, dtype: DTypeLike | None = None, edges: EdgeTreatment = EdgeTreatment.MISSING, missing_value=nan) ndarray[source]
- pygetm.domain.centers_to_supergrid_2d(source: ArrayLike, *, dtype: DTypeLike | None = None, edges_x: EdgeTreatment = EdgeTreatment.MISSING, edges_y: EdgeTreatment = EdgeTreatment.MISSING, missing_value=nan) ndarray[source]
- pygetm.domain.coriolis(lat: ArrayLike) ndarray[source]
Calculate Coriolis parameter f for the given latitude.
- Parameters:
lat – latitude (°North)
- Returns:
Coriolis parameter f (rad s-1)
- pygetm.domain.create_cartesian(x: ArrayLike, y: ArrayLike, *, interfaces: bool = False, central_lon: float | None = None, central_lat: float | None = None, **kwargs) Domain[source]
Create Cartesian domain from x and y coordinates.
- Parameters:
x – array with x coordinates (m). It can have shape
(nx,)or(ny, nx). Coordinates are interpreted to be positioned at cell interfaces if ifinterfaces=True, at cell centers otherwise.y – array with y coordinates (m). It can have shape
(ny,)or(ny, nx). Coordinates are interpreted to be positioned at cell interfaces if ifinterfaces=True, at cell centers otherwise.interfaces – coordinates are given at cell interfaces rather than cell centers.
central_lon – longitude of the center of the domain (°East). The center is
[x.min() + x.max()]/2, [y.min() + y.max()]/2.central_lat – latitude of the center of the domain (°North). The center is
[x.min() + x.max()]/2, [y.min() + y.max()]/2.**kwargs – additional arguments passed to
Domain
- pygetm.domain.create_spherical(lon: ArrayLike, lat: ArrayLike, *, interfaces: bool = False, **kwargs) Domain[source]
Create spherical domain from longitudes and latitudes.
- Parameters:
lon – array with longitude coordinates (°East). It can have shape
(nx,)or(ny, nx). Coordinates are interpreted to be positioned at cell interfaces if ifinterfaces=True, at cell centers otherwise.lat – array with latitude coordinates (°North). It can have shape
(ny,)or(ny, nx). Coordinates are interpreted to be positioned at cell interfaces if ifinterfaces=True, at cell centers otherwise.interfaces – coordinates are given at cell interfaces rather than cell centers.
**kwargs – additional arguments passed to
Domain
- pygetm.domain.create_spherical_at_resolution(minlon: float, maxlon: float, minlat: float, maxlat: float, resolution: float, **kwargs) Domain[source]
Create spherical domain encompassing the specified longitude range and latitude range and desired resolution in m.
- Parameters:
minlon – minimum longitude (°East)
maxlon – maximum longitude (°East)
minlat – minimum latitude (°North)
maxlat – maximum latitude (°North)
resolution – maximum grid cell length and width (m)
**kwargs – additional arguments passed to
Domain
- pygetm.domain.expand_2d(source: ArrayLike, *, dtype: DTypeLike | None = None, edges_x: EdgeTreatment = EdgeTreatment.MISSING, edges_y: EdgeTreatment = EdgeTreatment.MISSING, missing_value=nan) ndarray[source]
- pygetm.domain.from_xarray(ds: Dataset, **kwargs) Domain[source]
Create domain from
xarray.Dataset. This dataset must represent the arguments toDomainby variables (for 2D arrays) or attributes (for scalars) with the same names. Such a dataset is produced byDomain.to_xarray().- Parameters:
ds – dataset with domain data
**kwargs – additional arguments passed to
Domain. These override corresponding values in the dataset.
- Returns:
domain created from the dataset