optimus.source.transducers

Common functionality for transducer sources.

Module Contents

Classes

_Transducer

Functions

transducer_field(source, medium, field_locations[, ...])

Calculate the field emitted by a transducer source.

calc_greens_functions_in_observation_points_numba(...)

Calculate the pressure field and its gradient for a

calc_field_from_point_sources(locations_source, ...)

Calculate the pressure field and its gradient of a point source,

calc_field_from_point_sources_numpy(locations_source, ...)

Calculate the pressure field and its gradient of a point source,

calc_field_from_point_sources_mp_source_para(...)

Computes the pressure and normal pressure gradient at field locations for

calc_field_from_point_sources_mp_field_para(...)

Computes the pressure and normal pressure gradient at field locations for

chunk_size_index(locations_source, locations_observation)

Computes the chunk sizes and indices used to calculate the incident field using

break_in_chunks(number_of_locations, chunk_size)

Computes the chunk sizes and indices used to calculate the incident field using

optimus.source.transducers.transducer_field(source, medium, field_locations, normals=None, verbose=False)

Calculate the field emitted by a transducer source.

Parameters
sourceoptimus.source.Source

The type of acoustic source used.

mediumoptimus.material.Material

The propagating medium.

field_locationsnumpy.ndarray

An array of size (3,N) with the coordinates of the locations at which the incident field is evaluated.

normalsnumpy.ndarray

Array of size (3,N) with the coordinates of the unit normal vectors for evaluation of the pressure normal gradient on the surface of scatterers.

verboseboolean

Verbosity of output. Default: False

Returns
transducer_Transducer

An object with the attributes ‘pressure’ and ‘normal_pressure_gradient’.

class optimus.source.transducers._Transducer(source, medium, field_locations, normals, verbose)
generate_source_points()

Generate the source points of the transducer. The field emitted from any transducer is modelled by a collection of point sources, each with weighting for its amplitude.

Returns
source_locationsnumpy.ndarray

An array of size (3,N_sourcepoints) with the location of each point source.

source_weightsnumpy.ndarray

An array of size (N_sourcepoints,) with the weighting of each point source.

define_source_points_in_reference_piston(radius)

Define the source points for a reference piston element, that is, the source points on a rectangular grid, located in the plane z=0, centered at the origin, and inside a disk of the specified radius. The resolution of the points is determined by the specified number of point sources per wavelength. If zero points per wavelength is specified, return the center of the disk as the only source point. The surface area weighting is uniform.

Parameters
radiusfloat

The radius of piston.

Returns
locations_inside_transducernumpy.ndarray

An array of size (3,N_points) with the locations of the point source inside the reference element.

surface_area_weightingfloat

The surface area weighting associated to each point source.

define_source_points_in_reference_bowl()

Define the source points for a reference spherical section bowl transducer, that is, the source points on a spherical section bowl whose apex is in contact with the z=0 plane, at the global origin and whose axis is the Cartesian x-axis. The bowl is defined by its outer radius and radius of curvature. A circular aperture may be defined by specifying an inner radius. The resolution of the points is determined by the specified number of point sources per wavelength which must be strictly positive. The surface area weighting is uniform.

Returns
locations_inside_transducernumpy.ndarray

An array of size (3,N_points) with the locations of the point source inside the reference element.

surface_area_weightingfloat

The surface area weighting associated to each point source.

define_source_points_in_reference_array()

Define the source points for a reference spherical section array transducer, that is, the source points on a spherical section bowl whose apex is in contact with the z=0 plane, at the global origin and whose axis is the Cartesian x-axis. The array is defined by the location of the element centroids and the radius of the elements. The resolution of the points is determined by the specified number of point sources per wavelength which must be strictly positive. The surface area weighting is uniform.

Returns
locations_inside_transducernumpy.ndarray

An array of size (3,N_points) with the locations of the point source inside the reference array.

surface_area_weightingfloat

The surface area weighting associated to each point source.

transform_source_points(source_locations_on_reference_source)

Transform the source points from the reference transducer to the actual location of the transducer, as specified by the source axis and the source location.

Parameters
source_locations_on_reference_sourcenumpy.ndarray

An array of size (3,N_sourcepoints) with the locations of the source points on the reference element of the transducer type.

Returns
source_locations_transformednumpy.ndarray

An array of size (3,N_sourcepoints) with the locations of the source points on the transducer.

calc_pressure_field()

Calculate the pressure field and the normal gradient of the transducer, in a collection of 3D observation points.

Parameters
source_locationsnumpy.ndarray

An array of size (3,N_sourcepoints) with the coordinates of the locations of the point sources used to discretise the acoustic source.

source_weightsnumpy.ndarray

An array of size (N_sourcepoints,) with the weighting assigned to each point source.

field_locationsnumpy.ndarray

An array of size (3,N_observationpoints) with the coordinates of the locations at which the incident field is evaluated.

normalsnumpy.ndarray

An array of size (3,N_observationpoints) with the coordinates of the normal vectors for evaluation of the pressure normal gradient on the surface of scatterers.

Returns
pressure: numpy.ndarray

An array of size (N_observationpoints,) with the pressure in the observation points.

normal_pressure_gradient: numpy.ndarray

An array of size (3,N_observationpoints) with the normal gradient of the pressure in the observation points.

optimus.source.transducers.calc_greens_functions_in_observation_points_numba(locations_source, locations_observation, wavenumber, source_weights)

Calculate the pressure field and its gradient for a summation of point sources describing a transducer, according to the Rayleigh integral formula. Use Numba for acceleration and parallelisation.

Parameters
locations_sourcenumpy.ndarray

An array of size (3,N_sourcepoints) with the locations of the source points.

locations_observationnumpy.ndarray

An array of size (3,N_observationpoints) with the locations of the observation points.

wavenumbercomplex

The wavenumber of the wave field.

source_weightsnumpy.ndarray

An array of size (N_sourcepoints,) with the weights of each source element.

Returns
greens_function_in_observation_pointsnumpy.ndarray

An array of size (N_observationpoints,) with the Green’s function of the wave field at the observation points with contribution of all source locations.

greens_gradient_in_observation_pointsnumpy.ndarray

An array of size (3,N_observationpoints) with the gradient of Green’s function of the wave field at the observation points with contribution of all source locations. When normals is None, gradient related quantities are not evaluated.

optimus.source.transducers.calc_field_from_point_sources(locations_source, locations_observation, frequency, density, wavenumber, source_weights, verbose)

Calculate the pressure field and its gradient of a point source, according to the Rayleigh integral formula. Use Numba or Multiprocessing acceleration.

Parameters
locations_sourcenumpy.ndarray

An array of size (3,N_sourcepoints) with the locations of the source points.

locations_observationnumpy.ndarray

An array of size (3,N_observationpoints) with the locations of the observation points.

frequencyfloat

The frequency of the wave field.

densityfloat

The density of the propagating medium.

wavenumbercomplex

The wavenumber of the wave field.

source_weightsnumpy.ndarray

An array of size (N_sourcepoints,) with the weights of each source element.

verboseboolean

Display output.

Returns
pressurenumpy.ndarray

An array of size (N_observationpoints,) with the pressure of the wave field in the observation points.

gradientnumpy.ndarray

An array of size (3,N_observationpoints) with the gradient of the pressure field in the observation points.

optimus.source.transducers.calc_field_from_point_sources_numpy(locations_source, locations_observation, frequency, density, wavenumber, source_weights)

Calculate the pressure field and its gradient of a point source, according to the Rayleigh integral formula.

Parameters
locations_sourcenumpy.ndarray

An array of size (3,N_sourcepoints) with the locations of the source points.

locations_observationnumpy.ndarray

An array of size (3,N_observationpoints) with the locations of the observation points.

frequencyfloat

The frequency of the wave field.

densityfloat

The density of the propagating medium.

wavenumbercomplex

The wavenumber of the wave field.

source_weightsnumpy.ndarray

An array of size (N_sourcepoints,) with the weights of each source element.

Returns
pressurenumpy.ndarray

An array of size (N_observationpoints,) with the pressure of the wave field in the observation points.

gradientnumpy.ndarray

An array of size (3,N_observationpoints) with the gradient of the pressure field in the observation points.

optimus.source.transducers.calc_field_from_point_sources_mp_source_para(parallelisation_index, locations_source, locations_observation, frequency, density, wavenumber, source_weights, chunks_index_source, chunks_index_field, number_of_observation_locations)

Computes the pressure and normal pressure gradient at field locations for selected source and field positions based on output from chunk_size_index. Used to calculate the incident field using multiprocessing and when parallelising over source locations.

Parameters
parallelisation_indexinteger

Index corresponding to the parallel job.

locations_sourcenumpy.ndarray

An array of size (3,N_sourcepoints) with the locations of the source points.

locations_observationnumpy.ndarray

An array of size (3,N_observationpoints) with the locations of the observation points.

frequencyfloat

The frequency of the wave field.

densityfloat

The density of the propagating medium.

wavenumbercomplex

The wavenumber of the wave field.

source_weightsnumpy.ndarray

An array of size (N_sourcepoints,) with the weights of each source element.

chunks_index_sourcenumpy.ndarray

An array of integers containing the indices to infer the source location chunks.

chunks_index_fieldnumpy.ndarray

An array of integers containing the indices to infer the field location chunks.

number_of_observation_locationsinteger

The total number of observation points.

Returns
pressurenumpy.ndarray

An array of size (N_observationpoints,) with the pressure of the wave field in the observation points.

gradientnumpy.ndarray

An array of size (3,N_observationpoints) with the gradient of the pressure field in the observation points.

optimus.source.transducers.calc_field_from_point_sources_mp_field_para(parallelisation_index, locations_source, locations_observation, frequency, density, wavenumber, source_weights, chunks_index_source, chunks_index_field)

Computes the pressure and normal pressure gradient at field locations for selected source and field positions based on output from chunk_size_index. Used to calculate the incident field using multiprocessing and when parallelising over field locations.

Parameters
parallelisation_indexinteger

Index corresponding to the parallel job.

locations_sourcenumpy.ndarray

An array of size (3,N_sourcepoints) with the locations of the source points.

locations_observationnumpy.ndarray

An array of size (3,N_observationpoints) with the locations of the observation points.

frequencyfloat

The frequency of the wave field.

densityfloat

The density of the propagating medium.

wavenumbercomplex

The wavenumber of the wave field.

source_weightsnumpy.ndarray

An array of size (N_sourcepoints,) with the weights of each source element.

chunks_index_sourcenumpy.ndarray

An array of integers containing the indices to infer the source location chunks.

chunks_index_fieldnumpy.ndarray

An array of integers containing the indices to infer the field location chunks.

Returns
pressurenumpy.ndarray

An array of size (N_observationpoints,) with the pressure of the wave field in the observation points.

gradientnumpy.ndarray

An array of size (3,N_observationpoints) with the gradient of the pressure field in the observation points.

optimus.source.transducers.chunk_size_index(locations_source, locations_observation)

Computes the chunk sizes and indices used to calculate the incident field using multiprocessing and when parallelising over source or observation locations. The chunks are allocated based on the global_parameters.incident_field.mem_per_core parameter.

Parameters
locations_sourcenumpy.ndarray

An array of size (3,N_sourcepoints) with the locations of the source points.

locations_observationnumpy.ndarray

An array of size (3,N_observationpoints) with the locations of the observation points.

Returns
chunks_index_sourcenumpy.ndarray

An array of integers describing the source point indices corresponding to each chunk for parallelisation.

chunks_index_observationnumpy.ndarray

An array of integers describing the observation point indices corresponding to each chunk for parallelisation.

number_of_source_chunksinteger

The number of source chunks used for parallelisation over source locations.

number_of_observation_chunksinteger

The number of source chunks used for parallelisation over observation locations.

optimus.source.transducers.break_in_chunks(number_of_locations, chunk_size)

Computes the chunk sizes and indices used to calculate the incident field using multiprocessing and when parallelising over source or observation locations. The chunks are allocated based on the global_parameters.incident_field.mem_per_core parameter.

Parameters
number_of_locationsinteger

The total number of source or observer locations.

chunk_sizeinteger

The size of the chunks the source or observation points are broken down into.

Returns
number_of_chunksinteger

The number of chunks the source or observation points are broken into.

chunks_indexnumpy.ndarray

An array of size (N,) containing the indices of the limits of the source or observation point chunks.