frequency_differencing#
- echopype.mask.frequency_differencing(source_Sv: Dataset | str | Path, storage_options: dict | None = {}, freqABEq: str | None = None, chanABEq: str | None = None) DataArray#
Create a mask based on the differences of Sv values using a pair of frequencies. This method is often referred to as the “frequency-differencing” or “dB-differencing” method.
- Parameters:
- source_Sv: xr.Dataset or str or pathlib.Path
If a Dataset this value contains the Sv data to create a mask for, else it specifies the path to a zarr or netcdf file containing a Dataset. This input must correspond to a Dataset that has the coordinate
channeland variablesfrequency_nominalandSv.- storage_options: dict, optional
Any additional parameters for the storage backend, corresponding to the path provided for
source_Sv- freqABEq: string, optional
The frequency differencing criteria. Only one of
freqABandchanABshould be provided, and not both.- chanAB: string, optional
The frequency differencing criteria in terms of channel names where channel names in the criteria are enclosed in double quotes. Only one of
freqABandchanABshould be provided, and not both.
- Returns:
- xr.DataArray
A DataArray containing the mask for the Sv data. Regions satisfying the thresholding criteria are filled with
True, else the regions are filled withFalse.
- Raises:
- ValueError
If neither
freqABEqorchanABEqare given- ValueError
If both
freqABEqandchanABEqare given- TypeError
If any input is not of the correct type
- ValueError
If either
freqABEqorchanABEqare provided and the extractedfreqABorchanABdoes not contain 2 distinct elements- ValueError
If
freqABEqcontains values that are not contained infrequency_nominal- ValueError
If
chanABEqcontains values that not contained inchannel- ValueError
If
operatoris not one of the following:">", "<", "<=", ">=", "=="- ValueError
If the path provided for
source_Svis not a valid path- ValueError
If
freqABEqorchanABEqis provided and the Dataset produced bysource_Svdoes not contain the coordinatechanneland variablefrequency_nominal
Notes
This function computes the frequency differencing as follows:
Sv_freqA - Sv_freqB operator diff. Thus, ifoperator = "<"anddiff = "5"the following would be calculated:Sv_freqA - Sv_freqB < 5.Examples
Compute frequency-differencing mask using a mock Dataset and channel selection:
>>> n = 5 # set the number of ping times and range samples ... >>> # create mock Sv data >>> Sv_da = xr.DataArray(data=np.stack([np.arange(n**2).reshape(n,n), np.identity(n)]), ... coords={"channel": ['chan1', 'chan2'], ... "ping_time": np.arange(n), "range_sample":np.arange(n)}) ... >>> # obtain mock frequency_nominal data >>> freq_nom = xr.DataArray(data=np.array([1.0, 2.0]), ... coords={"channel": ['chan1', 'chan2']}) ... >>> # construct mock Sv Dataset >>> Sv_ds = xr.Dataset(data_vars={"Sv": Sv_da, "frequency_nominal": freq_nom}) ... >>> # compute frequency-differencing mask using channel names >>> echopype.mask.frequency_differencing(source_Sv=Sv_ds, storage_options={}, ... freqABEq=None, chanABEq = '"chan1" - "chan2">=10.0dB') <xarray.DataArray 'mask' (ping_time: 5, range_sample: 5)> array([[False, False, False, False, False], [False, False, False, False, False], [ True, True, True, True, True], [ True, True, True, True, True], [ True, True, True, True, True]]) Coordinates: * ping_time (ping_time) int64 0 1 2 3 4 * range_sample (range_sample) int64 0 1 2 3 4