sources#
Measured multichannel data management and simulation of acoustic sources.
Container for processing time data in |
|
Container to process and manage time-domain data with support for masking samples and channels. |
|
Define a fixed point source emitting a signal, intended for simulations. |
|
Define a fixed point source with dipole characteristics. |
|
Define a fixed spherical harmonic source emitting a signal. |
|
Define a fixed line source with a signal. |
|
Define a moving |
|
Define a moving point source with dipole characteristics. |
|
A moving |
|
Simulate uncorrelated white or pink noise signals at multiple channels. |
|
Combine signals from multiple sources by mixing their outputs. |
|
Blockwise convolution of a source signal with an impulse response (IR). |
|
|
Compute the spherical Hankel function of the first kind. |
|
Calculate the azimuthal and elevation angles between the microphones and the source. |
|
Calculate the spherical harmonic radiation pattern at microphone positions. |
- acoular.sources.spherical_hn1(n, z)#
Compute the spherical Hankel function of the first kind.
The spherical Hankel function of the first kind, \(h_n^{(1)}(z)\), is defined as
\[h_n^{(1)}(z) = j_n(z) + i \cdot y_n(z)\]with the complex unit \(i\), the spherical Bessel function of the first kind as
\[j_n(z) = \sqrt{\frac{\pi}{2z}} J_{n + 1/2}(z),\]and the spherical Bessel function of the second kind as
\[y_n(z) = \sqrt{\frac{\pi}{2z}} Y_{n + 1/2}(z),\]where \(Y_n\) is the Bessel function of the second kind.
- Parameters:
- Returns:
- complex or
numpy.ndarray Value of the spherical Hankel function of the first kind for the given order
nand argumentz. Ifzis array-like, an array of the same shape is returned.
- complex or
See also
scipy.special.spherical_jn()Computes the spherical Bessel function of the first kind.
scipy.special.spherical_yn()Computes the spherical Bessel function of the second kind.
Notes
The function relies on
scipy.special.spherical_jn()for the spherical Bessel function of the first kind andscipy.special.spherical_yn()for the spherical Bessel function of the second kind.The input
nmust be a non-negative integer; otherwise, the behavior is undefined.
Examples
>>> import acoular as ac >>> >>> ac.sources.spherical_hn1(0, 1.0) np.complex128(0.8414709848078965-0.5403023058681398j) >>> ac.sources.spherical_hn1(1, [1.0, 2.0]) array([0.30116868-1.38177329j, 0.43539777-0.350612j ])
- acoular.sources.get_radiation_angles(direction, mpos, sourceposition)#
Calculate the azimuthal and elevation angles between the microphones and the source.
The function computes the azimuth (
azi) and elevation (ele) angles between each microphone position and the source position, taking into account the orientation of the spherical harmonics provided by the parameterdirection.- Parameters:
- direction
numpy.ndarrayof shape(3,) Unit vector representing the spherical harmonic orientation. It should be a 3-element array corresponding to the
x,y, andzcomponents of the direction.- mpos
numpy.ndarrayof shape(3, N) Microphone positions in a 3D Cartesian coordinate system. The array should have 3 rows (the
x,yandzcoordinates) andNcolumns (one for each microphone).- sourceposition
numpy.ndarrayof shape(3,) Position of the source in a 3D Cartesian coordinate system. It should be a 3-element array corresponding to the
x,y, andzcoordinates of the source.
- direction
- Returns:
- azi
numpy.ndarrayof shape(N,) Azimuth angles in radians between the microphones and the source. The range of the values is \([0, 2\pi)\).
- ele
numpy.ndarrayof shape(N,) Elevation angles in radians between the microphones and the source. The range of the values is \([0, \pi]\).
- azi
See also
numpy.linalg.norm()Computes the norm of a vector.
numpy.arctan2Computes the arctangent of two variables, preserving quadrant information.
Notes
The function accounts for a coordinate system transformation where the
z-axis in Acoular corresponds to they-axis in spherical coordinates, and they-axis in Acoular corresponds to thez-axis in spherical coordinates.The elevation angle (
ele) is adjusted to the range \([0, \pi]\) by adding \(\pi/2\) after the initial calculation.
Examples
>>> import acoular as ac >>> import numpy as np >>> >>> direction = [1, 0, 0] >>> mpos = np.array([[1, 2], [0, 0], [0, 1]]) # Two microphones >>> sourceposition = [0, 0, 0] >>> azi, ele = ac.sources.get_radiation_angles(direction, mpos, sourceposition) >>> azi array([0. , 5.8195377]) >>> ele array([4.71238898, 4.71238898])
- acoular.sources.get_modes(lOrder, direction, mpos, sourceposition=None)#
Calculate the spherical harmonic radiation pattern at microphone positions.
This function computes the spherical harmonic radiation pattern values at each microphone position for a given maximum spherical harmonic order (
lOrder), orientation (direction), and optional source position (sourceposition).- Parameters:
- lOrder
int The maximum order of spherical harmonics to compute. The resulting modes will include all orders up to and including
lOrder.- direction
numpy.ndarrayof shape(3,) Unit vector representing the orientation of the spherical harmonics. Should contain the
x,y, andzcomponents of the direction.- mpos
numpy.ndarrayof shape(3, N) Microphone positions in a 3D Cartesian coordinate system. The array should have 3 rows (the
x,yandzcoordinates) andNcolumns (one for each microphone).- sourceposition
numpy.ndarrayof shape(3,), optional Position of the source in a 3D Cartesian coordinate system. If not provided, it defaults to the origin
[0, 0, 0].
- lOrder
- Returns:
numpy.ndarrayof shape(N, (lOrder+1) ** 2)Complex values representing the spherical harmonic radiation pattern at each microphone position (
Nmicrophones) for each spherical harmonic mode.
See also
get_radiation_angles()Computes azimuth and elevation angles between microphones and the source.
scipy.special.sph_harmComputes spherical harmonic values.
Notes
The azimuth (
azi) and elevation (ele) angles between the microphones and the source are calculated using theget_radiation_angles()function.Spherical harmonics (
sph_harm) are computed for each mode(l, m), wherelis the degree (ranging from0tolOrder) andmis the order (ranging from-lto+l).For negative orders (m < 0), the conjugate of the spherical harmonic is computed and scaled by the imaginary unit
1j.
Examples
>>> import acoular as ac >>> import numpy as np >>> >>> lOrder = 2 >>> direction = [0, 0, 1] # Orientation along z-axis >>> mpos = np.array([[1, -1], [1, -1], [0, 0]]) # Two microphones >>> sourcepos = [0, 0, 0] # Source at origin >>> >>> modes = ac.sources.get_modes(lOrder, direction, mpos, sourcepos) >>> modes.shape (2, 9)
- class acoular.sources.TimeSamples#
Bases:
SamplesGeneratorContainer for processing time data in
*.h5or NumPy array format.The
TimeSamplesclass provides functionality for loading, managing, and accessing time-domain data stored in HDF5 files or directly provided as a NumPy array. This data can be accessed iteratively through theresult()method, which returns chunks of the time data for further processing.See also
MaskedTimeSamplesExtends the functionality of class
TimeSamplesby enabling the definition of start and stop samples as well as the specification of invalid channels.
Notes
Examples
Data can be loaded from a HDF5 file as follows:
>>> from acoular import TimeSamples >>> file = <some_h5_file.h5> >>> ts = TimeSamples(file=file) >>> print(f'number of channels: {ts.num_channels}') number of channels: 56
Alternatively, the time data can be specified directly as a NumPy array. In this case, the
dataandsample_freqattributes must be set manually.>>> import numpy as np >>> data = np.random.rand(1000, 4) >>> ts = TimeSamples(data=data, sample_freq=51200)
Chunks of the time data can be accessed iteratively via the
result()generator. The last block will be shorter than the block size if the number of samples is not a multiple of the block size.>>> blocksize = 512 >>> generator = ts.result(num=blocksize) >>> for block in generator: ... print(block.shape) (512, 4) (488, 4)
- file = Union(None, File(filter=['*.h5'], exists=True), desc='name of data file')#
Full path to the
.h5file containing time-domain data.
- basename = Property(depends_on=['file'], desc='basename of data file')#
Basename of the
.h5file, set automatically from thefileattribute.
- num_channels = CInt(0, desc='number of input channels')#
Number of input channels in the time data, set automatically based on the
loaded dataorspecified array.
- num_samples = CInt(0, desc='number of samples')#
Total number of time-domain samples, set automatically based on the
loaded dataorspecified array.
- data = Any(transient=True, desc='the actual time data array')#
A 2D NumPy array containing the time-domain data, shape (
num_samples,num_channels).
- h5f = Instance(H5FileBase, transient=True)#
HDF5 file object.
- metadata = Dict(desc='metadata contained in .h5 file')#
Metadata loaded from the HDF5 file, if available.
- digest = Property(depends_on=['basename', '_datachecksum', 'sample_freq', 'num_channels', 'num_samples'])#
A unique identifier for the samples, based on its properties. (read-only)
- result(num=128)#
Generate blocks of time-domain data iteratively.
The
result()method is a Python generator that yields blocks of time-domain data of the specified size. Data is either read from an HDF5 file (iffileis set) or from a NumPy array (ifdatais directly provided).- Parameters:
- num
int, optional The size of each block to be yielded, representing the number of time-domain samples per block.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) representing a block of time-domain data. The last block may have fewer thannumsamples if the total number of samples is not a multiple ofnum.
- Raises:
OSErrorIf no samples are available (i.e.,
num_samplesis0).
Examples
Create a generator and access blocks of data:
>>> import numpy as np >>> from acoular.sources import TimeSamples >>> ts = TimeSamples(data=np.random.rand(1000, 4), sample_freq=51200) >>> generator = ts.result(num=256) >>> for block in generator: ... print(block.shape) (256, 4) (256, 4) (256, 4) (232, 4)
Note that the last block may have fewer that
numsamples.
- class acoular.sources.MaskedTimeSamples#
Bases:
TimeSamplesContainer to process and manage time-domain data with support for masking samples and channels.
The
MaskedTimeSamplesclass extends the functionality ofTimeSamplesby allowing the definition ofstartandstopindices for valid samples and by supporting invalidation of specific channels. This makes it suitable for use cases where only a subset of the data is of interest, such as analyzing specific time segments or excluding faulty sensor channels.See also
TimeSamplesThe parent class for managing unmasked time-domain data.
Notes
Channels specified in
invalid_channelsare excluded from processing and not included in the generator output.Examples
Data can be loaded from a HDF5 file and invalid channels can be specified as follows:
>>> from acoular import MaskedTimeSamples >>> file = <some_h5_file.h5> >>> ts = MaskedTimeSamples(file=file, invalid_channels=[0, 1]) >>> print(f'number of valid channels: {ts.num_channels}') number of valid channels: 54
Alternatively, the time data can be specified directly as a numpy array. In this case, the
dataandsample_freqattributes must be set manually.>>> from acoular import MaskedTimeSamples >>> import numpy as np >>> data = np.random.rand(1000, 4) >>> ts = MaskedTimeSamples(data=data, sample_freq=51200)
Chunks of the time data can be accessed iteratively via the
result()generator:>>> block_size = 512 >>> generator = ts.result(num=block_size) >>> for block in generator: ... print(block.shape) (512, 4) (488, 4)
- start = CInt(0, desc='start of valid samples')#
Index of the first sample to be considered valid. Default is
0.
- stop = Union(None, CInt, desc='stop of valid samples')#
Index of the last sample to be considered valid. If
None, all remaining samples from thestartindex onward are considered valid. Default isNone.
- invalid_channels = List(int, desc='list of invalid channels')#
List of channel indices to be excluded from processing. Default is
[].
- channels = Property(depends_on=['invalid_channels', 'num_channels_total'], desc='channel mask')#
A mask or index array representing valid channels. Automatically updated based on the
invalid_channelsandnum_channels_totalattributes.
- num_channels_total = CInt(0, desc='total number of input channels')#
Total number of input channels, including invalid channels. (read-only).
- num_samples_total = CInt(0, desc='total number of samples per channel')#
Total number of samples, including invalid samples. (read-only).
- num_channels = Property( …#
Number of valid input channels after excluding
invalid_channels. (read-only)
- num_samples = Property( …#
Number of valid time-domain samples, based on
startandstopindices. (read-only)
- digest = Property(depends_on=['basename', 'start', 'stop', 'invalid_channels', '_datachecksum'])#
A unique identifier for the samples, based on its properties. (read-only)
- result(num=128)#
Generate blocks of valid time-domain data iteratively.
The
result()method is a Python generator that yields blocks of valid time-domain data based on the specifiedstartandstopindices and the valid channels.- Parameters:
- num
int, optional The size of each block to be yielded, representing the number of time-domain samples per block. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) representing a block of valid time-domain data. The last block may have fewer thannumsamples if thenumber of valid samplesis not a multiple ofnum.
- Raises:
Examples
Access valid data in blocks:
>>> import numpy as np >>> from acoular.sources import MaskedTimeSamples >>> >>> data = np.random.rand(1000, 4) >>> ts = MaskedTimeSamples(data=data, start=100, stop=900) >>> >>> generator = ts.result(num=256) >>> for block in generator: ... print(block.shape) (256, 4) (256, 4) (256, 4) (32, 4)
Note that the last block may have fewer that
numsamples.
- class acoular.sources.PointSource#
Bases:
SamplesGeneratorDefine a fixed point source emitting a signal, intended for simulations.
The
PointSourceclass models a stationary sound source that generates a signal detected by microphones. It includes support for specifying the source’s location, handling signal behaviors for pre-padding, and integrating environmental effects on sound propagation. The output is being generated via theresult()generator.See also
SignalGeneratorFor defining custom emitted signals.
MicGeomFor specifying microphone geometries.
EnvironmentFor modeling sound propagation effects.
Notes
The signal is adjusted to account for the distances between the source and microphones.
The
prepaddingattribute allows control over how the signal behaves for time indices beforestart_t.Environmental effects such as sound speed are included through the
envattribute.
Examples
To define a point source emitting a signal at a specific location, we first programmatically set a microphone geomertry as in
MicGeom:>>> import numpy as np >>> >>> # Generate a (3,3) grid of points in the x-y plane >>> x = np.linspace(-1, 1, 3) # Generate 3 points for x, from -1 to 1 >>> y = np.linspace(-1, 1, 3) # Generate 3 points for y, from -1 to 1 >>> >>> # Create a meshgrid for 3D coordinates, with z=0 for all points >>> X, Y = np.meshgrid(x, y) >>> Z = np.zeros_like(X) # Set all z-values to 0 >>> >>> # Stack the coordinates into a single (3,9) array >>> points = np.vstack([X.ravel(), Y.ravel(), Z.ravel()]) >>> points array([[-1., 0., 1., -1., 0., 1., -1., 0., 1.], [-1., -1., -1., 0., 0., 0., 1., 1., 1.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0.]])
Now, to set the actual point source (
ps), we define a microphone geomerity (mg), using the positional data frompoints, and a sine generator (sg) with a total number of 6 samples.>>> from acoular import PointSource, SineGenerator, MicGeom >>> mg = MicGeom(pos_total=points) >>> sg = SineGenerator(freq=1000, sample_freq=51200, num_samples=6) >>> ps = PointSource(signal=sg, loc=(0.5, 0.5, 1.0), mics=mg)
We choose a blocksize of 4 and generate the output signal at the microphones in blocks:
>>> for block in ps.result(num=4): ... print(block.shape) (4, 9) (2, 9)
The first block has shape (4,9) for 4 samples and 9 microphones. The second block has shape (2,9), since of a total of 6 samples only 2 remained.
- signal = Instance(SignalGenerator)#
Instance of the
SignalGeneratorclass defining the emitted signal.
- loc = Tuple((0.0, 0.0, 1.0), desc='source location')#
Coordinates
(x, y, z)of the source in a left-oriented system. Default is(0.0, 0.0, 1.0).
- num_channels = Delegate('mics', 'num_mics')#
Number of output channels, automatically set based on the
microphone geometry.
- mics = Instance(MicGeom, desc='microphone geometry')#
MicGeomobject defining the positions of the microphones.
- env = 2)#
An
Environmentor derived object providing sound propagation details, such asspeed of sound in the medium. Default isEnvironment.
- start_t = Float(0.0, desc='signal start time')#
Start time of the signal in seconds. Default is
0.0.
- start = Float(0.0, desc='sample start time')#
Start time of data acquisition at the microphones in seconds. Default is
0.0.
- prepadding = Enum('loop', 'zeros', desc='Behaviour for negative time indices.')#
Behavior of the signal for negative time indices, i.e. if (
start<start_t):'loop': Repeat thesignalfrom its end.'zeros': Use zeros, recommended for deterministic signals.
Default is
'loop'.
- up = Int(16, desc='upsampling factor')#
Internal upsampling factor for finer signal resolution. Default is
16.
- num_samples = Delegate('signal')#
Total number of samples in the emitted signal, derived from the
signalgenerator.
- sample_freq = Delegate('signal')#
Sampling frequency of the signal, derived from the
signalgenerator.
- digest = Property( …#
A unique identifier for the current state of the source, based on its properties. (read-only)
- result(num=128)#
Generate output signal at microphones in blocks, incorporating propagation effects.
The
result()method provides a generator that yields blocks of the signal detected at microphones. The signal is adjusted for the distances between the source and microphones, as well as any environmental propagation effects.- Parameters:
- num
int, optional Number of samples per block to be yielded. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the signal detected at the microphones. The last block may have fewer samples ifnum_samplesis not a multiple ofnum.
- Raises:
ValueErrorIf the source and a microphone are located at the same position.
RuntimeErrorIf signal processing or propagation cannot be performed.
- class acoular.sources.SphericalHarmonicSource#
Bases:
PointSourceDefine a fixed spherical harmonic source emitting a signal.
The
SphericalHarmonicSourceclass models a stationary sound source that emits a signal with spatial properties represented by spherical harmonics. This source can simulate directionality and orientation in sound emission, making it suitable for advanced acoustic simulations.The output is being generated via the
result()generator.- lOrder = Int(0, desc='Order of spherical harmonic') # noqa: N815#
Order of the spherical harmonic representation. Default is
0.
- alpha = CArray(desc='coefficients of the (lOrder,) spherical harmonic mode')#
Coefficients of the spherical harmonic modes for the given
lOrder.
- direction = Tuple((1.0, 0.0, 0.0), desc='Spherical Harmonic orientation')#
Vector defining the orientation of the spherical harmonic source. Default is
(1.0, 0.0, 0.0).
- prepadding = Enum('loop', desc='Behaviour for negative time indices.')#
Behavior of the signal for negative time indices. Currently only supports loop. Default is
'loop'.
- transform(signals)#
Apply spherical harmonic transformation to input signals.
The
transform()method modifies the input signals using the spherical harmonic modes, taking into account the specified coefficients (alpha), order (lOrder), and source orientation (direction).- Parameters:
- signals
numpy.ndarray Input signal array of shape (
num_samples,num_channels).
- signals
- Returns:
numpy.ndarrayTransformed signal array of the same shape as
signals.
See also
get_modes()Method for computing spherical harmonic modes.
Notes
The spherical harmonic modes are computed using the
get_modes()function, which requires the microphone positions, source position, and source orientation.The transformation applies the spherical harmonic coefficients (
alpha) to the signal in the frequency domain.
- result(num=128)#
Generate output signal at microphones in blocks, incorporating propagation effects.
The
result()method provides a generator that yields blocks of the signal detected at microphones. The signal is adjusted for the distances between the source and microphones, as well as any environmental propagation effects.- Parameters:
- num
int, optional Number of samples per block to be yielded. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the signal detected at the microphones. The last block may have fewer samples ifnum_samplesis not a multiple ofnum.
- Raises:
IndexErrorIf no more samples are available from the signal source.
- class acoular.sources.MovingPointSource#
Bases:
PointSourceDefine a moving
point sourceemitting asignal.The
MovingPointSourceclass models a sound source that follows aspecified trajectorywhile emitting asignal. This allows for the simulation of dynamic acoustic scenarios, e.g. sources changing position over time such as vehicles in motion.See also
PointSourceFor modeling stationary point sources.
TrajectoryFor specifying source motion paths.
- conv_amp = Bool(False, desc='determines if convective amplification is considered')#
Determines whether convective amplification is considered. When
True, the amplitude of the signal is adjusted based on the relative motion between the source and microphones. Default isFalse.
- trajectory = Instance(Trajectory, desc='trajectory of the source')#
Instance of the
Trajectoryclass specifying the source’s motion. The trajectory defines the source’s position and velocity at any given time.
- prepadding = Enum('loop', desc='Behaviour for negative time indices.')#
Behavior of the signal for negative time indices. Currently only supports
'loop'. Default is'loop'.
- digest = Property( …#
A unique identifier for the current state of the source, based on its properties. (read-only)
- get_moving_direction(direction, time=0)#
Calculate the moving direction of the source along its trajectory.
This method computes the updated direction vector for the moving source, considering both translation along the
trajectoryand rotation defined by thereference vector. If thereference vectoris (0, 0, 0), only translation is applied. Otherwise, the method incorporates rotation into the calculation.- Parameters:
- direction
numpy.ndarray The initial orientation of the source, specified as a three-dimensional array.
- time
float, optional The time at which the
trajectoryposition and velocity are evaluated. Defaults to0.
- direction
- Returns:
numpy.ndarrayThe updated direction vector of the moving source after translation and, if applicable, rotation. The output is a three-dimensional array.
Notes
The method computes the translation direction vector based on the
trajectory’s velocity at the specified time.If the
reference vectoris non-zero, the method constructs a rotation matrix to compute the new source direction based on thetrajectory’s motion and thereference vector.The rotation matrix ensures that the new orientation adheres to the right-hand rule and remains orthogonal.
- result(num=128)#
Generate the output signal at microphones in blocks, accounting for source motion.
The
result()method provides a generator that yields blocks of the signal received at microphones. It incorporates thesource's trajectory, convective amplification (if enabled), and environmental propagation effects.- Parameters:
- num
int, optional Number of samples per block to be yielded. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the signal detected at the microphones. The last block may have fewer samples ifnum_samplesis not a multiple ofnum.
- Raises:
IndexErrorIf no more samples are available from the signal source.
Notes
The method iteratively solves for the emission times of the signal at each microphone using the Newton-Raphson method.
Convective amplification is applied if
conv_amp= True, modifying the signal’s amplitude based on the relative motion between the source and microphones.The signal’s emission time is calculated relative to the trajectory’s position and velocity at each step.
- class acoular.sources.PointSourceDipole#
Bases:
PointSourceDefine a fixed point source with dipole characteristics.
The
PointSourceDipoleclass simulates a fixed point source with dipole characteristics by superimposing two nearby inversely phased monopoles. This is particularly useful for acoustic simulations where dipole sources are required.The generated output is available via the
result()generator.See also
PointSourceFor modeling stationary point sources.
Notes
The dipole’s output is calculated as the superposition of two monopoles: one shifted forward and the other backward along the
directionvector, with inverse phases. This creates the characteristic dipole radiation pattern.- direction = Tuple((0.0, 0.0, 1.0), desc='dipole orientation and distance of the inversely phased monopoles')#
Vector defining the orientation of the dipole lobes and the distance between the inversely phased monopoles. The magnitude of the vector determines the monopoles’ separation:
distance = [lowest wavelength in spectrum] * [magnitude] * 1e-5
Use vectors with magnitudes on the order of
1.0or smaller for best results. Default is(0.0, 0.0, 1.0)(z-axis orientation).Note: Use vectors with order of magnitude around
1.0or less for good results.
- prepadding = Enum('loop', desc='Behaviour for negative time indices.')#
Behavior of the signal for negative time indices. Currently only supports
'loop'. Default is'loop'.
- digest = Property( …#
A unique identifier for the current state of the source, based on its properties. (read-only)
- result(num=128)#
Generate output signal at microphones in blocks.
- Parameters:
- num
int, optional Number of samples per block to yield. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the signal detected at the microphones. The last block may have fewer samples ifnum_samplesis not a multiple ofnum.
- Raises:
IndexErrorIf no more samples are available from the source.
Notes
If samples are needed for times earlier than the source’s
start_t, the signal is taken from the end of the signal array, effectively looping the signal for negative indices.
- class acoular.sources.MovingPointSourceDipole#
Bases:
PointSourceDipole,MovingPointSourceDefine a moving point source with dipole characteristics.
This class extends the functionalities of
PointSourceDipoleandMovingPointSourceto simulate a dipole source that moves along adefined trajectory. It incorporates both rotational and translational dynamics for the dipole lobes, allowing simulation of complex directional sound sources.- Key Features:
Combines dipole characteristics with source motion.
Supports rotation of the dipole directivity via the
rvecattribute.Calculates emission times using Newton-Raphson iteration.
See also
PointSourceDipoleFor stationary dipole sources.
MovingPointSourceFor moving point sources without dipole characteristics.
- digest = Property( …#
A unique identifier for the current state of the source, based on its properties. (read-only)
- rvec = CArray(dtype=float, shape=(3,), value=np.array((0, 0, 0)), desc='reference vector')#
A reference vector, perpendicular to the x and y-axis of moving source, defining the axis of rotation for the dipole directivity. If set to
(0, 0, 0), the dipole is only translated along thetrajectorywithout rotation. Default is(0, 0, 0).
- get_emission_time(t, direction)#
Calculate the emission time and related properties for a moving source.
- Parameters:
- t
numpy.ndarray The current receiving time at the microphones.
- direction
floatornumpy.ndarray Direction vector for the source’s dipole directivity.
- t
- Returns:
- tuple
A tuple containing:
- te
numpy.ndarray Emission times for each microphone.
- te
- rm
numpy.ndarray Distances from the source to each microphone.
- rm
- Mr
numpy.ndarray Radial Mach numbers for the source’s motion.
- Mr
- xs
numpy.ndarray Source coordinates at the calculated emission times.
- xs
Warning
Ensure that the maximum iteration count (
100) is sufficient for convergence in all scenarios, especially for high Mach numbers or long trajectories.Notes
The emission times are computed iteratively using the Newton-Raphson method. The iteration terminates when the time discrepancy (
eps) is below a threshold (epslim) or after 100 iterations.
- result(num=128)#
Generate the output signal at microphones in blocks.
- Parameters:
- num
int, optional Number of samples per block to yield. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the signal detected at the microphones. The last block may have fewer samples ifnum_samplesis not a multiple ofnum.
Notes
Radial Mach number adjustments are applied if
conv_ampis enabled.
- class acoular.sources.LineSource#
Bases:
PointSourceDefine a fixed line source with a signal.
The
LineSourceclass models a fixed line source composed of multiple monopole sources arranged along a specified direction. Each monopole can have its own source strength, and the coherence between them can be controlled.- Key Features:
Specify the
orientation,length, andnumberof monopoles in the line source.Control the
source strengthof individual monopoles.Support for
coherent or incoherentmonopole sources.
The output signals at microphones are generated block-wise using the
result()generator.See also
PointSourceFor modeling stationary point sources.
Notes
For incoherent sources, a unique seed is set for each monopole to generate independent signals.
- direction = Tuple((0.0, 0.0, 1.0), desc='Line orientation ')#
Vector to define the orientation of the line source. Default is
(0.0, 0.0, 1.0).
- length = Float(1, desc='length of the line source')#
Vector to define the length of the line source in meters. Default is
1.0.
- num_sources = Int(1)#
Number of monopole sources in the line source. Default is
1.
- source_strength = CArray(desc='coefficients of the source strength')#
Strength coefficients for each monopole source.
- coherence = Enum('coherent', 'incoherent', desc='coherence mode')#
Coherence mode for the monopoles (
'coherent'or'incoherent').
- digest = Property( …#
A unique identifier for the current state of the source, based on its properties. (read-only)
- result(num=128)#
Generate the output signal at microphones in blocks.
- Parameters:
- num
int, optional Number of samples per block to yield. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the signal detected at the microphones. The last block may have fewer samples ifnum_samplesis not a multiple ofnum.
- class acoular.sources.MovingLineSource#
Bases:
LineSource,MovingPointSourceA moving
line sourcewith an arbitrary signal.The
MovingLineSourceclass models aline sourcecomposed of multiple monopoles that move along atrajectory. It supportscoherent and incoherentsources and considers Doppler effects due to motion.- Key Features:
Specify the
trajectoryand rotation of theline source.Compute emission times considering motion and source
direction.Generate block-wise microphone output with moving source effects.
See also
LineSourceFor
line sourcesconsisting ofcoherent or incoherentmonopoles.MovingPointSourceFor moving point sources without dipole characteristics.
- digest = Property( …#
A unique identifier for the current state of the source, based on its properties. (read-only)
- rvec = CArray(dtype=float, shape=(3,), value=np.array((0, 0, 0)), desc='reference vector')#
A reference vector, perpendicular to the x and y-axis of moving source, defining the axis of rotation for the line source directivity. If set to
(0, 0, 0), the line source is only translated along thetrajectorywithout rotation. Default is(0, 0, 0).
- get_emission_time(t, direction)#
Calculate the emission time for a moving line source based on its trajectory.
This method computes the time at which sound waves are emitted from the line source at a specific point along its
trajectory. It also determines the distances from the source to each microphone and calculates the radial Mach number, which accounts for the Doppler effect due to the motion of the source.- Parameters:
- t
float The current receiving time at the microphones, specified in seconds.
- direction
numpy.ndarray The current direction vector of the line source, specified as a 3-element array representing the orientation of the line.
- t
- Returns:
- te
numpy.ndarray The computed emission times for each microphone, specified as an array of floats.
- rm
numpy.ndarray The distances from the line source to each microphone, represented as an array of absolute distances.
- Mr
numpy.ndarray The radial Mach number, which accounts for the Doppler effect, calculated for each microphone.
- xs
numpy.ndarray The position of the line source at the computed emission time, returned as a 3-element array.
- te
Notes
This method performs Newton-Raphson iteration to find the emission time where the sound wave from the source reaches the microphones.
The distance between the line source and microphones is computed using Euclidean geometry.
The radial Mach number (
Mr) is calculated using the velocity of the source and the speed of sound in the medium (c).The method iterates until the difference between the computed emission time and the current time is sufficiently small (within a defined threshold).
- result(num=128)#
Generate the output signal at microphones in blocks.
- Parameters:
- num
int, optional Number of samples per block to yield. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the signal detected at the microphones. The last block may have fewer samples ifnum_samplesis not a multiple ofnum.
Bases:
SamplesGeneratorSimulate uncorrelated white or pink noise signals at multiple channels.
The
UncorrelatedNoiseSourceclass generates noise signals (e.g., white or pink noise) independently at each channel. It supports a user-defined random seed for reproducibility and adapts the number of channels based on the provided microphone geometry. The output is generated block-by-block through theresult()generator.See also
SignalGeneratorFor defining noise types and properties.
MicGeomFor specifying microphone geometries.
Notes
The type of noise is defined by the
signalattribute, which must be an instance of aSignalGenerator-derived class that supports aseedparameter.Each channel generates independent noise, with optional pre-defined random seeds via the
seedattribute.If no seeds are provided, they are generated automatically based on the number of channels and the signal seed.
Examples
To simulate uncorrelated white noise at multiple channels:
>>> from acoular import UncorrelatedNoiseSource, WNoiseGenerator, MicGeom >>> import numpy as np >>> >>> # Define microphone geometry >>> mic_positions = np.array([[0, 0, 0], [1, 0, 0], [0, 1, 0]]).T # Three microphones >>> mics = MicGeom(pos_total=mic_positions) >>> >>> # Define white noise generator >>> noise_gen = WNoiseGenerator(sample_freq=51200, num_samples=1024, rms=1.0, seed=42) >>> >>> # Create the noise source >>> noise_source = UncorrelatedNoiseSource(signal=noise_gen, mics=mics) >>> >>> # Generate noise output block-by-block >>> for block in noise_source.result(num=256): ... print(block.shape) (256, 3) (256, 3) (256, 3) (256, 3)
The output blocks contain noise signals for each of the 3 channels. The number of blocks depends on the total number of samples and the block size.
Instance of a
NoiseGenerator-derived class. For example: -WNoiseGeneratorfor white noise. -PNoiseGeneratorfor pink noise.
Array of random seed values for generating uncorrelated noise at each channel. If left empty, seeds will be automatically generated as
np.arange(self.num_channels) + signal.seed. The size of the array must match thenumber of output channels.
Number of output channels, automatically determined by the number of microphones defined in the
micsattribute. Corresponds to the number of uncorrelated noise signals generated.
MicGeomobject specifying the positions of microphones. This attribute is used to define the microphone geometry and thenumber of channels.
Start time of the generated noise signal in seconds. Determines the time offset for the noise output relative to the start of data acquisition. Default is
0.0.
Start time of data acquisition at the microphones in seconds. This value determines when the generated noise begins relative to the acquisition process. Default is
0.0.
Total number of samples in the noise signal, derived from the
signalgenerator. This value determines the length of the output signal for all channels.
Sampling frequency of the generated noise signal in Hz, derived from the
signalgenerator. This value defines the temporal resolution of the noise output.
A unique identifier for the current state of the source, based on its properties. (read-only)
Generate uncorrelated noise signals at microphones in blocks.
The
result()method produces a Python generator that yields blocks of noise signals generated independently for each channel. This method supports customizable block sizes and ensures that the last block may have fewer samples if the total number of samples is not an exact multiple of the block size.- Parameters:
- num
int, optional Number of samples per block to be yielded. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing uncorrelated noise signals. The last block may be shorter if the total number of samples is not a multiple ofnum.
- Raises:
ValueErrorIf the shape of the
seedarray does not match the number of channels.
Notes
Each channel’s noise signal is generated using a unique random seed.
The type and characteristics of the noise are defined by the
signalattribute.
- class acoular.sources.SourceMixer#
Bases:
SamplesGeneratorCombine signals from multiple sources by mixing their outputs.
The
SourceMixerclass takes signals generated by multipleSamplesGeneratorinstances and combines them into a single mixed output. The signals are weighted (if weights are provided) and added block-by-block, supporting efficient streaming.See also
SamplesGeneratorBase class for signal generators.
Notes
All sources must have the same sampling frequency, number of channels, and number of samples for proper mixing.
The weights for the sources can be specified to control their relative contributions to the mixed output. If no weights are provided, all sources are equally weighted.
Examples
Mix a stationary point source emitting a sine signal with two pink noise emitting point sources circling it and white noise for each channel:
>>> import numpy as np >>> import acoular as ac >>> >>> # Generate positional microphone data for a 3x3 grid in the x-y plane at z=0 >>> mic_positions = [] >>> for i in range(3): ... for j in range(3): ... mic_positions.append([i - 1, j - 1, 0]) # Center the grid at the origin >>> >>> # Convert positions to the format required by MicGeom >>> mg = ac.MicGeom(pos_total=np.array(mic_positions).T) >>> >>> # Generate positional data for trajectories of two moving sources >>> # Trajectory 1: Circle in x-y plane at z=1 >>> args = 2 * np.pi * np.arange(10) / 10 # Discrete points around the circle >>> x = np.cos(args) >>> y = np.sin(args) >>> z = np.ones_like(x) # Constant height at z=1 >>> >>> locs1 = np.array([x, y, z]) >>> # Map time indices to positions for Trajectory 1 >>> points1 = {time: tuple(pos) for time, pos in enumerate(locs1.T)} >>> tr1 = ac.Trajectory(points=points1) >>> >>> # Trajectory 2: Same circle but with a 180-degree phase shift >>> locs2 = np.roll(locs1, 5, axis=1) # Shift the positions by half the circle >>> # Map time indices to positions for Trajectory 2 >>> points2 = {time: tuple(pos) for time, pos in enumerate(locs2.T)} >>> tr2 = ac.Trajectory(points=points2) >>> >>> # Create signal sources >>> # Pink noise sources with different RMS values and random seeds >>> pinkNoise1 = ac.PNoiseGenerator(sample_freq=51200, num_samples=1024, rms=1.0, seed=42) >>> pinkNoise2 = ac.PNoiseGenerator(sample_freq=51200, num_samples=1024, rms=0.5, seed=24) >>> >>> # Moving sources emitting pink noise along their respective trajectories >>> pinkSource1 = ac.MovingPointSource(trajectory=tr1, signal=pinkNoise1, mics=mg) >>> pinkSource2 = ac.MovingPointSource(trajectory=tr2, signal=pinkNoise2, mics=mg) >>> >>> # White noise source generating uncorrelated noise for each microphone channel >>> whiteNoise = ac.WNoiseGenerator(sample_freq=51200, num_samples=1024, rms=1.0, seed=73) >>> whiteSources = ac.UncorrelatedNoiseSource(signal=whiteNoise, mics=mg) >>> >>> # Stationary point source emitting a sine wave >>> sineSignal = ac.SineGenerator(freq=1200, sample_freq=51200, num_samples=1024) >>> sineSource = ac.PointSource(signal=sineSignal, loc=(0, 0, 1), mics=mg) >>> >>> # Combine all sources in a SourceMixer with specified weights >>> sources = [pinkSource1, pinkSource2, whiteSources, sineSource] >>> mixer = ac.SourceMixer(sources=sources, weights=[1.0, 1.0, 0.3, 2.0]) >>> >>> # Generate and process the mixed output block by block >>> for block in mixer.result(num=256): # Generate blocks of 256 samples ... print(block.shape) Pink noise filter depth set to maximum possible value of 10. Pink noise filter depth set to maximum possible value of 10. (256, 9) (256, 9) (256, 9) (256, 9)
The output contains blocks of mixed signals. Each block is a combination of the four signals, weighted according to the provided weights.
- sources = List(Instance(SamplesGenerator, ()))#
List of
SamplesGeneratorinstances to be mixed. Each source provides a signal that will be combined in the output. All sources must have the same sampling frequency, number of channels, and number of samples. The list must contain at least one source.
- sample_freq = Property(depends_on=['sdigest'])#
Sampling frequency of the mixed signal in Hz. Derived automatically from the first source in
sources. If no sources are provided, default is0.
- num_channels = Property(depends_on=['sdigest'])#
Number of channels in the mixed signal. Derived automatically from the first source in
sources. If no sources are provided, default is0.
- num_samples = Property(depends_on=['sdigest'])#
Total number of samples in the mixed signal. Derived automatically from the first source in
sources. If no sources are provided, default is0.
- weights = CArray(desc='channel weights')#
Array of amplitude weights for the sources. If not set, all sources are equally weighted. The size of the weights array must match the number of sources in
sources. For example, with two sources,weights = [1.0, 0.5]would mix the first source at full amplitude and the second source at half amplitude.
- sdigest = Str()#
Internal identifier for the combined state of all sources, used to track changes in the sources for reproducibility and caching.
- digest = Property(depends_on=['sdigest', 'weights'])#
A unique identifier for the current state of the source, based on the states of the sources and the weights. (read-only)
- validate_sources()#
Ensure that all sources are compatible for mixing.
This method checks that all sources in
sourceshave the same sampling frequency, number of channels, and number of samples. AValueErroris raised if any mismatch is detected.- Raises:
ValueErrorIf any source has incompatible attributes.
- result(num)#
Generate uncorrelated the mixed signal at microphones in blocks.
The
result()method combines signals from all sources block-by-block, applying the specified weights to each source. The output blocks contain the mixed signal for all channels.- Parameters:
- num
int Number of samples per block to be yielded.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the mixed signal. The last block may have fewer samples if the total number of samples is not a multiple ofnum.
- Raises:
ValueErrorIf the sources are not compatible for mixing.
- class acoular.sources.PointSourceConvolve#
Bases:
PointSourceBlockwise convolution of a source signal with an impulse response (IR).
The
PointSourceConvolveclass extendsPointSourceto simulate the effects of sound propagation through a room or acoustic environment by convolving the input signal with a specifiedconvolution kernel(the IR).The convolution is performed block-by-block to allow efficient streaming and processing of large signals.
See also
PointSourceBase class for point sources.
TimeConvolveClass used for performing time-domain convolution.
Notes
The input
convolution kernelmust be provided as a time-domain array.The second dimension of
kernelmust either be1(a single kernel applied to all channels) or match thenumber of channelsin the output.Convolution is performed using the
TimeConvolveclass.
Examples
Convolve a stationary sine wave source with a room impulse response (RIR):
>>> import numpy as np >>> import acoular as ac >>> >>> # Define microphone geometry: 4 microphones in a 2x2 grid at z=0 >>> mic_positions = np.array([[0, 0, 0], [1, 0, 0], [0, 1, 0], [1, 1, 0]]).T >>> mg = ac.MicGeom(pos_total=mic_positions) >>> >>> # Generate a sine wave signal >>> sine_signal = ac.SineGenerator(freq=1000, sample_freq=48000, num_samples=1000) >>> >>> # Define an impulse response kernel (example: 100-tap random kernel) >>> kernel = np.random.randn(100, 1) # One kernel for all channels >>> >>> # Create the convolving source >>> convolve_source = PointSourceConvolve( ... signal=sine_signal, ... loc=(0, 0, 1), # Source located at (0, 0, 1) ... kernel=kernel, ... mics=mg, ... ) >>> >>> # Generate the convolved signal block by block >>> for block in convolve_source.result(num=256): # 256 samples per block ... print(block.shape) (256, 4) (256, 4) (256, 4) (256, 4) (75, 4)
The last block has fewer samples.
- kernel = CArray(dtype=float, desc='Convolution kernel.')#
Convolution kernel in the time domain. The array must either have one column (a single kernel applied to all channels) or match the number of output channels in its second dimension.
- start_t = Enum(0.0, desc='signal start time')#
Start time of the signal in seconds. Default is
0.0.
- start = Enum(0.0, desc='sample start time')#
Start time of the data acquisition the the microphones in seconds. Default is
0.0.
- prepadding = Enum(None, desc='Behavior for negative time indices.')#
Behavior for negative time indices. Default is
None.
- up = Enum(None, desc='upsampling factor')#
Upsampling factor for internal use. Default is
None.
- digest = Property(depends_on=['mics.digest', 'signal.digest', 'loc', 'kernel'])#
Unique identifier for the current state of the source, based on microphone geometry, input signal, source location, and kernel. (read-only)
- result(num=128)#
Generate the convolved signal at microphones in blocks.
The
result()method produces blocks of the output signal by convolving the input signal with the specified kernel. Each block contains the signal for all output channels (microphones).- Parameters:
- num
int, optional The number of samples per block to yield. Default is
128.
- num
- Yields:
numpy.ndarrayA 2D array of shape (
num,num_channels) containing the convolved signal for all microphones. The last block may contain fewer samples if the total number of samples is not a multiple ofnum.
Notes
The input signal is expanded to match the number of microphones, if necessary.
Convolution is performed using the
TimeConvolveclass to ensure efficiency.