Tools

Useful tools for radar system analysis

This script requires that ‘numpy’ and ‘scipy’ be installed within the Python environment you are running this script in.

This file can be imported as a module and contains the following functions:

  • roc_pd - Calculate probability of detection (Pd) in receiver operating

    characteristic (ROC)

  • roc_snr - Calculate the minimal SNR for certain probability of

    detection (Pd) and probability of false alarm (Pfa) in receiver operating characteristic (ROC)

██████╗  █████╗ ██████╗  █████╗ ██████╗ ███████╗██╗███╗   ███╗██╗  ██╗
██╔══██╗██╔══██╗██╔══██╗██╔══██╗██╔══██╗██╔════╝██║████╗ ████║╚██╗██╔╝
██████╔╝███████║██║  ██║███████║██████╔╝███████╗██║██╔████╔██║ ╚███╔╝ 
██╔══██╗██╔══██║██║  ██║██╔══██║██╔══██╗╚════██║██║██║╚██╔╝██║ ██╔██╗ 
██║  ██║██║  ██║██████╔╝██║  ██║██║  ██║███████║██║██║ ╚═╝ ██║██╔╝ ██╗
╚═╝  ╚═╝╚═╝  ╚═╝╚═════╝ ╚═╝  ╚═╝╚═╝  ╚═╝╚══════╝╚═╝╚═╝     ╚═╝╚═╝  ╚═╝
radarsimpy.tools.marcumq(a, x, m=1)[source]

Calculates the generalized Marcum Q function.

The Marcum Q function is defined as:

Q_m(a, x) = 1 - F_ncx2(m * 2, a^2, x^2)

Parameters:
  • a (float) – Non-centrality parameter.

  • x (float) – Threshold value.

  • m (int) – Order of the function, positive integer (default is 1).

Returns:

Generalized Marcum Q function value.

Return type:

float

References:
radarsimpy.tools.pd_swerling0(npulses, snr, thred)[source]

Calculates the probability of detection (Pd) for Swerling 0 target model.

Parameters:
  • npulses (int) – Number of pulses.

  • snr (float) – Signal-to-noise ratio.

  • thred (float) – Detection threshold.

Returns:

Probability of detection (Pd).

Return type:

float

Notes:
  • For npulses <= 50, uses the Marcum Q function and modified Bessel functions.

  • For npulses > 50, employs an approximation based on statistical parameters.

References:
  • Swerling, P. (1953). Probability of Detection for Fluctuating Targets. IRE Transactions on Information Theory, 6(3), 269-308.

radarsimpy.tools.pd_swerling1(npulses, snr, thred)[source]

Calculates the probability of detection (Pd) for Swerling 1 target model.

Parameters:
  • npulses (int) – Number of pulses.

  • snr (float) – Signal-to-noise ratio.

  • thred (float) – Detection threshold.

Returns:

Probability of detection (Pd).

Return type:

float

Notes:
  • Swerling 1 assumes a target made up of many independent scatterers of roughly equal areas.

  • The RCS varies according to a chi-squared probability density function with two degrees

    of freedom (m = 1).

  • The radar cross section is constant from pulse-to-pulse but varies independently from

    scan to scan.

References:
  • Swerling, P. (1953). Probability of Detection for Fluctuating Targets. IRE Transactions on Information Theory, 6(3), 269-308.

radarsimpy.tools.pd_swerling2(npulses, snr, thred)[source]

Calculates the probability of detection (Pd) for Swerling 2 target model.

Parameters:
  • npulses (int) – Number of pulses.

  • snr (float) – Signal-to-noise ratio.

  • thred (float) – Detection threshold.

Returns:

Probability of detection (Pd).

Return type:

float

Notes:
  • Swerling 2 assumes a target made up of many independent scatterers of roughly equal areas.

  • The radar cross section (RCS) varies from pulse to pulse.

  • Statistics follow a chi-squared probability density function with two degrees of freedom.

References:
  • Swerling, P. (1953). Probability of Detection for Fluctuating Targets. IRE Transactions on Information Theory, 6(3), 269-308.

radarsimpy.tools.pd_swerling3(npulses, snr, thred)[source]

Calculates the probability of detection (Pd) for Swerling 3 target model.

Parameters:
  • npulses (int) – Number of pulses.

  • snr (float) – Signal-to-noise ratio.

  • thred (float) – Detection threshold.

Returns:

Probability of detection (Pd).

Return type:

float

Notes:
  • Swerling 3 assumes a target made up of one dominant isotropic reflector superimposed

    by several small reflectors.

  • The radar cross section (RCS) varies from pulse to pulse but remains constant within

    a single scan.

  • The statistical properties follow a density of probability based on the Chi-squared

    distribution with four degrees of freedom (m = 2).

References:
  • Swerling, P. (1953). Probability of Detection for Fluctuating Targets. IRE Transactions on Information Theory, 6(3), 269-308.

radarsimpy.tools.pd_swerling4(npulses, snr, thred)[source]

Calculates the probability of detection (Pd) for Swerling 4 target model.

Parameters:
  • npulses (int) – Number of pulses.

  • snr (float) – Signal-to-noise ratio.

  • thred (float) – Detection threshold.

Returns:

Probability of detection (Pd).

Return type:

float

Notes:
  • Swerling 4 assumes a target made up of one dominant isotropic reflector

    superimposed by several small reflectors.

  • The radar cross section (RCS) varies from pulse to pulse rather than from scan to scan.

  • The statistical properties follow a density of probability based on the Chi-squared

    distribution with four degrees of freedom (m = 2).

References:
  • Swerling, P. (1953). Probability of Detection for Fluctuating Targets. IRE Transactions on Information Theory, 6(3), 269-308.

radarsimpy.tools.roc_pd(pfa, snr, npulses=1, stype='Coherent')[source]

Calculate probability of detection (Pd) in receiver operating characteristic (ROC)

Parameters:
  • pfa (float or numpy.1darray) – Probability of false alarm (Pfa)

  • snr (float or numpy.1darray) – Signal to noise ratio in decibel (dB)

  • npulses (int) – Number of pulses for integration (default is 1)

  • stype (str) –

    Signal type (default is Coherent)

    • Coherent: Non-fluctuating coherent

    • Real: Non-fluctuating real signal

    • Swerling 0: Non-coherent Swerling 0, Non-fluctuating non-coherent

    • Swerling 1: Non-coherent Swerling 1

    • Swerling 2: Non-coherent Swerling 2

    • Swerling 3: Non-coherent Swerling 3

    • Swerling 4: Non-coherent Swerling 4

    • Swerling 5: Non-coherent Swerling 5, Non-fluctuating non-coherent

Returns:

probability of detection (Pd). if both pfa and snr are floats, pd is a float if pfa or snr is a 1-D array, pd is a 1-D array if both pfa and snr are 1-D arrays, pd is a 2-D array

Return type:

float or 1-D array or 2-D array

Reference

Mahafza, Bassem R. Radar systems analysis and design using MATLAB. Chapman and Hall/CRC, 2005.

radarsimpy.tools.roc_snr(pfa, pd, npulses=1, stype='Coherent')[source]

Calculate the minimal SNR for certain probability of detection (Pd) and probability of false alarm (Pfa) in receiver operating characteristic (ROC) with Secant method

Parameters:
  • pfa (float or numpy.1darray) – Probability of false alarm (Pfa)

  • pd (float or numpy.1darray) – Probability of detection (Pd)

  • npulses (int) – Number of pulses for integration (default is 1)

  • stype (str) –

    Signal type (default is Coherent)

    • Coherent : Non-fluctuating coherent

    • Real : Non-fluctuating real signal

    • Swerling 0 : Non-fluctuating non-coherent

    • Swerling 1 : Non-coherent Swerling 1

    • Swerling 2 : Non-coherent Swerling 2

    • Swerling 3 : Non-coherent Swerling 3

    • Swerling 4 : Non-coherent Swerling 4

    • Swerling 5 : Same as Swerling 0

Returns:

Minimal signal to noise ratio in decibel (dB) if both pfa and pd are floats, SNR is a float if pfa or pd is a 1-D array, SNR is a 1-D array if both pfa and pd are 1-D arrays, SNR is a 2-D array

Return type:

float or 1-D array or 2-D array

Reference

Secant method:

The x intercept of the secant line on the the Nth interval

\[m_n = a_n - f(a_n)*(b_n - a_n)/(f(b_n) - f(a_n))\]

The initial interval [a_0,b_0] is given by [a,b]. If f(m_n) == 0 for some intercept m_n then the function returns this solution. If all signs of values f(a_n), f(b_n) and f(m_n) are the same at any iterations, the secant method fails and return None.