dorado.particle_track

Lower-level methods to manage parameters and particle movement.

Particle class for managing the definition of particle attributes and parameters of the domain as well as iterative movement of the particles through the domain.

Project Homepage: https://github.com/passaH2O/dorado

class dorado.particle_track.Particles(params)[source]

Class for the particle(s) that is(are) going to be routed.

__init__(params)[source]

Check input parameters and assign default values where/if needed.

Methods require a class of parameters (dorado.particle_track.modelParams) to be passed to the Particles class. e.g. particle = Particles(modelParams)

Initialization tries to assign each value from the parameter class, otherwise an error is raised or default values are assigned when possible/sensible

__weakref__

list of weak references to the object (if defined)

clear_walk_data()[source]

Manually reset self.walk_data back to None.

generate_particles(Np_tracer, seed_xloc, seed_yloc, seed_time=0, method='random', previous_walk_data=None)[source]

Generate a set of particles in defined locations.

After a dorado.particle_track.Particles class has been initialized, this function can be called to create some particles within a given region. If particles are to be seeded in different groupings in different regions, this function can be called multiple times in succession prior to moving them via dorado.particle_track.run_iteration.

Walk data is stored within dorado.particle_track.Particles as Particles.walk_data so if this method is called in succession, even without the previous_walk_data flag, the new particle seed locations will be appended to the walk_data attribute of dorado.particle_track.Particles.

Inputs :

Np_tracerint

Number of particles to generate.

seed_xloclist

List of x-coordinates over which to initially distribute the particles.

seed_yloclist

List of y-coordinates over which to initially distribute the particles.

seed_timefloat, int, optional

Seed time to apply to all newly seeded particles. This can be useful if you are generating new particles and adding them to a pre-existing set of particles and you want to run them in a group to a “target-time”. The default value is 0, as new particles have not traveled for any time, but a savy user might have need to seed this value with something else.

methodstr, optional

Specify the type of particle generation you wish to use. Current options are ‘random’ and ‘exact’ for seeding particles either randomly across a set of points, or exactly at a set of defined locations. If unspecified, the default is ‘random’.

previous_walk_datadict, optional

Dictionary of all prior x locations, y locations, and travel times. This input parameter should only be used if ‘new’ or ‘independent’ walk data exists (e.g. data loaded from a .json file). Otherwise the walk data stored in self.walk_data is likely to contain the previous information. Order of indices is previous_walk_data[field][particle][iter], where e.g. [‘travel_times’][5][10] is the travel time of the 5th particle at the 10th iteration.

run_iteration(target_time=None, max_iter=10000.0)[source]

Run an iteration of the particle routing.

Runs an iteration of the particle routing. Returns at each step the particle’s locations and travel times.

Inputs :

target_timefloat, optional

The travel time (seconds) each particle should aim to have at end of this iteration. If left undefined, then just one iteration is run and the particles will be out of sync in time. Note that this loop will terminate before the target_time if the particle exceeds the maximum number of permitted iterations (max_iter)

max_iterint, optional

The maximum number of iterations a particle is allowed to take in order to try to reach a specified target_time. This is an optional parameter with a default value of 10,000.

Outputs :

all_walk_datadict

Dictionary of all x and y locations and travel times. Order of indices is previous_walk_data[field][particle][iter], where [‘travel_times’][5][10] is the travel time of the 5th particle at the 10th iteration

dorado.particle_track.coord2ind(coordinates, raster_origin, raster_size, cellsize)[source]

Convert geographical coordinates into raster index coordinates.

Assumes geographic coordinates are projected onto a Cartesian grid. Accepts coordinates in meters or decimal degrees.

Inputs :

coordinateslist

List [] of (x,y) pairs or tuples of coordinates to be converted from starting units (e.g. meters UTM) into raster index coordinates used in particle routing functions.

raster_origintuple

Tuple of the (x,y) raster origin in physical space, i.e. the coordinates of lower left corner. For rasters loaded from a GeoTIFF, lower left corner can be obtained using e.g. gdalinfo

raster_sizetuple

Tuple (L,W) of the raster dimensions, i.e. the output of numpy.shape(raster).

cellsizefloat or int

Length along one square cell face.

Outputs :

indslist

List [] of tuples (x,y) of raster index coordinates

dorado.particle_track.exposure_time(walk_data, region_of_interest, verbose=True)[source]

Measure exposure time distribution of particles in a specified region.

Function to measure the exposure time distribution (ETD) of particles to the specified region. For steady flows, the ETD is exactly equivalent to the residence time distribution. For unsteady flows, if particles make multiple excursions into the region, all of those times are counted.

Inputs :

walk_datadict

Output of a previous function call to run_iteration.

region_of_interestint array

Binary array the same size as input arrays in modelParams class with 1’s everywhere inside the region in which we want to measure exposure time, and 0’s everywhere else.

verbosebool, optional

Toggles whether or not messages print to the console. If False, nothing is output, if True, messages are output. Default value is True. Errors are always raised.

Outputs :

exposure_timeslist

List of exposure times to region of interest, listed in order of particle ID.

dorado.particle_track.gen_input_check(Np_tracer, seed_xloc, seed_yloc, seed_time, method)[source]

Check the inputs provided to generate_particles().

This function does input type checking and either succeeds or returns an error if the types provided cannot be converted to those that we expect in dorado.particle_track.generate_particles().

Inputs :

Np_tracerint

Number of particles to generate.

seed_xloclist

List of x-coordinates over which to initially distribute the particles

seed_yloclist

List of y-coordinates over which to initially distribute the particles

seed_timeint, float

Value to set as initial travel time for newly seeded particles.

methodstr, optional

Type of particle generation to use, either ‘random’ or ‘exact’.

Outputs :

Np_tracerint

Number of particles to generate.

seed_xloclist

List of x-coordinates over which to initially distribute the particles

seed_yloclist

List of y-coordinates over which to initially distribute the particles

dorado.particle_track.ind2coord(walk_data, raster_origin, raster_size, cellsize)[source]

Convert raster index coordinates into geographical coordinates.

Appends the walk_data dictionary from the output of run_iteration with additional fields ‘xcoord’ and ‘ycoord’ in projected geographic coordinate space. Locations align with cell centroids. Units of output coordinates match those of raster_origin and cellsize, can be meters or decimal degrees.

Inputs :

walk_datadict

Dictionary of all prior x locations, y locations, and travel times (the output of run_iteration)

raster_origintuple

Tuple of the (x,y) raster origin in physical space, i.e. the coordinates of lower left corner. For rasters loaded from a GeoTIFF, lower left corner can be obtained using e.g. gdalinfo

raster_sizetuple

Tuple (L,W) of the raster dimensions, i.e. the output of numpy.shape(raster).

cellsizefloat or int

Length along one square cell face.

Outputs :

walk_datadict

Same as the input walk_data dictionary, with the added ‘xcoord’ and ‘ycoord’ fields representing the particles geographic position at each iteration.

class dorado.particle_track.modelParams[source]

Parameter class with attributes and grid particles will be routed on.

The parameters class, dorado.particle_track.modelParams, must be populated with user-defined attributes of the grid the particles will be modeled on.

Required Parameters:

dxfloat

Length along one square cell face

depthnumpy.ndarray

Array of water depth values, if absent then the stage and topography arrays will be used to compute it

stagenumpy.ndarray

Array of water stage values, if absent then the depth and topography arrays will be used to compute it

qxnumpy.ndarray

Array of the x-component of flow discharge

qynumpy.ndarray

Array of the y-component of flow discharge

unumpy.ndarray

Array of the x-component of flow velocity

vnumpy.ndarray

Array of the y-component of flow velocity

Optional Parameters:

topographynumpy.ndarray

Array of cell elevation values

modelstr

Name of the hydrodynamic model input being used (e.g. ‘DeltaRCM’)

thetafloat

First of two weighting parameters for the weighted random walk. Default value is 1.0, higher values give higher weighting probabilities to cells with greater water depths

gammafloat

Second of two weighting parameters for the weighted random walk. Default value is 0.05. Gamma must be in the range [0,1]. Gamma == 1 means that the random walk weights are independent of the discharge values, and instead are based on the water surface gradient (the stage). Gamma == 0 means that the random walk weights are not dependent on the surface gradient, and instead are based on the inertial forces (the flow discharge).

diff_coefffloat

Diffusion/dispersion coefficient for use in travel time computation. If set to 0.0, flow is purely advection with no diffusion. Higher values lead to more spread in exit age distribution. Max diffusion time in any given step is 0.5*diff_coeff percent. Default value is 0.2 (i.e. max of 10%)

dry_depthfloat

Minimum depth for a cell to be considered wet, default value is 0.1m

cell_typenumpy.ndarray

Array of the different types of cells in the domain where 2 = land, 1 = channel, 0 = ocean, and -1 = edge. If not explicitly defined then the values are estimated based on the depth array and the defined dry_depth

steepest_descentbool

Toggle for routing based on a steepest descent rather than the weighted random walk. If True, then the highest weighted cells are used to route the particles. Default value is False.

verbosebool, optional

Toggles whether or not warnings and print output are output to the console. If False, nothing is output, if True, messages are output. Default value is True. Errors are always raised.

This list of expected parameter values can also be obtained by querying the class attributes with dir(modelParams), modelParams.__dict__, or vars(modelParams).

__init__()[source]

Create the expected variables for the modelParams class.

Variables are initialized as NoneType they need to be assigned by the user. Due to the wide variety of data formats produced by differeny hydrodynamic models, this is process is not automated and must be handled on a case-by-case basis.

__weakref__

list of weak references to the object (if defined)

dorado.particle_track.nourishment_area(walk_data, raster_size, sigma=0.7, clip=99.5)[source]

Determine the nourishment area of a particle injection

Function will measure the regions of the domain ‘fed’ by a seed location, as indicated by the history of particle travel locations in walk_data. Returns a heatmap raster, in which values indicate number of occasions each cell was occupied by a particle (after spatial filtering to reduce noise in the random walk, and normalizing by number of particles).

Inputs :

walk_datadict

Dictionary of all prior x locations, y locations, and travel times (the output of run_iteration)

raster_sizetuple

Tuple (L,W) of the domain dimensions, i.e. the output of numpy.shape(raster).

sigmafloat, optional

Degree of spatial smoothing of the area, implemented using a Gaussian kernal of the same sigma, via scipy.ndimage.gaussian_filter Default is light smoothing with sigma = 0.7 (to turn off smoothing, set sigma = 0)

clipfloat, optional

Percentile at which to truncate the distribution. Particles which get stuck can lead to errors at the high-extreme, so this parameter is used to normalize by a “near-max”. Default is the 99.5th percentile. To use true max, specify clip = 100

Outputs :

visit_freqnumpy.ndarray

Array of normalized particle visit frequency, with cells in the range [0, 1] representing the number of instances particles visited that cell. If sigma > 0, the array values include spatial filtering

dorado.particle_track.nourishment_time(walk_data, raster_size, sigma=0.7, clip=99.5)[source]

Determine the nourishment time of a particle injection

Function will measure the average length of time particles spend in each area of the domain for a given seed location, as indicated by the history of particle travel times in walk_data. Returns a heatmap raster, in which values indicate the average length of time particles tend to remain in each cell (after spatial filtering to reduce noise in the random walk).

Inputs :

walk_datadict

Dictionary of all prior x locations, y locations, and travel times (the output of run_iteration)

raster_sizetuple

Tuple (L,W) of the domain dimensions, i.e. the output of numpy.shape(raster).

sigmafloat, optional

Degree of spatial smoothing of the area, implemented using a Gaussian kernal of the same sigma, via scipy.ndimage.gaussian_filter Default is light smoothing with sigma = 0.7 (to turn off smoothing, set sigma = 0)

clipfloat, optional

Percentile at which to truncate the distribution. Particles which get stuck can lead to errors at the high-extreme, so this parameter is used to normalize by a “near-max”. Default is the 99.5th percentile. To use true max, specify clip = 100

Outputs :

mean_timenumpy.ndarray

Array of mean occupation time, with cell values representing the mean time particles spent in that cell. If sigma > 0, the array values include spatial filtering

dorado.particle_track.unstruct2grid(coordinates, quantity, cellsize, k_nearest_neighbors=3, boundary=None, crop=True)[source]

Convert unstructured model outputs into gridded arrays.

Interpolates model variables (e.g. depth, velocity) from an unstructured grid onto a Cartesian grid using inverse-distance-weighted interpolation. Assumes projected (i.e. “flat”) geographic coordinates. Accepts coordinates in meters or decimal degrees. Extent of output rasters are based on extent of coordinates. (Function modeled after ANUGA plot_utils code)

Inputs :

coordinateslist

List [] of (x,y) pairs or tuples of coordinates at which the interpolation quantities are located (e.g. centroids or vertices of an unstructured hydrodynamic model).

quantitylist

List [] of data to be interpolated with indices matching each (x,y) location given in coordinates. If quantity is depth, list would be formatted as [d1, d2, … , dn].

cellsizefloat or int

Length along one square cell face.

k_nearest_neighborsint, optional

Number of nearest neighbors to use in the interpolation. If k>1, inverse-distance-weighted interpolation is used.

boundarylist, optional

List [] of (x,y) coordinates used to delineate the boundary of interpolation. Points outside the polygon will be assigned as nan. Format needs to match requirements of matplotlib.path.Path()

cropbool, optional

If a boundary is specified, setting crop to True will eliminate any all-NaN borders from the interpolated rasters.

Outputs :

interp_funcfunction

Nearest-neighbor interpolation function for gridding additional quantities. Quicker to use this output function on additional variables (e.g. later time-steps of an unsteady model) than to make additional function calls to unstruct2grid. Function assumes data have the same coordinates. It is used as follows: “new_gridded_quantity = interp_func(new_quantity)”.

gridded_quantitynumpy.ndarray

Array of quantity after interpolation.