BELLHOP Python API Documentation¶
This is the Python API documentation for BELLHOP, an underwater acoustic propagation modeling toolbox.
API Reference¶
Core Functions¶
Underwater acoustic propagation modeling toolbox.
This toolbox uses the Bellhop acoustic propagation model. For this model to work, the conplete bellhop-hub package must be built and installed and bellhop.exe should be in your PATH.
- bellhop.bellhop.read_env2d(fname)[source]¶
Read a 2D underwater environment from a BELLHOP .env file.
This function parses a BELLHOP .env file and returns a Python data structure that is compatible with create_env2d(). This enables round-trip testing and compatibility between file-based and programmatic environment definitions.
- Parameters:
fname – path to .env file (with or without .env extension)
- Returns:
environment dictionary compatible with create_env2d()
The environment dictionary used in this code contains a large number of parameters, documented here to keep the code later more concise:
ENV parameters¶
- name:
environment title/name
- type:
‘2D’ (fixed for 2D environments)
- frequency:
acoustic frequency in Hz
- soundspeed:
sound speed profile (scalar for constant, array for depth-dependent)
- soundspeed_interp:
interpolation method (‘linear’, ‘spline’, ‘quadrilateral’)
- bottom_soundspeed:
bottom sediment sound speed in m/s
- bottom_soundspeed_shear:
bottom sediment sound speed in m/s
- bottom_density:
bottom sediment density in kg/m³
- bottom_absorption:
bottom sediment absorption in dB/wavelength
- bottom_absorption_shear:
bottom sediment absorption in dB/wavelength
- bottom_roughness:
bottom roughness RMS in meters
- surface:
surface altimetry profile (None if flat surface)
- surface_interp:
surface interpolation method (‘linear’, ‘curvilinear’)
- top_boundary_condition:
(‘vacuum’, ‘acousto-elastic’, ‘rigid’, ‘from-file’)
- volume_attenuation:
(‘none’, ‘thorp’, ‘francois-garrison’, ‘biological’)
- attenuation_units:
(‘nepers per meter’, ‘frequency dependent’, ‘dB per meter’, ‘frequency scaled dB per meter’, ‘dB per wavelength’, ‘quality factor’, ‘loss parameter’)
- tx_depth:
transmitter depth(s) in meters
- tx_directionality:
transmitter beam pattern (None if omnidirectional)
- rx_depth:
receiver depth(s) in meters
- rx_range:
receiver range(s) in meters
- depth:
maximum water depth in meters
- depth_interp:
bathymetry interpolation method (‘linear’, ‘curvilinear’)
- min_angle:
minimum beam angle in degrees
- max_angle:
maximum beam angle in degrees
- nbeams:
number of beams (0 for automatic)
- step_size:
(maximum) step size to trace rays in meters (0 for automatic)
- box_depth:
box extent to trace rays in meters (auto-calculated based on max depth data if not specified)
- box_range:
box extent to trace rays in meters (auto-calculated based on max receiver range if not specified)
- tx_type:
point (default) or line
- beam_type:
todo
- grid:
rectilinear or irregular
Supported ENV file formats:
Standard BELLHOP format with various boundary conditions
Constant or depth-dependent sound speed profiles
Compressed vector notation (e.g., “0.0 5000.0 /” for linearly spaced values)
Comments (lines with ! are handled correctly)
Different top/bottom boundary options (halfspace, file-based, etc.)
Unit conversions performed:
Receiver ranges: km → m
Bottom density: g/cm³ → kg/m³
All other units preserved as in ENV file
Examples:
>>> import bellhop as bh >>> env = bh.read_env2d('examples/Munk/MunkB_ray.env') >>> print(env['name']) 'Munk profile' >>> print(env['frequency']) 50.0
>>> # Use with existing functions >>> checked_env = bh.check_env2d(env) >>> rays = bh.compute_rays(env)
>>> # Round-trip compatibility >>> env_orig = bh.create_env2d(name="test", frequency=100) >>> # ... write to file via BELLHOP ... >>> env_read = bh.read_env2d("test.env") >>> assert env_read['frequency'] == env_orig['frequency']
Limitations:
External files (.ssp, .bty, .ati, .sbp) are noted but not automatically loaded
Some advanced BELLHOP features may not be fully supported
Assumes standard 2D BELLHOP format (not BELLHOP3D)
- bellhop.bellhop.read_ssp(fname)[source]¶
Read a 2D sound speed profile (.ssp) file used by BELLHOP.
This function reads BELLHOP’s .ssp files which contain range-dependent sound speed profiles. The file format is: - Line 1: Number of range profiles (NPROFILES) - Line 2: Range coordinates in km (space-separated) - Line 3+: Sound speed values, one line per depth point across all ranges
- Parameters:
fname – path to .ssp file (with or without .ssp extension)
- Returns:
for single-profile files: numpy array with [depth, soundspeed] pairs; for multi-profile files: pandas DataFrame with range-dependent sound speed data
Return format:
Single-profile files (1 range): Returns a 2D numpy array with [depth, soundspeed] pairs, compatible with create_env2d() soundspeed parameter.
Multi-profile files (>1 ranges): Returns a pandas DataFrame where:
Columns: Range coordinates (in meters, converted from km in file)
Index: Depth indices (0, 1, 2, … for each depth level in the file)
Values: Sound speeds (m/s)
This DataFrame can be directly assigned to create_env2d() soundspeed parameter for range-dependent acoustic modeling.
Note on depths: For multi-profile files, depth indices are used (0, 1, 2, …) since the actual depth coordinates come from the associated BELLHOP .env file. Users can modify the DataFrame index if actual depth values are known.
Examples:
>>> import bellhop as bh >>> # Single-profile file >>> ssp1 = bh.read_ssp("single_profile.ssp") # Returns numpy array >>> env = bh.create_env2d() >>> env["soundspeed"] = ssp1 >>> >>> # Multi-profile file >>> ssp2 = bh.read_ssp("tests/MunkB_geo_rot/MunkB_geo_rot.ssp") # Returns DataFrame >>> env = bh.create_env2d() >>> env["soundspeed"] = ssp2 # Range-dependent sound speed
File format example:
30 -50 -5 -1 -.8 -.75 -.6 -.4 -.2 0 0.2 0.4 0.6 0.8 1.0 1.2 1.4 1.6 1.8 2.0 2.2 2.4 2.6 2.8 3.0 3.2 3.4 3.6 3.8 4.0 10.0 1500 1500 1548.52 1530.29 1526.69 1517.78 1509.49 1504.30 1501.38 1500.14 1500.12 1501.02 1502.57 1504.62 1507.02 1509.69 1512.55 1515.56 1518.67 1521.85 1525.10 1528.38 1531.70 1535.04 1538.39 1541.76 1545.14 1548.52 1551.91 1551.91 1500 1500 1548.52 1530.29 1526.69 1517.78 1509.49 1504.30 1501.38 1500.14 1500.12 1501.02 1502.57 1504.62 1507.02 1509.69 1512.55 1515.56 1518.67 1521.85 1525.10 1528.38 1531.70 1535.04 1538.39 1541.76 1545.14 1548.52 1551.91 1551.91
- bellhop.bellhop.read_bty(fname)[source]¶
Read a bathymetry (.bty) file used by BELLHOP.
This function reads BELLHOP’s .bty files which define the bottom depth profile. The file format is: - Line 1: Interpolation type (‘L’ for linear, ‘C’ for curvilinear) - Line 2: Number of points - Line 3+: Range (km) and depth (m) pairs
- Parameters:
fname – path to .bty file (with or without .bty extension)
- Returns:
numpy array with [range, depth] pairs compatible with create_env2d()
The returned array can be assigned to env[“depth”] for range-dependent bathymetry.
Examples:
>>> import bellhop as bh >>> bty,bty_interp = bh.read_bty("tests/MunkB_geo_rot/MunkB_geo_rot.bty") >>> env = bh.create_env2d() >>> env["depth"] = bty >>> env["depth_interp"] = bty_interp >>> arrivals = bh.calculate_arrivals(env)
File format example:
'L' 5 0 3000 10 3000 20 500 30 3000 100 3000
- bellhop.bellhop.read_refl_coeff(fname)[source]¶
Read a reflection coefficient (.brc) file used by BELLHOP.
This function reads BELLHOP’s .brc files which define the reflection coefficient data. The file format is: - Line 1: Number of points - Line 2+: THETA(j) RMAG(j) RPHASE(j)
Where: - THETA(): Angle (degrees) - RMAG(): Magnitude of reflection coefficient - RPHASE(): Phase of reflection coefficient (degrees)
- Parameters:
fname – path to .brc/.trc file (extension required)
- Returns:
numpy array with [theta, rmag, rphase] triplets compatible with create_env2d()
The returned array can be assigned to env[“bottom_reflection_coefficient”] or env[“top_reflection_coefficient”] .
Example:
>>> import bellhop as bh >>> brc = bh.read_refl_coeff("tests/MunkB_geo_rot/MunkB_geo_rot.brc") >>> env = bh.create_env2d() >>> env["bottom_reflection_coefficient"] = brc >>> arrivals = bh.calculate_arrivals(env)
File format example:
3 0.0 1.00 180.0 45.0 0.95 175.0 90.0 0.90 170.0
- bellhop.bellhop.create_env2d(**kv)[source]¶
Create a new 2D underwater environment.
- Parameters:
**kv (dict) – Keyword arguments for environment configuration.
- Returns:
env – A new 2D underwater environment dictionary.
- Return type:
dict
Example
To see all the parameters available and their default values:
>>> import bellhop as bh >>> env = bh.create_env2d() >>> bh.print_env(env)
The environment parameters may be changed by passing keyword arguments or modified later using dictionary notation:
>>> import bellhop as bh >>> env = bh.create_env2d(depth=40, soundspeed=1540) >>> bh.print_env(env) >>> env['depth'] = 25 >>> env['bottom_soundspeed'] = 1800 >>> bh.print_env(env)
The default environment has a constant sound speed. A depth dependent sound speed profile be provided as a Nx2 array of (depth, sound speed):
>>> import bellhop as bh >>> env = bh.create_env2d(depth=20, >>>. soundspeed=[[0,1540], [5,1535], [10,1535], [20,1530]])
A range-and-depth dependent sound speed profile can be provided as a Pandas frame:
>>> import bellhop as bh >>> import pandas as pd >>> ssp2 = pd.DataFrame({ 0: [1540, 1530, 1532, 1533], # profile at 0 m range 100: [1540, 1535, 1530, 1533], # profile at 100 m range 200: [1530, 1520, 1522, 1525] }, # profile at 200 m range index=[0, 10, 20, 30]) # depths of the profile entries in m >>> env = bh.create_env2d(depth=20, soundspeed=ssp2)
The default environment has a constant water depth. A range dependent bathymetry can be provided as a Nx2 array of (range, water depth):
>>> import bellhop as bh >>> env = bh.create_env2d(depth=[[0,20], [300,10], [500,18], [1000,15]])
- bellhop.bellhop.check_env2d(env)[source]¶
Check the validity of a 2D underwater environment definition.
- Parameters:
env – environment definition
Exceptions are thrown with appropriate error messages if the environment is invalid.
>>> import bellhop as bh >>> env = bh.create_env2d() >>> check_env2d(env)
- bellhop.bellhop.print_env(env)[source]¶
Display the environment in a human readable form.
- Parameters:
env – environment definition
>>> import bellhop as bh >>> env = bh.create_env2d(depth=40, soundspeed=1540) >>> bh.print_env(env)
- bellhop.bellhop.plot_env(env, surface_color='dodgerblue', bottom_color='peru', tx_color='orangered', rx_color='midnightblue', rx_plot=None, **kwargs)[source]¶
Plots a visual representation of the environment.
- Parameters:
env – environment description
surface_color – color of the surface (see Bokeh colors)
bottom_color –
color of the bottom (see Bokeh colors)
tx_color –
color of transmitters (see Bokeh colors)
rx_color –
color of receviers (see Bokeh colors)
rx_plot – True to plot all receivers, False to not plot any receivers, None to automatically decide
Other keyword arguments applicable for bellhop.plot.plot() are also supported.
The surface, bottom, transmitters (marker: ‘*’) and receivers (marker: ‘o’) are plotted in the environment. If rx_plot is set to None and there are more than 2000 receivers, they are not plotted.
>>> import bellhop as bh >>> env = bh.create_env2d(depth=[[0, 40], [100, 30], [500, 35], [700, 20], [1000,45]]) >>> bh.plot_env(env)
- bellhop.bellhop.plot_ssp(env, **kwargs)[source]¶
Plots the sound speed profile.
- Parameters:
env – environment description
Other keyword arguments applicable for bellhop.plot.plot() are also supported.
If the sound speed profile is range-dependent, this function only plots the first profile.
>>> import bellhop as bh >>> env = bh.create_env2d(soundspeed=[[ 0, 1540], [10, 1530], [20, 1532], [25, 1533], [30, 1535]]) >>> bh.plot_ssp(env)
- bellhop.bellhop.compute_arrivals(env, model=None, debug=False, fname_base=None)[source]¶
Compute arrivals between each transmitter and receiver.
- Parameters:
env – environment definition
model – propagation model to use (None to auto-select)
debug – generate debug information for propagation model
fname_base – base file name for Bellhop working files, default (None), creates a temporary file
- Returns:
arrival times and coefficients for all transmitter-receiver combinations
>>> import bellhop as bh >>> env = bh.create_env2d() >>> arrivals = bh.compute_arrivals(env) >>> bh.plot_arrivals(arrivals)
- bellhop.bellhop.compute_eigenrays(env, tx_depth_ndx=0, rx_depth_ndx=0, rx_range_ndx=0, model=None, debug=False, fname_base=None)[source]¶
Compute eigenrays between a given transmitter and receiver.
- Parameters:
env – environment definition
tx_depth_ndx – transmitter depth index
rx_depth_ndx – receiver depth index
rx_range_ndx – receiver range index
model – propagation model to use (None to auto-select)
debug – generate debug information for propagation model
fname_base – base file name for Bellhop working files, default (None), creates a temporary file
- Returns:
eigenrays paths
>>> import bellhop as bh >>> env = bh.create_env2d() >>> rays = bh.compute_eigenrays(env) >>> bh.plot_rays(rays, width=1000)
- bellhop.bellhop.compute_rays(env, tx_depth_ndx=0, model=None, debug=False, fname_base=None)[source]¶
Compute rays from a given transmitter.
- Parameters:
env – environment definition
tx_depth_ndx – transmitter depth index
model – propagation model to use (None to auto-select)
debug – generate debug information for propagation model
fname_base – base file name for Bellhop working files, default (None), creates a temporary file
- Returns:
ray paths
>>> import bellhop as bh >>> env = bh.create_env2d() >>> rays = bh.compute_rays(env) >>> bh.plot_rays(rays, width=1000)
- bellhop.bellhop.compute_transmission_loss(env, tx_depth_ndx=0, mode=_Strings.coherent, model=None, tx_type='default', debug=False, fname_base=None)[source]¶
Compute transmission loss from a given transmitter to all receviers.
- Parameters:
env – environment definition
tx_depth_ndx – transmitter depth index
mode – coherent, incoherent or semicoherent
model – propagation model to use (None to auto-select)
tx_type – point or line
debug – generate debug information for propagation model
fname_base – base file name for Bellhop working files, default (None), creates a temporary file
- Returns:
complex transmission loss at each receiver depth and range
>>> import bellhop as bh >>> env = bh.create_env2d() >>> tloss = bh.compute_transmission_loss(env, mode=bh.incoherent) >>> bh.plot_transmission_loss(tloss, width=1000)
- bellhop.bellhop.arrivals_to_impulse_response(arrivals, fs, abs_time=False)[source]¶
Convert arrival times and coefficients to an impulse response.
- Parameters:
arrivals – arrivals times (s) and coefficients
fs – sampling rate (Hz)
abs_time – absolute time (True) or relative time (False)
- Returns:
impulse response
If abs_time is set to True, the impulse response is placed such that the zero time corresponds to the time of transmission of signal.
>>> import bellhop as bh >>> env = bh.create_env2d() >>> arrivals = bh.compute_arrivals(env) >>> ir = bh.arrivals_to_impulse_response(arrivals, fs=192000)
- bellhop.bellhop.plot_arrivals(arrivals, dB=False, color='blue', **kwargs)[source]¶
Plots the arrival times and amplitudes.
- Parameters:
arrivals – arrivals times (s) and coefficients
dB – True to plot in dB, False for linear scale
color –
line color (see Bokeh colors)
Other keyword arguments applicable for bellhop.plot.plot() are also supported.
>>> import bellhop as bh >>> env = bh.create_env2d() >>> arrivals = bh.compute_arrivals(env) >>> bh.plot_arrivals(arrivals)
- bellhop.bellhop.plot_rays(rays, env=None, invert_colors=False, **kwargs)[source]¶
Plots ray paths.
- Parameters:
rays – ray paths
env – environment definition
invert_colors – False to use black for high intensity rays, True to use white
If environment definition is provided, it is overlayed over this plot using default parameters for bellhop.plot_env().
Other keyword arguments applicable for bellhop.plot.plot() are also supported.
>>> import bellhop as bh >>> env = bh.create_env2d() >>> rays = bh.compute_eigenrays(env) >>> bh.plot_rays(rays, width=1000)
- bellhop.bellhop.plot_transmission_loss(tloss, env=None, **kwargs)[source]¶
Plots transmission loss.
- Parameters:
tloss – complex transmission loss
env – environment definition
If environment definition is provided, it is overlayed over this plot using default parameters for bellhop.plot_env().
Other keyword arguments applicable for bellhop.plot.image() are also supported.
>>> import bellhop as bh >>> import numpy as np >>> env = bh.create_env2d( rx_depth=np.arange(0, 25), rx_range=np.arange(0, 1000), min_angle=-45, max_angle=45 ) >>> tloss = bh.compute_transmission_loss(env) >>> bh.plot_transmission_loss(tloss, width=1000)
- bellhop.bellhop.pyplot_env(env, surface_color='dodgerblue', bottom_color='peru', tx_color='orangered', rx_color='midnightblue', rx_plot=None, **kwargs)[source]¶
Plots a visual representation of the environment with matplotlib.
- Parameters:
env – environment description
surface_color –
color of the surface (see Bokeh colors)
bottom_color –
color of the bottom (see Bokeh colors)
tx_color –
color of transmitters (see Bokeh colors)
rx_color –
color of receviers (see Bokeh colors)
rx_plot – True to plot all receivers, False to not plot any receivers, None to automatically decide
Other keyword arguments applicable for bellhop.plot.plot() are also supported.
The surface, bottom, transmitters (marker: ‘*’) and receivers (marker: ‘o’) are plotted in the environment. If rx_plot is set to None and there are more than 2000 receivers, they are not plotted.
>>> import bellhop as bh >>> env = bh.create_env2d(depth=[[0, 40], [100, 30], [500, 35], [700, 20], [1000,45]]) >>> bh.plot_env(env)
- bellhop.bellhop.pyplot_ssp(env, **kwargs)[source]¶
Plots the sound speed profile with matplotlib.
- Parameters:
env – environment description
Other keyword arguments applicable for bellhop.plot.plot() are also supported.
If the sound speed profile is range-dependent, this function only plots the first profile.
>>> import bellhop as bh >>> env = bh.create_env2d(soundspeed=[[ 0, 1540], [10, 1530], [20, 1532], [25, 1533], [30, 1535]]) >>> bh.plot_ssp(env)
- bellhop.bellhop.pyplot_arrivals(arrivals, dB=False, color='blue', **kwargs)[source]¶
Plots the arrival times and amplitudes with matplotlib.
- Parameters:
arrivals – arrivals times (s) and coefficients
dB – True to plot in dB, False for linear scale
color –
line color (see Bokeh colors)
Other keyword arguments applicable for bellhop.plot.plot() are also supported.
>>> import bellhop as bh >>> env = bh.create_env2d() >>> arrivals = bh.compute_arrivals(env) >>> bh.plot_arrivals(arrivals)
- bellhop.bellhop.pyplot_rays(rays, env=None, invert_colors=False, **kwargs)[source]¶
Plots ray paths with matplotlib
- Parameters:
rays – ray paths
env – environment definition
invert_colors – False to use black for high intensity rays, True to use white
If environment definition is provided, it is overlayed over this plot using default parameters for bellhop.plot_env().
Other keyword arguments applicable for bellhop.plot.plot() are also supported.
>>> import bellhop as bh >>> env = bh.create_env2d() >>> rays = bh.compute_eigenrays(env) >>> bh.plot_rays(rays, width=1000)
- bellhop.bellhop.pyplot_transmission_loss(tloss, env=None, **kwargs)[source]¶
Plots transmission loss with matplotlib.
- Parameters:
tloss – complex transmission loss
env – environment definition
If environment definition is provided, it is overlayed over this plot using default parameters for bellhop.plot_env().
Other keyword arguments applicable for bellhop.plot.image() are also supported.
>>> import bellhop as bh >>> import numpy as np >>> env = bh.create_env2d( rx_depth=np.arange(0, 25), rx_range=np.arange(0, 1000), min_angle=-45, max_angle=45 ) >>> tloss = bh.compute_transmission_loss(env) >>> bh.plot_transmission_loss(tloss, width=1000)
- bellhop.bellhop.models(env=None, task=None)[source]¶
List available models.
- Parameters:
env – environment to model
task – arrivals/eigenrays/rays/coherent/incoherent/semicoherent
- Returns:
list of models that can be used
>>> import bellhop as bh >>> bh.models() ['bellhop'] >>> env = bh.create_env2d() >>> bh.models(env, task=coherent) ['bellhop']
Plotting Utilities¶
Internal Constants (bellhop/constants.py)¶
- class bellhop.constants._Strings(*values)[source]¶
String definitions to avoid hard-coding magic strings in the source code
This helps prevent typos and permits autocomplete (if your editor is smart enough).
- linear = 'linear'¶
- spline = 'spline'¶
- curvilinear = 'curvilinear'¶
- quadrilateral = 'quadrilateral'¶
- pchip = 'pchip'¶
- hexahedral = 'hexahedral'¶
- arrivals = 'arrivals'¶
- nlinear = 'nlinear'¶
- eigenrays = 'eigenrays'¶
- rays = 'rays'¶
- coherent = 'coherent'¶
- incoherent = 'incoherent'¶
- semicoherent = 'semicoherent'¶
- vacuum = 'vacuum'¶
- acousto_elastic = 'acousto-elastic'¶
- rigid = 'rigid'¶
- from_file = 'from-file'¶
- flat = 'flat'¶
- line = 'line'¶
- point = 'point'¶
- rectilinear = 'rectilinear'¶
- irregular = 'irregular'¶
- thorp = 'thorp'¶
- francois_garrison = 'francois-garrison'¶
- biological = 'biological'¶
- class bellhop.constants._Maps[source]¶
Mappings from Bellhop single-char input file options to readable Python options
These are also defined with reverse mappings in the form:
>>> _Maps.interp["S"] "spline"
>>> _Maps.interp_rev["spline"] "S"
- interp = {' ': 'default', 'C': _Strings.linear, 'H': _Strings.hexahedral, 'N': _Strings.nlinear, 'P': _Strings.pchip, 'Q': _Strings.quadrilateral, 'S': _Strings.spline}¶
- bty_interp = {'C': _Strings.curvilinear, 'L': _Strings.linear}¶
- boundcond = {' ': 'default', 'A': _Strings.acousto_elastic, 'F': _Strings.from_file, 'R': _Strings.rigid, 'V': _Strings.vacuum}¶
- attunits = {' ': 'default', 'F': 'frequency dependent', 'L': 'loss parameter', 'M': 'dB per meter', 'N': 'nepers per meter', 'Q': 'quality factor', 'W': 'dB per wavelength', 'm': 'frequency scaled dB per meter'}¶
- volatt = {'': 'none', ' ': 'default', 'B': 'biological', 'F': 'francois-garrison', 'T': 'thorp'}¶
- bottom = {' ': 'default', '*': _Strings.from_file, '_': _Strings.flat, '~': _Strings.from_file}¶
- source = {' ': 'default', 'R': _Strings.point, 'X': _Strings.line}¶
- grid = {' ': 'default', 'I': _Strings.irregular, 'R': _Strings.rectilinear}¶
- beam = {' ': 'default', 'B': 'gaussian-cartesian', 'G': 'hat-cartesian', '^': 'hat-cartesian', 'b': 'gaussian-ray', 'g': 'hat-ray'}¶
- interp_rev = {'default': ' ', _Strings.hexahedral: 'H', _Strings.linear: 'C', _Strings.nlinear: 'N', _Strings.pchip: 'P', _Strings.quadrilateral: 'Q', _Strings.spline: 'S'}¶
- bty_interp_rev = {_Strings.curvilinear: 'C', _Strings.linear: 'L'}¶
- boundcond_rev = {_Strings.acousto_elastic: 'A', 'default': ' ', _Strings.from_file: 'F', _Strings.rigid: 'R', _Strings.vacuum: 'V'}¶
- attunits_rev = {'dB per meter': 'M', 'dB per wavelength': 'W', 'default': ' ', 'frequency dependent': 'F', 'frequency scaled dB per meter': 'm', 'loss parameter': 'L', 'nepers per meter': 'N', 'quality factor': 'Q'}¶
- volatt_rev = {'biological': 'B', 'default': ' ', 'francois-garrison': 'F', 'none': '', 'thorp': 'T'}¶
- bottom_rev = {'default': ' ', _Strings.flat: '_', _Strings.from_file: '*'}¶
- source_rev = {'default': ' ', _Strings.line: 'X', _Strings.point: 'R'}¶
- grid_rev = {'default': ' ', _Strings.irregular: 'I', _Strings.rectilinear: 'R'}¶
- beam_rev = {'default': ' ', 'gaussian-cartesian': 'B', 'gaussian-ray': 'b', 'hat-cartesian': '^', 'hat-ray': 'g'}¶