import radarsimpy
print("`RadarSimPy` used in this example is version: " + str(radarsimpy.__version__))
`RadarSimPy` used in this example is version: 15.2.0
FMCW Radar with Corner Reflector and CFAR Detection¶
An FMCW radar transmits a linear chirp; mixing the received echo with the transmitted signal produces a beat frequency proportional to target range:
$$f_{beat} = \frac{2BR}{cT_c}, \qquad \Delta R = \frac{c}{2B}$$
Multiple chirps allow Doppler extraction via phase progression across pulses:
$$v_r = \frac{f_d \lambda}{2} = \frac{f_d c}{2f_c}$$
CFAR (Constant False Alarm Rate) detection adapts the threshold to the local noise/clutter level, maintaining a fixed $P_{fa}$ across the range-Doppler map. OS-CFAR (Ordered Statistic) uses the $k$-th sorted sample from reference cells, making it robust to outliers and nearby interfering targets.
This Example¶
- Radar: 77 GHz FMCW, 100 MHz BW (ΔR = 1.5 m), 256 chirps
- Target: Corner reflector at 50 m, approaching at −5 m/s
- Processing: Range-Doppler FFT → 2D OS-CFAR detection
- Visualization: 3D range-Doppler surface overlaid with adaptive CFAR threshold
Create Radar Model¶
Import Required Modules¶
import numpy as np
import plotly.graph_objs as go
from plotly.subplots import make_subplots
from IPython.display import Image, display
from radarsimpy import Radar, Transmitter, Receiver
# 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)))
Transmitter¶
| Parameter | Value | Notes |
|---|---|---|
| Frequency sweep | 76.95–77.05 GHz | 100 MHz BW, ΔR = 1.5 m |
| Chirp duration | 80 µs | |
| TX power | 25 dBm | |
| PRP | 100 µs | PRF = 10 kHz |
| Pulses | 256 | Δv ≈ 0.076 m/s |
# Define transmitter antenna location
tx_channel = dict(
location=(0, 0, 0), # Position at origin
)
# Configure FMCW transmitter
tx = Transmitter(
f=[77e9 - 50e6, 77e9 + 50e6], # Frequency sweep: 76.95-77.05 GHz (100 MHz BW)
t=[0, 80e-6], # Chirp duration: 0-80 μs (explicit time range)
tx_power=25, # Transmit power: 25 dBm (316 mW)
prp=100e-6, # Pulse repetition period: 100 μs
pulses=256, # Number of chirps: 256
channels=[tx_channel], # Transmitter antenna configuration
)
Receiver¶
Monostatic configuration — co-located with transmitter.
| Parameter | Value |
|---|---|
| Sampling rate | 2 MHz |
| Noise figure | 8 dB |
| RF gain | 20 dB |
| Baseband gain | 30 dB |
| Load resistor | 500 Ω |
# Define receiver antenna location
rx_channel = dict(
location=(0, 0, 0), # Co-located with transmitter (monostatic)
)
# Configure radar receiver
rx = Receiver(
fs=2e6, # Sampling rate: 2 MHz
noise_figure=8, # Noise figure: 8 dB
rf_gain=20, # RF gain: 20 dB
load_resistor=500, # Load resistance: 500 Ω
baseband_gain=30, # Baseband gain: 30 dB
channels=[rx_channel], # Receiver antenna configuration
)
Create Radar System¶
Combine transmitter and receiver to form the complete FMCW radar.
# Create complete radar system
radar = Radar(transmitter=tx, receiver=rx)
Corner Reflector Target¶
Corner reflectors provide high, predictable RCS ($\sigma \propto a^4/\lambda^2$) with wide angular tolerance — ideal for radar calibration and algorithm testing.
| Parameter | Value |
|---|---|
| Model | Trihedral corner reflector (STL) |
| Range | 50 m |
| Speed | −5 m/s (approaching) |
| Expected Doppler | $f_d = 2v f_c/c \approx -2.57$ kHz |
# Configure corner reflector target
target_1 = {
"model": "../models/cr.stl", # Trihedral corner reflector 3D model
"unit": "m", # Model units in meters
"location": (50, 0, 0), # Position: 50m along x-axis
"speed": (-5, 0, 0), # Velocity: 5 m/s approaching (closing)
}
# Create target list for simulation
targets = [target_1]
Visualize Corner Reflector¶
import pymeshlab
ms = pymeshlab.MeshSet()
ms.load_new_mesh(target_1["model"])
t_mesh = ms.current_mesh()
v_matrix = np.array(t_mesh.vertex_matrix())
f_matrix = np.array(t_mesh.face_matrix())
fig = go.Figure()
fig.add_trace(
go.Mesh3d(
x=v_matrix[:, 0],
y=v_matrix[:, 1],
z=v_matrix[:, 2],
i=f_matrix[:, 0],
j=f_matrix[:, 1],
k=f_matrix[:, 2],
color="lightsteelblue",
flatshading=False,
lighting=dict(ambient=0.4, diffuse=0.9, specular=0.3, roughness=0.5),
lightposition=dict(x=1000, y=500, z=1000),
)
)
fig.update_layout(
scene=dict(
aspectmode="data",
xaxis=dict(visible=False),
yaxis=dict(visible=False),
zaxis=dict(visible=False),
camera=dict(eye=dict(x=-1.25, y=-1.25, z=1.25)),
),
height=500,
margin=dict(l=0, r=0, b=0, t=0),
)
show(fig)
Simulate Radar Returns¶
Ray tracing density 0.2 rays/λ² — higher than typical to capture the corner reflector's triple-bounce geometry accurately. Output: complex baseband [channels, pulses, samples].
# Import radar simulator and timing module
from radarsimpy.simulator import sim_radar
import time
# Start timing
tic = time.time()
# Simulate radar returns from corner reflector
# density=0.2: Higher ray density for accurate corner reflector RCS
data = sim_radar(radar, targets, density=0.2)
# Extract baseband I/Q signals and add system noise
baseband = data["baseband"] + data["noise"] # Complex samples (I + jQ)
# End timing
toc = time.time()
# Display execution time
print("Exec time:", toc - tic, "s")
Exec time: 0.24310517311096191 s
Range-Doppler Processing¶
range_doppler_fft performs both FFTs in a single call:
- Range FFT (fast-time) → beat frequency → range bins
- Doppler FFT (slow-time) → phase progression → velocity bins
Both use Chebyshev windows with 60 dB sidelobe suppression.
# Import signal processing modules
from scipy import signal
import radarsimpy.processing as proc
# Create Chebyshev windows for sidelobe suppression (60 dB)
range_window = signal.windows.chebwin(radar.sample_prop["samples_per_pulse"], at=60)
doppler_window = signal.windows.chebwin(
radar.radar_prop["transmitter"].waveform_prop["pulses"], at=60
)
# Perform combined range-Doppler FFT
# Input: baseband [channels, pulses, samples]
# Output: range_doppler [channels, Doppler_bins, range_bins]
range_doppler = proc.range_doppler_fft(baseband, rwin=range_window, dwin=doppler_window)
Visualize Range-Doppler Map¶
Strong peak expected at ~50 m range, −5 m/s velocity.
max_range = (
3e8
* radar.radar_prop["receiver"].bb_prop["fs"]
* radar.radar_prop["transmitter"].waveform_prop["pulse_length"]
/ radar.radar_prop["transmitter"].waveform_prop["bandwidth"]
/ 2
)
unambiguous_speed = (
3e8
/ radar.radar_prop["transmitter"].waveform_prop["prp"][0]
/ 77e9
/ 2
)
range_axis = np.linspace(
0, max_range, radar.sample_prop["samples_per_pulse"], endpoint=False
)
doppler_axis = np.linspace(
-unambiguous_speed,
0,
radar.radar_prop["transmitter"].waveform_prop["pulses"],
endpoint=False,
)
fig = go.Figure()
fig.add_trace(
go.Surface(
x=range_axis,
y=doppler_axis,
z=20 * np.log10(np.abs(range_doppler[0, :, :])),
colorscale="Rainbow",
)
)
fig.update_layout(
title="Range-Doppler Map",
height=600,
scene=dict(
xaxis=dict(title="Range (m)", range=[0, 100]),
yaxis=dict(title="Velocity (m/s)"),
zaxis=dict(title="Amplitude (dB)"),
),
margin=dict(l=0, r=0, b=60, t=100),
)
show(fig)
2D OS-CFAR Detection¶
OS-CFAR uses the $k$-th ordered sample from the reference window to estimate the noise floor, making it robust to interfering targets unlike CA-CFAR.
| Parameter | Value | Purpose |
|---|---|---|
| Guard cells | 2 | Prevent target energy leaking into reference |
| Trailing cells | 20 | Reference window for noise estimation |
| $P_{fa}$ | 10⁻⁴ | 1 false alarm per 10 000 cells |
| $k$ | 1500 | Ordered statistic rank |
| Detector | linear | Power-domain detection |
# Average range-Doppler map across channels (if multiple)
rdop_avg = np.mean(np.abs(range_doppler), axis=0)
# Apply 2D OS-CFAR detection
cfar = proc.cfar_os_2d(
rdop_avg, # Input: Range-Doppler magnitude map
guard=2, # Guard cells: 2 (protect against target spillover)
trailing=20, # Trailing cells: 20 (reference for threshold)
pfa=1e-4, # Probability of false alarm: 0.0001
k=1500, # Ordered statistic rank
detector="linear" # Detector type: linear (power detection)
)
# CFAR output is adaptive threshold surface matching input dimensions
Radar Returns vs. CFAR Threshold¶
- Rainbow surface: radar returns (signal + noise)
- Viridis surface (semi-transparent): adaptive CFAR threshold
Detection occurs where the signal surface exceeds the threshold surface. The CFAR threshold tracks the local noise floor and rises near the target due to the guard region.
fig = go.Figure()
fig.add_trace(
go.Surface(
x=range_axis,
y=doppler_axis,
z=20 * np.log10(rdop_avg),
colorscale="Rainbow",
name="Radar Returns",
showscale=True,
)
)
fig.add_trace(
go.Surface(
x=range_axis,
y=doppler_axis,
z=20 * np.log10(cfar),
colorscale="Viridis",
name="CFAR Threshold",
showscale=False,
opacity=0.7,
)
)
fig.update_layout(
title="Range-Doppler Map with 2D OS-CFAR Threshold",
height=600,
scene=dict(
xaxis=dict(title="Range (m)", range=[0, 100]),
yaxis=dict(title="Velocity (m/s)"),
zaxis=dict(title="Amplitude (dB)"),
camera=dict(eye=dict(x=1.5, y=-1.5, z=1.2)),
),
margin=dict(l=0, r=0, b=60, t=100),
)
show(fig)
Summary¶
- Combined range-Doppler FFT with Chebyshev windowing (60 dB) resolves the corner reflector at 50 m / −5 m/s.
- 2D OS-CFAR adapts the detection threshold to local noise, maintaining $P_{fa} = 10^{-4}$ across the map.
- OS-CFAR's $k$-th ordered statistic makes it more robust than CA-CFAR when multiple targets or clutter edges are present.
- The detection margin (signal peak minus threshold) directly indicates target SNR.
Things to Try¶
| Experiment | Parameter to change | Observable effect |
|---|---|---|
| CA-CFAR comparison | Replace with cfar_ca_2d |
Less robust near interfering targets |
| Lower false-alarm rate | pfa=1e-6 |
Higher threshold, may miss weak targets |
| More guard cells | guard=5 |
Better target isolation, wider blind zone |
| Multiple targets | Add 2nd reflector at different range/velocity | Test masking / OS-CFAR robustness |
| Wider bandwidth | f=[77e9-250e6, 77e9+250e6] |
Finer range resolution (0.3 m) |
| More chirps | pulses=512 |
Finer Doppler resolution (0.038 m/s) |