Noise handling

The lisainstrument.noisy package provides noise definitions and noise generation tools

The definition and generation of noise is separated. The modules noise_defs and noise_defs_lisa provide a collection of basic and specialized noise types, via a unified protocol allowing code to work with noises of all types. The module noise_gen_numpy provides a simple noise generator based on numpy arrays. The streams package also contains noise generators in the streams.noise and streams.noise_alt modules. All generators use the unified noise definitions.

The noise definition interface also provides a method to analytically compute the PSD of the sampled noise as produced by the generators. This is not needed for the simulator, except for testing, but may be useful on its own.

General purpose noise types

The noisy.noise_defs module contains classes for parameters of each noise type and its PSD.

The module module does not generate noise, but provides the definitions used in the modules dedicated to noise generation. See noise_gen_numpy for a numpy-based noise stream generator, and streams.noise_alt for a generator based on the streams package.

The noise definitions can be divided in two categories: elementary noises that the noise generator can generate directly, and noises that first have to be expressed in terms of elementary noises. All noise classes have a method elementary() that returns an equivalent elementary noise. For example, the non-elementary NoiseDefRed describes filtered white noise that can be represented in terms of the elementary noises NoiseDefFiltered and NoiseDefWhite. In addition, the elementary method can apply optimizations. This is used to convert noises with zero amplitude into an NoiseDefZero, allowing noise generators to optimize away zero-amplitude noises.

All noises have a sample rate and a name, stored in the base class NoiseDefBase. The name is used when deriving random seeds from a master seed. Naming noises differently makes them uncorrelated. Noise classes used within a non-elementary noise instance are generated with a seed that also incorporates the names of all parent classes. For example, several noise classes create a NoiseDefWhite instance named “white” within their elementary() method. Yet those will be generated with different seeds, unless the names of the parent class also agree.

By design, the names of noises and the master seed of the generator fully determine the seed. The order of instance creation does not matter, neither the amount of instances. There is also no random element. The noise parameters can also enter the seed generation via optimizations. For example if a sum involving only one noise with non-zero amplitude is replaced by the nonzero component.

The intended usage is to name each noise definition uniquely (with the exception of noises used internally within other noises) to obtain uncorrelated noises when passing the definition to the generator. Computing correlated noises is not within the scope of this module. Using the same name for two noise definitions can lead to unintended results and should be avoided.

The noise classes provide a method discrete_psd_twosided() which computes the two-sided PSD of the discrete noise samples that would be generated (defined in the limit of infinite signal length). Those methods do not generate noises but rely on analytically computed PSDs. Note that the discrete PSD is a continuous function of frequency, but differs from the continuous signal PSDs that the various noises aim to approximate. In particular, the discrete PSD is periodic in frequency, with period given by the sample rate.

Instances of noise definitions can be added together. Generating noise samples from the resulting noise definition will result in a noise with PSD that is the sum of the individual PSDs. The noise is uncorrelated with the noises that would be generated for the individual noise definitions, since the seed of the sum contains the name of the sum, which is computed from the names of the summed noises. Noise definitions can also be multiplied by constant scalars. This results in a noise definition with ASD scaled by the factor. As with sums, generating noise from a scaled noise defintion is uncorrelated with noise generated for the unscaled definition.

If arithmetic combinations of the same noise realization are needed, one should instead first generate the individual noises, and then directly operate on the resulting arrays with noise samples.

The noises described by the definitions in this module refer to quantities with different units as defined in the documentation of each noise class. Some noise definitions do not refer to specific physical quantities, but represent general purpose spectral shapes. On their own, those general purpose noises are formally dimensionless, but dimensions will be re-interpreted when used in the context of implementing other dimensionful models.

Frequency parameters such as sample rates are all in [Hz]. Also the PSD methods employ frequencies in units [Hz]. Thus, PSDs of dimensionless quantities are all in units [1/Hz] and ASDs in units [1/sqrt(Hz)]. For noises of a dimensionful quantity with unit [u], the ASD unit is [u/sqrt(Hz)] and the PSD unit [u^2/Hz].

All noises in this module follow a Gaussian distribution for any individual noise sample, with zero mean.

class lisainstrument.noisy.noise_defs.NoiseDefAlpha(asd_onesided_at_1hz: float, f_min_hz: float, f_max_hz: float, alpha: float, f_sample: float, name: str)[source]

Defines colored noise with PSD approximating \(1/f^\alpha\) power law.

In the limit of frequencies small compared to the sample rate, the PSD scales according to a power law \(f^{-\alpha}\), but only within a given frequency range, outside of which the PSD approaches a constant. The exact discrete PSD can be obtained from discrete_psd_twosided(). The power law exponent is restricted to the range \(0.01 \le \alpha \le 2\).

For compatibility with previous noise generation code, the amplitude is specified in a very roundabout way: The alpha noise is defined as white noise filtered by a chain of IIR filters, each of which has a frequency response approximating the desired global power law within a small frequency segment, while approaching constant values above and below the segment. The amplitude of the noise obtained after applying all filters is scaled with a factor based on the last segment only, such that the power law model that the last filter would approximate by itself has the specified amplitude at a fiducial frequency of 1 Hz.

This definition is not very useful, in particular since the choice of segments is computed internally and thus hidden from the user. However, the exact discrete PSD can be obtained from the discrete_psd_twosided() method.

This is a general-purpose noise formally defined for dimensionless quantities. The ASD has units [1/sqrt(Hz)].

property alpha: float

Exponent of the power-law falloff for the PSD

property asd_onesided_at_1hz: float

The incomprehensibly defined amplitude [1/sqrt(Hz)]

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

The alpha noise is equivalent to white noise filtered by a chain of suitable IIR filters.

property f_max_hz: float

Frequency above which PSD approaches another constant [Hz]

property f_min_hz: float

Frequency below which PSD approaches a constant [Hz]

class lisainstrument.noisy.noise_defs.NoiseDefBase(f_sample: float, name: str)[source]

Base class for all noise definitions

This stores the sample rate and a name. The name will be used while generating noise from the definition to derive a seed for the random number generator used. It is up to the user to specify unique names to ensure uncorrelated noises. Note that defining correlations between noises is not in the scope of this noise definition class.

The sample rate is used when generating noises from the definition, and to compute the corresponding discrete PSD.

All noises definitions must implement the discrete_psd_twosided() method, providing the exact discrete PSD of the generated noise.

The class defines summation and multiplication operators. Summation of two noise definitions results in a noise definition for the sum of two uncorrelated noises. Multiplication of a noise definition with a scalar value results in a noise definition with ASD scaled by that factor.

abstract discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

abstract elementary() NoiseDefBase[source]

Return equivalent definition in terms of elementary noises

property f_sample: float

Sampling rate at which noise will be generated [Hz]

property name: str

Noise name used in generation of random seeds

class lisainstrument.noisy.noise_defs.NoiseDefCumSum(base: NoiseDefBase, name: str)[source]

Defines a noise obtained by taking the cumulative sum of another noise

Noise samples will be generated by computing the cumulative sum of noise samples generated from the base noise definition. The result is also scaled.

The discrete PSD in the limit of low frequencies is the base PSD divided by frequency squared.

Warning: depending on the base noise, this might not lead to a valid stationary noise with defined PSD.

property base: NoiseDefBase

The base noise for which to compute the cumulative sum

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

Note: we compute the PSD by using the frequncy response of an IIR filter equivalent to the cumulative sum.

elementary() NoiseDefBase[source]

Return equivalent definition in terms of elementary noises

class lisainstrument.noisy.noise_defs.NoiseDefFalloff6(asd_onesided_at_1hz: float, f_min_hz: float, f_sample: float, name: str)[source]

Defines noise with PSD approximating \(f^{-6}\) falloff

Below a cutoff frequency, the PSD falloff changes to a powerlaw \(f^{-4}\). NoiseDefFalloff6 is equivalent to NoiseDefCumSum based on NoiseDefInfrared. The amplitude parameter asd_onesided_at_1hz and the frequency cutoff f_min_Hz are passed on to NoiseDefInfrared.

Warning: it is currently unclear if NoiseDefInfrared is a valid base for NoiseDefCumSum. This definition might not result in stationary noise with well-defined PSD. In any case, the low-frequency divergence can lead to huge dynamic range even for finite length noise signals. Estimates of the PSD using Welch method might also become biased and won’t converge with increasing signal length.

property asd_onesided_at_1hz: float

Noise amplitude [1/sqrt(Hz)], see class documentation

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

property f_min_hz: float

Low frequency cutoff [Hz]

class lisainstrument.noisy.noise_defs.NoiseDefFalloff8(asd_onesided_at_1hz: float, f_min_hz: float, f_sample: float, name: str)[source]

Defines noise with PSD approximating \(f^{-8}\) falloff

Below a cutoff frequency, the PSD falloff changes to a powerlaw \(f^{-6}\). NoiseDefFalloff8 is equivalent to NoiseDefCumSum based on NoiseDefFalloff6. The amplitude parameter asd_onesided_at_1hz and the frequency cutoff f_min_Hz are passed on to NoiseDefFalloff6.

Warning: it is currently unclear if NoiseDefFalloff6 is a valid noise or at least a valid base for NoiseDefCumSum. This definition might not result in stationary noise with well-defined PSD. In any case, the low-frequency divergence can lead to huge dynamic range even for finite length noise signals. Estimates of the PSD using Welch method might also become biased and won’t converge with increasing signal length.

property asd_onesided_at_1hz: float

Noise amplitude [1/sqrt(Hz)], see class documentation

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

property f_min_hz: float

Low frequency cutoff [Hz]

class lisainstrument.noisy.noise_defs.NoiseDefFiltered(base: NoiseDefBase, iirfilters: list[DefFilterIIR], burn_in_size: int, name: str)[source]

Defines a noise obtained by applying one or more IIR filters to another noise

The discrete PSD of this noise is the one of the base noise times the squared frequency response of each IIR filter.

When generating finite amount of samples for filtered noise using infinite response filters, all samples are in principle affected by the impact of the left boundary. Depending on the filter coefficients, the impact decreases to insignificance after a certain amount of samples, denoted here as the burn-in size. This amount needs to be specified by the user together with the filter. The noise generator will generate the same amount of extra samples at the beginning, which are then discarded.

Note: not all IIR filters can be used to generate stationary noise. It is up to the user to ensure the filters are meaningful.

property base: NoiseDefBase

The definition of the base noise to be filtered

property burn_in_size: int

Number of extra samples that should be discarded to avoid edge effects

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent definition in terms of elementary noises

property iirfilters: list[DefFilterIIR]

The list of IIR filters to apply to the base noise

class lisainstrument.noisy.noise_defs.NoiseDefGradient(base: NoiseDefBase, name: str)[source]

Defines a noise obtained by finite differencing of another noise

Noise samples will be generated by applying a second order central finite difference operator to noise samples generated from the base noise definition. The result is also scaled. The discrete PSD in the limit of low frequencies is the base PSD times the frequency squared.

Note: when approaching the Nyquist frequency, the discrete PSD deviates from this simple scaling. The exact PSD can be obtained using the discrete_psd_twosided() method.

property base: NoiseDefBase

The base noise to which the finite difference should be applied

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

The formula used was obtained by direct calculation from the PSD definition applied to central finite differences, in the limit of infinitly long signals.

elementary() NoiseDefBase[source]

Return equivalent definition in terms of elementary noises

class lisainstrument.noisy.noise_defs.NoiseDefInfraRed(asd_onesided_at_1hz: float, f_min_hz: float, f_sample: float, name: str)[source]

Defines infrared noise with ASD approximating \(f^{-2}\) falloff

Below a cutoff frequency, the ASD falloff changes to a powerlaw \(f^{-1}\). NoiseDefInfraRed is equivalent to NoiseDefCumSum based on NoiseDefRed. The amplitude parameter asd_onesided_at_1hz and the frequency cutoff f_min_Hz are passed on to NoiseDefRed.

property asd_onesided_at_1hz: float

Noise amplitude [1/sqrt(Hz)], see class documentation

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

property f_min_hz: float

Low frequency cutoff [Hz]

class lisainstrument.noisy.noise_defs.NoiseDefPink(asd_onesided_at_1hz: float, f_min_hz: float, f_sample: float, name: str)[source]

Pink noise with PSD falloff \(1/f\)

This is just an alias for alpha noise with alpha=1 and f_max_hz = f_sample / 2 see NoiseDefAlpha for details.

class lisainstrument.noisy.noise_defs.NoiseDefPowerLaw(asd_onesided_at_1hz: float, asd_exp: float, f_min_hz: float, f_sample: float, name: str)[source]

Defines noise with ASD approximating \(f^\gamma\) power law.

The exact discrete noise PSD depends on the power law exponent, based on which different noise models are chosen as follows. For \(\gamma = -1\), NoiseDefPowerLaw is equivalent to NoiseDefRed. In the range \(-1 < \gamma < 0\), it is equivalent to NoiseDefAlpha with \(\alpha = -2 \gamma\). For \(\gamma = 0\), it is equivalent to NoiseDefWhite. Finally, for \(\gamma > 0\), NoiseDefPowerLaw is equivalent to NoiseDefGradient based on another NoiseDefPowerLaw definition with exponent \(\gamma - 1\).

For the case where NoiseDefAlpha is used, the upper frequency cutoff (f_max_hz) is set to the Nyquist frequency.

For the cases where NoiseDefAlpha or NoiseDefRed are used, the lower frequency cutoff parameter of these is set to the f_min_hz parameter of the NoiseDefPowerLaw. For the cases where NoiseDefWhite is used, this parameter is ignored.

The amplitude is specified by the asd_onesided_at_1hz parameter, which is used for the asd_onesided_at_1hz parameter of NoiseDefAlpha or NoiseDefRed, and the asd_onesided_const parameter of NoiseDefWhite.

Note that due to the different models used for different parameter ranges, the PSD does not depend on the exponent \(\gamma\) in a continuous fashion.

property asd_exp: float

The PSD power law exponent \(\gamma\).

property asd_onesided_at_1hz: float

Noise amplitude [1/sqrt(Hz)], see class documentation

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

property f_min_hz: float

Lower cutoff frequency (see above) [Hz]

class lisainstrument.noisy.noise_defs.NoiseDefRed(asd_onesided_at_1hz: float, f_min_hz: float, f_sample: float, name: str)[source]

Defines red (aka Brownian) noise with ASD approximating \(f^{-1}\) falloff.

In the limit of frequencies small compared to the sample rate, the PSD scales with \(f^{-2}\) (the ASD with \(f^{-1}\)), but only until a cutoff frequency, below which the PSD approaches a constant. The cutoff frequency needs to be strictly positive. The exact discrete PSD can be obtained from discrete_psd_twosided().

The noise amplitude is specified in terms of the amplitude of the above idealized scaling law at a fiducal frequency of 1 Hz. Note however that the actual discrete ASD of the noise samples at the same frequency differs from this value, depending on the sample rate. The two agree only in the limit where the sample rate is much bigger than the fiducial frequency of 1 Hz.

This is a general-purpose noise formally defined for dimensionless quantities. The ASD has units [1/sqrt(Hz)].

property asd_onesided_at_1hz: float

Noise amplitude [1/sqrt(Hz)], see class documentation

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

The red noise is equivalent to white noise filtered by a particular IIR filter

property f_min_hz: float

Cutoff frequency below which PSD approaches a constant [Hz]

class lisainstrument.noisy.noise_defs.NoiseDefRedAndViolet(asd_onesided_at_1hz: float, f_min_hz: float, f_knee_hz: float, f_sample: float, name: str)[source]

Defines a sum of red and violet noise components

This is used as base for the definitions of several LISA instrument noises.

This noise approximates a model PSD given by

\[S^\text{ffd}_\mathrm{bl}(f) = A^2 \qty[ f^2 + \frac{f_\mathrm{knee}^4}{f^2} ] = A^2 f^2 + \qty(A f_\mathrm{knee}^2)^2 f^{-2}.\]

The \(f^2\) component is implemented using NoiseDefViolet, and the \(f^{-2}\) using NoiseDefRed. Note that the latter has a low frequency cutoff. The exact resulting discrete PSD can be obtained using the discrete_psd_twosided() method.

property asd_onesided_at_1hz: float

ASD amplitude coefficient \(A\) [1/sqrt(Hz)]

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

property f_knee_hz: float

Transition frequency to red part [Hz]

property f_min_hz: float

Low frequency cutoff [Hz]

class lisainstrument.noisy.noise_defs.NoiseDefScaled(base: NoiseDefBase, factor: float, name: str)[source]

Defines a noise obtained by scaling another noise by a constant

Generating this noise is equivalent (in the statistical sense) to generating the base noise and multiply the samples by the scaling factor. Note however that a different seed value will by used, such that generating noise samples from the base definition and the scaled definition will yield uncorrelated noises.

The discrete PSD is the original one times the scaling factor squared, the ASD is the original one time the scaling factor. The ASD of the scaled noise definition depends pn the unit of the base ASD and the unit of the scaling factor. The latter is up to the user and might be dimensionless or not.

property base: NoiseDefBase

The definition of the noise to be scaled

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent definition in terms of elementary noises

property factor: float

The scale factor to apply to the noise samples

class lisainstrument.noisy.noise_defs.NoiseDefSum(components: list[NoiseDefBase], name: str)[source]

Defines a noise that is the sum of two or more uncorrelated noises.

This describes a noise definition with a PSD given by the sum of the individual noise definitions. Such a PSD describes the PSD of a noise obtained by generating uncorrelated noises with the PSDs of the individual zero-mean noise definitions and adding the noise samples. All valid noise definitions described in this module have zero mean.

The unit of the summed noise ASD is the same as the individual ASD units, which have to be all the same. The latter is not enforced, it is up to the user to make sure not to add noises with mismatched units.

Warning: when one ore more of the added noises are ill-defined, e.g. NoiseDefCumSum with an ill-suited base noise, the mean value is not zero or even defined, due to low-frequency divergence. Summing might thus worsen such problems.

property components: list[NoiseDefBase]

The definitions of the noises to be added

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent definition in terms of elementary noises

class lisainstrument.noisy.noise_defs.NoiseDefViolet(asd_onesided_at_1hz: float, f_sample: float, name: str)[source]

Defines violet noise with ASD approximating \(\propto f\)

The amplitude is specified such that in the limit of infinite sample rate, the discrete ASD at a fiducial frequency of 1 Hz is the given value. Note that the actual discrete PSD deviates from the ideal (continous) model PSD. The exact PSD can be obtained using the discrete_psd_twosided() method.

property asd_onesided_at_1hz: float

Noise amplitude [1/sqrt(Hz)], see class documentation

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

NoiseDefViolet is equivalent to NoiseDefGradient based on NoiseDefWhite.

class lisainstrument.noisy.noise_defs.NoiseDefWhite(asd_onesided_const: float, f_sample: float, name: str)[source]

Defines a Gaussian white noise

The one-sided ASD of this noise is a constant \(A\).

This is a general-purpose noise formally defined for dimensionless quantities. The ASD has units [1/sqrt(Hz)].

Individual samples of noise series generated for this definition are uncorrelated with each other and follow a Gaussian normal distribution with zero mean and standard deviation of \(\sigma = \sqrt{f_s S}\), where \(f_s\) is the sample rate and \(S=A^2/2\) is the constant two-sided PSD.

property asd_onesided_const: float

The one-sided constant ASD [1/sqrt(Hz)]

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent definition in terms of elementary noises

property psd_twosided_const: float

The two-sided constant PSD [1/Hz]

sample_distribution(x: ndarray) ndarray[source]

Evaluate the zero-mean Gaussian distribution for single noise samples

property sample_rms: float

The RMS for the Gaussian distribution of a single noise sample

class lisainstrument.noisy.noise_defs.NoiseDefZero(f_sample: float = 1.0, name: str = 'ZERO')[source]

Defines zero noise

The one-sided ASD of this noise is zero. Noise generated for this definition will simply be samples that are all zero. However, the purpose of this noise definition is to simplify logistics, not to generate zeros.

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent definition in terms of elementary noises

lisainstrument.noisy.noise_defs.asd_onesided_from_psd_twosided(psd2)[source]

Compute the one-sided amplitude spectral density (ASD) from a two-sided PSD

Parameters:

psd2 – Two-sided PSD [1/Hz]

Returns:

One-sided ASD [1/sqrt(Hz)]

lisainstrument.noisy.noise_defs.common_sample_rate(noises: list[NoiseDefBase]) float[source]

Return common sample rate of noises or raise ValueError

Parameters:

noises – list of noise definitions with identical sample rate

Returns:

Common sample rate

lisainstrument.noisy.noise_defs.iir_freq_resp(iir_filt: DefFilterIIR, f: ArrayLike, f_samp: float) np.ndarray[source]

Compute the frequency response of a IIR filter

Parameters:
  • iir_filt – The filter definition containing the coefficients

  • f – Array with frequencies for which to compute response [Hz]

  • f_samp – Sampling rate [Hz]

Returns:

Array with frequency response

Note: in principle the frequency units are arbitrary (but the same for f and f_samp) but the general convention in this module is that frequencies are in Hz.

lisainstrument.noisy.noise_defs.psd_onesided_from_asd_onesided(asd1)[source]

Compute the one-sided PSD from the one-sided ASD

Parameters:

asd1 – One-sided ASD [1/sqrt(Hz)]

Returns:

One-sided PSD [1/Hz]

lisainstrument.noisy.noise_defs.psd_onesided_from_psd_twosided(psd2)[source]

Compute the one-sided power spectral density (PSD) from a two-sided PSD

Parameters:

psd2 – Two-sided PSD [1/Hz]

Returns:

One-sided PSD [1/Hz]

lisainstrument.noisy.noise_defs.psd_twosided_from_asd_onesided(asd1)[source]

Compute two-sided PSD from one-sided ASD

Parameters:

asd1 – One-sided ASD [1/sqrt(Hz)]

Returns:

Two-sided PSD [1/Hz]

lisainstrument.noisy.noise_defs.psd_twosided_from_psd_onesided(psd1)[source]

Compute the two-sided PSD from a one-sided PSD

Parameters:

psd1 – One-sided PSD [1/Hz]

Returns:

Two-sided PSD [1/Hz]

LISA noise types

This module collects all physical noise models for the LISA Instrument

The noise definitions are based on the general-purpose noise definitions in lisainstrument.noise_defs module. This only provides definitions, to generate noise samples, see module noisy.noise_gen_numpy and stream.noise_alt

class lisainstrument.noisy.noise_defs_lisa.LaserNoiseShape(value)[source]

Enum class listing the choices for the spectral shape in NoiseDefLaser

class lisainstrument.noisy.noise_defs_lisa.NoiseDefAngularJitter(angle_asd_onesided: float, f_min_hz: float, f_knee_hz: float, f_sample: float, name: str)[source]

Defines LISA instrument jitter for one angular degree of freedom

This defines the angular jitter in terms of angular velocity, approximating a model PSD given by

\[ \begin{align}\begin{aligned}S^\text{rad/s}_\mathrm{jitter}(f) &= (2 \pi A)^2 \qty[ f^2 + \frac{f_\mathrm{knee}^4}{f^2} ]\\&= (2 \pi A)^2 f^2 + (2 \pi A f_\mathrm{knee}^2)^2 f^{-2}.\end{aligned}\end{align} \]

Generated noise samples have unit [Rad/s], and the ASD has units [(Rad/s)/sqrt(Hz)].

Note that the amplitude coefficient is given in terms of the jitter of the angle, not the angular velocity. The power spectral density in angle is given by

\[S^\text{rad}_\mathrm{jitter}(f) = A^2 \qty[ 1 + \qty(\frac{f_\mathrm{knee}}{f})^4 ],\]

which is converted to angular velocity by mutliplying by \((2 \pi f)^2\),

The noise is implemented using NoiseDefRedAndViolet as a base, see class documentation for further details.

property angle_asd_onesided: float

Amplitude coefficient \(A\) for ASD of angle [Rad/sqrt(Hz)]

Defines LISA instrument backlink noise

This noise approximates a model PSD for fractional frequency deviations given by

\[S^\text{ffd}_\mathrm{bl}(f) = \qty(\frac{2 \pi A}{c})^2 \qty[ f^2 + \frac{f_\mathrm{knee}^4}{f^2} ] = \qty(\frac{2 \pi A}{c})^2 f^2 + \qty(\frac{2 \pi A f_\mathrm{knee}^2}{c})^2 f^{-2}.\]

The exact discrete PSD can be obtained using the discrete_psd_twosided() method.

Note that the amplitude coefficient refers to the displacement amplitude and not the fractional frequency amplitude. The displacemnt amplitude is also the parameter used to create NoiseDefBacklink instances. The amplitude coefficient in terms of fractional frequency deviation will be available as a property asd_onesided_at_1hz.

The PSD in terms of the displacement is given by

\[S^\text{m}_\mathrm{bl}(f) = A^2 \qty[ 1 + \qty(\frac{f_\mathrm{knee}}{f})^4 ].\]

Multiplying by \((2 \pi f / c)^2\) yields the fractional frequency deviations PSD.

Because this is a optical pathlength noise expressed as fractional frequency deviation, it should be multiplied by the beam frequency to obtain the beam frequency fluctuations.

The noise is implemented using NoiseDefRedAndViolet as a base, see class documentation for further details.

property displacement_asd_onesided: float

ASD amplitude of displacement \(A\) [m/sqrt(Hz)]

class lisainstrument.noisy.noise_defs_lisa.NoiseDefClock(asd_onesided_at_1hz: float, f_sample: float, name: str, f_min_hz: float = 1e-05)[source]

LISA instrument clock noise

The power spectral density in fractional frequency deviations is a pink noise,

\[S^\text{ffd}_q(f) = A^2 f^{-1}.\]

Clock noise saturates below 1E-5 Hz, as the low-frequency part is modeled by deterministing clock drifts.

This model is just an alias for pink noise with f_min_hz=1e-5, see NoiseDefPink for details.

class lisainstrument.noisy.noise_defs_lisa.NoiseDefDWS(angle_asd_onesided: float, f_sample: float, name: str)[source]

Defines LISA instrument DWS measurement noise.

This defines the DWS measurement noise in terms of angular velocity, approximating a model PSD given by

\[S^\text{rad/s}_\mathrm{dws}(f) = (2 \pi A)^2 f^2.\]

Generated noise samples have unit [Rad/s], and the ASD has units [(Rad/s) / sqrt(Hz)].

Note that the amplitude coefficient \(A\) refers to the ASD of angles, not angular velocity. The PSD in angle is given by

\[S^\text{rad}_\mathrm{dws}(f) = A^2,\]

which is converted to the angular velocity PSD by mutliplying by \((2 \pi f)^2\),

Technical note: this noise class is implemented as a subclass of NoiseDefViolet, which is a dimensionless general-purpose noise. In the context of NoiseDefLongitudinalJitter, the NoiseDefViolet noise is to be considered dimensionful with same units as NoiseDefDWS. The amplitude parameter asd_onesided_at_1hz of NoiseDefViolet is set to \((2 \pi A)\) and has units of [Rad/sqrt(Hz)].

property angle_asd_onesided: float

The amplitude coefficient \(A\) for ASD of angle [Rad/sqrt(Hz)]

class lisainstrument.noisy.noise_defs_lisa.NoiseDefLaser(asd_onesided_const: float, shape: LaserNoiseShape | str, f_min_hz: float, f_sample: float, name: str, f_knee_hz: float = 0.002)[source]

Defines LISA instrument laser noise

This is a white noise with an infrared relaxation towards low frequencies, approximating the continuous PSD model

\[S^\text{Hz}_p(f) = A^2 \qty[ 1 + \qty(\frac{f_\mathrm{knee}}{f})^4 ] = A^2 + A^2 \frac{f_\mathrm{knee}^4}{f^4}.\]

The low-frequency part (infrared relaxation) can be disabled, in which case the PSD model becomes

\[S_p(f) = A^2.\]

The laser noise is defined as noise of the frequency. The ASD therefore has units of [Hz/sqrt(Hz)].

property asd_onesided_const: float

ASD amplitude coefficient \(A\) [1/sqrt(Hz)]

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

property f_knee_hz: float

Transition frequency to infrared part [Hz]

property f_min_hz: float

Low frequency cutoff [Hz]

property shape: LaserNoiseShape

Which spectral shape option to use

class lisainstrument.noisy.noise_defs_lisa.NoiseDefLongitudinalJitter(displacement_asd_onesided: float, f_sample: float, name: str)[source]

Defines LISA instrument MOSA longitudinal jitter noise

This computes the jitter in terms of the velocity jitter given by the continous one-sided PSD

\[S^\text{vel}_\mathrm{jitter}(f) = (2 \pi A)^2 f^2.\]

This noise is not dimensionless. Generated noise samples have unit of velocity, and the PSD has units of [(m/s)^2 / Hz)].

The corresponding onse-sided PSD in terms of displacement is given by white noise

\[S^\text{disp}_\mathrm{jitter}(f) = A^2,\]

which is converted to velocities by multiplying by \((2 \pi f)^2\).

Note we specify the amplitude of the velocity jitter in terms of the coefficient \(A\) that is is the amplitude of the one-sided ASD of the displacement jitter. The unit of \(A\) is thus [m/sqrt(Hz)].

Technical note: this noise class is implemented as a subclass of NoiseDefViolet, which is a dimensionless general-purpose noise. In the context of NoiseDefLongitudinalJitter, the NoiseDefViolet noise is to be considered dimensionful with same units as NoiseDefLongitudinalJitter. The amplitude parameter asd_onesided_at_1hz of NoiseDefViolet is set to \((2 \pi A)\)

property displacement_asd_onesided: float

The amplitude coefficient \(A\) [m/sqrt(Hz)]

class lisainstrument.noisy.noise_defs_lisa.NoiseDefMocTimeCorrelation(asd_onesided_const: float, f_sample: float, name: str)[source]

Defines LISA instrument MOC time correlation noise

High-level noise model for the uncertainty we have in computing the MOC time correlation (or time couples), i.e., the equivalent TCB times for the equally-sampled TPS timestamps.

Assumed to be a white noise in timing,

\[S^\text{s}_\mathrm{moc}(f) = A^2.\]

The generated noise has units of [s], and the ASD units of [s/sqrt(Hz)].

The noise is implemented as an alias for the formally dimensionless NoiseDefWhite. However, the latter is to be interpreted here as noise for a dimensionful quantity with unit [s], such that the ASD amplitude has unit [s/sqrt(Hz)].

class lisainstrument.noisy.noise_defs_lisa.NoiseDefModulation(asd_onesided_at_1hz: float, f_min_hz: float, f_sample: float, name: str)[source]

Defines LISA instrument modulation noise

The PSD model for the fractional frequency deviations reads

\[S^\text{ffd}_M(f) = A^2 f^{2/3}.\]

Note that these are fractional frequency deviations wrt. the GHz modulation frequencies.

The model is approximated using NoiseDefPowerLaw noise with asd_exp = 1/3. See NoiseDefPowerLaw for details.

class lisainstrument.noisy.noise_defs_lisa.NoiseDefOMS(displacement_asd_onesided: float, f_min_hz: float, f_knee_hz: float, f_sample: float, name: str)[source]

Defines LISA instrument optical metrology system (OMS) noise allocation

This noise approximates a model PSD for fractional frequency deviations given by

\[S^\text{ffd}_\mathrm{bl}(f) = \qty(\frac{2 \pi A}{c})^2 \qty[ f^2 + \frac{f_\mathrm{knee}^4}{f^2} ] = \qty(\frac{2 \pi A}{c})^2 f^2 + \qty(\frac{2 \pi A f_\mathrm{knee}^2}{c})^2 f^{-2}.\]

Note that the amplitude coefficient refers to the displacement amplitude and not the fractional frequency amplitude. The corresponding PSD for the displacement is given by

\[S^\text{m}_\mathrm{bl}(f) = A^2 \qty[ 1 + \qty(\frac{f_\mathrm{knee}}{f})^4 ].\]

Multiplying by \((2 \pi f / c)^2\) yields the fractional frequency deviations PSD.

Note that the level of this noise depends on the interferometer and the type of beatnote.

Warning

This corresponds to the overall allocation for the OMS noise from the Performance Model. It is a collection of different noises, some of which are duplicates of standalone noises we already implement in the simulation (e.g., backlink noise).

This noise is implemented as an alias for NoiseDefBacklink which is currently using the exact same PSD model.

class lisainstrument.noisy.noise_defs_lisa.NoiseDefRanging(asd_onesided_const: float, f_sample: float, name: str)[source]

Defines LISA instrument stochastic ranging noise

This models timing jitter as white noise

\[S^ ext{s}_R(f) = A^2.\]

The noise is implemented as an alias for NoiseDefWhite. However, the latter is to be interpreted as noise for a dimensionful quantity with unit [s], such that the ASD amplitude has unit [s/sqrt(Hz)].

class lisainstrument.noisy.noise_defs_lisa.NoiseDefTestMass(accel_asd_onesided: float, shape: TestMassNoiseShape | str, f_knee_hz: float, f_break_hz: float, f_relax_hz: float, f_min_hz: float, f_sample: float, name: str)[source]

Defines LISA instrument testmass noise

This models the noise caused by test-mass acceleration, more precisely, it describes the correponding noise of the testmass velocity. Generated noise samples have unit [m/s], and the ASD has unit [(m/s)/sqrt(Hz)].

There are two different choices for the spectral shape, specified via the TestMassNoiseShape enum class.

The first model, named “original”, models the PSD of the test mass acceleration as

\[S^\text{acc}_\delta(f) = A^2 \qty[ 1 + \qty(\frac{f_\mathrm{knee}}{f})^2] \qty[ 1 + \qty(\frac{f}{f_\mathrm{break}})^4].\]

Multiplying by \(1 / (2 \pi f)^2\) yields the PSD of the velocity,

\[ \begin{align}\begin{aligned}S^\text{vel}_\delta(f) &= \qty(\frac{A}{2 \pi})^2 \qty[ f^{-2} + \frac{f_\mathrm{knee}^2}{f^4} + \frac{f^2}{f_\mathrm{break}^4} + \frac{f_\mathrm{knee}^2}{f_\mathrm{break}^4} ]\\&= \qty(\frac{A f_\mathrm{knee}}{2 \pi})^2 f^{-4} + \qty(\frac{A}{2 \pi})^2 f^{-2} + \qty(\frac{A f_\mathrm{knee}}{2 \pi f_\mathrm{break}^2})^2 + \qty(\frac{A}{2 \pi f_\mathrm{break}^2})^2 f^2,\end{aligned}\end{align} \]

Note that we parametrize the velocity ASD amplitude using the acceleration amplitude coefficient \(A\), which has units [(m/s^2)/sqrt(Hz)].

The PSD model corresponds to the incoherent sum of an infrared, a red, a white, and a violet noise. Noise generated for NoiseDefTestMass is equivalent to the incoherent sum of noises generated for the NoiseDefInfrared, NoiseDefRed, NoiseDefWhite, and NoiseDefViolet noise definitions. The resulting actual PSD for NoiseDefTestMass approximates the above model PSD. Note however that NoiseDefInfrared and NoiseDefRed have an additional low frequency cutoff parameter f_min_hz, which has to be specified for NoiseDefTestMass as well. Below the cutoff, the PSD deviates from the model above. See NoiseDefRed for details. The resulting exact discrete PSD of noise generated for NoiseDefTestMass can be computed using the discrete_psd_twosided() method.

The second PSD model, named “lowfreq-relax”, multiplies the “original” model PSD by a factor modifying the low-frequency behavior. The PSD of the acceleration becomes

\[S^\text{acc}_\delta(f) = \ldots \times \qty[ 1 + \qty(\frac{f_\mathrm{relax}}{f})^4 ].\]

This corresponds to additional terms in the PSD of the testmass velocity noise

\[ \begin{align}\begin{aligned}S\text{vel}_\delta(f) &= \ldots \times \qty[ 1 + \qty(\frac{f_\mathrm{relax}}{f})^4 ]\\&= \ldots + \qty(\frac{A f_\mathrm{knee} f_\mathrm{relax}^2}{2 \pi})^2 f^{-8} + \qty(\frac{A f_\mathrm{relax}^2}{2 \pi})^2 f^{-6}\\&\qquad + \qty(\frac{A f_\mathrm{knee} f_\mathrm{relax}^2}{2 \pi f_\mathrm{break}^2})^2 f^{-4} + \qty(\frac{A f_\mathrm{relax}^2}{2 \pi f_\mathrm{break}^2})^2 f^{-2}.\end{aligned}\end{align} \]

Warning

This PSD strongly diverges for low frequencies, adding terms with falloff \(f^{-6}\) and \(f^{-8}\). The latter are represented by the noise definitions NoiseDefFalloff6 and NoiseDefFalloff8. Although both have a low frequency cutoff parameter f_min_hz, they still diverge below that cutoff. See NoiseDefFalloff6 for details. For both, it is unclear if they represent a valid PSD for stationary noise. In any case, including low frequencies causes huge dynamic range of the PSD. This makes it difficult to measure the PSD using the Welch method. It may also result in unphysically large test mass displacements (for which the PSD picks up another factor \(f^{-2}\)). Results obtained with the “lowfreq-relax” model should be treated with care.

property accel_asd_onesided: float

Acceleration ASD amplitude coefficient \(A\) [(m/s^2)/sqrt(Hz)]

discrete_psd_twosided(f: ndarray) ndarray[source]

Compute twosided PSD of discrete sampled data

elementary() NoiseDefBase[source]

Return equivalent noise definition using elementary noises

property f_break_hz: float

Transition frequency to TODO [Hz]

property f_knee_hz: float

Transition frequency to TODO [Hz]

property f_min_hz: float

Low frequency cutoff [Hz]

property f_relax_hz: float

Transition frequency to TODO [Hz]

property shape: TestMassNoiseShape

Which spectral shape option to use

class lisainstrument.noisy.noise_defs_lisa.TestMassNoiseShape(value)[source]

Enum class listing the choices for the spectral shape in NoiseDefTestMass

Simulator noise handling

The simulator collects the parameters of all noises via a dedicated module. This interface is also exposed to users, via the Instrument.get_noise_defs() method, allowing to obtain each noise type in a form that can be used to directly compute the PSD or generate noise samples.

The instru.instru_noises module handles the configuration of all noise parameters in the simulator

The InstruNoisesConfig class collects all noise related parameters. The InstruNoises class provides noise definitions in terms of NoiseDefBase, which can be used to generate noises.

class lisainstrument.instru.instru_noises.InstruNoises(config: InstruNoisesConfig, seed: int)[source]

Class providing the definition for any noise in the instrument simulator

The purpose is to provide a unified view of all noises. Instead of working with heterogenuous parametrizations of the different noises, InstruNoise instances provide a method defining each noise via the NoiseDefBase protocol. The latter can be used directly to generate noise with any of the noise generators in the lisainstrument.noisy or lisainstrument.streams packages. Besides the noise definitions, this also stores a master seed value.

modified(**args) InstruNoises[source]

Return another instance with same seed but selected parameters modified

Definition of backlink noise used for given MOSA

noise_def_clock(sc: SatID) NoiseDefBase[source]

Definition of clock noise used for given spacecraft

noise_def_dws_eta(mosa: MosaID) NoiseDefBase[source]

Definition of DWS eta angle noise used for given MOSA

noise_def_dws_phi(mosa: MosaID) NoiseDefBase[source]

Definition of DWS phi angle noise used for given MOSA

noise_def_laser(mosa: MosaID) NoiseDefBase[source]

Definition of laser noise used for given MOSA

noise_def_moc_time_correlation(sc: SatID) NoiseDefBase[source]

Definition of MOC time correlation noise used for given spacecraft

noise_def_modulation(mosa: MosaID) NoiseDefBase[source]

Definition of modulation noise used for given MOSA

noise_def_mosa_jitter_etas(mosa: MosaID) NoiseDefBase[source]

Definition of eta-angle noise used for given MOSA

noise_def_mosa_jitter_phis(mosa: MosaID) NoiseDefBase[source]

Definition of phi-angle noise used for given MOSA

noise_def_mosa_jitter_xs(mosa: MosaID) NoiseDefBase[source]

Definition of longitudinal noise used for given MOSA

noise_def_oms_ref_carrier(mosa: MosaID) NoiseDefBase[source]

Definition of OMS RFI carrier noise used for given MOSA

noise_def_oms_ref_usb(mosa: MosaID) NoiseDefBase[source]

Definition of OMS RFI usb noise used for given MOSA

noise_def_oms_sci_carrier(mosa: MosaID) NoiseDefBase[source]

Definition of OMS ISI carrier noise used for given MOSA

noise_def_oms_sci_usb(mosa: MosaID) NoiseDefBase[source]

Definition of OMS ISI usb noise used for given MOSA

noise_def_oms_tmi_carrier(mosa: MosaID) NoiseDefBase[source]

Definition of OMS TMI carrier noise used for given MOSA

noise_def_oms_tmi_usb(mosa: MosaID) NoiseDefBase[source]

Definition of OMS TMI usb noise used for given MOSA

noise_def_ranging(mosa: MosaID) NoiseDefBase[source]

Definition of ranging noise used for given MOSA

noise_def_sc_jitter_etas(sc: SatID) NoiseDefBase[source]

Definition of spacecraft eta-angle noise used for given spacecraft

noise_def_sc_jitter_phis(sc: SatID) NoiseDefBase[source]

Definition of spacecraft phi-angle noise used for given spacecraft

noise_def_sc_jitter_thetas(sc: SatID) NoiseDefBase[source]

Definition of spacecraft theta-angle noise used for given spacecraft

noise_def_testmass(mosa: MosaID) NoiseDefBase[source]

Definition of testmass noise used for given MOSA

property seed: int

Master seed value

class lisainstrument.instru.instru_noises.InstruNoisesConfig(physics_fs: float, telemetry_fs: float, noises_fmin: float, laser_asds: dict[MosaID, float], laser_shape: LaserNoiseShape, modulation_asds: dict[MosaID, float], ranging_asds: dict[MosaID, float], backlink_asds: dict[MosaID, float], backlink_fknees: dict[MosaID, float], testmass_asds: dict[MosaID, float], testmass_fknees: dict[MosaID, float], testmass_fbreak: dict[MosaID, float], testmass_frelax: dict[MosaID, float], testmass_shape: TestMassNoiseShape, oms_sci_carrier_asds: dict[MosaID, float], oms_sci_usb_asds: dict[MosaID, float], oms_tmi_carrier_asds: dict[MosaID, float], oms_tmi_usb_asds: dict[MosaID, float], oms_ref_carrier_asds: dict[MosaID, float], oms_ref_usb_asds: dict[MosaID, float], oms_fknees: dict[MosaID, float], mosa_jitter_x_asds: dict[MosaID, float], mosa_jitter_phi_asds: dict[MosaID, float], mosa_jitter_eta_asds: dict[MosaID, float], mosa_jitter_phi_fknees: dict[MosaID, float], mosa_jitter_eta_fknees: dict[MosaID, float], dws_asds: dict[MosaID, float], moc_time_correlation_asds: dict[SatID, float], clock_asds: dict[SatID, float], clock_fmin: float, sc_jitter_phi_asds: dict[SatID, float], sc_jitter_eta_asds: dict[SatID, float], sc_jitter_theta_asds: dict[SatID, float], sc_jitter_phi_fknees: dict[SatID, float], sc_jitter_eta_fknees: dict[SatID, float], sc_jitter_theta_fknees: dict[SatID, float])[source]

Immutable data class collecting all parameters needed to generate the noises

Simple noise generation

Besides internal use within the simulator, the noisy package also provides a user-facing interface for generating noise samples for any of the noise types.

The lisainstrument.noisy.noise_gen_numpy module provides numpy-based noise generation

Given a noise definition of type NoiseDefBase, create a NoiseGen instance, which provides a corresponding noise generation function. This noise generator has internal state and can be called repeatedly to obtain chunks of noise samples.

Internally, the code is organized around stateless generators that require the state as explicit parameters. This is suitable for use in the streams package. There is a class for each elementary noise definition type, implementing the noise generation. Those classes implement a protocol NoiseGenStateless. By elementary noises we refer to noises that can be part of noise definitions returned by calls to the elementary() method of any noise.

class lisainstrument.noisy.noise_gen_numpy.NoiseGen(noisedef: NoiseDefBase, seed: int | str)[source]

Noise generator with internal state for repeated calls

This is the user-friendly wrapper of the stateless NoiseGenStateless noise generator.

The sequence produced by the generator is fully determined by the combination of the seed value, the name of the noise, the type of the noise, and to some extend on the noise parameters. The intended use case is for generating uncorrelated noises, providing a master seed to all generated noises, each of which has a definition with a unique name. When generating two or more realizations of the same noise, create a new noise definition with a unique name for each.

We stress that different parameters alone, with everything else fixed, do not necessarily lead to uncorrelated noises. This only happens when the amplitude of a noise or a internally generated noise component becomes zero, triggering some optimizations.

generate(size: int) ndarray[source]

Generate an array with noise samples.

Concatenating the results of repeated calls is guaranteed to be the same as the result of a single call with the combined size.

Parameters:

size – number of noise samples to generate

Returns:

Numpy 1D array with nosie samples