Passive Acoustics Module

The passive acoustics module provides a set of functions for analyzing and visualizing passive acoustic monitoring data deployed in water bodies. This package reads in raw .wav files and conducts basic acoustics analysis and visualization.

To start using the module, import it directly from MHKiT: from mhkit import acoustics. The analysis functions are available directly from the main import, while the I/O and graphics submodules are available from acoustics.io and acoustics.graphics, respectively. The base functions are intended to be used on top of the I/O submodule, and include functionality to calibrate data, create spectral densities, sound pressure levels, and time or band aggregate spectral data.

This module contains key functions for passive acoustics analysis, designed to process and analyze sound pressure data from .wav files in the frequency and time domains. The functions herein build on each other, with a structured flow that facilitates the calculation of sound pressure levels, spectral densities, and banded averages, based on input audio data.

The following functionality is provided:

  1. Frequency Validation and Warning:

    • _fmax_warning: Ensures specified maximum frequency does not exceed the Nyquist frequency, adjusting if necessary to avoid aliasing.

  2. Shallow Water Cutoff Frequency:

    • minimum_frequency: Calculates the minimum frequency cutoff based on water depth and the speed of sound in water and seabed materials.

  3. Spectral Density Calculations:

    • sound_pressure_spectral_density: Computes the mean square sound pressure spectral density using FFT binning with Hanning windowing and 50% overlap.

  4. Calibration:

    • apply_calibration: Applies calibration adjustments to the spectral density data using a sensitivity curve, filling missing values as specified.

  5. Spectral Density Level Calculation:

    • sound_pressure_spectral_density_level: Converts mean square spectral density values to sound pressure spectral density levels in dB.

  6. Spectral Density Aggregation:

    • band_aggregate: Aggregates spectral density data into fractional octave bands using specified statistical methods (e.g., median, mean).

    • time_aggregate: Aggregates spectral density data into specified time windows using similar statistical methods.

  7. Sound Pressure Level Calculation:

    • sound_pressure_level: Computes the overall sound pressure level within a frequency band from mean square spectral density.

  8. Frequency-Banded Sound Pressure Level:

    • _band_sound_pressure_level: Helper function for calculating sound pressure levels over specified frequency bandwidths.

    • third_octave_sound_pressure_level and decidecade_sound_pressure_level: Compute sound pressure levels across third-octave and decidecade bands, respectively.

minimum_frequency

Estimate the shallow water cutoff frequency based on the speed of sound in the water column and the speed of sound in the seabed material (generally ranges from 1450 - 1800 m/s)

sound_pressure_spectral_density

Calculates the mean square sound pressure spectral density from audio samples split into FFTs with a specified bin length in seconds, using Hanning windowing with 50% overlap.

apply_calibration

Applies custom calibration to spectral density values.

sound_pressure_spectral_density_level

Calculates the sound pressure spectral density level from the mean square sound pressure spectral density.

band_aggregate

Reorganizes spectral density level frequency tensor into fractional octave bands and applies a function to them.

time_aggregate

Reorganizes spectral density level frequency tensor into time windows and applies a function to them.

sound_pressure_level

Calculates the sound pressure level in a specified frequency band from the mean square sound pressure spectral density.

third_octave_sound_pressure_level

Calculates the sound pressure level in third octave bands directly from the mean square sound pressure spectral density.

decidecade_sound_pressure_level

Calculates the sound pressure level in decidecade bands directly from the mean square sound pressure spectral density.

IO

This submodule provides input/output functions for passive acoustics data, focusing on hydrophone recordings stored in WAV files. The main functionality includes reading and processing hydrophone data from various manufacturers and exporting audio files for easy playback and analysis.

Supported Hydrophone Models

  • SoundTrap (Ocean Instruments)

  • icListen (Ocean Sonics)

Functions Overview

  1. Data Reading:

    • read_hydrophone: Main function to read a WAV file from a hydrophone and convert it to either a voltage or pressure time series, depending on the availability of sensitivity data.

    • read_soundtrap: Wrapper for reading Ocean Instruments SoundTrap hydrophone files, automatically using appropriate metadata.

    • read_iclisten: Wrapper for reading Ocean Sonics icListen hydrophone files, including metadata processing to apply hydrophone sensitivity for direct sound pressure calculation.

  2. Audio Export:

    • export_audio: Converts processed sound pressure data back into a WAV file format, with optional gain adjustment to improve playback quality.

  3. Data Extraction:

    • _read_wav_metadata: Extracts metadata from a WAV file, including bit depth and other header information.

    • _calculate_voltage_and_time: Converts raw WAV data into voltage values and generates a time index based on the sampling frequency.

read_hydrophone

Read .wav file from a hydrophone.

read_soundtrap

Read .wav file from an Ocean Instruments SoundTrap hydrophone.

read_iclisten

Read .wav file from an Ocean Sonics icListen "Smart" hydrophone.

export_audio

Creates human-scaled audio file from underwater recording.

mhkit.acoustics.io.read_hydrophone(filename: str | Path, peak_voltage: int | float, sensitivity: int | float | None = None, gain: int | float = 0, start_time: str = '2024-01-01T00:00:00') DataArray[source]

Read .wav file from a hydrophone. Returns voltage timeseries if sensitivity not provided, returns pressure timeseries if it is provided.

Parameters:
  • filename (str or pathlib.Path) – Input filename

  • peak_voltage (int or float) – Peak voltage supplied to the analog to digital converter (ADC) in V. (Or 1/2 of the peak to peak voltage).

  • sensitivity (int or float) – Hydrophone calibration sensitivity in dB re 1 V/uPa. Should be negative. Default: None.

  • gain (int or float) – Amplifier gain in dB re 1 V/uPa. Default 0.

  • start_time (str) – Start time in the format yyyy-mm-ddTHH:MM:SS

Returns:

out (numpy.array) – Sound pressure [Pa] or Voltage [V] indexed by time[s]

mhkit.acoustics.io.read_soundtrap(filename: str, sensitivity: int | float | None = None, gain: int | float = 0) DataArray[source]

Read .wav file from an Ocean Instruments SoundTrap hydrophone. Returns voltage timeseries if sensitivity not provided, returns pressure timeseries if it is provided.

Parameters:
  • filename (str) – Input filename.

  • sensitivity (int or float, optional) – Hydrophone calibration sensitivity in dB re 1 V/μPa. Should be negative. Default is None.

  • gain (int or float) – Amplifier gain in dB re 1 V/μPa. Default is 0.

Returns:

out (xarray.DataArray) – Sound pressure [Pa] or Voltage [V] indexed by time[s].

mhkit.acoustics.io.read_iclisten(filename: str, sensitivity: int | float | None = None, use_metadata: bool = True) DataArray[source]

Read .wav file from an Ocean Sonics icListen “Smart” hydrophone. Returns voltage timeseries if sensitivity not provided, returns pressure timeseries if it is provided.

Parameters:
  • filename (str) – Input filename.

  • sensitivity (int or float, optional) – Hydrophone calibration sensitivity in dB re 1 V/μPa. Should be negative. Default is None.

  • use_metadata (bool) – If True and sensitivity is None, applies sensitivity value stored in the .wav file’s LIST block. If False and sensitivity is None, a sensitivity value isn’t applied.

Returns:

out (xarray.DataArray) – Sound pressure [Pa] or Voltage [V] indexed by time[s].

mhkit.acoustics.io.export_audio(filename: str, pressure: DataArray, gain: int | float = 1) None[source]

Creates human-scaled audio file from underwater recording.

Parameters:
  • filename (str) – Output filename for the WAV file (without extension).

  • pressure (xarray.DataArray) –

    Sound pressure data with attributes:
    • ’values’ (numpy.ndarray): Pressure data array.

    • ’sensitivity’ (int or float): Sensitivity of the hydrophone in dB.

    • ’fs’ (int or float): Sampling frequency in Hz.

  • gain (int or float, optional) – Gain to multiply the original time series by. Default is 1.

Returns:

None

Graphics

This submodule provides essential plotting functions for visualizing passive acoustics data. The functions allow for customizable plotting of sound pressure spectral density levels across time and frequency dimensions.

Each plotting function leverages the flexibility of Matplotlib, allowing for passthrough of Matplotlib keyword arguments via **kwargs, making it easy to modify plot aspects such as color, scale, and label formatting.

Key Functions

  1. plot_spectrogram:

    • Generates a spectrogram plot from sound pressure spectral density level data, with a logarithmic frequency scale by default for improved readability of acoustic data.

  2. plot_spectra:

    • Produces a spectral density plot with a log-transformed x-axis, allowing for clear visualization of spectral density across frequency bands.

plot_spectrogram

Plots the spectrogram of the sound pressure spectral density level.

plot_spectra

Plots spectral density.

mhkit.acoustics.graphics.plot_spectrogram(spsdl: DataArray, fmin: int = 10, fmax: int = 100000, fig: figure = None, ax: Axes = None, **kwargs) Tuple[figure, Axes][source]

Plots the spectrogram of the sound pressure spectral density level.

Parameters:
  • spsdl (xarray DataArray (time, freq)) – Mean square sound pressure spectral density level in dB rel 1 uPa^2/Hz

  • fmin (int) – Lower frequency band limit (lower limit of the hydrophone). Default: 10 Hz

  • fmax (int) – Upper frequency band limit (Nyquist frequency). Default: 100000 Hz

  • fig (matplotlib.pyplot.figure) – Figure handle to plot on

  • ax (matplotlib.pyplot.axis) – Figure axis containing plot objects

  • kwargs (dict) – Dictionary of matplotlib function keyword arguments

Returns:

  • fig (matplotlib.pyplot.figure) – Figure handle of plot

  • ax (matplotlib.pyplot.Axes) – Figure plot axis

mhkit.acoustics.graphics.plot_spectra(spsdl: DataArray, fmin: int = 10, fmax: int = 100000, fig: figure = None, ax: Axes = None, **kwargs) Tuple[figure, Axes][source]

Plots spectral density. X axis is log-transformed.

Parameters:
  • spsdl (xarray DataArray (time, freq)) – Mean square sound pressure spectral density level in dB rel 1 uPa^2/Hz

  • fmin (int) – Lower frequency band limit (lower limit of the hydrophone). Default: 10 Hz

  • fmax (int) – Upper frequency band limit (Nyquist frequency). Default: 100000 Hz

  • fig (matplotlib.pyplot.figure) – Figure handle to plot on

  • ax (matplotlib.pyplot.Axes) – Figure axis containing plot objects

  • kwargs (dict) – Dictionary of matplotlib function keyword arguments

Returns:

  • fig (matplotlib.pyplot.figure) – Figure handle of plot

  • ax (matplotlib.pyplot.Axes) – Figure plot axis