API reference
This is the API reference for the echoSMs package.
Each type of model is contained in a separate Python class (with name ending in Model
), but with common calling signatures across all model classes, as defined in ScatterModelBase
. There are also classes to provide ready access to the benchmark models and reference model definitions. There are also utility functions.
ScatterModelBase
Bases: ABC
Base class for a class that provides a scattering model.
All scattering models should inherit from this class, have a name that ends with 'Model', and provide initialisation and calculate_ts_single() functions.
Initialise.
Attributes:
Name | Type | Description |
---|---|---|
long_name |
str
|
The long name of the model. |
short_name |
str
|
A short version of the model's long name, typically an acronym. |
analytical_type |
str
|
Whether the model implements an |
boundary_types |
list[str]
|
The types of boundary conditions that the model provides, e.g., 'fixed rigid', 'pressure release', 'fluid filled' |
shapes |
list[str]
|
The target shapes that the model can represent. |
max_ka |
float
|
An approximate maximum ka value that will result in accurate target strength results. Note that ka is often not the only parameter that determines the accuracy of the model (e.g., aspect ratio and incident angle can also affect the accuracy). |
no_expand_parameters |
list[str]
|
The model parameters that are not expanded into Pandas DataFrame columns or
Xarray DataArray coordinates. They will instead end up as a dict in the DataFrame or
DataArray |
Source code in src/echosms/scattermodelbase.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single()
abstractmethod
validate_parameters(p)
abstractmethod
Validate the model parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
p
|
dict
|
Dict containing the model parameters. |
required |
Raises:
Type | Description |
---|---|
ValueError
|
If any of the model parameters are invalid. |
KeyError
|
If any required model parameters are not present. |
Source code in src/echosms/scattermodelbase.py
DCMModel
Bases: ScatterModelBase
Modal series deformed cylinder model (DCM).
This class contains methods to calculate acoustic scatter from finite straight cylinders with various boundary conditions.
Source code in src/echosms/dcmmodel.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single(medium_c, medium_rho, a, b, theta, f, boundary_type, target_c=None, target_rho=None, validate_parameters=True, **kwargs)
Calculate the scatter from a finite cylinder using the modal series deformed cylinder model.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
medium_c
|
float
|
Sound speed in the fluid medium surrounding the target [m/s]. |
required |
medium_rho
|
float
|
Density of the fluid medium surrounding the target [kg/m³]. |
required |
a
|
float
|
Radius of the cylinderical target [m]. |
required |
b
|
float
|
Length of the cylinderical target [m]. |
required |
theta
|
float
|
Pitch angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
f
|
float
|
Frequency to calculate the scattering at [Hz]. |
required |
boundary_type
|
str
|
The model type. Supported model types are given in the |
required |
target_c
|
float
|
Sound speed in the fluid inside the sphere [m/s].
Only required for |
None
|
target_rho
|
float
|
Density of the fluid inside the sphere [kg/m³].
Only required for |
None
|
validate_parameters
|
bool
|
Whether to validate the model parameters. |
True
|
Returns:
Type | Description |
---|---|
float
|
The target strength (re 1 m²) of the target [dB]. |
Notes
The class implements the code in Section B.1 of Jech et al. (2015).
References
Jech, J.M., Horne, J.K., Chu, D., Demer, D.A., Francis, D.T.I., Gorska, N., Jones, B., Lavery, A.C., Stanton, T.K., Macaulay, G.J., Reeder, D.B., Sawada, K., 2015. Comparisons among ten models of acoustic backscattering used in aquatic ecosystem research. Journal of the Acoustical Society of America 138, 3742–3764. https://doi.org/10.1121/1.4937607
Source code in src/echosms/dcmmodel.py
42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 |
|
validate_parameters(params)
Validate the model parameters.
See here for calling details.
Source code in src/echosms/dcmmodel.py
DWBA models
There are several models that use the distorted-wave Born approximation, documented below:
DWBA
Bases: ScatterModelBase
Distorted-wave Born approximation scattering model.
Note
The DWBA model is not yet functional.
Source code in src/echosms/dwbamodel.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single(theta, phi, f, target_rho, target_c, validate_parameters=True)
Distorted-wave Born approximation scattering model.
Implements the distorted-wave Born approximation model for calculating the acoustic backscatter from weakly scattering bodies.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
theta
|
float
|
Pitch angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
phi
|
float
|
Roll angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
f
|
float
|
Frequency to run the model at [Hz] |
required |
target_rho
|
iterable[float]
|
Densities of each material. Must have at least the same number of entries as unique
integers in |
required |
target_c
|
iterable[float]
|
Sound speed of each material. Must have at least the same number of entries as unique
integers in |
required |
validate_parameters
|
bool
|
Whether to validate the model parameters. |
True
|
Returns:
Type | Description |
---|---|
float
|
The target strength (re 1 m²) [dB] of the target. |
Notes
This class implements the method presented in Chu et al. (1993).
References
Chu, D., Foote, K. G., & Stanton, T. K. (1993). Further analysis of target strength measurements of Antarctic krill at 38 and 120 kHz: Comparison with deformed cylinder model and inference or orientation distribution. The Journal of the Acoustical Society of America, 93(5), 2985–2988. https://doi.org/10.1121/1.405818
Source code in src/echosms/dwbamodel.py
validate_parameters(params)
PT-DWBA
Bases: ScatterModelBase
Phase-tracking distorted-wave Born approximation scattering model.
Source code in src/echosms/ptdwbamodel.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single(volume, theta, phi, f, voxel_size, rho, c, validate_parameters=True, **kwargs)
Phase-tracking distorted-wave Born approximation scattering model.
Implements the phase-tracking distorted-wave Born approximation model for calculating the acoustic backscatter from weakly scattering bodies.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
volume
|
Numpy ndarray[int]
|
The object to be modelled as a 3D volume of voxels. Array contents should be 0
for the surrounding medium, then increasing by 1 for each additional material
type (i.e., 1, 2, 3, etc).
Increasing axes indices correspond to increasing x, y, and z values. |
required |
theta
|
float
|
Pitch angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
phi
|
float
|
Roll angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
f
|
float
|
Frequency to run the model at [Hz] |
required |
voxel_size
|
iterable[float]
|
The size of the voxels in |
required |
rho
|
iterable[float]
|
Densities of each material. Must have at least the same number of entries as unique
integers in |
required |
c
|
iterable[float]
|
Sound speed of each material. Must have at least the same number of entries as unique
integers in |
required |
validate_parameters
|
bool
|
Whether to validate the model parameters. |
True
|
Returns:
Type | Description |
---|---|
float
|
The target strength (re 1 m²) [dB] of the target. |
Notes
This class implements the method presented in Jones et. al. (2009). The code is based closely on the Matlab code in Jones (2006).
References
Jones, B. A. (2006). Acoustic scattering of broadband echolocation signals from prey of Blainville's beaked whales: Modeling and analysis. Master of Science, Massachusetts Institute of Technology. https://doi.org/10.1575/1912/1283
Jones, B. A., Lavery, A. C., & Stanton, T. K. (2009). Use of the distorted wave Born approximation to predict scattering by inhomogeneous objects: Application to squid. The Journal of the Acoustical Society of America, 125(1), 73-88. https://doi.org/10.1121/1.3021298
Source code in src/echosms/ptdwbamodel.py
31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 |
|
validate_parameters(params)
SDWBA
Bases: ScatterModelBase
Stochastic distorted-wave Born approximation scattering model.
Note
The SDWBA model is not yet functional.
Source code in src/echosms/sdwbamodel.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single(theta, phi, f, target_rho, target_c, validate_parameters=True)
Stochastic distorted-wave Born approximation scattering model.
Implements the stochastic distorted-wave Born approximation model for calculating the acoustic backscatter from weakly scattering bodies.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
theta
|
float
|
Pitch angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
phi
|
float
|
Roll angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
f
|
float
|
Frequency to run the model at [Hz] |
required |
target_rho
|
iterable[float]
|
Densities of each material. Must have at least the same number of entries as unique
integers in |
required |
target_c
|
iterable[float]
|
Sound speed of each material. Must have at least the same number of entries as unique
integers in |
required |
validate_parameters
|
bool
|
Whether to validate the model parameters. |
True
|
Returns:
Type | Description |
---|---|
float
|
The target strength (re 1 m²) [dB] of the target. |
Notes
This class implements the method presented in Demer & Conti (2003), Demer & Conti (2004), and Conti & Demer (2006).
References
Demer, D. A., & Conti, S. G. (2003). Reconciling theoretical versus empirical target strengths of krill: Effects of phase variability on the distorted-wave Born approximation. ICES Journal of Marine Science, 60, 429-434. https://doi.org/10.1016/S1054-3139(03)00002-X
Demer, D. A., & Conti, S. G. (2004). Reconciling theoretical versus empirical target strengths of krill: Effects of phase variability on the distorted-wave Born approximation. ICES Journal of Marine Science, 61(1), 157-158. https://doi.org/10.1016/j.icesjms.2003.12.003
Conti, S. G., & Demer, D. A. (2006). Improved parameterization of the SDWBA for estimating krill target strength. ICES Journal of Marine Science, 63(5), 928-935. https://doi.org/10.1016/j.icesjms.2006.02.007
Source code in src/echosms/sdwbamodel.py
validate_parameters(params)
ESModel
Bases: ScatterModelBase
Elastic sphere (ES) scattering model.
This class calculates acoustic backscatter from elastic spheres.
Source code in src/echosms/esmodel.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single(medium_c, medium_rho, a, f, target_longitudinal_c, target_transverse_c, target_rho, validate_parameters=True, **kwargs)
Calculate the backscatter from an elastic sphere for one set of parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
medium_c
|
float
|
Sound speed in the fluid medium surrounding the sphere [m/s]. |
required |
medium_rho
|
float
|
Density of the fluid medium surrounding the sphere [kg/m³]. |
required |
a
|
float
|
Radius of the sphere [m]. |
required |
f
|
float
|
Frequency to calculate the scattering at [Hz]. |
required |
target_longitudinal_c
|
float
|
Longitudinal sound speed in the material inside the sphere [m/s]. |
required |
target_transverse_c
|
float
|
Transverse sound speed in the material inside the sphere [m/s]. |
required |
target_rho
|
float
|
Density of the material inside the sphere [kg/m³]. |
required |
validate_parameters
|
bool
|
Whether to validate the model parameters. |
True
|
Returns:
Type | Description |
---|---|
float
|
The target strength (re 1 m²) of the sphere [dB]. |
Notes
The class implements the code in MacLennan (1981).
References
MacLennan, D. N. (1981). The Theory of Solid Spheres as Sonar Calibration Targets. Scottish Fisheries Research Report Number 22. Department of Agriculture and Fisheries for Scotland.
Source code in src/echosms/esmodel.py
validate_parameters(params)
Validate the model parameters.
See here for calling details.
Source code in src/echosms/esmodel.py
KAModel
Bases: ScatterModelBase
Kirchhoff approximation (KA) scattering model.
This class calculates acoustic scatter from arbitrary surfaces.
Source code in src/echosms/kamodel.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single(medium_c, theta, phi, f, mesh, boundary_type, validate_parameters=True, **kwargs)
Calculate the scatter using the ka model for one set of parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
medium_c
|
float
|
Sound speed in the fluid medium surrounding the target [m/s]. |
required |
theta
|
float
|
Pitch angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
phi
|
float
|
Roll angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
f
|
float
|
Frequency to calculate the scattering at [Hz]. |
required |
mesh
|
Any
|
The triangular mesh that defines the scattering surface. This parameter must provide attributes with names of:
A suitable library for creating and manipulating triangular meshes is trimesh. |
required |
boundary_type
|
str
|
The boundary type. Supported types are given in the |
required |
validate_parameters
|
bool
|
Whether to validate the model parameters. |
True
|
Returns:
Type | Description |
---|---|
float
|
The target strength (re 1 m²) of the target [dB]. |
Notes
The class implements the code in Foote (1985).
References
Foote, K. G. (1985). Rather-high-frequency sound scattering of swimbladdered fish. The Journal of the Acoustical Society of America, 78(2), 688–700. https://doi.org/10.1121/1.392438
Source code in src/echosms/kamodel.py
validate_parameters(params)
MSSModel
Bases: ScatterModelBase
Modal series solution (MSS) scattering model.
This class calculates acoustic scatter from spheres and shells with various
boundary conditions, as listed in the boundary_types
class attribute.
Source code in src/echosms/mssmodel.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single(medium_c, medium_rho, a, f, boundary_type, target_c=None, target_rho=None, shell_c=None, shell_rho=None, shell_thickness=None, validate_parameters=True, **kwargs)
Calculate the scatter using the mss model for one set of parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
medium_c
|
float
|
Sound speed in the fluid medium surrounding the target [m/s]. |
required |
medium_rho
|
float
|
Density of the fluid medium surrounding the target [kg/m³]. |
required |
a
|
float
|
Radius of the spherical target [m]. |
required |
f
|
float
|
Frequency to calculate the scattering at [Hz]. |
required |
boundary_type
|
str
|
The boundary type. Supported types are given in the |
required |
target_c
|
float
|
Sound speed in the fluid inside the sphere [m/s].
Only required for |
None
|
target_rho
|
float
|
Density of the fluid inside the sphere [kg/m³].
Only required for |
None
|
shell_c
|
float
|
Sound speed in the spherical shell [m/s].
Only required for |
None
|
shell_rho
|
float
|
Density in the spherical shell [kg/m³].
Only required for |
None
|
shell_thickness
|
float
|
Thickness of the spherical shell [m]. This value is subtracted from |
None
|
validate_parameters
|
bool
|
Whether to validate the model parameters. |
True
|
Returns:
Type | Description |
---|---|
float
|
The target strength (re 1 m²) of the target [dB]. |
Notes
The class implements the code in Section A.1 of Jech et al. (2015).
References
Jech, J.M., Horne, J.K., Chu, D., Demer, D.A., Francis, D.T.I., Gorska, N., Jones, B., Lavery, A.C., Stanton, T.K., Macaulay, G.J., Reeder, D.B., Sawada, K., 2015. Comparisons among ten models of acoustic backscattering used in aquatic ecosystem research. Journal of the Acoustical Society of America 138, 3742–3764. https://doi.org/10.1121/1.4937607
Source code in src/echosms/mssmodel.py
50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 |
|
validate_parameters(params)
Validate the model parameters.
See here for calling details.
Source code in src/echosms/mssmodel.py
PSMSModel
Bases: ScatterModelBase
Prolate spheroidal modal series (PSMS) scattering model.
Note
The fluid filled boundary type implementation is currently only accurate for weakly scattering interiors. Support for strongly scattering (e.g., gas-filled) will come later.
Source code in src/echosms/psmsmodel.py
calculate_ts(data, expand=False, inplace=False, multiprocess=False)
Calculate the target strength (TS) for many parameters.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
Pandas DataFrame, Xarray DataArray or dict
|
Requirements for the different input data types are:
|
required |
multiprocess
|
bool
|
Split the ts calculation across CPU cores. Multiprocessing is currently provided by
mapply with little customisation. For more
sophisticated uses it may be preferred to use a multiprocessing package of your choice
directly on the |
False
|
expand
|
bool
|
Only applicable if |
False
|
inplace
|
bool
|
Only applicable if |
False
|
Returns:
Type | Description |
---|---|
None, list[float], Series, or DataFrame
|
The return type and value are determined by the type of the input variable (
|
Source code in src/echosms/scattermodelbase.py
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 |
|
calculate_ts_single(medium_c, medium_rho, a, b, theta, f, boundary_type, target_c=None, target_rho=None, validate_parameters=True)
Prolate spheroid modal series (PSMS) solution model.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
medium_c
|
float
|
Sound speed in the fluid medium surrounding the target [m/s]. |
required |
medium_rho
|
float
|
Density of the fluid medium surrounding the target [kg/m³]. |
required |
a
|
float
|
Prolate spheroid major axis radius [m]. |
required |
b
|
float
|
Prolate spheroid minor axis radius [m]. |
required |
theta
|
float
|
Pitch angle to calculate the scattering as per the echoSMs coordinate system [°]. |
required |
f
|
float
|
Frequency to calculate the scattering at [Hz]. |
required |
boundary_type
|
str
|
The model type. Supported model types are given in the |
required |
target_c
|
float
|
Sound speed in the fluid inside the target [m/s].
Only required for |
None
|
target_rho
|
float
|
Density of the fluid inside the target [kg/m³].
Only required for |
None
|
validate_parameters
|
bool
|
Whether to validate the input parameters. |
True
|
Returns:
Type | Description |
---|---|
float
|
The target strength (re 1 m²) of the target [dB]. |
Notes
The backscattered target strength of a pressure release or fluid-filled prolate spheroid is calculated using the PSMS method of Furusawa (1988) and corrections in Furusawa et al. (1994).
References
Furusawa, M. (1988). "Prolate spheroidal models for predicting general trends of fish target strength," J. Acoust. Soc. Jpn. 9, 13-24. Furusawa, M., Miyanohana, Y., Ariji, M., and Sawada, Y. (1994). “Prediction of krill target strength by liquid prolate spheroid model,” Fish. Sci., 60, 261-265.
Source code in src/echosms/psmsmodel.py
43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 |
|
validate_parameters(params)
Validate the model parameters.
See here for calling details.
Source code in src/echosms/psmsmodel.py
ReferenceModels
Provide access to reference scattering model parameters.
Reference models are the models and parameters defined in Jech et al. (2015). The parameters are stored in a TOML-formatted file in the echoSMs repository and this class provides easy access to the data in that file. Additional reference models may be defined in the future and added to the TOML file (for example, entries have been added for all known calibration sphere sizes).
Attributes:
Name | Type | Description |
---|---|---|
definitions |
dict
|
A dict representation of the |
Raises:
Type | Description |
---|---|
TOMLDecodeError
|
If the |
KeyError
|
If the |
References
Jech, J.M., Horne, J.K., Chu, D., Demer, D.A., Francis, D.T.I., Gorska, N., Jones, B., Lavery, A.C., Stanton, T.K., Macaulay, G.J., Reeder, D.B., Sawada, K., 2015. Comparisons among ten models of acoustic backscattering used in aquatic ecosystem research. Journal of the Acoustical Society of America 138, 3742–3764. https://doi.org/10.1121/1.4937607
Source code in src/echosms/referencemodels.py
names()
Names of all model definitions.
Returns:
Type | Description |
---|---|
iterable of str
|
All model names in the |
parameters(name)
Model parameters for a particular model.
Model parameters are a subset of the model specification where the metadata items have been removed.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
name
|
str
|
The name of a model in the |
required |
Returns:
Type | Description |
---|---|
dict
|
The model parameters for the requested model or an empty set if no model with that name. |
Source code in src/echosms/referencemodels.py
specification(name)
Model definitions for a particular model.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
name
|
str
|
The name of a model in the |
required |
Returns:
Type | Description |
---|---|
dict
|
The model definitions for the requested model or an empty set if no model with that name. |
Source code in src/echosms/referencemodels.py
BenchmarkData
Convenient interface to the benchmark dataset.
This dataset contains the TS results from Jech et al. (2015).
Attributes:
Name | Type | Description |
---|---|---|
angle_dataset |
Pandas DataFrame
|
The angle dataset from the benchmark model runs. |
freq_dataset |
Pandas DataFrame
|
The frequency dataset from the benchmark model runs. |
Notes
The column names in the source benchmark files have been changed to be the same as those used in the ReferenceModels model definitions.
References
Jech, J.M., Horne, J.K., Chu, D., Demer, D.A., Francis, D.T.I., Gorska, N., Jones, B., Lavery, A.C., Stanton, T.K., Macaulay, G.J., Reeder, D.B., Sawada, K., 2015. Comparisons among ten models of acoustic backscattering used in aquatic ecosystem research. Journal of the Acoustical Society of America 138, 3742-3764. https://doi.org/10.1121/1.4937607
Source code in src/echosms/benchmarkdata.py
Utilities
Miscellaneous utility functions.
Neumann(m)
Neumann number.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
m
|
int
|
The input integer. |
required |
Returns:
Type | Description |
---|---|
int
|
The Neumann number. |
as_dataarray(params, no_expand=[])
Convert model parameters from dict form to a Xarray DataArray.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
params
|
dict
|
The model parameters. |
required |
no_expand
|
list
|
Key values of the non-expandable model parameters in |
[]
|
Returns:
Type | Description |
---|---|
DataArray
|
Returns a multi-dimensional DataArray generated from the Cartesian product of all
expandable items in the input dict. Non-expandable items are added to the DataArray
attrs property. Expandable items are those that can be sensibly expanded into
DataArray coordinates. Not all models have non-expandable items.
The array is named |
Source code in src/echosms/utils.py
as_dataframe(params, no_expand=[])
Convert model parameters from dict form to a Pandas DataFrame.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
params
|
dict
|
The model parameters. |
required |
no_expand
|
list
|
Key values of the non-expandable model parameters in |
[]
|
Returns:
Type | Description |
---|---|
DataFrame
|
Returns a Pandas DataFrame generated from the Cartesian product of all expandable items in the input dict. DataFrame column names are obtained from the dict keys. Non-expandable items are added to the DataFrame attrs property. Expandable items are those that can be sensibly expanded into DataFrame columns. Not all models have non-expandable items. |
Source code in src/echosms/utils.py
as_dict(params)
Convert model parameters from DataFrame or DataArray to dict.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
params
|
dict | DataFrame | DataArray
|
The model parameters |
required |
Raises:
Type | Description |
---|---|
TypeError:
|
If the input data type is not supported. |
Returns:
Type | Description |
---|---|
dict
|
A dict containing the model parameters. |
Source code in src/echosms/utils.py
h1(n, z, derivative=False)
Spherical Hankel function of the first kind or its' derivative.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
n
|
int
|
Order (n ≥ 0). |
required |
z
|
float
|
Argument of the Hankel function. |
required |
derivative
|
if True, the value of the derivative (rather than the function itself) is returned. |
False
|
Returns:
Type | Description |
---|---|
complex
|
Value of the spherical Hankel function |
Raises:
Type | Description |
---|---|
ValueError
|
For negative n values. |
Notes
The value of the Hankel function is calculated from spherical Bessel functions [1].
The derivative is computed from spherical Hankel functions [2].
Source code in src/echosms/utils.py
pro_ang1(m, n, c, eta, norm=False)
Prolate spheroidal angular function of the first kind and derivative.
Calculates the prolate spheroidal angular function of the first kind and its'
derivative with respect to eta
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
m
|
int
|
The order parameter (≥ 0) |
required |
n
|
int
|
The degree parameter (≥ |
required |
c
|
float
|
The size parameter. |
required |
eta
|
float
|
The angular coordinate, η, where |η| ≤ 1. |
required |
norm
|
If |
False
|
Returns:
Type | Description |
---|---|
tuple[float, float]
|
The value of the prolate spheroidal angular function and its' derivative. |
Notes
This method uses the prolate spheroidal wave function code for non complex
arguments (van Buren & Boisvert, 2002, and van Buren & Boisvert, 2024), available on
github. This code is in Fortran90
and was interfaced to Python using numpy.f2py
then wrapped with the current method to
provide a similar calling convention as the Scipy function of the same name.
References
Van Buren, A. L., & Boisvert, J. E. (2002). Accurate calculation of prolate spheroidal radial functions of the first kind and their first derivatives. Quarterly of Applied Mathematics, 60(3), 589-599. https://doi.org/10.1090/qam/1914443
Van Buren, A. L., & Boisvert, J. E. (2004). Improved Calculation of Prolate Spheroidal Radial Functions of the Second Kind and Their First Derivatives. Quarterly of Applied Mathematics, 62(3), 493-507. https://doi.org/10.1090/qam/2086042
Source code in src/echosms/utils.py
pro_rad1(m, n, c, xi)
Prolate spheroidal radial function of the first kind and derivative.
Calculates the prolate spheroidal radial function of the first kind and its'
derivative with respect to xi
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
m
|
int
|
The order parameter (≥ 0). |
required |
n
|
int
|
The degree parameter (≥ |
required |
c
|
float
|
The size parameter. |
required |
xi
|
float
|
The radial coordinate, ξ, where ξ ≥ 1. |
required |
Returns:
Type | Description |
---|---|
tuple[float, float]
|
The value of the prolate spheroidal radial function and its' derivative. |
Notes
This method uses the prolate spheroidal wave function code for non complex
arguments (van Buren & Boisvert, 2002, and van Buren & Boisvert, 2024), available on
github. This code is in Fortran90
and was interfaced to Python using numpy.f2py
then wrapped with the current method to
provide a similar calling convention as the Scipy function of the same name.
References
Van Buren, A. L., & Boisvert, J. E. (2002). Accurate calculation of prolate spheroidal radial functions of the first kind and their first derivatives. Quarterly of Applied Mathematics, 60(3), 589-599. https://doi.org/10.1090/qam/1914443
Van Buren, A. L., & Boisvert, J. E. (2004). Improved Calculation of Prolate Spheroidal Radial Functions of the Second Kind and Their First Derivatives. Quarterly of Applied Mathematics, 62(3), 493-507. https://doi.org/10.1090/qam/2086042
Source code in src/echosms/utils.py
pro_rad2(m, n, c, xi)
Prolate spheroidal radial function of the second kind and derivative.
Calculates the prolate spheroidal radial function of the second kind and its'
derivative with respect to xi
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
m
|
int
|
The order parameter (≥ 0). |
required |
n
|
int
|
The degree parameter (≥ |
required |
c
|
float
|
The size parameter. |
required |
xi
|
float
|
The radial coordinate, ξ, where ξ ≥ 1. |
required |
Returns:
Type | Description |
---|---|
tuple[float, float]
|
The value of the prolate spheroidal radial function and its' derivative. |
Notes
This method uses the prolate spheroidal wave function code for non complex
arguments (van Buren & Boisvert, 2002, and van Buren & Boisvert, 2024), available on
github. This code is in Fortran90
and was interfaced to Python using numpy.f2py
then wrapped with the current method to
provide a similar calling convention as the Scipy function of the same name.
References
Van Buren, A. L., & Boisvert, J. E. (2002). Accurate calculation of prolate spheroidal radial functions of the first kind and their first derivatives. Quarterly of Applied Mathematics, 60(3), 589-599. https://doi.org/10.1090/qam/1914443
Van Buren, A. L., & Boisvert, J. E. (2004). Improved Calculation of Prolate Spheroidal Radial Functions of the Second Kind and Their First Derivatives. Quarterly of Applied Mathematics, 62(3), 493-507. https://doi.org/10.1090/qam/2086042
Source code in src/echosms/utils.py
spherical_jnpp(n, z)
Second derivative of the spherical Bessel function.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
n
|
int
|
Order (n ≥ 0) |
required |
z
|
float
|
Argument of the Bessel function. |
required |
Returns:
Type | Description |
---|---|
float
|
The second derivative of the spherical Bessel function. |
Source code in src/echosms/utils.py
split_dict(d, s)
Split a dict into two dicts based on a list of keys.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
d
|
dict
|
Dict to be split. |
required |
s
|
list
|
List of dict keys to use for splitting |
required |
Returns:
Type | Description |
---|---|
tuple(dict, dict)
|
The |
Source code in src/echosms/utils.py
wavelength(c, f)
Calculate the acoustic wavelength.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
c
|
float
|
Sound speed [m/s] |
required |
f
|
float
|
Frequency [Hz] |
required |
Returns:
Type | Description |
---|---|
float
|
The acoustic wavelength [m]. |
wavenumber(c, f)
Calculate the acoustic wavenumber.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
c
|
float
|
Sound speed [m/s] |
required |
f
|
float
|
Frequency [Hz] |
required |
Returns:
Type | Description |
---|---|
float
|
The acoustic wavenumber [m⁻¹]. |