grids#

Implement support for multidimensional grids and integration sectors.

Inheritance diagram of acoular.grids.Grid, acoular.grids.RectGrid, acoular.grids.RectGrid3D, acoular.grids.ImportGrid, acoular.grids.LineGrid, acoular.grids.MergeGrid

Inheritance diagram of acoular.grids.Sector

Grid

Abstract base class for grid geometries.

RectGrid

Provides a 2D Cartesian grid for beamforming results.

RectGrid3D

Provide a cartesian 3D grid for the beamforming results.

ImportGrid

Load a 3D grid from an XML file.

LineGrid

Define a 3D grid for a line geometry.

MergeGrid

Base class for merging multiple grid geometries.

Sector

Abstract base class for all sector types.

SingleSector

Base class for single sector types.

RectSector

Class for defining a rectangular sector.

RectSector3D

Class for defining a cuboid sector.

CircSector

Class for defining a circular sector.

PolySector

Class for defining a polygon sector.

ConvexSector

Class for defining a convex hull sector.

MultiSector

Class for defining a sector consisting of multiple sectors.

Polygon(x, y)

Create an object representing a general polygon in a 2D plane.

in_hull(p, hull[, border, tol])

Test if points in p are in hull, in- or excluding the border.

acoular.grids.in_hull(p, hull, border=True, tol=0)#

Test if points in p are in hull, in- or excluding the border.

Parameters:
pnumpy.ndarray of floats, shape (N, K)

Coordinates of N points in K dimensions.

hullnumpy.ndarray of floats, shape (M, K), or Delaunay object

Coordinates of M points in K dimensions for which Delaunay triangulation will be computed.

borderbool, optional

Points in p on the border of hull will be kept in the return if True. If False, only points inside hull will be kept. Default is True.

tolfloat, optional

Tolerance allowed in the inside-triangle check. Default is 0.

Returns:
numpy.ndarray of bools

An array of boolean values indicating which points in p are inside the hull, same shape as p. Each entry is True if the corresponding point is inside the hull (or on the border, if border=True), and False otherwise.

Notes

This function uses Delaunay triangulation to determine if a point is inside the convex hull, which is efficient and robust for arbitrary shapes in higher-dimensional spaces.

Examples

>>> from acoular.grids import in_hull
>>> import numpy as np
>>> from scipy.spatial import Delaunay
>>> points = np.array([[0, 0], [1, 0], [0, 1], [1, 1]])
>>> hull = Delaunay(points)
>>> p = np.array([[0.5, 0.5], [2, 2]])
>>> in_hull(p, hull)
array([ True, False])
class acoular.grids.Polygon(x, y)#

Bases: object

Create an object representing a general polygon in a 2D plane.

This class allows defining a polygon by specifying the coordinates of its vertices and provides methods for checking whether a set of points lies inside the polygon, or if a point is closer to a side or vertex of the polygon.

Parameters:
xarray_like

Array of x-coordinates of the vertices that define the polygon. These coordinates should form a closed shape (i.e., the last point should be the same as the first point).

yarray_like

Array of y-coordinates of the vertices that define the polygon. These coordinates should correspond to the x-coordinates, forming a closed shape.

Attributes:
xnumpy.ndarray

Array of x-coordinates of the polygon vertices.

ynumpy.ndarray

Array of y-coordinates of the polygon vertices.

is_inside(xpoint, ypoint, smalld=1e-12)#

Check if a point or set of points are inside the polygon.

Parameters:
xpointfloat or array_like

Array of x-coordinates of the points to be tested.

ypointfloat or array_like

Array of y-coordinates of the points to be tested.

smalldfloat, optional

Tolerance used for floating point comparisons when checking if a point is exactly on a polygon’s edge. The default value is 1e-12.

Returns:
float or array_like

The distance from the point to the nearest point on the polygon. The values returned have the following meanings: - mindst < 0: Point is outside the polygon. - mindst = 0: Point is on an edge of the polygon. - mindst > 0: Point is inside the polygon.

Notes

The method uses an improved algorithm based on Nordbeck and Rydstedt for determining whether a point is inside a polygon [16].

class acoular.grids.Grid#

Bases: ABCHasStrictTraits

Abstract base class for grid geometries.

This class defines a common interface for all grid geometries and provides tools to query grid properties and related data. It is intended to serve as a base class for specialized grid implementations and should not be instantiated directly as it lacks concrete functionality.

Unit of length

The source code is agnostic to the unit of length. The positions’ coordinates are assumed to be in meters. This is consistent with the standard Environment class which uses the speed of sound at 20°C at sea level under standard atmosphere pressure in m/s. If the positions’ coordinates are provided in a unit other than meter, it is advisable to change the c attribute to match the given unit.

size = Property(desc='overall number of grid points')#

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

shape = Property(desc='grid shape as tuple')#

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

pos = Property(desc='x, y, z positions of grid points')#

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

digest = Property#

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

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.

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>
class acoular.grids.RectGrid#

Bases: Grid

Provides a 2D Cartesian grid for beamforming results.

This grid is composed of square or nearly square cells and lies on a plane perpendicular to the z-axis. It is defined by the lower and upper x- and y-limits and a constant z-coordinate.

x_min = Float(-1.0, desc='minimum  x-value')#

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

x_max = Float(1.0, desc='maximum  x-value')#

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

y_min = Float(-1.0, desc='minimum  y-value')#

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

y_max = Float(1.0, desc='maximum  y-value')#

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

z = Float(1.0, desc='position on z-axis')#

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

increment = Float(0.1, desc='step size')#

The side length of each cell. Default is 0.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)

extent = Property(desc='grid extent as (x_min, x_max, y_min, y_max)')#

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

digest = Property( #

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

extend()#

Return the grid’s extension in matplotlib.pyplot.imshow compatible form.

Returns:
tuple of floats

(x_min, x_max, y_min, y_max) representing the grid’s extent.

Notes

This method is deprecated. Use the extent property instead.

index(x, y)#

Find the indices of a grid point near a given coordinate.

Parameters:
xfloat

The x coordinate of interest.

yfloat

The y coordinate of interest.

Returns:
tuple of ints

Indices corresponding to the nearest grid point.

Raises:
ValueError

If the coordinates are outside the grid boundaries.

indices(*r)#

Find the indices of a subdomain in the grid.

Supports rectangular, circular, and polygonal subdomains.

Parameters:
rtuple of floats
Defines the subdomain shape and dimensions:
  • If 3 values are provided: center (x1, y1) and radius r2 define a circle.

  • If 4 values are provided: corners (x1, y1) and (x2, y2) define a rectangle.

  • If more than 4 values are provided: vertices (xn, yn) define a polygon.

Returns:
tuple

A 2-tuple of indices or slices corresponding to the subdomain.

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.0.

z_max = Float(1.0, desc='maximum  z-value')#

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

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 of floats with respective increments in x-, y-, and z-direction. Default is 0.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, 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))
class acoular.grids.ImportGrid#

Bases: Grid

Load a 3D grid from an XML file.

This class is used to import a 3D grid defined in an XML file. The grid’s positions and subgrid names are parsed and stored for further processing.

file = Union(None, File(filter=['*.xml'], exists=True), desc='name of the xml file to import')#

Name of the .xml-file from which to read the data.

subgrids = CArray(desc='names of subgrids for each point')#

Names of subgrids for each point. This is an optional property, typically used when grids are divided into named subregions.

digest = Property(depends_on=['_gpos'])#

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

class acoular.grids.LineGrid#

Bases: Grid

Define a 3D grid for a line geometry.

The LineGrid class represents a grid where points are arranged linearly in 3D space. The grid is defined by a starting location (loc), a direction vector (direction), a total length (length), and the number of points (num_points) along the line.

Notes

  • The distance between points is length / ( num_points - 1).

  • The direction vector is normalized to ensure consistency.

Examples

Create a line grid with 5 points along the x-axis, starting at (0, 0, 0), with a length of 4 meters:

>>> import acoular as ac
>>> grid = ac.LineGrid()
>>> grid.loc = (0.0, 0.0, 0.0)
>>> grid.direction = (1.0, 0.0, 0.0)
>>> grid.length = 4
>>> grid.num_points = 5
>>> grid.pos
array([[0., 1., 2., 3., 4.],
       [0., 0., 0., 0., 0.],
       [0., 0., 0., 0., 0.]])
loc = Tuple((0.0, 0.0, 0.0))#

Starting point of the grid in 3D space. Default is (0.0, 0.0, 0.0).

direction = Tuple((1.0, 0.0, 0.0), desc='Line orientation ')#

A vector defining the orientation of the line in 3D space. Default is (1.0, 0.0, 0.0).

length = Float(1, desc='length of the line source')#

Total length of the line. Default is 1.0.

num_points = Int(1, desc='length of the line source')#

Number of grid points along the line. Default is 1.

size = Property(desc='overall number of grid points')#

The total number of grid points. Automatically updated when other grid-defining attributes are set. (read-only)

pos = Property(desc='x, y, z positions of grid points')#

A (3, size) array containing the x, y, and z positions of the grid points. (read-only) All positions’ coordinates are in meters by default (see here).

digest = Property( #

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

class acoular.grids.MergeGrid#

Bases: Grid

Base class for merging multiple grid geometries.

The MergeGrid class allows the combination of multiple grid geometries into a single unified grid. Each input grid is assigned a subdomain in the resulting grid, and all properties, such as positions and identifiers, are appropriately merged.

Notes

  • The merged grid eliminates duplicate points based on their positions.

  • Each subgrid retains its original grid properties, such as digest and size.

Examples

Merging two simple grids:

>>> import acoular as ac
>>> grid1 = ac.LineGrid(loc=(0, 0, 0), direction=(1, 0, 0), length=1, num_points=3)
>>> grid2 = ac.LineGrid(loc=(0, 0, 0), direction=(0, 1, 0), length=1, num_points=3)
>>> merged_grid = ac.MergeGrid()
>>> merged_grid.grids = [grid1, grid2]
>>> merged_grid.size
5
>>> merged_grid.pos
array([[0. , 0. , 0. , 0.5, 1. ],
       [0. , 0.5, 1. , 0. , 0. ],
       [0. , 0. , 0. , 0. , 0. ]])
grids = List(desc='list of grids')#

A list of Grid objects to be merged. Each grid is treated as a subdomain in the resulting merged grid.

grid_digest = Str(desc='digest of the merged grids')#

A list of unique digests for each grid being merged. (read-only)

subgrids = Property(desc='names of subgrids for each point')#

Names of subgrids corresponding to each point in the merged grid. (read-only)

digest = Property(depends_on=['grids', 'grid_digest'])#

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

class acoular.grids.Sector#

Bases: ABCHasStrictTraits

Abstract base class for all sector types.

The Sector class defines the common interface for all sector implementations. It serves as the base class for creating diverse sector geometries, each capable of determining whether specific grid points fall within its bounds.

When used directly, this class represents a sector encompassing the entire grid, meaning all positions are considered valid.

Notes

This class is designed to be subclassed. Derived classes should override the contains() method to implement specific sector geometries (e.g., circular, rectangular, or polygon shapes).

Examples

Load example data and set different Sectors for integration in the sector integration example.

contains(pos)#

Check whether the given coordinates lie within the sector’s bounds.

This method determines if each column of the input array pos corresponds to a point that falls within the sector. For this base class, all points are considered within the sector.

Parameters:
posnumpy.ndarray

A 2D array with shape (3, N), where N is the number of grid points. Each column represents the x, y, z coordinates of a grid point.

Returns:
numpy.ndarray of bools

A 1D array of length N, where each entry indicates whether the corresponding column in pos lies within the sector’s bounds.

Examples

>>> import numpy as np
>>> import acoular as ac
>>> sector = ac.Sector()
>>> positions = np.array([[0, 1], [0, 0], [0, 0]])  # Two grid points
>>> sector.contains(positions)
array([ True,  True])
class acoular.grids.SingleSector#

Bases: Sector

Base class for single sector types.

Defines the common interface for all single sector classes. This class can serve as a base for various single sector implementations. When used directly, it defines a sector that encompasses the whole grid. It includes attributes for handling border inclusion, tolerance for sector borders, and default behavior when no points are inside the sector.

Examples

Load example data and set diffrent Sectors for intergration in the sector integration example.

include_border = Bool(True, desc='include points on the border')#

If True, grid points lying on the sector border are included in the sector. Default is True.

abs_tol = Float(1e-12, desc='absolute tolerance for sector border')#

The absolute tolerance to apply when determining if a grid point lies on the sector border. Default is 1e-12.

default_nearest = Bool(True, desc='``contains`` method return nearest grid point to center if none inside sector')#

If True, the contains method (as in RectSector.contains(), RectSector3D.contains(), CircSector.contains(), and PolySector.contains()) returns the nearest grid point if no grid points are inside the sector. Default is True.

class acoular.grids.RectSector#

Bases: SingleSector

Class for defining a rectangular sector.

Defines a rectangular sector either for 2D grids (rectangle in the XY-plane) or for 3D grids (rectangular cylindrical sector parallel to the z-axis). The sector is bounded by the specified x_min, x_max, y_min, and y_max positions, defining the lower and upper bounds of the rectangle along the x and y axes.

Examples

Load example data and set diffrent Sectors for intergration in the sector integration example.

x_min = Float(-1.0, desc='minimum x position of the rectangle')#

The minimum x position of the rectangle. Default is -1.0.

x_max = Float(1.0, desc='maximum x position of the rectangle')#

The maximum x position of the rectangle. Default is 1.0.

y_min = Float(-1.0, desc='minimum y position of the rectangle')#

The minimum y position of the rectangle. Default is -1.0.

y_max = Float(1.0, desc='maximum y position of the rectangle')#

The maximum y position of the rectangle. Default is 1.0.

contains(pos)#

Check if the coordinates in a given array lie within the rectangular sector.

If no coordinate is inside, the nearest one to the rectangle center is returned if default_nearest is True.

Parameters:
posarray of floats

A (3, N) array containing the positions of N grid points.

Returns:
numpy.ndarray of bools

An array of shape (N,) indicating which of the given positions lie within the given sector.

Examples

>>> import acoular as ac
>>> grid = ac.RectGrid(increment=2)
>>> sec = ac.RectSector(x_min=0, y_min=0)
>>> sec.contains(grid.pos)
array([False, False, False,  True])
class acoular.grids.RectSector3D#

Bases: RectSector

Class for defining a cuboid sector.

This class extends the RectSector class to define a cuboid sector, which can be used for 3D grids. The cuboid sector is defined by its bounds along the x, y, and z axes.

Examples

Load example data and set diffrent Sectors for intergration in the sector integration example.

z_min = Float(-1.0, desc='minimum z position of the cuboid')#

The lower z position of the cuboid. Default is -1.0.

z_max = Float(1.0, desc='maximum z position of the cuboid')#

The upper z position of the cuboid. Default is 1.0.

contains(pos)#

Check if the coordinates in a given array lie within the cuboid sector.

The method checks if the points in the provided position array are within the cuboid defined by the bounds along the x, y, and z axes. If no point is inside the sector, and if default_nearest is True, the nearest point to the center of the cuboid is returned.

Parameters:
posarray of floats

A (3, N) array containing the positions of N grid points, where each point is represented by its x, y, and z coordinates.

Returns:
numpy.ndarray of bools

A boolean array of shape shape (N,) indicating which of the given positions lie within the cuboid sector. True if the grid point is inside the cuboid, otherwise False.

Examples

>>> import acoular as ac
>>> grid = ac.RectGrid3D(increment=2)
>>> sec = ac.RectSector3D(x_min=0, y_min=0, z_min=0)
>>> sec.contains(grid.pos)
array([False, False, False, False, False, False, False,  True])
class acoular.grids.CircSector#

Bases: SingleSector

Class for defining a circular sector.

Defines a circular sector, which can be used for both 2D grids (as a circle in the XY-plane) or for 3D grids (as a cylindrical sector parallel to the z-axis). The sector is defined by its center position (x, y) and its radius r.

Examples

Load example data and set diffrent Sectors for intergration in the sector integration example.

x = Float(0.0, desc='x position of the circle center')#

The x position of the circle center. Default is 0.0.

y = Float(0.0, desc='y position of the circle center')#

The y position of the circle center. Default is 0.0.

r = Float(1.0, desc='radius of the circle')#

Radius of the circle. Default is 1.0.

contains(pos)#

Check if the coordinates in a given array lie within the circular sector.

The method calculates the squared distance of each point from the center of the circle and checks if it lies within the sector, considering the sector’s radius r. If no point is inside and default_nearest is True, the nearest point outside the sector will be returned.

Parameters:
posarray of floats

A (3, N) array containing the positions of N grid points, where each point is represented by its x, y, and z coordinates.

Returns:
numpy.ndarray of bools

A boolean array of shape shape (N,) indicating which of the given positions lie within the circular sector. True if the grid point is inside the circular sector, otherwise False.

Examples

>>> import acoular as ac
>>> grid = ac.RectGrid(increment=1)
>>> grid.pos
array([[-1., -1., -1.,  0.,  0.,  0.,  1.,  1.,  1.],
       [-1.,  0.,  1., -1.,  0.,  1., -1.,  0.,  1.],
       [ 1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.]])
>>> sec = ac.CircSector(x=1, y=1, r=0.5)
>>> sec.contains(grid.pos)
array([False, False, False, False, False, False, False, False,  True])
class acoular.grids.PolySector#

Bases: SingleSector

Class for defining a polygon sector.

Inherits from SingleSector and provides functionality to define a polygonal sector on a 2D grid.

Notes

The polygon is specified by the Polygon class.

Examples

Load example data and set diffrent Sectors for intergration in the sector integration example.

edges = List(Float)#

List of coordinates representing the polygon’s vertices. The coordinates must define a closed polygon like x1, y1, x2, y2, ... xn, yn.

contains(pos)#

Check if the coordinates in a given array lie within the polygon sector.

If no coordinate is inside, the nearest one to the rectangle center is returned if default_nearest is True.

Parameters:
posarray of floats

A (3, N) array containing the positions of N grid points where each point is represented by its x, y, and z coordinates.

Returns:
numpy.ndarray of bools

A boolean array of shape (N,) indicating which of the given positions lie within the polygon sector. True if the grid point is inside the polygon, otherwise False.

Examples

>>> import acoular as ac
>>> grid = ac.RectGrid(increment=1)
>>> grid.pos
array([[-1., -1., -1.,  0.,  0.,  0.,  1.,  1.,  1.],
       [-1.,  0.,  1., -1.,  0.,  1., -1.,  0.,  1.],
       [ 1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.]])
>>> sec = ac.PolySector(edges=[0, 0, 1, 0, 1, 1, 0, 1])
>>> sec.contains(grid.pos)
array([False, False, False, False,  True,  True, False,  True,  True])
class acoular.grids.ConvexSector#

Bases: SingleSector

Class for defining a convex hull sector.

This class defines a convex hull sector for 2D grids. The sector is created using a list of edge coordinates edges which represent the vertices of a polygon. The convex hull is the smallest convex shape that contains all the given vertices.

Examples

Load example data and set diffrent Sectors for intergration in the sector integration example.

edges = List(Float)#

List of edge coordinates that define the convex hull. The coordinates must define a closed polygon that forms the convex hull like x1, y1, x2, y2, … xn, yn.

contains(pos)#

Check if the coordinates in a given array lie within the convex sector.

If no coordinate is inside, the nearest one to the rectangle center is returned if default_nearest is True.

Parameters:
posarray of floats

Array containing the positions of N grid points, shape (3, N).

Returns:
numpy.ndarray of bools

An array of shape (N,) indicating which of the given positions lie within the given sector.

Examples

>>> import acoular as ac
>>> grid = ac.RectGrid(increment=1)
>>> grid.pos
array([[-1., -1., -1.,  0.,  0.,  0.,  1.,  1.,  1.],
       [-1.,  0.,  1., -1.,  0.,  1., -1.,  0.,  1.],
       [ 1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.]])
>>> sec = ac.ConvexSector(edges=[0, 0, 1, 0, 1, 1, 0, 1])
>>> sec.contains(grid.pos)
array([False, False, False, False,  True,  True, False,  True,  True])
class acoular.grids.MultiSector#

Bases: Sector

Class for defining a sector consisting of multiple sectors.

This class allows the combination of several individual sectors into one.

Examples

Load example data and set diffrent Sectors for intergration in the sector integration example.

sectors = List(Instance(Sector))#

List of Sector objects to be mixed, each defining a different sector.

contains(pos)#

Check if the coordinates in a given array lie within any of the sub-sectors.

This method iterates over the list of sectors, checking if each point in the given position array lies within any of the defined sectors.

Parameters:
posarray of floats

A (3, N) array containing the positions of N grid points, where each point is represented by its x, y, and z coordinates.

Returns:
numpy.ndarray of bools

A boolean array of shape (N,) indicating which of the given positions lie within any of the defined sectors. True if the grid point is inside the circular sector, False if otherwise.

Examples

>>> import acoular as ac
>>> grid = ac.RectGrid(increment=1)
>>> grid.pos
array([[-1., -1., -1.,  0.,  0.,  0.,  1.,  1.,  1.],
       [-1.,  0.,  1., -1.,  0.,  1., -1.,  0.,  1.],
       [ 1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.,  1.]])
>>> sec1 = ac.RectSector(x_min=0, y_min=0)
>>> sec2 = ac.CircSector(x=1, y=1, r=0.5)
>>> multi_sec = ac.MultiSector(sectors=[sec1, sec2])
>>> multi_sec.contains(grid.pos)
array([False, False, False, False,  True,  True, False,  True,  True])