opensoundscape.localization package
Submodules
opensoundscape.localization.audiomoth_sync module
utilities for synchronizing files recorded with audiomoth GPS firmware
Audio synchronization is a preprocessing step before time-difference-of-arrival based acoustic localization. These tools use the .CSV metadata files generated by audiomoth GPS firmware to generate audio files starting at an exact timestamp and having the exact desired sampling rate by interpolating from the original samples and using GPS timestamps.
Example: ```python audio = Audio.from_file(‘/path/to/20240801_103000.WAV’) pps_table = pps_table = pd.read_csv(‘/path/to/20240801_103000.CSV’, index_col=0)
# create correspondence between GPS timestamps and WAV file sample positions sample_timestamp_table = associate_pps_samples_timestamps(csv_file)
# Resample the audio second-by-second using the GPS timestamps to achieve nominal samping rate resampled_audio = correct_sample_rate(
audio, sample_timestamp_table, desired_sr=48000
)
# save sample_timestamp_table.to_csv(‘samples_and_timestamps.csv’) resampled_audio.save(‘resampled_audio.wav’) ```
See also:` tutorials/audiomoth_sync.py` for an example of syncing an entire dataset
NOTE: This module is in _beta_. Report issues to the opensoundscape GitHub repository.
- opensoundscape.localization.audiomoth_sync.associate_pps_samples_timestamps(pps_table, expected_sr=48000, sr_tol=50, cpu_clock_counter_col='TIMER_COUNT') pandas.DataFrame [source]
Using the CSV File from Audiomoth GPS firmware, creates a table aligning GPS_TIME_SEC values with the corresponding PPS and SAMPLES_WRITTEN Values
This is the information we need to interpolate audio files to achieve “synchronized” files with the exact desired start time and sampling rate throughout the file.
The interpretation of the CSV metadata file can be confusing, because it integrates information from (a) the Audiomoth internal state (samples written and internal clock); (b) a PPS signal, a voltage impulse received from the GPS once per second which is used to provide accurate timing information, i.e. the voltage impulse’s timing corresponds to the timestamp and other metadata that the GPS will subsequently send in an RMC packet; and (c) RMC data, the actual data from the GPS (timestamp, position, positional uncertainty) which arrives some time after the PPS signal (e.g 0.4 seconds). Confusingly, the rows of the CSV metadata file do not contain RMC data that corresponds to the PPS data in the same row. Instead, the firmware writes a line to the csv as soon as it receives a pps impulse (meaning that it hasn’t received the associated RMC data yet). The reason for this is that the RMC packet could be dropped and never arrive. The row includes the current number of samples written to the WAV file at the time of receiving the PPS signal (SAMPLES_WRITTEN), and the internal AudioMoth clock time (AUDIOMOTH_TIME). It also includes includes the _most recent previous_ RMC data (LAST_RMC_GPS_TIME) and the internal AudioMoth clock time for when the RMC packet was recieved (LAST_RMC_AUDIOMOTH_TIME). In a typical scenario, where PPS and RMC signals are arriving every second, this means that LAST_RMC_GPS_TIME[N+1] gives the accurate timestamp at which TOTAL_SAMPLES[N] had been recorded. However, PPS and/or RMC signals might not always be received every second.
This function performs several checks to try to determine whether a value in the LAST_RMC_GPS_TIME column of row N+1 corresponds to the PPS from row N. If so, row N’s SAMPLES_WRITTEN value gives the number of samples written at the time of the PPS, whose accurate time stamp is given by the LAST_RMC_GPS_TIME of row N+1.
If any of the validation checks fail, we discard the RMC data from that row. As long as there is eventually a good pairing of RMC data (LAST_RMC_GPS_TIME) with a PPS, we will know the number of samples written and the exact time interval during which they were recorded. This allows us to interpolate the audio signal to achieve the nominal sampling rate.
We also check for buffer overflows, which are recorded in the metadata if they occur. A buffer overflow occurs if the memory buffer of samples arriving from the ADC fills up and is cleared before the samples are written to the WAV file on the SD card. If a buffer overflow occurs, the WAV file will be missing samples, leading to misalignment of all samples after the buffer overflow compared to their expected timing. They might occur if the SD card is slow or too many files are being written simultaneously. This function intentionally fails if buffer overlfows are detected. However, it may be possible to correct for buffer overflows by detecting when they occur and inserting zero-valued samples for the missing samples so that subsequent samples occur at the correct time.
- Parameters:
pps_table – dataframe containing the .CSV metadata file written by AudioMoth GPS firmware Set first column as index, e.g. pps_table = pd.read_csv(path, index_col=0). (one CSV file has PPS data for a one WAV audio file)
expected_sr – expected sample rate of the input audio
sr_tol – tolerance for observed vs expected sample rate
cpu_clock_counter_col – name of the column containing the CPU clock counter - default is “TIMER_COUNT” for the public Audiomoth GPS sync firmware - specify “COUNTER” for older firmware such as AudioMothGPSDeploy_1_0_8_Hardware_1_1
- Returns:
a DataFrame with rows in which GPS_TIME_SEC is matched with the corresponding PPS, SAMPLES_WRITTEN, and TIMER_COUNT values
columns: [“PPS”,”SAMPLES_WRITTEN”, “GPS_TIME_SEC”, “GPS_TIME_STR”, “CPU_CLOCK_COUNTER”]
- opensoundscape.localization.audiomoth_sync.correct_sample_rate(audio, sample_timestamp_df, desired_sr, expected_sr=48000, interpolation_method='nearest', sr_warning_tolerance=100)[source]
Resample audio file to a nominal sample rate using PPS metadata
We will resample the audio piece by piece, taking the samples recorded between each sequential pair of accurate GPS timestamps and resampling to the nominal sampling rate
Error/fluxuation in sampling rate is corrected for segements of audio between each RMC datapoint For example, if pps pings are received every second and desired_sr = 48000 the first second is corrected to be exactly 48000, the second second is corrected to be exactly 48000 samples. These are concatenated to form the entire audio file.
- Parameters:
audio (-) – Audio objecte
sample_timestamp_df (-) –
a dataframe of samples written to the wav file and the corresponding RMC time, typically created with associate_pps_samples_timestamps()
must contain 2 required columns: SAMPLES_WRITTEN and RMC_TIME
desired_sr (-) – the desired sample_rate for the files.
expected_sr (-) – the expected sample rate of the input; gives warning if input sr differs by more than sr_warning_tolerance Hz on any interval
interpolation_method (-) –
Can be e.g. “linear”, “quadratic”, “nearest”. See scipy.interpolate.interp1d docs.
nearest is slightly faster than linear
sr_warning_tolerance (-) – gives warning if input sr differs by more than sr_warning_tolerance Hz on any interval
- Returns:
Audio object with desired sample rate and updated .metadata[‘recording_start_time’]
- Effects:
writes out a .wav file with corrected sample_rate
- opensoundscape.localization.audiomoth_sync.parse_RMC_time_to_datetime(timestamp: str) datetime [source]
Takes in the GPS-audiomoth formatted timestamp and returns datetime.datetime localized with UTC timezone. Assumes the input time string is given in UTC.
- Parameters:
timestamp – string in the format used by our GPS-audiomoths ‘%Y-%m-%dT%H:%M:%S.%f’ e.g. “2022-05-14T10:30:00.725”
Returns: datetime object localized to UTC timezone
- opensoundscape.localization.audiomoth_sync.parse_RMC_time_to_seconds(timestamp: str) float [source]
Takes in the GPS-audiomoth formatted timestamp and returns unix epoch time in seconds.
- Parameters:
timestamp – string in the format used by audiomoth GPS firmware ‘%Y-%m-%dT%H:%M:%S.%f’ e.g. “2022-05-14T10:30:00.725”
Returns: time from unix-epoch in seconds (to 10^-6s resolution) as a float
opensoundscape.localization.localization_algorithms module
Tools for localizing audio events from synchronized recording arrays
- opensoundscape.localization.localization_algorithms.calc_speed_of_sound(temperature=20)[source]
Calculate speed of sound in air, in meters per second
Calculate speed of sound for a given temperature in Celsius (Humidity has a negligible effect on speed of sound and so this functionality is not implemented)
- Parameters:
temperature – ambient air temperature in Celsius
- Returns:
the speed of sound in air in meters per second
- opensoundscape.localization.localization_algorithms.gillette_localize(receiver_locations, arrival_times, speed_of_sound)[source]
Uses the Gillette and Silverman [1] localization algorithm to localize a sound event from a set of TDOAs. :param receiver_locations: a list of [x,y] or [x,y,z] locations for each receiver
locations should be in meters, e.g., the UTM coordinate system.
- Parameters:
arrival_times – a list of TDOA times (arrival times) for each receiver The times should be in seconds.
speed_of_sound – speed of sound in m/s
- Returns:
a tuple of (x,y,z) coordinates of the sound source
- Return type:
coords
Algorithm from: [1] M. D. Gillette and H. F. Silverman, “A Linear Closed-Form Algorithm for Source Localization From Time-Differences of Arrival,” IEEE Signal Processing Letters
- opensoundscape.localization.localization_algorithms.least_squares_localize(receiver_locations, arrival_times, speed_of_sound)[source]
Use a least squares optimizer to perform TDOA localization on a sound event, based on a set of TDOA times and receiver locations. :param receiver_locations: a list of [x,y,z] locations for each receiver
locations should be in meters, e.g., the UTM coordinate system.
- Parameters:
arrival_times – a list of TDOA times (onset times) for each recorder The first receiver is the reference receiver, with arrival time 0.
sound (speed of) – speed of sound in m/s
- Returns:
The solution (x,y,z) in meters. In the same number of dimensions as the receiver locations.
- opensoundscape.localization.localization_algorithms.localize(receiver_locations, tdoas, algorithm, speed_of_sound)[source]
Perform TDOA localization on a sound event. If there are not enough receivers to localize the event, return None. :param receiver_locations: a list of [x,y,z] locations for each receiver
locations should be in meters, e.g., the UTM coordinate system.
- Parameters:
tdoas – a list of TDOA times (onset times) for each recorder The times should be in seconds.
speed_of_sound – speed of sound in m/s
algorithm – the algorithm to use for localization Options: ‘soundfinder’, ‘gillette’, ‘least_squares’. See the documentation for the functions soundfinder_localize, gillette_localize, and least_squares_localize for more detail on how each algorithm works.
- Returns:
The estimated source location in meters, in the same number of dimensions as the receiver locations.
- opensoundscape.localization.localization_algorithms.lorentz_ip(u, v=None)[source]
Compute Lorentz inner product of two vectors
For vectors u and v, the Lorentz inner product for 3-dimensional case is defined as
u[0]*v[0] + u[1]*v[1] + u[2]*v[2] - u[3]*v[3]
Or, for 2-dimensional case as
u[0]*v[0] + u[1]*v[1] - u[2]*v[2]
- Parameters:
u – vector with shape either (3,) or (4,)
v – vector with same shape as x1; if None (default), sets v = u
- Returns:
value of Lorentz IP
- Return type:
float
- opensoundscape.localization.localization_algorithms.soundfinder_localize(receiver_locations, arrival_times, speed_of_sound, invert_alg='gps', center=True, pseudo=True)[source]
Use the soundfinder algorithm to perform TDOA localization on a sound event Localize a sound event given relative arrival times at multiple receivers. This function implements a localization algorithm from the equations described in [1]. Localization can be performed in a global coordinate system in meters (i.e., UTM), or relative to recorder locations in meters.
This implementation follows [2] with corresponding variable names.
- Parameters:
receiver_locations – a list of [x,y,z] locations for each receiver locations should be in meters, e.g., the UTM coordinate system.
arrival_times – a list of TDOA times (onset times) for each recorder The times should be in seconds.
sound (speed of) – speed of sound in m/s
invert_alg – what inversion algorithm to use (only ‘gps’ is implemented)
center – whether to center recorders before computing localization result. Computes localization relative to centered plot, then translates solution back to original recorder locations. (For behavior of original Sound Finder, use True)
pseudo – whether to use the pseudorange error (True) or sum of squares discrepancy (False) to pick the solution to return (For behavior of original Sound Finder, use False. However, in initial tests, pseudorange error appears to perform better.)
- Returns:
The solution (x,y,z) in meters.
[1] Wilson, David R., Matthew Battiston, John Brzustowski, and Daniel J. Mennill. “Sound Finder: A New Software Approach for Localizing Animals Recorded with a Microphone Array.” Bioacoustics 23, no. 2 (May 4, 2014): 99–112. https://doi.org/10.1080/09524622.2013.827588.
[2] Global locationing Systems handout, 2002 http://web.archive.org/web/20110719232148/http://www.macalester.edu/~halverson/math36/GPS.pdf
- opensoundscape.localization.localization_algorithms.travel_time(source, receiver, speed_of_sound)[source]
Calculate time required for sound to travel from a souce to a receiver
- Parameters:
source – cartesian location [x,y] or [x,y,z] of sound source, in meters
receiver – cartesian location [x,y] or [x,y,z] of sound receiver, in meters
speed_of_sound – speed of sound in m/s
- Returns:
time in seconds for sound to travel from source to receiver
opensoundscape.localization.position_estimate module
- class opensoundscape.localization.position_estimate.PositionEstimate(location_estimate, class_name=None, receiver_files=None, receiver_locations=None, tdoas=None, cc_maxs=None, start_timestamp=None, receiver_start_time_offsets=None, duration=None, distance_residuals=None)[source]
Bases:
object
class created by localization algorithms to store estimated location of a sound source
Information about the sound event: - location_estimate: 2 or 3 element array of floats, estimated location of the sound source - class_name: string, class name of the sound source - start_timestamp: timestamp of the start of the event - duration: duration of the event in seconds
Also contains information about the receivers used for localization, and intermediate outputs of the localization algorithm: - receiver_files: list of file paths to audio files used for localization - receiver_start_time_offsets: list of floats, time from start of audio to start of event
for each receiver
receiver_locations: list of receiver locations
tdoas: list of time differences of arrival computed with cross correlation
cc_maxs: list of cross correlation maxima
- Parameters:
estimate (location) – 2 or 3 element array of floats, estimated location of the sound source
args (see SpatialEvent() for other)
- classmethod from_dict(dictionary)[source]
Recover PositionEstimate from dictionary, eg loaded from json
- load_aligned_audio_segments(start_offset=0, end_offset=0)[source]
Load audio segments from each receiver at start times offset by tdoas
This is useful for checking for correct alignment: the sound source should be aligned across all receivers if the tdoas (time differences of arrival) are correct.
Note: requires self.receiver_start_time_offsets to determine starting position for each file
- Parameters:
start_offset – float, amount of time before start of event to include in audio segment
end_offset – float, amount of time after end of event to include in audio segment
- Returns:
list of Audio objects, one for each receiver
Example:
`python low_f = 1000 # Hz high_f = 5000 # Hz audio_segments = example.load_aligned_audio_segments() specs = [Spectrogram.from_audio(a).bandpass(low_f,high_f) for a in audio_segments] plt.pcolormesh(np.vstack([s.spectrogram for s in specs]),cmap='Greys') `
- property residual_rms
- opensoundscape.localization.position_estimate.df_to_positions(df)[source]
convert a pd DataFrame to list of PositionEstimate objects
works best if df was created using events_to_df
- Parameters:
df – pd.DataFrame with columns for each attribute of a PositionEstimate object
- Returns:
list of PositionEstimate objects
See also: events_to_df, PositionEstimate.from_dict, PositionEstimate.to_dict
- opensoundscape.localization.position_estimate.positions_to_df(list_of_events)[source]
convert a list of PositionEstimate objects to pd DataFrame
- Parameters:
list_of_events – list of PositionEstimate objects
- Returns:
pd.DataFrame with columns for each attribute of a PositionEstimate object
See also: df_to_events, PositionEstimate.from_dict, PositionEstimate.to_dict
opensoundscape.localization.spatial_event module
- class opensoundscape.localization.spatial_event.SpatialEvent(receiver_files, receiver_locations, max_delay, min_n_receivers=3, receiver_start_time_offsets=None, start_timestamp=None, duration=None, class_name=None, bandpass_range=None, cc_threshold=0, cc_filter=None, speed_of_sound=343)[source]
Bases:
object
Class that estimates the location of a single sound event
Uses receiver locations and time-of-arrival of sounds to estimate sound source location
- estimate_location(localization_algorithm='gillette', use_stored_tdoas=True)[source]
Estimate spatial location of this event.
This method first estimates the time delays (TDOAS) using cross-correlation, then estimates the location from those TDOAS. Localization is performed in 2d or 3d according to the dimensions of self.receiver_locations (x,y) or (x,y,z) Note: if self.tdoas or self.receiver_locations is None, first calls self._estimate_delays() to estimate the time delays.
If you want to change some parameters of the localization (e.g. try a different cc_threshold, or bandpass_range), you can set the appropriate attribute (e.g. self.cc_threshold = 0.01) before calling self.estimate_location().
- Parameters:
localization_algorithm (-) – algorithm to use for estimating the location of a sound event from the locations and time delays of a set of detections. Options are ‘gillette’ or ‘soundfinder’. Default is ‘gillette’.
use_stored_tdoas (-) –
if True, uses the tdoas stored in self.tdoas to estimate the location.
If False, first calls self._estimate_delays() to estimate the tdoas. default: True
- Returns:
PositionEstimate object with .location_estimate and other attributes - if localization is not successful, .location_estimate attribute of returned object is
None
- opensoundscape.localization.spatial_event.calculate_tdoa_residuals(receiver_locations, tdoas, location_estimate, speed_of_sound)[source]
Calculate the residual distances of the TDOA localization algorithm
The residual represents the discrepancy between (difference in distance of each reciever to estimated location) and (observed tdoa), and has units of meters. Residuals are calculated as follows:
- expected = calculated time difference of arrival between reference and
another receiver, based on the locations of the receivers and estimated event location
observed = observed tdoas provided to localization algorithm
residual time = expected - observed (in seconds)
residual distance = speed of sound * residual time (in meters)
- Parameters:
receiver_location – The list of coordinates (in m) of each receiver, as [x,y] for 2d or or [x,y,z] for 3d.
tdoas – List of time delays of arival for the sound at each receiver, relative to the first receiver in the list (tdoas[0] should be 0)
location_estimate – The estimated location of the sound, as (x,y) or (x,y,z) in meters
speed_of_sound – The speed of sound in m/s
- Returns:
np.array containing the residuals in units of meters, one per receiver
- opensoundscape.localization.spatial_event.df_to_events(df)[source]
convert a pd DataFrame to list of SpatialEvent objects
works best if df was created using events_to_df
- Parameters:
df – pd.DataFrame with columns for each attribute of a SpatialEvent object
- Returns:
list of SpatialEvent objects
See also: events_to_df, SpatialEvent.from_dict, SpatialEvent.to_dict
- opensoundscape.localization.spatial_event.events_to_df(list_of_events)[source]
convert a list of SpatialEvent objects to pd DataFrame
- Parameters:
list_of_events – list of SpatialEvent objects
- Returns:
pd.DataFrame with columns for each attribute of a SpatialEvent object
See also: df_to_events, SpatialEvent.from_dict, SpatialEvent.to_dict
opensoundscape.localization.synchronized_recorder_array module
- class opensoundscape.localization.synchronized_recorder_array.GetStartTimestamp[source]
Bases:
object
- class opensoundscape.localization.synchronized_recorder_array.SynchronizedRecorderArray(file_coords, speed_of_sound=343)[source]
Bases:
object
Class with utilities for localizing sound events from array of recorders
- localize_detections()[source]
Attempt to localize a sound event for each detection of each class. First, creates candidate events with: create_candidate_events()
Create SpatialEvent objects for all simultaneous, spatially clustered detections of a class
Then, attempts to localize each candidate event via time delay of arrival information: For each candidate event:
calculate relative time of arrival with generalized cross correlation (event.estimate_delays())
- if enough cross correlation values exceed a threshold, attempt to localize the event
using the time delays and spatial locations of each receiver with event.estimate_location()
- if the residual distance rms value is below a cutoff threshold, consider the event
to be successfully localized
- check_files_missing_coordinates(detections)[source]
Check that all files in detections have coordinates in file_coords :returns:
a list of files that are in detections but not in file_coords
- create_candidate_events(detections, min_n_receivers, max_receiver_dist, cc_threshold, bandpass_ranges, cc_filter, max_delay=None)[source]
Takes the detections dictionary and groups detections that are within max_receiver_dist of each other.
- Parameters:
detections –
a dictionary of detections, with multi-index (file,start_time,end_time,start_timestamp), and one column per class with 0/1 values for non-detection/detection where start_timestamp contains timzeone-aware datetime.datetime objects corresponding to start_time - If start_timestamp index not present, attempts to automatically determine the start_timestamps file from audio file metadata (‘recording_start_time’ field). This will generally succeed with files generated by AudioMoth or other devices for which metadata parsing is supported and recording start time is given in the metadata.
Example of adding start_timestamp to multi-index of detections df manually: ```python # assume we have a csv with columns file, start_time, end_time, class1, … # e.g., the output of and opensoundscape CNN prediction workflow detections = pd.read_csv(‘my_opso_detections.csv’,index_col=[0,1,2]) # let’s assume all files started recording at the same time for this example: import pytz import datetime tz = pytz.timezone(‘America/New_York’) detections = detections.reset_index(drop=False) detections[‘start_timestamp’]=[
datetime.datetime(2024,1,1,10,0,0).astimezonoe(tz)+ datetime.timedelta(detections.at[i,’start_time’]) for i in detections.index
] # add the start_timestamp column into the multi-index detections = detections.reset_index(drop=False).set_index(
[“file”, “start_time”, “end_time”, “start_timestamp”]
min_n_receivers – if fewer nearby receivers have a simultaneous detection, do not create candidate event
max_receiver_dist – the maximum distance between recorders to consider a detection as a single event
bandpass_ranges – dictionary of bandpass ranges for each class, or a single range for all classes, or None - if None, does not bandpass audio before cross-correlation - if a single range [low_f,high_f], uses this bandpass range for SpatialEvents of all classes - if a dictionary, must have keys for each class in detections.columns
max_delay – the maximum delay (in seconds) to consider between receivers for a single event if None, defaults to max_receiver_dist / self.speed_of_sound
- Returns:
a list of SpatialEvent objects to attempt to localize
- localize_detections(detections, max_receiver_dist, localization_algorithm='gillette', max_delay=None, min_n_receivers=3, cc_threshold=0, cc_filter='phat', bandpass_ranges=None, residual_threshold=numpy.inf, return_unlocalized=False, num_workers=1)[source]
Attempt to localize positions of sound sources from clip-level detections across multiple recorders
Wraps self.create_candidate_events() and self.localize_events()
Algorithm
The user provides a table of class detections from each recorder with timestamps. The object’s self.file_coords dataframe contains a table listing the spatial location of the recorder for each unique audio file in the table of detections. The audio recordings must be synchronized such that timestamps from each recording correspond to the exact same real-world time.
Localization of sound events proceeds in four steps:
Grouping of detections into candidate events (self.create_candidate_events()):
Simultaneous and spatially clustered detections of a class are selected as targets for localization of a single real-world sound event.
For each detection of a species, the grouping algorithm treats the reciever with the detection as a “reference receiver”, then selects all detections of the species at the same time and within max_receiver_dist of the reference reciever (the “surrounding detections”). This selected group of simulatneous, spatially-clustered detections of a class beomes one “candidate event” for subsequent localization.
If the number of recorders in the candidate event is fewer than min_n_receivers, the candidate event is discarded.
This step creates a highly redundant set of candidate events to localize, because each detection is treated separately with its recorder as the ‘reference recorder’. Thus, the localized events created by this algorithm may contain multiple instances representing the same real-world sound event.
Estimate time delays with cross correlation:
For each candidate event, the time delay between the reference reciever’s detection and the surrounding recorders’ detections is estimated through generalized cross correlation.
If bandpass_ranges are provided, cross correlation is performed on audio that has been bandpassed to class-specific low and high frequencies.
If the max value of the cross correlation is below cc_threshold, the corresponding time delay is discarded and not used during localization. This provides a way of filtering out undesired time delays that do not correspond to two recordings of the same sound event.
If the number of estimated time delays in the candidate event is fewer than min_n_receivers after filtering by cross correlation threshold, the candidate event is discarded.
Estimate locations
The location of the event is estimated based on the locations and time delays of each detection.
location estimation from the locations and time delays at a set of receivers is performed using one of two algorithms, described in localization_algorithm below.
Filter by spatial residual error
The residual errors represent descrepencies between (a) time of arrival of the event at a reciever and (b) distance from reciever to estimated location.
Estimated locations are discarded if the root mean squared spatial residual is greater than residual_rms_threshold
- param detections:
a dictionary of detections, with multi-index (file,start_time,end_time), and one column per class with 0/1 values for non-detection/detection The times in the index imply the same real world time across all files: eg 0 seconds assumes that the audio files all started at the same time, not on different dates/times
- param max_receiver_dist:
float (meters) Radius around a recorder in which to use other recorders for localizing an event. Simultaneous detections at receivers within this distance (meters) of a receiver with a detection will be used to attempt to localize the event.
- param max_delay:
float, optional Maximum absolute value of time delay estimated during cross correlation of two signals For instance, 0.2 means that the maximal cross-correlation in the range of delays between -0.2 to 0.2 seconds will be used to estimate the time delay. if None (default), the max delay is set to max_receiver_dist / self.speed_of_sound
- param min_n_receivers:
int Minimum number of receivers that must detect an event for it to be localized [default: 3]
- param localization_algorithm:
str, optional algorithm to use for estimating the location of a sound event from the locations and time delays of a set of detections. [Default: ‘gillette’] Options:
‘gillette’: linear closed-form algorithm of Gillette and Silverman 2008 [1]
‘soundfinder’: GPS location algorithm of Wilson et al. 2014 [2]
- param cc_threshold:
float, optional Threshold for cross correlation: if the max value of the cross correlation is below this value, the corresponding time delay is discarded and not used during localization. Default of 0 does not discard any delays.
- param cc_filter:
str, optional Filter to use for generalized cross correlation. See signalprocessing.gcc function for options. Default is “phat”.
- param bandpass_ranges:
dict, None, or 2-element list/tuple, optional. Any of: - dict: {“class name”: [low_f, high_f]} for audio bandpass filtering during cross correlation with keys for each class - list/tuple: [low_f,high_f] for all classes - None [Default]: does not bandpass audio.
Note: Bandpassing audio to the frequency range of the relevant sound is recommended for best cross correlation results.
- param residual_threshold:
discard localized events if the root mean squared residual of the
- param TDOAs exceeds this value:
np.inf does not filter out any
- type TDOAs exceeds this value:
a distance in meters
- param events by residual]:
- param return_unlocalized:
bool, optional. If True, returns a second value, the list of PositionEstimates that either failed to localize or was filtered out, for example because too few receivers had detections, or too few receivers passed the cc_threshold, or the TDOA residuals were too high.
- param num_workers:
int, optional. Number of workers to use for parallelization. Default is 1 (no parallelization, no multiprocessing). If > 1, uses joblib.Parallel to parallelize
- returns:
- If return_unlocalized is False,
returns a list of localized positons, each of which is a PositionEstimate object with a .location_estimate attribute (= None if localization failed)
- If return_unlocalized is True, returns 2 lists:
list of localized positions, list of un-localized positions
[1] M. D. Gillette and H. F. Silverman, “A Linear Closed-Form Algorithm for Source Localization From Time-Differences of Arrival,” IEEE Signal Processing Letters
[2] Wilson, David R., Matthew Battiston, John Brzustowski, and Daniel J. Mennill. “Sound Finder: A New Software Approach for Localizing Animals Recorded with a Microphone Array.” Bioacoustics 23, no. 2 (May 4, 2014): 99–112. https://doi.org/10.1080/09524622.2013.827588.
- make_nearby_files_dict(r_max)[source]
create dictinoary listing nearby files for each file
pre-generate a dictionary listing all close files for each audio file dictionary will have a key for each audio file, and value listing all other receivers within r_max of that receiver
eg {ARU_0.mp3: [ARU_1.mp3, ARU_2.mp3…], ARU_1… }
Note: could manually create this dictionary to only list _simulataneous_ nearby files if the detection dataframe contains files from different times
- The returned dictionary is used in create_candidate_events as a look-up table for
recordings nearby a detection in any given file
- Parameters:
r_max – maximum distance from each recorder in which to include other recorders in the list of ‘nearby recorders’, in meters
- Returns:
dictionary with keys for each file and values = list of nearby recordings