API Documentation¶
Environment Specification¶
Ocean Environment Specification Generic enironment specification for ocean sound propagation.
- class pygenray.environment.OceanEnvironment2D(sound_speed=None, bathymetry=None, lat=35, flat_earth_transform=True, verbose=False)¶
Ocean Environment Specification (2D).
An ocean acoustic environment specification built around xarray scientific computing environment. Flat earth transformation is computed for single latitude.
- Parameters:
sound_speed (xr.DataArray) – 1D or 2D array with coordinate dimensions (depth,) or (depth,range.Order of dimensions is arbitrary. Units of sound speed should be in [m/s] and units of coordinates should be in [m]. Default is range-independent Munk profile for range of 100km.
bathymetry (xr.DataArray) – 1D array of bottom depth with coordinate dimension (range,). Units of bathymetry should be in [m]. Range grid does not have to be aligned with sound speed range grid. Default is flat bottom at 5000m depth.
lat (float) – latitude in degrees. Used for flat earth transformation. Deffault is 35 degrees.
flat_earth_transform (bool) – whether to transform sound speed and bathymetry to flat earth.
verbose (bool) – whether to print initialization progress messages. Default is False.
- sound_speed¶
2D array of sound speed with dimensions (range, depth) with units [m]
- Type:
xr.DataArray
- bathymetry¶
1D array of bathymetry with dimension (range,) in units [m]
- Type:
xr.DataArray
- dcdz¶
2D array of dc/dz with dimensions (range, depth) in units [m/s/m]
- Type:
xr.DataArray
- bottom_angle¶
1D array of bottom slope angles in degrees at each range. The angle is computed as the arctangent of the gradient of bathymetry with respect to range. The gradient is computed usimg numpy’s gradient function (2nd order)
- Type:
np.ndarray
- bottom_angle_call¶
function that returns the bottom slope angle in degrees at a given range. The angle is computed as the arctangent of the gradient of bathymetry with respect to range. The function takes a single argument, range in meters, and returns the angle in degrees.
- Type:
callable
- flat_earth_transform(lat)¶
Earth flattening transformation, change depths and sound speeds so that spherical shell can be done as an x-z slice (using WGS-84).
Single latitude is used. If latitude changes sufficiently over track, consider using flat_earth_rd() ``range dependenant’’ version instead.
- Parameters:
lat (float) – latitude in degrees. Positive is north, negative is south.
- flat_earth_transform_rd()¶
Earth flattening transformation, change depths and sound speeds so that spherical shell can be done as an x-z slice (using WGS-84).
earth flattening is computed for each range / lat value independently.
- plot(**kwargs)¶
plot 2D slice of environment
kwargs are passed to matplotlib.pyplot.pcolormesh through xarray
Launch Rays¶
- pygenray.launch_rays.shoot_ray(source_depth, source_range, launch_angle, receiver_range, num_range_save, environment, rtol=1e-09, terminate_backwards=True, debug=True, flatearth=True)¶
Integrate rays for given environment and launch angles. Different launch angle initial conditions are mapped to available CPUS.
- Parameters:
source_depth (
float) – array of source depths (meters)source_range (
float) – array of source ranges (meters)launch_angle (
float) – array of source angles (degrees), should be 1D with shape (k,)receiver_range (
float) – receiver range (meters)num_range_save (
int) – The number of range values to save the ray state at. This value is unrelated to the numerical integration. The ray state value that is at the end bounds of the range integration is saved exactly. All other values are interpolated to a range grid with num_range_save points between the source and receiver range.environment (
OceanEnvironment2D) – OceanEnvironment object specifying sound speed and bathymetry.rtol (float) – relative tolerance for the ODE solver, default is 1e-6
terminate_backwards (
bool) – whether to terminate ray if it bounces backwards)debug (
bool) – whether to print debug information, default is Falseflatearth (
bool) – whether to transform environment to flat earth coordinates. Default is True.
- Returns:
ray – pr.Ray object
- Return type:
pr.Ray
- pygenray.launch_rays.shoot_rays(source_depth, source_range, launch_angles, receiver_range, num_range_save, environment, rtol=1e-09, terminate_backwards=True, n_processes=None, debug=True, flatearth=True)¶
Integrate rays for given environment and launch angles. Different launch angle initial conditions are mapped to available CPUS.
- Parameters:
source_depth (
float) – array of source depths (meters)source_range (
float) – array of source ranges (meters)launch_angles (
array) – array of source angles (degrees)receiver_range (
float) – receiver range (meters)num_range_save (
int) – The number of range values to save the ray state at. This value is unrelated to the numerical integration. The ray state value that is at the end bounds of the range integration is saved exactly. All other values are interpolated to a range grid with num_range_save points between the source and receiver range.environment (
OceanEnvironment2D) – OceanEnvironment object specifying sound speed and bathymetry.rtol (float) – relative tolerance for the ODE solver, default is 1e-9
terminate_backwards (
bool) – whether to terminate ray if it bounces backwardsn_processes (
int) – number of processes to use, Default of None (mp.cpu_count)debug (
bool) – whether to print debug information, default is Falseflatearth (
bool) – whether to transform environment to flat earth coordinates. Default is True.
- Returns:
ray (np.array) – 2D array of ray state at each x_eval point, shape (4, n_eval), where n_eval is the number of evaluation points
n_bott (int) – number of bottom bounces
n_surf (int) – number of surface bounces
Eigen Ray Solving¶
Tools and methods for calculating eigenrays for specifed receiver depths.
- pygenray.eigenrays.find_eigenrays(rays, receiver_depths, source_depth, source_range, receiver_range, num_range_save, environment, ztol=1, max_iter=20, num_workers=None, **kwargs)¶
Given an initial ray fan, find eigenrays with [regula falsi](https://en.wikipedia.org/wiki/Regula_falsi#The_regula_falsi_(false_position)_method) method of root finding.
- Parameters:
rays (pr.RayFan) – RayFan object containing sweep of rays to be used for finding eigenrays. Can be computed with pr.shoot_rays().
receiver_depths (array like) – one dimensional array, or list containing receiver depths
source_depth (float) – source depth in meters
source_range (float) – source range in meters
receiver_range (float) – receiver range in meters
num_range_save (int) – number of range values to save the ray state at
environment (pr.OceanEnvironment2D) – OceanEnvironment2D object containing environment parameters for ray tracing.
ztol (float, optional) – depth tolerance for eigenrays, by default 1 m
max_iter (int, optional) – maximum number of root finding iterations, default 20
num_workers (int, optional) – number of workers for parallel processing, by default None (uses all available cores, i.e. mp.cpu_count())
kwargs (keyword arguments) – additional keyword arguments passed to pr.shoot_ray
- Returns:
erays – dictionary of eigen rays. Key values are values in receiver_depths.
- Return type:
Integration Processes¶
Ray equations and the methods needed to solve them. These functions are the primary bottle neck for the numerical integration and are optimized for speed using Numba just in time compilation. The ray equations are derived from the Hamiltonian formulation for ray theory [Colosi2016].
where \(t\) is the travel time, \(z\) is the depth, and \(p\) is the ray parameter \((\frac{sin(\theta)}{c})\), and range, \(x\) is the independant variable.
References
Colosi, J. A. (2016). Sound Propagation through the Stochastic Ocean, Cambridge University Press, 443 pages.
- pygenray.integration_processes.bilinear_interp(x, y, x_grid, y_grid, values)¶
Perform bilinear interpolation on a 2D grid.
Fast, purely functional bilinear interpolation for scattered points on a regular 2D grid using Numba JIT compilation for performance.
- Parameters:
x (float) – The x-coordinate at which to interpolate.
y (float) – The y-coordinate at which to interpolate.
x_grid (array_like) – 1-D array of x-coordinates of the grid points, must be sorted in ascending order.
y_grid (array_like) – 1-D array of y-coordinates of the grid points, must be sorted in ascending order.
values (array_like) – 2-D array of shape (len(x_grid), len(y_grid)) containing the values at each grid point.
- Returns:
The interpolated value at point (x, y).
- Return type:
Notes
This function uses bilinear interpolation, which linearly interpolates first in one dimension, then in the other. The interpolation is performed using the four nearest grid points surrounding the query point.
If the query point lies outside the grid bounds, it is clamped to the nearest edge of the grid before interpolation.
The function is compiled with Numba’s JIT compiler for improved performance.
Examples
>>> import numpy as np >>> x_grid = np.array([0.0, 1.0, 2.0]) >>> y_grid = np.array([0.0, 1.0, 2.0]) >>> values = np.array([[1, 2, 3], [4, 5, 6], [7, 8, 9]]) >>> result = bilinear_interp(0.5, 0.5, x_grid, y_grid, values) >>> print(result) # Should be 3.0
- pygenray.integration_processes.bottom_bounce(x, y, cin, cpin, rin, zin, depths, depth_ranges)¶
Bottom event: only trigger when approaching bottom from above
- pygenray.integration_processes.derivsrd(x, y, cin, cpin, rin, zin, depths, depth_ranges)¶
Compute the differential equations for ray propagation. The ray equations are derived from the Hamiltonian formulation for ray theory [Colosi2016a], which consist of three coupled ODEs with range as the independant varibale, given by equations (4), (5), and (6).
\[y = \left [ t, z, p \right ]^T\]where \(t\) is the travel time, \(z\) is the depth, and \(p\) is the ray parameter \((\frac{sin(\theta)}{c})\), and range, \(x\) is the independant variable.
(4)¶\[\begin{split}\frac{dT}{dx} = \frac{1}{c\sqrt{1-c^2 \ p_z^2}} \\\end{split}\](5)¶\[\begin{split}\frac{dz}{dx} = \frac{c \ p_z}{ \sqrt{1-c^2 \ p_z^2}} \\\end{split}\](6)¶\[\begin{split}\frac{dp_z}{dx} = -\frac{1}{c^2}\frac{1}{\sqrt{1-c^2 \ p_z^2}}\frac{\partial c}{\partial z} \\\end{split}\]- Parameters:
x (
float) – horizontal range (meters)y (
array) – ray variables, [travel time, depth, ray parameter (sin(θ)/c)]cin (
array) – 2D array of sound speed valuescpin (
array) – 2D array of dc/dzrin (
array) – range coordinate for c arrayszin (
array) – depth coordinate for c arraysdepths (
array) – array of bathymetry valuesdepth_ranges (
array) – array of bathymetry value ranges. Does not have to match rin grid.
- Returns:
dydx – derivative of ray variables with respect to horizontal range, [dT/dx, dz/dx, dp/dx]
- Return type:
array
References
[Colosi2016a]Colosi, J. A. (2016). Sound Propagation through the Stochastic Ocean, Cambridge University Press, 443 pages.
- pygenray.integration_processes.linear_interp(x, xin, yin)¶
Perform linear interpolation on a 1D grid.
Fast, purely functional linear interpolation for scattered points on a regular 1D grid using Numba JIT compilation for performance.
- Parameters:
x_interp (float) – The x-coordinate at which to interpolate.
xin (array_like) – 1-D array of x-coordinates of the grid points, must be sorted in ascending order.
yin (array_like) – 1-D array of shape (len(x_grid),) containing the values at each grid point.
- Returns:
y_interp – The interpolated value at point x.
- Return type:
Notes
This function uses linear interpolation between the two nearest grid points surrounding the query point.
If the query point lies outside the grid bounds, it is clamped to the nearest edge of the grid before interpolation.
The function is compiled with Numba’s JIT compiler for improved performance.
Examples
>>> import numpy as np >>> x_grid = np.array([0.0, 1.0, 2.0]) >>> values = np.array([1.0, 4.0, 7.0]) >>> result = linear_interp(0.5, x_grid, values) >>> print(result) # Should be 2.5
- pygenray.integration_processes.ray_angle(x, y, cin, rin, zin)¶
calculate angle of ray for specific ray state
- Parameters:
x (
float) – horizontal range (meters)y (
array) – ray variables, [travel time, depth, ray parameter (sin(θ)/c)]cin (
array) – 2D array of sound speed valuesrin (
array) – range coordinate for c arrayszin (
array) – depth coordinate for c arrays
- Returns:
theta (float) – angle of ray (degrees)
c (float) – sound speed at ray state (m/s)
- pygenray.integration_processes.ray_bounding_box_event(x, y, cin, cpin, rin, zin, depths, depth_ranges)¶
Ray Bounding Box Event - trigger when ray position goes outside of bounding box. Bounding box is defined as the box where sound speed is defined.
- Returns:
bbox – True if ray is outside bounding box, False otherwise
- Return type:
- pygenray.integration_processes.surface_bounce(x, y, cin, cpin, rin, zin, depths, depth_ranges)¶
Surface event: only trigger when approaching surface from below
- pygenray.integration_processes.vertical_ray(x, y, cin, cpin, rin, zin, depths, depth_ranges)¶
Vertical ray event: trigger when ray is vertical (θ = 90 degrees)
Ray Objects¶
- class pygenray.ray_objects.EigenRays(receiver_depths, eigenray_dict, environment, num_eigenrays, num_eigenrays_found, failed_eray_theta_brackets)¶
EigenRays Object - python object that store all parameters associated with eigen rays for given receiver depths
- Parameters:
receiver_depths (list) – List of receiver depths for which eigen rays are computed.
eigenray_dict (dict) – dictionary of eigen rays. Key values are indices of receiver depths, and values are lists of pr.Ray objects.
environment (pr.OceanEnvironment2D) – OceanEnvironment2D environment used for ray tracing.
num_eigenrays (dict) – Total number of eigen rays from the RayFan. (i.e. number of zero crossings of (z-rd) and launch angle)
num_eigenrays_found (dict) – Total number of eigen rays found for each receiver depth.
failed_eray_theta_brackets (dict) – Dictionary of failed eigen ray theta brackets. Keys are receiver depth indices, and values are lists of tuples (theta1, theta2) that bracket an eigenray from the Ray Fan, but for which an eigen ray wasn’t found for the given ztol and iteration limit.
- receiver_depths¶
List of receiver depths for which eigen rays are computed. This is used to index the eigen rays.
- Type:
- rs¶
dictionary of eigenray ranges. keys are range depth indices. values are arrays of shape (M,N), where M is number of eigen rays and N is number of range steps
- Type:
- ts¶
dictionary of eigenray times. keys are range depth indices. values are arrays of shape (M,N), where M is number of eigen rays and N is number of range steps
- Type:
- zs¶
dictionary of eigenray depths. keys are range depth indices. values are arrays of shape (M,N), where M is number of eigen rays and N is number of range steps
- Type:
- ps¶
dictionary of eigenray ray parameters (sin(θ)/c). keys are range depth indices. values are arrays of shape (M,N), where M is number of eigen rays and N is number of range steps
- Type:
- received_angles¶
dictionary of eigenray launch angles. keys are range depth indices. values are arrays of shape (M,), where M is number of eigen rays
- Type:
- launch_angles¶
dictionary of eigenray launch angles. keys are range depth indices. values are arrays of shape (M,), where M is number of eigen rays
- Type:
- n_bottom¶
dictionary of number of bottom reflections for eigen rays. keys are range depth indices. values are arrays of shape (M,), where M is number of eigen rays
- Type:
- n_surface¶
dictionary of number of surface reflections for eigen rays. keys are range depth indices. values are arrays of shape (M,), where M is number of eigen rays
- Type:
- ray_id¶
Ray ID string with boundary indicator.
- Type:
string
- plot(ridxs=[0], **kwargs)¶
Plot all eigen rays in time-depth space
- Parameters:
ridxs (list) – list of receiver depth indices to plot. Default is [0], which plots the first receiver depth.
- plot_ducted(**kwargs)¶
Plot all eigen rays that don’t interact with boundaries
- class pygenray.ray_objects.Ray(r, y, n_bottom, n_surface, launch_angle=None, source_depth=None)¶
Single Ray Object - python object that store all parameters associated with a single ray.
- Parameters:
- plot(**kwargs)¶
Plot ray in time-depth space
- class pygenray.ray_objects.RayFan(Rays)¶
RayFan Object - python object that store all parameters associated with a ray fan.
- Parameters:
Rays (
list)
- compute_rayids()¶
compute the ray IDs for each ray in ray fan.
- plot_ray_fan(**kwargs)¶
plot ray fan
- plot_time_front(include_lines=False, range_idx=-1, add_colorbar=True, ray_id=False, **kwargs)¶
plot time front. Key word arguments are passed to plt.scatter.
- Parameters:
include_lines (bool) – if True, lines between received rays on time plots are plotted
range_idx (int) – index of the range to plot the time front for
add_colorbar (bool) – if True, a colorbar is added to the plot, Default True
ray_id (bool) – if true, than the colors of the time-front arrivals are the ray IDs
Multi-processing functions¶
mostly for internal use
- class pygenray.ray_objects.EigenRays(receiver_depths, eigenray_dict, environment, num_eigenrays, num_eigenrays_found, failed_eray_theta_brackets)¶
EigenRays Object - python object that store all parameters associated with eigen rays for given receiver depths
- Parameters:
receiver_depths (list) – List of receiver depths for which eigen rays are computed.
eigenray_dict (dict) – dictionary of eigen rays. Key values are indices of receiver depths, and values are lists of pr.Ray objects.
environment (pr.OceanEnvironment2D) – OceanEnvironment2D environment used for ray tracing.
num_eigenrays (dict) – Total number of eigen rays from the RayFan. (i.e. number of zero crossings of (z-rd) and launch angle)
num_eigenrays_found (dict) – Total number of eigen rays found for each receiver depth.
failed_eray_theta_brackets (dict) – Dictionary of failed eigen ray theta brackets. Keys are receiver depth indices, and values are lists of tuples (theta1, theta2) that bracket an eigenray from the Ray Fan, but for which an eigen ray wasn’t found for the given ztol and iteration limit.
- receiver_depths¶
List of receiver depths for which eigen rays are computed. This is used to index the eigen rays.
- Type:
- rs¶
dictionary of eigenray ranges. keys are range depth indices. values are arrays of shape (M,N), where M is number of eigen rays and N is number of range steps
- Type:
- ts¶
dictionary of eigenray times. keys are range depth indices. values are arrays of shape (M,N), where M is number of eigen rays and N is number of range steps
- Type:
- zs¶
dictionary of eigenray depths. keys are range depth indices. values are arrays of shape (M,N), where M is number of eigen rays and N is number of range steps
- Type:
- ps¶
dictionary of eigenray ray parameters (sin(θ)/c). keys are range depth indices. values are arrays of shape (M,N), where M is number of eigen rays and N is number of range steps
- Type:
- received_angles¶
dictionary of eigenray launch angles. keys are range depth indices. values are arrays of shape (M,), where M is number of eigen rays
- Type:
- launch_angles¶
dictionary of eigenray launch angles. keys are range depth indices. values are arrays of shape (M,), where M is number of eigen rays
- Type:
- n_bottom¶
dictionary of number of bottom reflections for eigen rays. keys are range depth indices. values are arrays of shape (M,), where M is number of eigen rays
- Type:
- n_surface¶
dictionary of number of surface reflections for eigen rays. keys are range depth indices. values are arrays of shape (M,), where M is number of eigen rays
- Type:
- ray_id¶
Ray ID string with boundary indicator.
- Type:
string
- plot(ridxs=[0], **kwargs)¶
Plot all eigen rays in time-depth space
- Parameters:
ridxs (list) – list of receiver depth indices to plot. Default is [0], which plots the first receiver depth.
- plot_ducted(**kwargs)¶
Plot all eigen rays that don’t interact with boundaries
- class pygenray.ray_objects.Ray(r, y, n_bottom, n_surface, launch_angle=None, source_depth=None)¶
Single Ray Object - python object that store all parameters associated with a single ray.
- Parameters:
- plot(**kwargs)¶
Plot ray in time-depth space
- class pygenray.ray_objects.RayFan(Rays)¶
RayFan Object - python object that store all parameters associated with a ray fan.
- Parameters:
Rays (
list)
- compute_rayids()¶
compute the ray IDs for each ray in ray fan.
- plot_ray_fan(**kwargs)¶
plot ray fan
- plot_time_front(include_lines=False, range_idx=-1, add_colorbar=True, ray_id=False, **kwargs)¶
plot time front. Key word arguments are passed to plt.scatter.
- Parameters:
include_lines (bool) – if True, lines between received rays on time plots are plotted
range_idx (int) – index of the range to plot the time front for
add_colorbar (bool) – if True, a colorbar is added to the plot, Default True
ray_id (bool) – if true, than the colors of the time-front arrivals are the ray IDs