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 = Float(-1.0, desc='minimum z-value')¶
The lower z-limit that defines the grid. Default is
-1
.
- z_max = Float(1.0, desc='maximum z-value')¶
The upper z-limit that defines the grid. Default is
1
.
- nxsteps = Property(desc='number of grid points along x-axis')¶
Number of grid points along x-axis. (read-only)
- nysteps = Property(desc='number of grid points along y-axis')¶
Number of grid points along y-axis. (read-only)
- nzsteps = Property(desc='number of grid points along z-axis')¶
Number of grid points along z-axis. (read-only)
- increment = Property(desc='step size')¶
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
offloats
with respective increments in x-, y-, and z-direction. Default is0.1
.
- digest = Property( …¶
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, z
float
The coordinates for which the indices is queried.
- x, y, z
- Returns:
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, z2
float
A box-shaped sector is assumed that is given by two corners
(x1,y1,z1)
and(x2,y2,z2)
.
- x1, y1, z1, x2, y2, z2
- Returns:
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))