The Observations Module

Functions for handling observations of SFR package output.

sfrmaker.observations.add_observations(sfrdata, data, flowline_routing=None, obstype=None, sfrlines_shapefile=None, rno_column_in_sfrlines='rno', x_location_column=None, y_location_column=None, line_id_column=None, rno_column=None, obstype_column=None, obsname_column='site_no')[source]

Add SFR observations to the observations DataFrame attribute of an sfrdata instance. Observations can by located on the SFR network by specifying reach number directly (rno_column), by x, y location (x_column_in_data and y_column in data), or by specifying the source hydrography lines that they are located on (line_id_column).

Parameters:
sfrdatasfrmaker.SFRData instance

SFRData instance with reach_data table attribute. To add observations from x, y coordinates, the reach_data table must have a geometry column with LineStrings representing each reach, or an sfrlines_shapefile is required. Reach numbers are assumed to be in an ‘rno’ column.

dataDataFrame, path to csv file, or list of DataFrames or file paths

Table with information on the observation sites to be located. Must have either reach numbers (rno_column), line_ids (line_id_column), or x and y locations (x_column_in_data and y_column_in_data).

obstypestr or list-like (optional)

Type(s) of observation to record, for MODFLOW-6 (default ‘downstream-flow’; see MODFLOW-6 IO documentation for more details). Alternatively, observation types can be specified by row in data, using the obstype_column_in_data argument.

x_location_columnstr (optional)

Column in data with site x-coordinates (in same CRS as SFR network).

y_location_columnstr (optional)

Column in data with site y-coordinates (in same CRS as SFR network).

sfrlines_shapefilestr (optional)

Shapefile version of SFRdata.reach_data. Only needed if SFRdata.reach_data doesn’t have LineString geometries for the reaches.

rno_column_in_sfrlinesstr (optional)

Column in sfrlines with reach numbers for matching lines with reaches in sfrdata, or reach numbers assigned to observation sites. (default ‘rno’)

line_id_columnstr

Column in data matching observation sites to line_ids in the source hydrography data.

rno_columnstr

Column in data matching observation sites to reach numbers in the SFR network.

flowline_routingdict

Optional dictionary of routing for source hydrography. Only needed if locating by line_id, and SFR network is a subset of the full source hydrography (i.e. some lines were dropped in the creation of the SFR packge, or if the sites are inflow points corresponding to lines outside of the model perimeter). In this case, observation points referenced to line_ids that are missing from the SFR network are placed at the first reach corresponding to the next downstream line_id that is represented in the SFR network.

obstype_columnstr (optional)

Column in data with MODFLOW-6 observation types. For adding observations of different types. If obstype and obstype_column_in_data are none, the default of ‘downstream-flow’ will be used.

obsname_columnstr

Column in data with unique identifier (e.g. site number or name) for observation sites.

Notes

Sites located by line_id (source hydrography) will be assigned to the last reach in the segment corresponding to the line_id. Locating by x, y or reach number is more accurate.

sfrmaker.observations.get_closest_reach(x, y, sfrlines, rno_column='rno')[source]

Get the SFR reach number closest to a point feature.

Parameters:
xscalar or list of scalars

x-coordinate(s) of point feature(s)

yscalar or list of scalars

y-coordinate(s) or point feature(s)

sfrlines: dataframe

DataFrame containing a geometry column with SFR line arcs, and a column rno_column with unique numbers for each reach.

rno_column: str

Column with unique number for each reach. default “rno”

thresholdnumeric

Distance threshold (in CRS units). Only return reaches within this distance.

Returns:
rnoint or list of ints

Reach numbers for reaches closest to each location defined by x, y.

distance: int or list of ints

Distance from each site in x, y to the corresponding closest SFR reach.

sfrmaker.observations.locate_sites(site_data, reach_data, active_area_shapefile=None, x_column_in_data=None, y_column_in_data=None, reach_id_col='rno', ireach_col='ireach', iseg_col='iseg', site_number_col='site_no', keep_columns=None, perimeter_buffer=1000, distance_threshold=1000)[source]

Get SFR reach locations corresponding to x, y points (e.g. measurement site locations).

Parameters:
site_data: (Geo)DataFrame or shapefile

DataFrame or shapefile with point locations and attribute data for stream flow observation sites. Point locations can be specified in a DataFrame by either x_column_in_data and y_column_in_data, or a ‘geometry’ column of shapely points. If shapefiles or GeoDataFrames are provided for both site_data and reach_data, they can be in any CRS, but both must have valid CRS references (.prj file or .crs attribute); otherwise site_data and reach_data are assumed to be in the same CRS.

reach_data: (Geo)DataFrame or shapefile

SFRData.reach_data DataFrame, or shapefile equivalent with line-arcs representing all segments and/or reaches. If shapefiles or GeoDataFrames are provided for both site_data and reach_data, they can be in any CRS, but both must have valid CRS references (.prj file or .crs attribute); otherwise site_data and reach_data are assumed to be in the same CRS.

active_area_shapefile: ESRI shapefile or shapely polygon (optional)

Shapefile or polygon, in same CRS as sfr_lines_shapefile, defining areal extent (perimeter) of SFR network.

x_column_in_datastr (optional)

Column in data with site x-coordinates (in same CRS as SFR network).

y_column_in_datastr (optional)

Column in data with site y-coordinates (in same CRS as SFR network).

reach_id_col: str

Column with 1-based unique number for each stream line-arc. default “rno”

ireach_col: str

Column with 1-based reach number within segment (i.e. in the modflow-2005 SFR data model), by default, ‘ireach’

iseg_col: str

Column with 1-based reach segment number (i.e. in the modflow-2005 SFR data model), by default, ‘iseg’

site_number_colstr

Name of column in sites_shapefile with number identifying each site to be located. default “site_no”

keep_columns: list of strings

List of columns in sites_shapefile to retain when writing output_csv_file and output_shape_file.

perimeter_bufferscalar

Exclude flows within this distance of perimeter defined by active_area. For example, a value of 1000 would mean that sites must be at least 1 km inside of the active area perimeter to be included.

distance_thresholdscalar

Only consider sites within this distance of a stream line-arc.

Returns:
locsDataFrame
sfrmaker.observations.write_gage_package(location_data, gage_package_filename=None, gage_namfile_entries_file=None, model=None, sitename_col='site_no', gage_package_unit=25, start_gage_unit=200)[source]
Parameters:
location_datapandas.DataFrame

Table of observation locations. Must have columns: ‘iseg’: segment number ‘ireach’: reach number sitename_col: specified by sitename_col argument (default ‘site_no’)

gage_package_filenamestr or pathlike

Name for the gage package file, which will be written to the model workspace (model.model_ws) of the supplied flopy model instance. This must also be manually added to the MODFLOW .nam file, for example:

GAGE          25 br_trans.gage

Or, a new namefile can be written from the supplied flopy model instance.

gage_namfile_entries_filestr or pathlike

Namefile entries for the gage package output files will be written to this file. The contents of this file must then be manually copied and pasted into the MODFLOW .nam file. By default, None, in which case this file is named after gage_package_filename.

modelflopy model instance

Used for writing the gage package input file via flopy.

sitename_colstr

Unique name or number for each gage site. By default, ‘site_no’

gage_package_unitint

MODFLOW unit number to assign to gage package. By default, 25

start_gage_unitint

Modflow unit numbers for each gage package output file will be assigned sequentially starting at this number. By default, 200

Returns:
sfrmaker.observations.write_mf6_sfr_obsfile(observation_locations, filename, sfr_output_filename, digits=5, print_input=True)[source]

Write MODFLOW-6 observation input for the SFR package.

Parameters:
observation_locationsDataFrame or str (filepath to csv file)

line_id : unique identifier for observation locations

filenamestr

File path for MODFLOW-6 SFR observation input file

sfr_output_filenamestr

File path that SFR observation output file

digitsint

the number of significant digits with which simulated values are written to the output file.

print_inputbool

keyword to indicate that the list of observation information will be written to the listing file immediately after it is read.

Returns:
writes filename