mfsetup.grid module

Code for creating and working with regular (structured) grids. Focus is on the 2D representation of the grid in the cartesian plane. For methods involving layering (in the vertical dimension), see the discretization module.

class mfsetup.grid.MFsetupGrid(delc, delr, top=None, botm=None, idomain=None, laycbd=None, lenuni=None, binary_grid_file=None, epsg=None, proj_str=None, prj=None, wkt=None, crs=None, xoff=0.0, yoff=0.0, xul=None, yul=None, angrot=0.0)[source]

Bases: StructuredGrid

Class representing a structured grid. Extends flopy.discretization.StructuredGrid to facilitate gis operations in a projected (real-word) coordinate reference system (CRS).

Parameters:
delcndarray

1D numpy array of grid spacing along a column (len nrow), in CRS units.

delrndarray

1D numpy array of grid spacing along a row (len ncol), in CRS units.

topndarray

2D numpy array of model top elevations

botmndarray

3D numpy array of model bottom elevations

idomainndarray

3D numpy array of model idomain values

laycbdndarray

(Modflow 2005 and earlier style models only): LAYCBD—is a flag, with one value for each model layer, that indicates whether or not a layer has a Quasi-3D confining bed below it. 0 indicates no confining bed, and not zero indicates a confining bed. LAYCBD for the bottom layer must be 0.

lenuniint, optional

MODFLOW length units variable. See the Online Guide to MODFLOW

epsgint, optional

EPSG code for the model CRS

proj_strstr, optional

PROJ string for model CRS. In general, a spatial reference ID (such as an EPSG code) or Well-Known Text (WKT) string is prefered over a PROJ string (see References)

prjstr, optional

Filepath for ESRI projection file (containing wkt) describing model CRS

wktstr, optional

Well-known text string describing model CRS.

crsobj, optional

A Python int, dict, str, or pyproj.crs.CRS instance passed to pyproj.crs.CRS.from_user_input() Can be any of:

  • PROJ string

  • Dictionary of PROJ parameters

  • PROJ keyword arguments for parameters

  • JSON string with PROJ parameters

  • CRS WKT string

  • An authority string [i.e. ‘epsg:4326’]

  • An EPSG integer code [i.e. 4326]

  • A tuple of (“auth_name”: “auth_code”) [i.e (‘epsg’, ‘4326’)]

  • An object with a to_wkt method.

  • A pyproj.crs.CRS class

xoff, yofffloat, float, optional

Model grid offset (location of lower left corner), by default 0.0, 0.0

xul, yulfloat, float, optional

Model grid offset (location of upper left corner), by default 0.0, 0.0

angrotfloat, optional

Rotation of the model grid, in degrees counter-clockwise about the lower left corner. Non-zero rotation values require input of xoff, yoff (xul, yul not supported). By default 0.0

References

https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems

property bbox

Shapely polygon bounding box of the model grid.

property botm
property bounds

Grid bounding box in order used by shapely.

property crs

pyproj.crs.CRS instance describing the coordinate reference system for the model grid.

property dataframe

Pandas DataFrame of grid cell polygons with i, j locations.

get_dataframe(layers=True)[source]

Get a pandas DataFrame of grid cell polygons with i, j locations.

Parameters:
layersbool

If True, return a row for each k, i, j location and a ‘k’ column; if False, only return i, j locations with no ‘k’ column. By default, True

Returns:
layersDataFrame

Pandas Dataframe with k, i, j and geometry column with a shapely polygon representation of each model cell.

get_intercell_connections(binary_grid_file=None)[source]

_summary_

Parameters:
binary_grid_filestr or pathlike

MODFLOW 6 binary grid file

Returns:
dfDataFrame

Intercell connections, with the following columns:

n

from zero-based node number

kn

from zero-based layer

in

from zero-based row

jn

from zero-based column

m

to zero-based node number

kn

to zero-based layer

in

to zero-based row

jn

to zero-based column

Raises:
ValueError

_description_

get_vertices(i, j)[source]

Get vertices for a single cell or sequence if i, j locations.

property intercell_connections

Pandas DataFrame of flow connections between grid cells.

property length_multiplier
property length_units
property polygons

Vertices for grid cell polygons.

property proj_str
property rotation
property size
property top
property transform

Rasterio Affine object (same as transform attribute of rasters).

property vertices

Vertices for grid cell polygons.

property wkt
write_bbox_shapefile(filename='grid_bbox.shp')[source]
write_shapefile(filename='grid.shp')[source]

Write a shapefile of the grid with just the row and column attributes.

property xul
property yul
mfsetup.grid.get_cellface_midpoint(grid, k, i, j, direction)[source]

Return the midpoint of vertical cell face within a structured grid. For example, the midpoint for the right cell face is halfway between the upper and lower right corners of the cell, halfway between the top and bottom edges.

mfsetup.grid.get_crs(crs=None, epsg=None, prj=None, wkt=None, proj_str=None)[source]

Get a pyproj CRS instance from various CRS representations.

mfsetup.grid.get_crs_length_units(crs)[source]
mfsetup.grid.get_grid_bounding_box(modelgrid)[source]

Get bounding box of potentially rotated modelgrid as a shapely Polygon object.

Parameters:
modelgridflopy.discretization.StructuredGrid instance
mfsetup.grid.get_ij(grid, x, y, local=False)[source]

Return the row and column of a point or sequence of points in real-world coordinates.

Parameters:
gridflopy.discretization.StructuredGrid instance
xscalar or sequence of x coordinates
yscalar or sequence of y coordinates
local: bool (optional)

If True, x and y are in local coordinates (defaults to False)

Returns:
irow or sequence of rows (zero-based)
jcolumn or sequence of columns (zero-based)
mfsetup.grid.get_intercell_connections(binary_grid_file)[source]

Get all of the connections between cells in a MODFLOW 6 structured grid.

Parameters:
binary_grid_filestr or pathlike

MODFLOW 6 binary grid file

Returns:
dfDataFrame

Intercell connections, with the following columns:

n

from zero-based node number

kn

from zero-based layer

in

from zero-based row

jn

from zero-based column

m

to zero-based node number

kn

to zero-based layer

in

to zero-based row

jn

to zero-based column

qidx

index position of flow in cell budget file

Raises:
ValueError

_description_

mfsetup.grid.get_kij_from_node3d(node3d, nrow, ncol)[source]

For a consecutive cell number in row-major order (row, column, layer), get the zero-based row, column position.

mfsetup.grid.get_nearest_point_on_grid(x, y, transform=None, xul=None, yul=None, dx=None, dy=None, rotation=0.0, offset='center', op=None)[source]
Parameters:
xfloat

x-coordinate of point

yfloat

y-coordinate of point

transformAffine instance, optional

Affine object instance describing grid

xulfloat

x-coordinate of upper left corner of the grid

yulfloat

y-coordinate of upper left corner of the grid

dxfloat

grid spacing in the x-direction (along rows)

dyfloat

grid spacing in the y-direction (along columns)

rotationfloat

grid rotation about the upper left corner, in degrees clockwise from the x-axis

offsetstr, {‘center’, ‘edge’}

Whether the point on the grid represents a cell center or corner (edge). This argument is only used if xul, yul, dx, dy and rotation are supplied. If an Affine transform instance is supplied, it is assumed to already incorporate the offset.

opfunction, optional

Function to convert fractional pixels to whole numbers (np.round, np.floor, np.ceiling). Defaults to np.round if offset == ‘center’; otherwise defaults to np.floor.

Returns:
x_nearest, y_nearestfloat

Coordinates of nearest grid cell center.

mfsetup.grid.get_point_on_national_hydrogeologic_grid(x, y, offset='edge', **kwargs)[source]

Given an x, y location representing the upper left corner of a model grid, return the upper left corner of the cell in the National Hydrogeologic Grid that contains it.

mfsetup.grid.get_transform(modelgrid)[source]

Get a rasterio Affine object from a Flopy modelgrid (same as transform attribute of rasters).

mfsetup.grid.rasterize(feature, grid, id_column=None, include_ids=None, exclude_ids=None, names_column=None, crs=None, **kwargs)[source]

Rasterize a feature onto the model grid, using the rasterio.features.rasterize method. Features are intersected if they contain the cell center.

Parameters:
featurestr (shapefile path), list of shapely objects,

or dataframe with geometry column

id_columnstr

Column with unique integer identifying each feature; values from this column will be assigned to the output raster.

include_idssequence

Subset of IDs in id_column to include

exclude_idssequence

Subset of IDs in id_column to exclude

names_columnstr, optional

By default, the IDs in id_column, or sequential integers are returned. This option allows another column of strings to be specified (i.e. feature names); in which case an array of the strings will be returned.

gridgrid.StructuredGrid instance
crsobj

A Python int, dict, str, or pyproj.crs.CRS instance passed to pyproj.crs.CRS.from_user_input() Can be any of:

  • PROJ string

  • Dictionary of PROJ parameters

  • PROJ keyword arguments for parameters

  • JSON string with PROJ parameters

  • CRS WKT string

  • An authority string [i.e. ‘epsg:4326’]

  • An EPSG integer code [i.e. 4326]

  • A tuple of (“auth_name”: “auth_code”) [i.e (‘epsg’, ‘4326’)]

  • An object with a to_wkt method.

  • A pyproj.crs.CRS class

**kwargskeyword arguments to rasterio.features.rasterize()

https://rasterio.readthedocs.io/en/stable/api/rasterio.features.html

Returns:
2D numpy array with intersected values
mfsetup.grid.setup_structured_grid(xoff=None, yoff=None, xul=None, yul=None, nrow=None, ncol=None, nlay=None, dxy=None, delr=None, delc=None, top=None, botm=None, rotation=0.0, parent_model=None, snap_to_parent=True, snap_to_NHG=False, features=None, features_shapefile=None, id_column=None, include_ids=None, buffer=1000, crs=None, epsg=None, prj=None, wkt=None, model_length_units=None, grid_file='grid.json', bbox_shapefile=None, **kwargs)[source]

_summary_

Parameters:
xoff_type_, optional

_description_, by default None

yoff_type_, optional

_description_, by default None

xul_type_, optional

_description_, by default None

yul_type_, optional

_description_, by default None

nrow_type_, optional

_description_, by default None

ncol_type_, optional

_description_, by default None

nlay_type_, optional

_description_, by default None

dxy_type_, optional

Specified uniform row/column spacing, in model grid (coordinate reference system) units, by default None

delrscalar or sequence, optional

Column spacing along a row, in model grid (coordinate reference system) units, by default None

delcscalar or sequence, optional

Row spacing along a column, in model grid (coordinate reference system) units, by default None

top_type_, optional

_description_, by default None

botm_type_, optional

_description_, by default None

rotation_type_, optional

_description_, by default 0.

parent_model_type_, optional

_description_, by default None

snap_to_parentbool, optional

_description_, by default True

snap_to_NHGbool, optional

_description_, by default False

features_type_, optional

_description_, by default None

features_shapefile_type_, optional

_description_, by default None

id_column_type_, optional

_description_, by default None

include_ids_type_, optional

_description_, by default None

bufferint, optional

_description_, by default 1000

crs_type_, optional

_description_, by default None

epsg_type_, optional

_description_, by default None

prj_type_, optional

_description_, by default None

wkt_type_, optional

_description_, by default None

model_length_units_type_, optional

_description_, by default None

grid_filestr, optional

_description_, by default ‘grid.json’

bbox_shapefile_type_, optional

_description_, by default None

Returns:
_type_

_description_

Raises:
ValueError

_description_

ValueError

_description_

ValueError

_description_

ValueError

_description_

ValueError

_description_

ValueError

_description_

mfsetup.grid.snap_to_cell_corner(x, y, modelgrid, corner='nearest')[source]

Move an x, y location to the nearest cell corner on a rectilinear modelgrid.

Parameters:
xfloat

x coordinate in coordinate reference system of modelgrid.

y_type_

y coordinate in coordinate reference system of modelgrid.

modelgridFlopy StructuredGrid instance
cornerstr, optional

‘upper left’, ‘lower right’ or ‘nearest’, by default ‘nearest’

Returns:
x_corner, y_corner

x, y location of cell corner in coordinate reference system of modelgrid.

Raises:
ValueError

If x, y are outside of the model domain, or if an invalid cell corner is specified.

mfsetup.grid.write_bbox_shapefile(modelgrid, outshp)[source]