Underwater acoustic propagation modeling

Underwater acoustic propagation modeling toolbox.

This toolbox currently uses the Bellhop acoustic propagation model. For this model to work, the acoustic toolbox must be installed on your computer and bellhop.exe should be in your PATH.

arlpy.uwapm.arrivals_to_impulse_response(arrivals, fs, abs_time=False)

Convert arrival times and coefficients to an impulse response.

Parameters:
  • arrivals – arrivals times (s) and coefficients
  • fs – sampling rate (Hz)
  • abs_time – absolute time (True) or relative time (False)
Returns:

impulse response

If abs_time is set to True, the impulse response is placed such that the zero time corresponds to the time of transmission of signal.

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> arrivals = pm.compute_arrivals(env)
>>> ir = pm.arrivals_to_impulse_response(arrivals, fs=192000)
arlpy.uwapm.check_env2d(env)

Check the validity of a 2D underwater environment definition.

Parameters:env – environment definition

Exceptions are thrown with appropriate error messages if the environment is invalid.

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> check_env2d(env)
arlpy.uwapm.compute_arrivals(env, model=None, debug=False)

Compute arrivals between each transmitter and receiver.

Parameters:
  • env – environment definition
  • model – propagation model to use (None to auto-select)
  • debug – generate debug information for propagation model
Returns:

arrival times and coefficients for all transmitter-receiver combinations

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> arrivals = pm.compute_arrivals(env)
>>> pm.plot_arrivals(arrivals)
arlpy.uwapm.compute_eigenrays(env, tx_depth_ndx=0, rx_depth_ndx=0, rx_range_ndx=0, model=None, debug=False)

Compute eigenrays between a given transmitter and receiver.

Parameters:
  • env – environment definition
  • tx_depth_ndx – transmitter depth index
  • rx_depth_ndx – receiver depth index
  • rx_range_ndx – receiver range index
  • model – propagation model to use (None to auto-select)
  • debug – generate debug information for propagation model
Returns:

eigenrays paths

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> rays = pm.compute_eigenrays(env)
>>> pm.plot_rays(rays, width=1000)
arlpy.uwapm.compute_rays(env, tx_depth_ndx=0, model=None, debug=False)

Compute rays from a given transmitter.

Parameters:
  • env – environment definition
  • tx_depth_ndx – transmitter depth index
  • model – propagation model to use (None to auto-select)
  • debug – generate debug information for propagation model
Returns:

ray paths

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> rays = pm.compute_rays(env)
>>> pm.plot_rays(rays, width=1000)
arlpy.uwapm.compute_transmission_loss(env, tx_depth_ndx=0, mode='coherent', model=None, debug=False)

Compute transmission loss from a given transmitter to all receviers.

Parameters:
  • env – environment definition
  • tx_depth_ndx – transmitter depth index
  • mode – coherent, incoherent or semicoherent
  • model – propagation model to use (None to auto-select)
  • debug – generate debug information for propagation model
Returns:

complex transmission loss at each receiver depth and range

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> tloss = pm.compute_transmission_loss(env, mode=pm.incoherent)
>>> pm.plot_transmission_loss(tloss, width=1000)
arlpy.uwapm.create_env2d(**kv)

Create a new 2D underwater environment.

A basic environment is created with default values. To see all the parameters available and their default values:

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> pm.print_env(env)

The environment parameters may be changed by passing keyword arguments or modified later using a dictionary notation:

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d(depth=40, soundspeed=1540)
>>> pm.print_env(env)
>>> env['depth'] = 25
>>> env['bottom_soundspeed'] = 1800
>>> pm.print_env(env)

The default environment has a constant sound speed. A depth dependent sound speed profile be provided as a Nx2 array of (depth, sound speed):

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d(depth=20, soundspeed=[[0,1540], [5,1535], [10,1535], [20,1530]])

A range-and-depth dependent sound speed profile can be provided as a Pandas frame:

>>> import arlpy.uwapm as pm
>>> import pandas as pd
>>> ssp2 = pd.DataFrame({
          0: [1540, 1530, 1532, 1533],     # profile at 0 m range
        100: [1540, 1535, 1530, 1533],     # profile at 100 m range
        200: [1530, 1520, 1522, 1525] },   # profile at 200 m range
        index=[0, 10, 20, 30])             # depths of the profile entries in m
>>> env = pm.create_env2d(depth=20, soundspeed=ssp2)

The default environment has a constant water depth. A range dependent bathymetry can be provided as a Nx2 array of (range, water depth):

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d(depth=[[0,20], [300,10], [500,18], [1000,15]])
arlpy.uwapm.models(env=None, task=None)

List available models.

Parameters:
  • env – environment to model
  • task – arrivals/eigenrays/rays/coherent/incoherent/semicoherent
Returns:

list of models that can be used

>>> import arlpy.uwapm as pm
>>> pm.models()
['bellhop']
>>> env = pm.create_env2d()
>>> pm.models(env, task=coherent)
['bellhop']
arlpy.uwapm.plot_arrivals(arrivals, dB=False, color='blue', **kwargs)

Plots the arrival times and amplitudes.

Parameters:
  • arrivals – arrivals times (s) and coefficients
  • dB – True to plot in dB, False for linear scale
  • color – line color (see Bokeh colors)

Other keyword arguments applicable for arlpy.plot.plot() are also supported.

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> arrivals = pm.compute_arrivals(env)
>>> pm.plot_arrivals(arrivals)
arlpy.uwapm.plot_env(env, surface_color='dodgerblue', bottom_color='peru', tx_color='orangered', rx_color='midnightblue', rx_plot=None, **kwargs)

Plots a visual representation of the environment.

Parameters:
  • env – environment description
  • surface_color

    color of the surface (see Bokeh colors)

  • bottom_color

    color of the bottom (see Bokeh colors)

  • tx_color

    color of transmitters (see Bokeh colors)

  • rx_color

    color of receviers (see Bokeh colors)

  • rx_plot – True to plot all receivers, False to not plot any receivers, None to automatically decide

Other keyword arguments applicable for arlpy.plot.plot() are also supported.

The surface, bottom, transmitters (marker: ‘*’) and receivers (marker: ‘o’) are plotted in the environment. If rx_plot is set to None and there are more than 2000 receivers, they are not plotted.

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d(depth=[[0, 40], [100, 30], [500, 35], [700, 20], [1000,45]])
>>> pm.plot_env(env)
arlpy.uwapm.plot_rays(rays, env=None, invert_colors=False, **kwargs)

Plots ray paths.

Parameters:
  • rays – ray paths
  • env – environment definition
  • invert_colors – False to use black for high intensity rays, True to use white

If environment definition is provided, it is overlayed over this plot using default parameters for arlpy.uwapm.plot_env().

Other keyword arguments applicable for arlpy.plot.plot() are also supported.

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d()
>>> rays = pm.compute_eigenrays(env)
>>> pm.plot_rays(rays, width=1000)
arlpy.uwapm.plot_ssp(env, **kwargs)

Plots the sound speed profile.

Parameters:env – environment description

Other keyword arguments applicable for arlpy.plot.plot() are also supported.

If the sound speed profile is range-dependent, this function only plots the first profile.

>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d(soundspeed=[[ 0, 1540], [10, 1530], [20, 1532], [25, 1533], [30, 1535]])
>>> pm.plot_ssp(env)
arlpy.uwapm.plot_transmission_loss(tloss, env=None, **kwargs)

Plots transmission loss.

Parameters:
  • tloss – complex transmission loss
  • env – environment definition

If environment definition is provided, it is overlayed over this plot using default parameters for arlpy.uwapm.plot_env().

Other keyword arguments applicable for arlpy.plot.image() are also supported.

>>> import arlpy.uwapm as pm
>>> import numpy as np
>>> env = pm.create_env2d(
        rx_depth=np.arange(0, 25),
        rx_range=np.arange(0, 1000),
        min_angle=-45,
        max_angle=45
    )
>>> tloss = pm.compute_transmission_loss(env)
>>> pm.plot_transmission_loss(tloss, width=1000)
arlpy.uwapm.print_env(env)

Display the environment in a human readable form.

Parameters:env – environment definition
>>> import arlpy.uwapm as pm
>>> env = pm.create_env2d(depth=40, soundspeed=1540)
>>> pm.print_env(env)