import radarsimpy
print("`RadarSimPy` used in this example is version: " + str(radarsimpy.__version__))
`RadarSimPy` used in this example is version: 15.2.0
Receiver Operating Characteristic (ROC) Analysis¶
ROC analysis quantifies the trade-off between detection probability $P_d$ and false alarm probability $P_{fa}$ as a function of SNR. RadarSimPy provides two tools:
roc_pd(pfa, snr, N, model)— given SNR, compute $P_d$roc_snr(pfa, pd, N, model)— given $P_d$, compute required SNR
Swerling target models supported:
| Model | RCS distribution | Fluctuation rate |
|---|---|---|
| Swerling 1 | Exponential (Rayleigh amp.) | Scan-to-scan |
| Swerling 2 | Exponential | Pulse-to-pulse |
| Swerling 3 | Chi-squared (4 DOF) | Scan-to-scan |
| Swerling 4 | Chi-squared (4 DOF) | Pulse-to-pulse |
| Swerling 5 | Constant RCS | None |
| Coherent | Ideal coherent integration | — |
| Real | Based on actual distributions | — |
N is the number of non-coherent integration pulses — more pulses improve detection but increase processing time and reduce Doppler ambiguity range.
Example 1: $P_d$ vs SNR — Swerling 3, N = 64¶
$P_{fa} \in \{10^{-4},\, 10^{-5}\}$, SNR from −10 to 10 dB.
from radarsimpy.tools import roc_pd, roc_snr
import numpy as np
import plotly.graph_objs as go
from IPython.display import Image, display
# Set to True for interactive plots; False renders a static JPEG (e.g. for HTML export)
INTERACTIVE = False
def show(fig):
if INTERACTIVE:
fig.show()
else:
display(Image(fig.to_image(format="jpg", scale=2)))
snr = np.arange(-10, 10, 0.1)
pfa = np.array([1e-4, 1e-5])
pd = roc_pd(pfa, snr, 64, "Swerling 3")
fig = go.Figure()
for i, p in enumerate(pfa):
fig.add_trace(go.Scatter(x=snr, y=pd[i, :], name=f"Pfa = {p}"))
fig.update_layout(
title="Swerling 3, N = 64",
yaxis=dict(title="Probability of detection Pd"),
xaxis=dict(title="SNR (dB)"),
height=600,
)
show(fig)
Example 2: Effect of Integration — Swerling 3, N = 20¶
Reducing from N = 64 to N = 20 lowers integration gain, shifting detection curves to the right (higher SNR required for the same $P_d$).
pd = roc_pd(pfa, snr, 20, "Swerling 3")
fig = go.Figure()
for i, p in enumerate(pfa):
fig.add_trace(go.Scatter(x=snr, y=pd[i, :], name=f"Pfa = {p}"))
fig.update_layout(
title="Swerling 3, N = 20",
yaxis=dict(title="Probability of detection Pd"),
xaxis=dict(title="SNR (dB)"),
)
show(fig)
Example 3: All Swerling Models — $P_{fa} = 10^{-6}$, N = 5¶
Coherent integration is the best-case bound. Swerling 1/2 are hardest to detect due to deep RCS fades; Swerling 5 (constant RCS) performs closest to the coherent bound.
SNR range: −10 to 30 dB.
snr = np.arange(-10, 30, 0.1)
pfa = 1e-6
N = 5
pd_sw1 = roc_pd(pfa, snr, N, "Swerling 1")
pd_sw2 = roc_pd(pfa, snr, N, "Swerling 2")
pd_sw3 = roc_pd(pfa, snr, N, "Swerling 3")
pd_sw4 = roc_pd(pfa, snr, N, "Swerling 4")
pd_sw5 = roc_pd(pfa, snr, N, "Swerling 5")
pd_coherent = roc_pd(pfa, snr, N, "Coherent")
pd_real = roc_pd(pfa, snr, N, "Real")
fig = go.Figure()
for name, pd_vals in [
("Swerling 1", pd_sw1),
("Swerling 2", pd_sw2),
("Swerling 3", pd_sw3),
("Swerling 4", pd_sw4),
("Swerling 5", pd_sw5),
("Coherent", pd_coherent),
("Real", pd_real),
]:
fig.add_trace(go.Scatter(x=snr, y=pd_vals, name=name))
fig.update_layout(
title=f"Pfa = {pfa}, N = {N}",
yaxis=dict(title="Probability of detection Pd"),
xaxis=dict(title="SNR (dB)"),
)
show(fig)
Example 4: Required SNR vs $P_d$ — $P_{fa} = 10^{-5}$, N = 10¶
Inverse calculation: given a detection requirement, how much SNR is needed? Useful for radar link budget design.
pfa = 1e-5
pd = np.arange(0.1, 0.9, 0.01)
N = 10
snr_real = roc_snr(1e-4, pd, N, "Real")
snr_sw1 = roc_snr(pfa, pd, N, "Swerling 1")
snr_sw2 = roc_snr(pfa, pd, N, "Swerling 2")
snr_sw3 = roc_snr(pfa, pd, N, "Swerling 3")
snr_sw4 = roc_snr(pfa, pd, N, "Swerling 4")
snr_sw5 = roc_snr(pfa, pd, N, "Swerling 5")
snr_coherent = roc_snr(pfa, pd, N, "Coherent")
snr_real = roc_snr(pfa, pd, N, "Real")
fig = go.Figure()
for name, snr_vals in [
("Swerling 1", snr_sw1),
("Swerling 2", snr_sw2),
("Swerling 3", snr_sw3),
("Swerling 4", snr_sw4),
("Swerling 5", snr_sw5),
("Coherent", snr_coherent),
("Real", snr_real),
]:
fig.add_trace(go.Scatter(x=pd, y=snr_vals, name=name))
fig.update_layout(
title=f"Pfa = {pfa}, N = {N}",
xaxis=dict(title="Probability of detection Pd"),
yaxis=dict(title="Minimal SNR (dB)"),
)
show(fig)
Generate SNR Lookup Table (Optional)¶
Uncomment cell below to pre-compute required SNR for N = 1…256 at $P_d = 0.5$, $P_{fa} = 10^{-4}$ and save to snr.csv. Useful for embedding in real-time systems that need fast lookups instead of calling roc_snr() at runtime.
# snr_sw1 = np.zeros(256)
# snr_sw2 = np.zeros(256)
# snr_sw3 = np.zeros(256)
# snr_sw4 = np.zeros(256)
# snr_sw5 = np.zeros(256)
# snr_ci = np.zeros(256)
# pd = 0.5
# pfa = 1e-4
# for N in range(1, 257):
# snr_sw1[N-1] = roc_snr(pfa, pd, N, 'Swerling 1')
# snr_sw2[N-1] = roc_snr(pfa, pd, N, 'Swerling 2')
# snr_sw3[N-1] = roc_snr(pfa, pd, N, 'Swerling 3')
# snr_sw4[N-1] = roc_snr(pfa, pd, N, 'Swerling 4')
# snr_sw5[N-1] = roc_snr(pfa, pd, N, 'Swerling 5')
# snr_ci[N-1] = roc_snr(pfa, pd, N, 'Coherent')
# mat = np.zeros((256, 6))
# mat[:, 0] = snr_sw1
# mat[:, 1] = snr_sw2
# mat[:, 2] = snr_sw3
# mat[:, 3] = snr_sw4
# mat[:, 4] = snr_sw5
# mat[:, 5] = snr_ci
# np.savetxt("../data/snr.csv", mat, delimiter=",")
Summary¶
roc_pd()androc_snr()let you work the detection problem in both directions.- Swerling 1/2 (Rayleigh-distributed RCS) require the most SNR; Swerling 5 and Coherent require the least.
- Increasing N improves detection at the cost of update rate and Doppler ambiguity.
Things to Try¶
| Experiment | Change | Observable effect |
|---|---|---|
| Integration sweep | N = 1, 10, 50, 100 | Observe diminishing returns |
| Stricter false alarm | $P_{fa} = 10^{-9}$ | SNR requirement rises significantly |
| Fix $P_d = 0.95$ | Use roc_snr for all models |
Quantify worst- vs best-case SNR gap |
| Generate design curve | Sweep N, fixed $P_d$/$P_{fa}$ | Integration gain vs N plot |