Coverage for python/bellhop/plot.py: 100%

127 statements  

« prev     ^ index     » next       coverage.py v7.11.0, created at 2025-10-22 12:04 +0000

1############################################################################## 

2# 

3# Copyright (c) 2025-, Will Robertson 

4# Copyright (c) 2018-2025, Mandar Chitre 

5# 

6# This file was originally part of arlpy, released under Simplified BSD License. 

7# It has been relicensed in this repository to be compatible with the Bellhop licence (GPL). 

8# 

9############################################################################## 

10 

11"""Plotting functions for the underwater acoustic propagation modeling toolbox. 

12""" 

13 

14from typing import Any, Optional 

15from sys import float_info as _fi 

16 

17import numpy as _np 

18import scipy.interpolate as _interp 

19import pandas as _pd 

20 

21import matplotlib.pyplot as _pyplt 

22import matplotlib.colors as _mplc 

23 

24from .environment import Environment 

25from .constants import _Strings 

26from .plotutils import figure as figure 

27 

28import bellhop.plotutils as _plt 

29 

30def plot_env(env: Environment, 

31 surface_color: str = 'dodgerblue', 

32 bottom_color: str = 'peru', 

33 source_color: str = 'orangered', 

34 receiver_color: str = 'midnightblue', 

35 receiver_plot: Optional[bool] = None, 

36 **kwargs: Any 

37 ) -> None: 

38 """Plots a visual representation of the environment. 

39 

40 Parameters 

41 ---------- 

42 env : dict 

43 Environment description 

44 surface_color : str, default='dodgerblue' 

45 Color of the surface (see `Bokeh colors <https://bokeh.pydata.org/en/latest/docs/reference/colors.html>`_) 

46 bottom_color : str, default='peru' 

47 Color of the bottom (see `Bokeh colors <https://bokeh.pydata.org/en/latest/docs/reference/colors.html>`_) 

48 source_color : str, default='orangered' 

49 Color of transmitters (see `Bokeh colors <https://bokeh.pydata.org/en/latest/docs/reference/colors.html>`_) 

50 receiver_color : str, default='midnightblue' 

51 Color of receivers (see `Bokeh colors <https://bokeh.pydata.org/en/latest/docs/reference/colors.html>`_) 

52 receiver_plot : bool, optional 

53 True to plot all receivers, False to not plot any receivers, None to automatically decide 

54 **kwargs 

55 Other keyword arguments applicable for `bellhop.plot.plot()` are also supported 

56 

57 Notes 

58 ----- 

59 The surface, bottom, transmitters (marker: '*') and receivers (marker: 'o') 

60 are plotted in the environment. If `receiver_plot` is set to None and there are 

61 more than 2000 receivers, they are not plotted. 

62 

63 Examples 

64 -------- 

65 >>> import bellhop as bh 

66 >>> env = bh.create_env(depth=[[0, 40], [100, 30], [500, 35], [700, 20], [1000,45]]) 

67 >>> bh.plot_env(env) 

68 """ 

69 

70 if env is not None: 

71 env.check() 

72 

73 min_x = 0.0 

74 max_x = float(_np.max(env['receiver_range'])) 

75 if max_x-min_x > 10000: 

76 divisor = 1000.0 

77 min_x /= divisor 

78 max_x /= divisor 

79 xlabel = 'Range (km)' 

80 else: 

81 divisor = 1.0 

82 xlabel = 'Range (m)' 

83 if env['surface'] is None: 

84 min_y = 0 

85 else: 

86 min_y = _np.min(env['surface'][:,1]) 

87 if _np.size(env['depth']) > 1: 

88 max_y = _np.max(env['depth'][:,1]) 

89 else: 

90 max_y = env['depth'] 

91 mgn_x = 0.01*(max_x-min_x) 

92 mgn_y = 0.1*(max_y-min_y) 

93 

94 oh = _plt.hold() 

95 if env['surface'] is None: 

96 xx = [min_x, max_x] 

97 yy = [0, 0] 

98 else: 

99 # linear and curvilinear options use the same altimetry, just with different normals 

100 s = env['surface'] 

101 xx = s[:,0]/divisor 

102 yy = -s[:,1] 

103 _plt.plot(xx, yy, xlabel=xlabel, ylabel='Depth (m)', xlim=(min_x-mgn_x, max_x+mgn_x), ylim=(-max_y-mgn_y, -min_y+mgn_y), color=surface_color, **kwargs) 

104 

105 if _np.size(env['depth']) == 1: 

106 xx = [min_x, max_x] 

107 yy = [-env['depth'], -env['depth']] 

108 else: 

109 # linear and curvilinear options use the same bathymetry, just with different normals 

110 s = env['depth'] 

111 xx = s[:,0]/divisor 

112 yy = -s[:,1] 

113 _plt.plot(xx, yy, color=bottom_color) 

114 

115 txd = env['source_depth'] 

116 _plt.plot([0]*_np.size(txd), -txd, marker='*', style='solid', color=source_color) 

117 

118 if receiver_plot is None: 

119 receiver_plot = _np.size(env['receiver_depth'])*_np.size(env['receiver_range']) < 2000 

120 if receiver_plot: 

121 rxr = env['receiver_range'] 

122 if _np.size(rxr) == 1: 

123 rxr = [rxr] 

124 for r in _np.array(rxr): 

125 rxd = env['receiver_depth'] 

126 _plt.plot([r/divisor]*_np.size(rxd), -rxd, marker='o', style='solid', color=receiver_color) 

127 

128 _plt.hold(oh if oh is not None else False) 

129 

130def plot_ssp(env: Environment, **kwargs: Any) -> None: 

131 """Plots the sound speed profile. 

132 

133 Parameters 

134 ---------- 

135 env : dict 

136 Environment description 

137 **kwargs 

138 Other keyword arguments applicable for `bellhop.plot.plot()` are also supported 

139 

140 Notes 

141 ----- 

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

143 

144 Examples 

145 -------- 

146 >>> import bellhop as bh 

147 >>> env = bh.create_env(soundspeed=[[ 0, 1540], [10, 1530], [20, 1532], [25, 1533], [30, 1535]]) 

148 >>> bh.plot_ssp(env) 

149 """ 

150 

151 if env is not None: 

152 env.check() 

153 

154 oh = _plt.hold() 

155 svp = env['soundspeed'] 

156 if isinstance(svp, _pd.DataFrame): 

157 svp = _np.hstack((_np.array([svp.index]).T, _np.asarray(svp))) 

158 if env['soundspeed_interp'] == _Strings.spline: 

159 ynew = _np.linspace(_np.min(svp[:,0]), _np.max(svp[:,0]), 100) 

160 tck = _interp.splrep(svp[:,0], svp[:,1], s=0) 

161 xnew = _interp.splev(ynew, tck, der=0) 

162 _plt.plot(xnew, -ynew, xlabel='Soundspeed (m/s)', ylabel='Depth (m)', hold=True, **kwargs) 

163 _plt.scatter(svp[:,1], -svp[:,0], **kwargs) 

164 else: 

165 for rr in range(1,svp.shape[1]): 

166 _plt.plot(svp[:,rr], -svp[:,0], xlabel='Soundspeed (m/s)', ylabel='Depth (m)', legend=f'Range {rr}', **kwargs) 

167 _plt.hold(oh if oh is not None else False) 

168 

169def plot_arrivals(arrivals: Any, dB: bool = False, color: str = 'blue', **kwargs: Any) -> None: 

170 """Plots the arrival times and amplitudes. 

171 

172 Parameters 

173 ---------- 

174 arrivals : pandas.DataFrame 

175 Arrivals times (s) and coefficients 

176 dB : bool, default=False 

177 True to plot in dB, False for linear scale 

178 color : str, default='blue' 

179 Line color (see `Bokeh colors <https://bokeh.pydata.org/en/latest/docs/reference/colors.html>`_) 

180 **kwargs 

181 Other keyword arguments applicable for `bellhop.plot.plot()` are also supported 

182 

183 Examples 

184 -------- 

185 >>> import bellhop as bh 

186 >>> env = bh.create_env() 

187 >>> arrivals = bh.compute_arrivals(env) 

188 >>> bh.plot_arrivals(arrivals) 

189 """ 

190 

191 t0 = min(arrivals.time_of_arrival) 

192 t1 = max(arrivals.time_of_arrival) 

193 oh = _plt.hold() 

194 if dB: 

195 min_y = 20*_np.log10(_np.max(_np.abs(arrivals.arrival_amplitude)))-60 

196 ylabel = 'Amplitude (dB)' 

197 else: 

198 ylabel = 'Amplitude' 

199 min_y = 0 

200 _plt.plot([t0, t1], [min_y, min_y], xlabel='Arrival time (s)', ylabel=ylabel, color=color, **kwargs) 

201 for _, row in arrivals.iterrows(): 

202 t = row.time_of_arrival.real 

203 y = _np.abs(row.arrival_amplitude) 

204 if dB: 

205 y = max(20*_np.log10(_fi.epsilon+y), min_y) 

206 _plt.plot([t, t], [min_y, y], color=color, **kwargs) 

207 _plt.hold(oh if oh is not None else False) 

208 

209def plot_rays(rays: Any, env: Optional[Environment] = None, invert_colors: bool = False, **kwargs: Any) -> None: 

210 """Plots ray paths. 

211 

212 Parameters 

213 ---------- 

214 rays : pandas.DataFrame 

215 Ray paths 

216 env : dict, optional 

217 Environment definition 

218 invert_colors : bool, default=False 

219 False to use black for high intensity rays, True to use white 

220 **kwargs 

221 Other keyword arguments applicable for `bellhop.plot.plot()` are also supported 

222 

223 Notes 

224 ----- 

225 If environment definition is provided, it is overlayed over this plot using default 

226 parameters for `bellhop.plot_env()`. 

227 

228 Examples 

229 -------- 

230 >>> import bellhop as bh 

231 >>> env = bh.create_env() 

232 >>> rays = bh.compute_eigenrays(env) 

233 >>> bh.plot_rays(rays, width=1000) 

234 """ 

235 rays = rays.sort_values('bottom_bounces', ascending=False) 

236 

237 # some edge cases to worry about here: rays.bottom_bounces could be all zeros? 

238 max_amp = _np.max(_np.abs(rays.bottom_bounces)) if len(rays.bottom_bounces) > 0 else 0.0 

239 max_amp = max_amp or 1.0 

240 

241 divisor = 1 

242 xlabel = 'Range (m)' 

243 r = [] 

244 for _, row in rays.iterrows(): 

245 r += list(row.ray[:,0]) 

246 if max(r)-min(r) > 10000: 

247 divisor = 1000 

248 xlabel = 'Range (km)' 

249 

250 oh = _plt.hold() 

251 for _, row in rays.iterrows(): 

252 c = float(_np.abs(row.bottom_bounces) / max_amp) 

253 if invert_colors: 

254 c = 1.0 - c 

255 cmap = _pyplt.get_cmap("gray") 

256 col_str = _mplc.to_hex(cmap(c)) 

257 _plt.plot(row.ray[:,0]/divisor, -row.ray[:,1], color=col_str, xlabel=xlabel, ylabel='Depth (m)', **kwargs) 

258 if env is not None: 

259 plot_env(env,title=None) 

260 _plt.hold(oh if oh is not None else False) 

261 

262def plot_transmission_loss(tloss: Any, env: Optional[Environment] = None, **kwargs: Any) -> None: 

263 """Plots transmission loss. 

264 

265 Parameters 

266 ---------- 

267 tloss : pandas.DataFrame 

268 Complex transmission loss 

269 env : dict, optional 

270 Environment definition 

271 **kwargs 

272 Other keyword arguments applicable for `bellhop.plot.image()` are also supported 

273 

274 Notes 

275 ----- 

276 If environment definition is provided, it is overlayed over this plot using default 

277 parameters for `bellhop.plot_env()`. 

278 

279 Examples 

280 -------- 

281 >>> import bellhop as bh 

282 >>> import numpy as np 

283 >>> env = bh.create_env( 

284 receiver_depth=np.arange(0, 25), 

285 receiver_range=np.arange(0, 1000), 

286 beam_angle_min=-45, 

287 beam_angle_max=45 

288 ) 

289 >>> tloss = bh.compute_transmission_loss(env) 

290 >>> bh.plot_transmission_loss(tloss, width=1000) 

291 """ 

292 

293 xr = (min(tloss.columns), max(tloss.columns)) 

294 yr = (-max(tloss.index), -min(tloss.index)) 

295 xlabel = 'Range (m)' 

296 if xr[1]-xr[0] > 10000: 

297 xr = (min(tloss.columns)/1000, max(tloss.columns)/1000) 

298 xlabel = 'Range (km)' 

299 oh = _plt.hold() 

300 _plt.image(20*_np.log10(_fi.epsilon+_np.abs(_np.flipud(_np.array(tloss)))), x=xr, y=yr, xlabel=xlabel, ylabel='Depth (m)', xlim=xr, ylim=yr, **kwargs) 

301 if env is not None: 

302 plot_env(env, receiver_plot=False, title=None) 

303 _plt.hold(oh if oh is not None else False) 

304 

305 

306### Export module names for auto-importing in __init__.py 

307 

308__all__ = [ 

309 name for name in globals() if not name.startswith("_") # ignore private names 

310]