RectGrid3D#

class acoular.grids.RectGrid3D

Bases: RectGrid

Provide a cartesian 3D grid for the beamforming results.

The grid has cubic or nearly cubic cells. It is defined by lower and upper x-, y- and z-limits.

z_min

The lower z-limit that defines the grid. Default is -1.0.

z_max

The upper z-limit that defines the grid. Default is 1.0.

nxsteps

Number of grid points along x-axis. (read-only)

nysteps

Number of grid points along y-axis. (read-only)

nzsteps

Number of grid points along z-axis. (read-only)

increment

The cell side length for the grid. This can either be a scalar (same increments in all 3 dimensions) or a (3,) array of tuple of floats with respective increments in x-, y-, and z-direction. Default is 0.1.

digest

A unique identifier for the grid, based on its properties. (read-only)

index(x, y, z)

Return the indices for a grid point near a certain coordinate.

This can be used to query results or coordinates at or near a certain coordinate. Raises an exception if the given coordinate is outside the grid.

Parameters:
x, y, zfloat

The coordinates for which the indices is queried.

Returns:
3-tuple of ints

The indices that give the grid point nearest to the given x, y, z coordinates from an array with the same shape as the grid.

Examples

Check which of the points in a simple 8-point rectangular grid is closest to the point (0.5, 0.5, 1.0).

>>> import acoular as ac
>>>
>>> grid = ac.RectGrid3D()
>>> grid.increment = 2
>>> grid.pos
array([[-1., -1., -1., -1.,  1.,  1.,  1.,  1.],
       [-1., -1.,  1.,  1., -1., -1.,  1.,  1.],
       [-1.,  1., -1.,  1., -1.,  1., -1.,  1.]])
>>>
>>> grid.index(0.5, 0.5, 1.0)
(1, 1, 1)
indices(x1, y1, z1, x2, y2, z2)

Return the indices for a subdomain in the grid.

Allows box-shaped subdomains. This can be used to mask or to query results from a certain sector or subdomain.

Parameters:
x1, y1, z1, x2, y2, z2float

A box-shaped sector is assumed that is given by two corners (x1,y1,z1) and (x2,y2,z2).

Returns:
3-tuple of numpy.s_ objects

The indices that can be used to mask/select the grid subdomain from an array with the same shape as the grid.

Examples

Get the indices of the grid points of a simple 27-point rectangular grid which are located inside the first octant.

>>> import acoular as ac
>>>
>>> grid = ac.RectGrid3D()
>>> grid.increment = 1
>>> grid.pos
array([[-1., -1., -1., -1., -1., -1., -1., -1., -1.,  0.,  0.,  0.,  0.,
         0.,  0.,  0.,  0.,  0.,  1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.,
         1.],
       [-1., -1., -1.,  0.,  0.,  0.,  1.,  1.,  1., -1., -1., -1.,  0.,
         0.,  0.,  1.,  1.,  1., -1., -1., -1.,  0.,  0.,  0.,  1.,  1.,
         1.],
       [-1.,  0.,  1., -1.,  0.,  1., -1.,  0.,  1., -1.,  0.,  1., -1.,
         0.,  1., -1.,  0.,  1., -1.,  0.,  1., -1.,  0.,  1., -1.,  0.,
         1.]])
>>>
>>> grid.indices(0, 0, 0, 1, 1, 1)
(slice(1, 3, None), slice(1, 3, None), slice(1, 3, None))
export_gpos(filename)

Export the grid positions to an XML file.

This method generates an XML file containing the positions of all grid points. Each point is represented by a <pos> element with Name, x, y, and z attributes. The generated XML is formatted to match the structure required for importing into the ImportGrid class.

Parameters:
filenamestr

The path to the file to which the grid positions will be written. The file extension must be .xml.

Raises:
OSError

If the file cannot be written due to permissions issues or invalid file paths.

Notes

  • The file will be saved in UTF-8 encoding.

  • The Name attribute for each point is set as "Point {i+1}", where i is the index of the grid point.

  • If subgrids are defined, they will be included as the subgrid attribute.

Examples

Export a grid with 100 points to an XML file:

>>> import acoular as ac
>>> import numpy as np
>>> grid = ac.ImportGrid()
>>> # Create some grid points
>>> points = np.arange(9).reshape(3, 3)
>>> grid.pos = points
>>> grid.export_gpos('grid_points.xml')

The generated grid_points.xml file will look like this:

<?xml version="1.1" encoding="utf-8"?><Grid name="grid_points">
  <pos Name="Point 1" x="0" y="1" z="2"/>
  <pos Name="Point 2" x="3" y="4" z="5"/>
  <pos Name="Point 3" x="6" y="7" z="8"/>
</Grid>
subdomain(sector)

Return the indices for a subdomain in the grid.

Allows arbitrary subdomains of type Sector.

Parameters:
sectorSector object

Sector describing the subdomain.

Returns:
tuple

A 2-tuple of arrays of integers or numpy.s_ objects that can be used to mask or select the specified subdomain from a grid-shaped array.

Notes

The numpy.where() method is used to determine the the indices.

x_min

The lower x-limit that defines the grid. Default is -1.

x_max

The upper x-limit that defines the grid. Default is 1.

y_min

The lower y-limit that defines the grid. Default is -1.

y_max

The upper y-limit that defines the grid. Default is 1.

z

The constant z-coordinate of the grid plane. Default is 1.0.

extent

The grid’s extension in matplotlib.pyplot.imshow compatible form. (read-only)

size

The total number of grid points. This property is automatically calculated based on other defining attributes of the grid. (read-only)

shape

The shape of the grid, represented as a tuple. Primarily useful for Cartesian grids. (read-only)

pos

The grid positions represented as a (3, size) array of floats. (read-only) All positions’ coordinates are in meters by default (see here).