optimus.source.array
Array sources.
Module Contents
Classes
Functions
|
Create an array source consisting of circular piston elements distributed on a spherical section bowl. |
- optimus.source.array.create_array(frequency, element_radius, velocity=1.0, source_axis=(1, 0, 0), number_of_point_sources_per_wavelength=6, location=(0, 0, 0), centroid_locations=None, centroid_locations_filename=None)
Create an array source consisting of circular piston elements distributed on a spherical section bowl.
- Parameters
- frequencyfloat
The frequency of the acoustic field.
- element_radiusfloat
The radius of elements which lie on the spherical section bowl.
- velocitycomplex, numpy.ndarray complex
Array of size (N,) with complex values for the normal velocities of the array elements. If one value is specified, this will be repeated for all array elements. Default : 1 m/s
- source_axistuple float
The direction vector of the axis of the bowl. Default: positive x direction
- number_of_point_sources_per_wavelengthinteger
The number of point sources per wavelength used to discretise each piston. Default: 6
- locationtuple float
The location of the centroid of the bowl. Default: global origin
- centroid_locationsnumpy.ndarray
An array of size (3, N) with the locations of the centroids of the piston elements. These must be specified in a local coordinate sytem where the axis of the transducer is the Cartesian positive z-axis and the focus of the transducer is (0,0,0).
- centroid_locations_filenamestr
Path and filename containing the centroid locations data. The file extension has to be “.dat”.
- class optimus.source.array._Array(frequency, element_radius, velocity, source_axis, number_of_point_sources_per_wavelength, location, centroid_locations, centroid_locations_filename)
Bases:
optimus.source.common.Source
- _calc_centroid_locations(centroid_locations, centroid_locations_filename)
Calculates centroid locations.
- Parameters
- centroid_locationsnumpy.ndarray
An array of size (3, N) with the locations of the centroids of the piston elements.
- centroid_locations_filenamestr
Path and filename containing the centroid locations data. The file extension has to be “.dat”.
- Returns
- centroid_locationsnumpy.ndarray
An array of size (3, N) with the locations of the centroids of the piston elements.
- _calc_radius_of_curvature(centroid_locations)
Calculates the radius of curvature of the array transducer from centroid locations.
- Parameters
- centroid_locationsnumpy.ndarray
An array of size (3, N) with the locations of the centroids of the piston elements.
- Returns
- radius of curvaturefloat
The radius of curvature of the array.
- pressure_field(medium, locations)
Calculate the pressure field in the specified locations.
- Parameters
- mediumoptimus.material.Material
The propagating medium.
- locationsnumpy.ndarray
An array of size (3,N) with the locations on which to evaluate the pressure field.
- Returns
- pressurenp.ndarray
An array of size (N,) with the pressure in the locations.
- normal_pressure_gradient(locations, normals, medium)
Calculate the normal gradient of the pressure field in the specified locations.
- Parameters
- mediumoptimus.material.Material
The propagating medium.
- locationsnumpy.ndarray
An array of size (3,N) with the locations on which to evaluate the pressure field.
- normalsnumpy.ndarray
An array of size (3,N) with the unit normal vectors at the locations on which to evaluate the pressure field.
- Returns
- gradientnumpy.ndarray
An array of size (3,N) with the normal gradient of the pressure in the locations.
- pressure_field_and_normal_gradient(medium, locations, normals)
Calculate the pressure field and the normal gradient of the pressure field in the specified locations.
- Parameters
- mediumoptimus.material.Material
The propagating medium.
- locationsnumpy.ndarray
An array of size (3,N) with the locations on which to evaluate the pressure field.
- normalsnumpy.ndarray
An array of size (3,N) with the unit normal vectors at the locations on which to evaluate the pressure field.
- Returns
- pressurenumpy.ndarray
An array of size (N,) with the pressure in the locations.
- gradientnumpy.ndarray
An array of size (3,N) with the normal gradient of the pressure in the locations.
- calc_surface_traces(medium, space_dirichlet=None, space_neumann=None, dirichlet_trace=True, neumann_trace=True)
Calculate the surface traces of the source field on the mesh.
- Parameters
- mediumoptimus.material.Material
The propagating medium.
- space_dirichlet, space_neumannbempp.api.FunctionSpace
The discrete spaces on the surface grid.
- dirichlet_trace, neumann_tracebool
Calculate the Dirichlet or Neumann trace of the field.
- Returns
- tracebempp.api.GridFunctions
The surface traces.