environments#

Implements acoustic environments with and without flow.

Inheritance diagram of acoular.environments.Environment

Inheritance diagram of acoular.environments.FlowField

Environment

A simple acoustic environment without flow.

UniformFlowEnvironment

An acoustic environment with uniform flow.

GeneralFlowEnvironment

An acoustic environment with a generic flow field.

FlowField

An abstract base class for a spatial flow field.

OpenJet

Analytical approximation of the flow field of an open jet.

RotatingFlow

Analytical approximation of a rotating flow field with velocity component in z-direction.

SlotJet

Analytical approximation of the flow field of a slot jet.

dist_mat(gpos, mpos)

Compute distance matrix.

cylToCart(x[, Q])

Return cartesian coordinate representation of an input array in cylindrical coordinate.

cartToCyl(x[, Q])

Return cylindrical coordinate representation of an input array in cartesian coordinate.

spiral_sphere(N[, Om, b])

Generate unit vectors equally distributed over a sphere or a portion of it.

acoular.environments.dist_mat(gpos, mpos)#

Compute distance matrix. (accelerated with numba).

Given an (3, N) array of the locations of points in the beamforming map grid in 3D cartesian coordinates and (3, M) array of the locations of microphones in 3D cartesian coordinates, the (N, M) matrix of the distances between each microphone and each point in the beamforming map grip.

Parameters:
gposnumpy.ndarray of floats

The locations of N points in the beamforming map grid in 3D cartesian coordinates, shape (3, N).

mposnumpy.ndarray of floats

The locations of M microphones in 3D cartesian coordinates, shape (3, M).

Returns:
numpy.ndarray of floats

Matrix of the distances between each microphone and each point in the beamforming map grid, shape (N, M).

acoular.environments.cartToCyl(x, Q=None)#

Return cylindrical coordinate representation of an input array in cartesian coordinate.

Return the cylindrical coordinate representation of an input position which was before transformed into a modified cartesian coordinate, which has flow into positive z direction.

Parameters:
xnumpy.ndarray of floats

Cartesian coordinates of N points, shape (3, N).

Qnumpy.ndarray of floats, optional

Orthogonal transformation matrix, shape (3, 3). If provided, the positional vectors are transformed via new_x = Q * x, before transforming those modified coordinates into cylindrical ones. Default is the identity matrix.

Returns:
numpy.ndarray of floats

Cylindrical representation of given N points in cartesian coodrinates as an array of shape (3, N) with new coordinates \((\phi, r, z)\).

acoular.environments.cylToCart(x, Q=None)#

Return cartesian coordinate representation of an input array in cylindrical coordinate.

Return the cartesian coordinate representation of a input position which was before transformed into a cylindrical coordinate, which has flow into positive z direction.

Parameters:
xnumpy.ndarray of floats

Cylindrical coordinates of N points, shape (3, N).

Qnumpy.ndarray of floats, optional

Orthogonal transformation matrix, shape (3, 3). If provided, the positional vectors are transformed via new_x = Q * x before transforming those modified coordinates into cartesian ones. Default is the identity matrix.

Returns:
numpy.ndarray of floats

Cartesian representation of given N points in cylindrical coodrinates as an array of shape (3, N) with coodinates \((x, y, z)\).

class acoular.environments.Environment#

Bases: HasStrictTraits

A simple acoustic environment without flow.

This class models an acoustic environment where the propagation of sound is considered to occur in a homogeneous medium without any flow effects (e.g., wind). It provides functionality for computing the travel time or distances between grid point locations and microphone locations.

Notes

  • The dist_mat() function is used internally to compute the pairwise distances between grid points and microphone positions efficiently.

  • This class assumes a static, homogeneous environment without accounting for factors like temperature gradients, humidity, or atmospheric turbulence.

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

A unique identifier based on the environment properties. (read-only)

c = Float(343.0, desc='speed of sound')#

The speed of sound in the environment. Default is 343.0, which corresponds to the approximate speed of sound at 20°C in dry air at sea level, if the unit is m/s.

roi = Union(None, CArray)#

The region of interest (ROI) for calculations. (Not needed for most types of environment.) Default is None.

class acoular.environments.UniformFlowEnvironment#

Bases: Environment

An acoustic environment with uniform flow.

This class models an acoustic environment where sound propagates in a medium with uniform flow. It extends the Environment class to account for the effects of flow on sound propagation, such as changes in travel times and distances due to advection by the flow field.

The flow is assumed to be uniform and steady, characterized by its Mach number (ma) and direction (fdv).

Notes

The effective distance is adjusted by solving a flow-dependent relationship that accounts for the cosine of the angle between the flow direction and the propagation path.

ma = Float(0.0, desc='flow mach number')#

The Mach number of the flow, defined as the ratio of the flow velocity to the speed of sound. Default is 0.0, which corresponds to no flow.

fdv = CArray(dtype=np.float64, shape=(3,), value=np.array((1.0, 0, 0)), desc='flow direction')#

A unit vector specifying the direction of the flow in 3D Cartesian coordinates. Default is (1.0, 0, 0), which corresponds to flow in the x-direction.

digest = Property( #

A unique identifier based on the environment properties. (read-only)

class acoular.environments.FlowField#

Bases: ABCHasStrictTraits

An abstract base class for a spatial flow field.

digest = Property#

A unique identifier based on the field properties. (read-only)

abstractmethod v(xx)#

Provide the flow field as a function of the location.

Parameters:
xxnumpy.ndarray of floats

Location in the fluid for which to provide the data, shape (3,).

class acoular.environments.SlotJet#

Bases: FlowField

Analytical approximation of the flow field of a slot jet.

This class provides an analytical model of the velocity field for a slot jet, based on the description in [1]. It describes the flow field originating from a slot nozzle, including the jet core and shear layer, and calculates the velocity field and its Jacobian matrix at a given location.

Notes

  • The slot jet is divided into two regions: the jet core and the shear layer. The model distinguishes between these regions based on the non-dimensionalized distance from the slot exit plane.

  • The flow field is aligned with the direction of the flow vector, while the plane vector helps define the orientation of the slot.

v0 = Float(0.0, desc='exit velocity')#

Exit velocity at the slot jet origin (nozzle). Default is 0.0.

origin = CArray(dtype=np.float64, shape=(3,), value=np.array((0.0, 0.0, 0.0)), desc='center of nozzle')#

The location of the slot nozzle center. Default is (0.0, 0.0, 0.0).

flow = CArray(dtype=np.float64, shape=(3,), value=np.array((1.0, 0.0, 0.0)), desc='flow direction')#

Unit vector representing the flow direction. Default is (1.0, 0.0, 0.0).

plane = CArray(dtype=np.float64, shape=(3,), value=np.array((0.0, 1.0, 0.0)), desc='slot center line direction')#

Unit vector parallel to the slot center plane, used to define the slot orientation. Default is (0.0, 1.0, 0.0).

B = Float(0.2, desc='nozzle diameter')#

Width of the slot (slot diameter). Default is 0.2.

l = Float(5.2, desc='flow establishment length')  # noqa: E741#

Non-dimensional length of the zone of flow establishment (jet core length). Default is 5.2.

digest = Property( #

A unique identifier based on the field properties. (read-only)

v(xx)#

Compute the velocity field and its Jacobian matrix at a given location.

This method provides the velocity vector and its Jacobian matrix at the given location xx in the fluid. The velocity is computed for the flow component in the direction of the slot jet, while the entrainment components are assumed to be zero.

Parameters:
xxnumpy.ndarray of floats

The 3D Cartesian coordinates of the location in the fluid where the velocity field is to be computed, shape (3,).

Returns:
velocity_vectornumpy.ndarray of floats

The velocity vector at the given location, shape (3,).

jacobian_matrixnumpy.ndarray of floats

The Jacobian matrix of the velocity vector field at the given location, shape (3,).

Notes

  • The velocity field is computed using a local coordinate system aligned with the flow direction and the slot orientation.

  • The velocity profile depends on whether the point lies within the jet core region or in the shear layer. In the jet core, the velocity is constant, while in the shear layer, it decays following a Gaussian distribution.

  • The Jacobian matrix provides the partial derivatives of the velocity vector components with respect to the spatial coordinates.

class acoular.environments.OpenJet#

Bases: FlowField

Analytical approximation of the flow field of an open jet.

This class provides a simplified analytical model of the velocity field for an open jet, based on the description in [1]. It calculates the velocity vector and its Jacobian matrix at a given location in the fluid domain, assuming flow in the x-direction only.

Notes

  • This is not a fully generic implementation, and is limited to flow in the x-direction only. No other directions are possible at the moment and flow components in the other direction are zero.

  • The flow field transitions from the jet core to the shear layer, with velocity decay modeled using a Gaussian profile in the shear layer.

Examples

>>> import acoular as ac
>>> import numpy as np
>>>
>>> jet = ac.OpenJet(v0=10.0, D=0.4, l=6.2)
>>> velocity, jacobian = jet.v(np.array((1.0, 0.1, 0.1)))
>>> velocity
array([9.62413564, 0.        , 0.        ])
>>> jacobian
array([[ -1.92660591, -23.25619062, -23.25619062],
       [  0.        ,   0.        ,   0.        ],
       [  0.        ,   0.        ,   0.        ]])
v0 = Float(0.0, desc='exit velocity')#

Exit velocity at the jet origin (nozzle). Default is 0.0.

origin = CArray(dtype=np.float64, shape=(3,), value=np.array((0.0, 0.0, 0.0)), desc='center of nozzle')#

The location of the nozzle center. Default is (0.0, 0.0, 0.0).

D = Float(0.2, desc='nozzle diameter')#

Diameter of the nozzle. Default is 0.2.

l = Float(6.2, desc='flow establishment length')  # noqa: E741#

Non-dimensional length of the zone of flow establishment (jet core length). Default is 6.2. [1]

digest = Property( #

A unique identifier based on the field properties. (read-only)

v(xx)#

Compute the velocity field and its Jacobian matrix at a given location.

This method calculates the velocity vector and its Jacobian matrix at the specified location xx in the fluid domain. The velocity is modeled only for the x-component of the flow, while the y- and z-components are assumed to be zero.

Parameters:
xxnumpy.ndarray of floats

The 3D Cartesian coordinates of the location in the fluid where the velocity field is to be computed, shape (3,).

Returns:
velocity_vectornumpy.ndarray of floats

The velocity vector at the specified location xx, shape (3,).

jacobian_matrixnumpy.ndarray of floats

The Jacobian matrix of the velocity vector field at the specified location xx, shape (3, 3).

Notes

  • The velocity field is determined based on whether the location is within the jet core region or in the shear layer. Within the jet core, the velocity is constant and equal to v0. In the shear layer, the velocity decays following a Gaussian distribution.

  • The Jacobian matrix provides the partial derivatives of the velocity components with respect to the spatial coordinates.

  • If the radial distance r from the jet axis is zero, the derivatives with respect to y and z are set to zero to avoid division by zero.

class acoular.environments.RotatingFlow#

Bases: FlowField

Analytical approximation of a rotating flow field with velocity component in z-direction.

This class provides an analytical model for a fluid flow field with a rigid-body-like rotation about the z-axis. The flow combines rotational motion in the x-y plane and a constant velocity component in the z-direction.

Notes

  • The rotation is assumed to be about the z-axis. The velocity components in the x-y plane are determined by the angular velocity omega, while the z-component is constant and set by v0.

  • The angular velocity omega is computed as: omega = 2 * np.pi * rps, with the rps given in revolutions per second (i.e. Hz).

Examples

>>> import acoular as ac
>>> import numpy as np
>>>
>>> flow = RotatingFlow(rps=1, v0=1.0)
>>> velocity, jacobian = flow.v(np.array((1.0, 1.0, 0.0)))
>>> velocity
array([-6.28318531,  6.28318531,  1.        ])
>>> jacobian
array([[ 0.        , -6.28318531,  0.        ],
       [ 6.28318531,  0.        ,  0.        ],
       [ 0.        ,  0.        ,  0.        ]])
rps = 2, #

Rotational speed in revolutions per second. Negative values indicate clockwise rigid-body-like rotation of the flow. Default is 0.0.

v0 = Float(0.0, desc='flow velocity')#

Constant flow velocity in the z-direction. Default is 0.0.

origin = CArray(dtype=np.float64, shape=(3,), value=np.array((0.0, 0.0, 0.0)), desc='center of rotation')#

The location of the center of rotation. Default is (0.0, 0.0, 0.0).

digest = Property( #

A unique identifier based on the field properties. (read-only)

omega = Property( #

Angular velocity (in radians per second) of the rotation. This is a derived property based on rps.

v(xx)#

Compute the rotating flow field and its Jacobian matrix at a given location.

This method calculates the velocity vector and its Jacobian matrix at the specified location xx in the fluid domain. The flow field consists of rotational components in the x-y plane and a constant velocity component in the z-direction.

Parameters:
xxnumpy.ndarray of floats

The 3D Cartesian coordinates of the location in the fluid where the velocity field is to be computed, shape (3,).

Returns:
velocity_vectornumpy.ndarray of floats
The velocity vector at the specified location xx, shape (3,). The components are:
  • U: Velocity in the x-direction (dependent on y-coordinate and omega).

  • V: Velocity in the y-direction (dependent on x-coordinate and omega).

  • W: Constant velocity in the z-direction (set by v0).

jacobian_matrixnumpy.ndarray of floats

The Jacobian matrix of the velocity vector field at the specified location xx. The matrix contains partial derivatives of each velocity component with respect to the spatial coordinates \((x, y, z)\), shape (3, 3).

Notes

The Jacobian matrix is constant for this flow field and represents the linear relationship between the velocity components and spatial coordinates in the x-y plane.

acoular.environments.spiral_sphere(N, Om=None, b=None)#

Generate unit vectors equally distributed over a sphere or a portion of it.

Internal helper function for the raycasting that returns an array of unit vectors of shape (N, 3) giving equally distributed directions on a part of sphere given by the center direction b and the solid angle Om.

The function uses spherical coordinates to distribute the points, the converts them to Cartesian coordinates. It also applies a transformation to reflect the points about a plane so that the direction defined by the vector b points toward the center of the sphere.

Parameters:
Nint

The number of points to generate on the sphere.

Omfloat, optional

The solid angle in steradians to cover on the sphere. Default is 2 * np.pi, which corresponds to a hemisphere. Smaller values result in covering a smaller portion of the hemisphere.

bnumpy.ndarray of floats, optional

A 3D unit vector specifying the desired center direction of the distribution. Points are mirrored such that this vector points toward the center of the sphere. Default is [0, 0, 1], which corresponds to the z-axis.

Returns:
numpy.ndarray of floats

An array of unit vectors representing points on the sphere, shape (3, N). Each column corresponds to a 3D Cartesian coordinate of a point.

Notes

  • The points are initially distributed using a spiral pattern in spherical coordinates. This ensures an approximately equal spacing between points over the specified portion of the sphere.

  • If a vector b is provided, the function mirrors the distribution using a Householder reflection so that b points toward the center.

  • The function avoids generating singularities at the poles by adjusting the spiral distribution formula.

Examples

Generate 100 points over a hemisphere:

>>> from acoular.environments import spiral_sphere
>>> points = spiral_sphere(100)
>>> points.shape
(3, 100)

Generate 50 points over half a hemisphere with the z-axis as the center direction:

>>> import numpy as np
>>> points = spiral_sphere(50, Om=np.pi, b=np.array((0, 0, 1)))
>>> points.shape
(3, 50)

Generate 200 points with a different direction vector:

>>> points = spiral_sphere(200, b=np.array((1, 0, 0)))
>>> points.shape
(3, 200)
class acoular.environments.GeneralFlowEnvironment#

Bases: Environment

An acoustic environment with a generic flow field.

This class provides the facilities to calculate the travel time (distances) between grid point locations and microphone locations in a generic flow field with non-uniform velocities that depend on the location. The algorithm for the calculation uses a ray-tracing approach that bases on rays cast from every microphone position in multiple directions and traced backwards in time. The result is interpolated within a tetrahedal grid spanned between these rays.

See also

scipy.interpolate.LinearNDInterpolator

Piecewise linear interpolator in N > 1 dimensions.

Examples

>>> import numpy as np
>>> import acoular as ac
>>>
>>> # Instantiate the flow field
>>> flow_field = ac.OpenJet(v0=10.0, D=0.4, l=3.121)
>>>
>>> # Create an instance of GeneralFlowEnvironment
>>> environment = ac.GeneralFlowEnvironment(
...     ff=flow_field,  # Use the custom flow field
...     N=300,  # Number of rays
...     Om=np.pi,  # Maximum solid angle
... )
ff = Instance(FlowField, desc='flow field')#

The flow field object describing the velocity field, which must be an instance of FlowField.

N = Int(200, desc='number of rays per Om')#

The number of rays used per solid angle \(\Omega\). Defaults to 200.

Om = Float(np.pi, desc='maximum solid angle')#

The maximum solid angle (in steradians) used in the ray-tracing algorithm. Default is numpy.pi.

digest = Property( #

A unique identifier based on the environment properties. (read-only)

idict = Dict#

A dictionary for storing precomputed interpolators to optimize repeated calculations. (internal use)

get_interpolator(roi, x0)#

Generate an interpolator for ray travel times based on a region of interest.

This method computes the ray trajectories starting from a given microphone position (x0) through a region of interest (roi). The rays’ paths are integrated numerically using a system of differential equations, and the resulting points are used to construct a convex hull. A linear interpolator is then created to estimate travel times for arbitrary points within the region.

Parameters:
roinumpy.ndarray of floats

Array representing the region of interest (ROI), where each column corresponds to a point in the 3D space \((x, y, z)\), shape (3, M).

x0numpy.ndarray of floats

Array representing the location of the microphone in 3D Cartesian coordinates, shape (3,).

Returns:
scipy.interpolate.LinearNDInterpolator object

A linear interpolator object for estimating travel times for 3D positions within the computed ray trajectories.