End-to-End CLI Example
This notebook demonstrates the full SERVAL pipeline from sky model construction to visibility generation using the command-line interface (CLI):
- Build a point-source sky model and save it to a Zarr store.
- Run
serval gencache fixed-skyto pre-compute the integrator cache. - Construct voltage beams and derive power beams, saving them to Zarr stores.
- Run
serval genvis from-voltage-beamsto generate visibilities (power beams derived on-the-fly). - Run
serval genvis from-power-beamsto generate visibilities. - Inspect the power-beam visibilities.
- Compare both visibility outputs and verify they are numerically identical.
import subprocess
import sys
from pathlib import Path
import astropy.coordinates as coords
import astropy.units as units
import matplotlib as mpl
import matplotlib.pyplot as plt
import numpy as np
import zarr
from astropy.stats import gaussian_fwhm_to_sigma
from astropy.time import Time
from serval.containers import TIRSPowerBeam, SkyModel, TIRSVoltageBeam, Visibilities
from serval.core import PowerBeamVisProjector, baseline_bandlimits
plt.style.use("./example_style.mplstyle")
mpl.rcParams["text.usetex"] = False
_ = np.seterr(divide="ignore")
Configuration
All parameters are defined here. We use a small number of frequency channels and modest harmonic bandlimits so the notebook runs quickly.
# HIRAX site coordinates
latitude_str = "-30d41m47.0s"
longitude_str = "21d34m20.0s"
latitude_rad = coords.Latitude(latitude_str).rad
longitude_rad = coords.Longitude(longitude_str).rad
# Frequency channels
frequencies_MHz = np.array([400.0, 500.0, 600.0, 700.0, 800.0])
# Harmonic bandlimits for voltage and power beams, very small for this example.
voltage_beam_lmax = 15
power_beam_lmax = 2 * voltage_beam_lmax
power_beam_mmax = power_beam_lmax
# HIRAX-like baselines in ENU (metres)
baselines_enu = [
(6.5, 0.0, 0.0), # East-West
(0.0, 8.5, 0.0), # North-South
(6.0, 8.5, 0.0),
]
# Pointing: zenith at nominal HIRAX position
pointing_altitude = np.pi / 2
pointing_azimuth = np.pi
# Dish diameter for the Gaussian beam model
D_dish = 6 * units.m
# Output data directory
DATA_DIR = Path("data/e2e_cli_example")
DATA_DIR.mkdir(parents=True, exist_ok=True)
# Baseline tag helper (matches CLI default --baseline-tag-ndigits=3)
def bl_tag(e, n, u, ndigits=3):
def fmt(x):
return f"{x:.{ndigits}f}".replace(".", "p")
return f"bl_E_{fmt(e)}m_N_{fmt(n)}m_U_{fmt(u)}m"
bl_tags = [bl_tag(*b) for b in baselines_enu]
print("Baseline tags:", bl_tags)
Baseline tags: ['bl_E_6p500m_N_0p000m_U_0p000m', 'bl_E_0p000m_N_8p500m_U_0p000m', 'bl_E_6p000m_N_8p500m_U_0p000m']
Step 1: Create a Point Source Sky Model
We place three fictitious point sources at ICRS positions with declinations near \(-30°\), which is the HIRAX latitude, so each source will transit close to zenith.
The sky model must have lmax at least as large as
\(\ell_\mathrm{beam} + \ell_\mathrm{baseline}\), where \(\ell_\mathrm{baseline}\)
is determined by the longest baseline at the highest frequency. We compute this
programmatically so the notebook is self-consistent.
# Determine the required sky lmax: power_beam_lmax + max(baseline_lmax)
max_baseline_lmax = max(
baseline_bandlimits(bl, frequencies_MHz.max(), latitude_rad, longitude_rad)[0]
for bl in baselines_enu
)
sky_lmax = power_beam_lmax + max_baseline_lmax
print(f"baseline_lmax (longest baseline at max freq): {max_baseline_lmax}")
print(f"sky_lmax required for integrator cache: {sky_lmax}")
baseline_lmax (longest baseline at max freq): 227
sky_lmax required for integrator cache: 257
# Three point sources near dec = -30 deg (close to HIRAX zenith)
src_coords = coords.SkyCoord(
ra = [ 45.0, 90.0, 135.0] * units.deg,
dec = [-30.0, -35.0, -25.0] * units.deg,
frame="icrs",
)
# Flat-spectrum, 1 Jy each
flux_Jy = np.ones((len(frequencies_MHz), len(src_coords)))
sky_model = SkyModel.from_point_source_catalog(
coords=src_coords,
flux_Jy=flux_Jy,
frequencies_MHz=frequencies_MHz,
lmax=sky_lmax,
coordinate_basis="ICRS",
n_jobs=1,
)
sky_zarr = DATA_DIR / "sky_model.zarr"
sky_model.to_zarr_store(sky_zarr)
print(f"Sky model saved to {sky_zarr}")
print(f"Sky alm shape: {sky_model.alm.shape} (n_freq, lmax+1, 2*lmax+1)")
Sky model saved to data/e2e_cli_example/sky_model.zarr
Sky alm shape: (5, 258, 515) (n_freq, lmax+1, 2*lmax+1)
Step 2: Generate the Integrator Cache (CLI)
We write a TOML configuration file and invoke
serval gencache fixed-sky as a subprocess. The command pre-integrates the sky and
baseline terms for each baseline and stores the result in a Zarr cache.
import tomlkit
cache_zarr = DATA_DIR / "integrator_cache.zarr"
gencache_toml = DATA_DIR / "gencache_config.toml"
gencache_cfg = tomlkit.document()
gencache_cfg.add("store_location", str(cache_zarr))
gencache_cfg.add("sky_model", str(sky_zarr))
gencache_cfg.add("latitude", latitude_str)
gencache_cfg.add("longitude", longitude_str)
gencache_cfg.add(
"baselines_enu",
[[f"{e}m", f"{n}m", f"{u}m"] for e, n, u in baselines_enu],
)
gencache_cfg.add("power_beam_lmax", power_beam_lmax)
gencache_cfg.add(
"frequencies",
[f"{f}MHz" for f in frequencies_MHz],
)
gencache_cfg.add("memory_per_task", "2 GiB")
gencache_cfg.add("sky_pol", "I")
gencache_toml.write_text(tomlkit.dumps(gencache_cfg))
print(tomlkit.dumps(gencache_cfg))
store_location = "data/e2e_cli_example/integrator_cache.zarr"
sky_model = "data/e2e_cli_example/sky_model.zarr"
latitude = "-30d41m47.0s"
longitude = "21d34m20.0s"
baselines_enu = [["6.5m", "0.0m", "0.0m"], ["0.0m", "8.5m", "0.0m"], ["6.0m", "8.5m", "0.0m"]]
power_beam_lmax = 30
frequencies = ["400.0MHz", "500.0MHz", "600.0MHz", "700.0MHz", "800.0MHz"]
memory_per_task = "2 GiB"
sky_pol = "I"
result = subprocess.run(
[
"serval", "--dist-mode", "serial",
"gencache", "fixed-sky",
"--config", str(gencache_toml),
],
capture_output=True,
text=True,
)
print(result.stdout)
if result.returncode != 0:
print("STDERR:", result.stderr)
raise RuntimeError(f"gencache failed (return code {result.returncode})")
print("Integrator cache generated successfully.")
Welcome to:
____ _____ ______ ___ _
/ ___|| ____| _ \ \ / / \ | |
\___ \| _| | |_) \ \ / / _ \ | |
___) | |___| _ < \ V / ___ \| |___
|____/|_____|_| \_\ \_/_/ \_\_____| v1.0.1b0
----------------------------------------------------------------
Distribution mode: serial (1 task).
----------------------------------------------------------------
Generating a fixed-sky integration cache...
----------------------------------------------------------------
Parsed, converted and merged parameters:
----------------------------------------------------------------
store_location = "data/e2e_cli_example/integrator_cache.zarr"
sky_model = "data/e2e_cli_example/sky_model.zarr"
baselines_enu = [
[6.5, 0.0, 0.0], [0.0, 8.5, 0.0], [6.0, 8.5, 0.0],
]
latitude = -0.5357530545837149
longitude = 0.3765063047496657
power_beam_lmax = 30
frequencies = [
400.0, 500.0, 600.0,
700.0, 800.0,
]
memory_per_task = 2.0
max_absms_per_chunk = 128
pointed_beam_mmax = 40
pointing_altitude = 1.5707963267948966
pointing_azimuth = 3.141592653589793
pointing_boresight = 0.0
pointing_contract = false
optimized_cache = true
start_chunk = 0
sky_pol = "I"
baseline_tag_ndigits = 3
----------------------------------------------------------------
Task Distribution Plan:
----------------------------------------------------------------
Split into 13 chunks across 1 task(s).
Finest possible chunking requires at least 0.08 GiB per chunk.
Maximium sky lmax used: 257
Chunk Statistics:
ms/chunk mem/chunk [GB] chunks/task
mean 51.23 1.83 13.00
std 25.68 0.35 NaN
min 24.00 0.75 13.00
25% 31.00 1.93 13.00
50% 42.00 1.97 13.00
75% 67.00 1.98 13.00
max 108.00 1.99 13.00
----------------------------------------------------------------
Per Polarisation & Baseline Statistics:
chunk frac. [%] chunks ms mem [GB]
Sky pol Baseline ENU [m]
I (0.0, 8.5, 0.0) 30.77 4 223 7.76
(6.0, 8.5, 0.0) 46.15 6 258 11.36
(6.5, 0.0, 0.0) 23.08 3 185 4.72
----------------------------------------------------------------
Plan generated in 0.05s.
[2026-03-15 08:28:42] [Task: 1/1] Task initialized with OMP_NUM_THREADS=unset
[2026-03-15 08:28:42] [Task: 1/1] Starting on chunk index 0 (1/13 for this task).
[2026-03-15 08:28:42] [Task: 1/1] Starting pol: I | baseline E: 6.5 m, N: 0.0 m, U: 0.0 m | Sky lmax: 184, baseline lmax: 154.
[2026-03-15 08:28:42] [Task: 1/1] Working on 83 m's. |m| ∈ [0->42).
[2026-03-15 08:28:42] [Task: 1/1] Computing 0.21G Gaunts (1.69 GB) took: 0.43s.
[2026-03-15 08:28:43] [Task: 1/1] Batch integrator: 0.43 s (at ~4.94 GFLOPS) | IO (0.013 GB): 0.09s
[2026-03-15 08:28:43] [Task: 1/1] Finished global chunk index 0 (local: 1/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:43] [Task: 1/1] Starting on chunk index 1 (2/13 for this task).
[2026-03-15 08:28:43] [Task: 1/1] Working on 134 m's. |m| ∈ [42->109).
[2026-03-15 08:28:43] [Task: 1/1] Computing 0.20G Gaunts (1.57 GB) took: 0.49s.
[2026-03-15 08:28:44] [Task: 1/1] Batch integrator: 0.44 s (at ~4.49 GFLOPS) | IO (0.020 GB): 0.13s
[2026-03-15 08:28:44] [Task: 1/1] Finished global chunk index 1 (local: 2/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:44] [Task: 1/1] Starting on chunk index 2 (3/13 for this task).
[2026-03-15 08:28:44] [Task: 1/1] Working on 152 m's. |m| ∈ [109->185).
[2026-03-15 08:28:44] [Task: 1/1] Computing 0.03G Gaunts (0.27 GB) took: 0.19s.
[2026-03-15 08:28:44] [Task: 1/1] Batch integrator: 0.12 s (at ~2.83 GFLOPS) | IO (0.023 GB): 0.16s
[2026-03-15 08:28:44] [Task: 1/1] Finished global chunk index 2 (local: 3/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:44] [Task: 1/1] Starting on chunk index 3 (4/13 for this task).
[2026-03-15 08:28:44] [Task: 1/1] Starting pol: I | baseline E: 0.0 m, N: 8.5 m, U: 0.0 m | Sky lmax: 222, baseline lmax: 192.
[2026-03-15 08:28:44] [Task: 1/1] Working on 59 m's. |m| ∈ [0->30).
[2026-03-15 08:28:45] [Task: 1/1] Computing 0.20G Gaunts (1.62 GB) took: 0.43s.
[2026-03-15 08:28:45] [Task: 1/1] Batch integrator: 0.44 s (at ~4.59 GFLOPS) | IO (0.009 GB): 0.06s
[2026-03-15 08:28:45] [Task: 1/1] Finished global chunk index 3 (local: 4/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:45] [Task: 1/1] Starting on chunk index 4 (5/13 for this task).
[2026-03-15 08:28:45] [Task: 1/1] Working on 72 m's. |m| ∈ [30->66).
[2026-03-15 08:28:46] [Task: 1/1] Computing 0.20G Gaunts (1.61 GB) took: 0.46s.
[2026-03-15 08:28:46] [Task: 1/1] Batch integrator: 0.42 s (at ~4.79 GFLOPS) | IO (0.011 GB): 0.07s
[2026-03-15 08:28:46] [Task: 1/1] Finished global chunk index 4 (local: 5/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:46] [Task: 1/1] Starting on chunk index 5 (6/13 for this task).
[2026-03-15 08:28:46] [Task: 1/1] Working on 98 m's. |m| ∈ [66->115).
[2026-03-15 08:28:47] [Task: 1/1] Computing 0.19G Gaunts (1.52 GB) took: 0.46s.
[2026-03-15 08:28:47] [Task: 1/1] Batch integrator: 0.44 s (at ~4.35 GFLOPS) | IO (0.015 GB): 0.09s
[2026-03-15 08:28:47] [Task: 1/1] Finished global chunk index 5 (local: 6/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:47] [Task: 1/1] Starting on chunk index 6 (7/13 for this task).
[2026-03-15 08:28:47] [Task: 1/1] Working on 216 m's. |m| ∈ [115->223).
[2026-03-15 08:28:48] [Task: 1/1] Computing 0.10G Gaunts (0.83 GB) took: 0.43s.
[2026-03-15 08:28:48] [Task: 1/1] Batch integrator: 0.27 s (at ~3.81 GFLOPS) | IO (0.033 GB): 0.21s
[2026-03-15 08:28:48] [Task: 1/1] Finished global chunk index 6 (local: 7/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:48] [Task: 1/1] Starting on chunk index 7 (8/13 for this task).
[2026-03-15 08:28:48] [Task: 1/1] Starting pol: I | baseline E: 6.0 m, N: 8.5 m, U: 0.0 m | Sky lmax: 257, baseline lmax: 227.
[2026-03-15 08:28:48] [Task: 1/1] Working on 47 m's. |m| ∈ [0->24).
[2026-03-15 08:28:49] [Task: 1/1] Computing 0.20G Gaunts (1.58 GB) took: 0.40s.
[2026-03-15 08:28:49] [Task: 1/1] Batch integrator: 0.44 s (at ~4.53 GFLOPS) | IO (0.007 GB): 0.04s
[2026-03-15 08:28:49] [Task: 1/1] Finished global chunk index 7 (local: 8/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:49] [Task: 1/1] Starting on chunk index 8 (9/13 for this task).
[2026-03-15 08:28:49] [Task: 1/1] Working on 54 m's. |m| ∈ [24->51).
[2026-03-15 08:28:50] [Task: 1/1] Computing 0.20G Gaunts (1.61 GB) took: 0.42s.
[2026-03-15 08:28:50] [Task: 1/1] Batch integrator: 0.44 s (at ~4.56 GFLOPS) | IO (0.008 GB): 0.05s
[2026-03-15 08:28:50] [Task: 1/1] Finished global chunk index 8 (local: 9/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:50] [Task: 1/1] Starting on chunk index 9 (10/13 for this task).
[2026-03-15 08:28:50] [Task: 1/1] Working on 62 m's. |m| ∈ [51->82).
[2026-03-15 08:28:50] [Task: 1/1] Computing 0.19G Gaunts (1.56 GB) took: 0.42s.
[2026-03-15 08:28:51] [Task: 1/1] Batch integrator: 0.43 s (at ~4.49 GFLOPS) | IO (0.009 GB): 0.09s
[2026-03-15 08:28:51] [Task: 1/1] Finished global chunk index 9 (local: 10/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:51] [Task: 1/1] Starting on chunk index 10 (11/13 for this task).
[2026-03-15 08:28:51] [Task: 1/1] Working on 76 m's. |m| ∈ [82->120).
[2026-03-15 08:28:51] [Task: 1/1] Computing 0.18G Gaunts (1.48 GB) took: 0.45s.
[2026-03-15 08:28:52] [Task: 1/1] Batch integrator: 0.41 s (at ~4.56 GFLOPS) | IO (0.011 GB): 0.08s
[2026-03-15 08:28:52] [Task: 1/1] Finished global chunk index 10 (local: 11/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:52] [Task: 1/1] Starting on chunk index 11 (12/13 for this task).
[2026-03-15 08:28:52] [Task: 1/1] Working on 106 m's. |m| ∈ [120->173).
[2026-03-15 08:28:52] [Task: 1/1] Computing 0.16G Gaunts (1.28 GB) took: 0.42s.
[2026-03-15 08:28:53] [Task: 1/1] Batch integrator: 0.38 s (at ~4.24 GFLOPS) | IO (0.016 GB): 0.12s
[2026-03-15 08:28:53] [Task: 1/1] Finished global chunk index 11 (local: 12/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:53] [Task: 1/1] Starting on chunk index 12 (13/13 for this task).
[2026-03-15 08:28:53] [Task: 1/1] Working on 170 m's. |m| ∈ [173->258).
[2026-03-15 08:28:53] [Task: 1/1] Computing 0.05G Gaunts (0.39 GB) took: 0.24s.
[2026-03-15 08:28:53] [Task: 1/1] Batch integrator: 0.16 s (at ~3.05 GFLOPS) | IO (0.026 GB): 0.15s
[2026-03-15 08:28:53] [Task: 1/1] Finished global chunk index 12 (local: 13/13) | MAXRSS: 2.01 GB | Est. time remaining: 00h00m.
Finished at 2026-03-15 08:28:53. Total time: 00h00m.
Integrator cache generated successfully.
Step 3: Generate Voltage and Power Beams
We model a single HIRAX dish with a Gaussian primary beam whose FWHM scales as \(\lambda / D\). X and Y voltage beams are created in the TIRS frame and normalised so that the power integral equals unity. We then form the auto-correlation power beam for each polarisation pair (XX, YY) and save them to a Zarr store.
def cocr_gaussian(freq, theta, phi, D=D_dish):
"""Ludwig-3 Gaussian beam model: returns (co-pol, cross-pol) arrays."""
wl = (freq * units.MHz).to("m", equivalencies=units.spectral())
sigma = (gaussian_fwhm_to_sigma * (wl / D)).to("").value
co_pol = np.exp(-(theta ** 2) / (4 * sigma ** 2))
cr_pol = 1e-3 * co_pol
return np.stack([co_pol, cr_pol], axis=0)
vbeam_x = TIRSVoltageBeam.from_ludwig3(
cocr_gaussian,
lmax=voltage_beam_lmax,
frequencies_MHz=frequencies_MHz,
polarisation="X",
latitude=latitude_rad,
longitude=longitude_rad,
altitude=pointing_altitude,
azimuth=pointing_azimuth,
).normalise("power_integral")
vbeam_y = TIRSVoltageBeam.from_ludwig3(
cocr_gaussian,
lmax=voltage_beam_lmax,
frequencies_MHz=frequencies_MHz,
polarisation="Y",
latitude=latitude_rad,
longitude=longitude_rad,
altitude=pointing_altitude,
azimuth=pointing_azimuth,
).normalise("power_integral")
print(f"Voltage beam alm shape: {vbeam_x.alm.shape} (n_pol_components, n_freq, lmax+1, 2*lmax+1)")
Voltage beam alm shape: (2, 5, 16, 31) (n_pol_components, n_freq, lmax+1, 2*lmax+1)
vbeam_zarr = DATA_DIR / "voltage_beams.zarr"
vbeam_x.to_zarr_store(vbeam_zarr, group_path="dish0_X")
vbeam_y.to_zarr_store(vbeam_zarr, group_path="dish0_Y")
vbeam_x.to_zarr_store(vbeam_zarr, group_path="dish1_X")
vbeam_y.to_zarr_store(vbeam_zarr, group_path="dish1_Y")
print(f"Voltage beams saved to {vbeam_zarr}")
Voltage beams saved to data/e2e_cli_example/voltage_beams.zarr
pbeam_zarr = DATA_DIR / "power_beams.zarr"
beam_pairs = [
("dish0_X_dish1_X", vbeam_x, vbeam_x),
("dish0_Y_dish1_Y", vbeam_y, vbeam_y),
]
for group_name, va, vb in beam_pairs:
pb = TIRSPowerBeam.from_voltage_beams(
va, vb,
sky_pol="I",
power_beam_lmax=power_beam_lmax,
power_beam_mmax=power_beam_mmax,
)
pb.to_zarr_store(pbeam_zarr, group_path=group_name)
print(f" Saved power beam '{group_name}', alm shape: {pb.alm.shape}")
print(f"Power beams saved to {pbeam_zarr}")
Saved power beam 'dish0_X_dish1_X', alm shape: (5, 31, 61)
Saved power beam 'dish0_Y_dish1_Y', alm shape: (5, 31, 61)
Power beams saved to data/e2e_cli_example/power_beams.zarr
Step 4: Generate Visibilities from Voltage Beams (CLI)
serval genvis from-voltage-beams accepts a voltage-beam store directly and derives the
power beams on-the-fly from pairs of TIRSVoltageBeam groups. The pol_baseline_mapping
entries are 3-tuples [bl_id, vbeam_i_group, vbeam_j_group] instead of the 2-tuples used
by from-power-beams.
We first save the voltage beams to their own Zarr store, then write a TOML config and invoke the CLI.
vis_vb_zarr = DATA_DIR / "visibilities_from_vbeams.zarr"
genvis_vb_toml = DATA_DIR / "genvis_vb_config.toml"
# 3-tuple mapping: [bl_id, vbeam_i_group, vbeam_j_group]
vb_pol_bl_mapping = [
[f"I/{tag}", "dish0_X", "dish1_X"]
for tag in bl_tags
] + [
[f"I/{tag}", "dish0_Y", "dish1_Y"]
for tag in bl_tags
]
genvis_vb_cfg = tomlkit.document()
genvis_vb_cfg.add("store_location", str(vis_vb_zarr))
genvis_vb_cfg.add("integrator_cache", str(cache_zarr))
genvis_vb_cfg.add("voltage_beams", str(vbeam_zarr))
genvis_vb_cfg.add("pol_baseline_mapping", vb_pol_bl_mapping)
genvis_vb_cfg.add("memory_per_task", "2 GiB")
genvis_vb_toml.write_text(tomlkit.dumps(genvis_vb_cfg))
print(tomlkit.dumps(genvis_vb_cfg))
store_location = "data/e2e_cli_example/visibilities_from_vbeams.zarr"
integrator_cache = "data/e2e_cli_example/integrator_cache.zarr"
voltage_beams = "data/e2e_cli_example/voltage_beams.zarr"
pol_baseline_mapping = [["I/bl_E_6p500m_N_0p000m_U_0p000m", "dish0_X", "dish1_X"], ["I/bl_E_0p000m_N_8p500m_U_0p000m", "dish0_X", "dish1_X"], ["I/bl_E_6p000m_N_8p500m_U_0p000m", "dish0_X", "dish1_X"], ["I/bl_E_6p500m_N_0p000m_U_0p000m", "dish0_Y", "dish1_Y"], ["I/bl_E_0p000m_N_8p500m_U_0p000m", "dish0_Y", "dish1_Y"], ["I/bl_E_6p000m_N_8p500m_U_0p000m", "dish0_Y", "dish1_Y"]]
memory_per_task = "2 GiB"
result = subprocess.run(
[
"serval", "--dist-mode", "serial",
"genvis", "from-voltage-beams",
"--config", str(genvis_vb_toml),
],
capture_output=True,
text=True,
)
print(result.stdout)
if result.returncode != 0:
print("STDERR:", result.stderr)
raise RuntimeError(f"genvis from-voltage-beams failed (return code {result.returncode})")
print("Visibilities from voltage beams generated successfully.")
Welcome to:
____ _____ ______ ___ _
/ ___|| ____| _ \ \ / / \ | |
\___ \| _| | |_) \ \ / / _ \ | |
___) | |___| _ < \ V / ___ \| |___
|____/|_____|_| \_\ \_/_/ \_\_____| v1.0.1b0
----------------------------------------------------------------
Distribution mode: serial (1 task).
----------------------------------------------------------------
Generating visibilities from voltage beams...
----------------------------------------------------------------
Parsed, converted and merged parameters:
----------------------------------------------------------------
store_location = "data/e2e_cli_example/visibilities_from_vbeams.zarr"
integrator_cache = "data/e2e_cli_example/integrator_cache.zarr"
voltage_beams = "data/e2e_cli_example/voltage_beams.zarr"
pol_baseline_mapping = [
["I/bl_E_6p500m_N_0p000m_U_0p000m", "dish0_X", "dish1_X"],
["I/bl_E_0p000m_N_8p500m_U_0p000m", "dish0_X", "dish1_X"],
["I/bl_E_6p000m_N_8p500m_U_0p000m", "dish0_X", "dish1_X"],
["I/bl_E_6p500m_N_0p000m_U_0p000m", "dish0_Y", "dish1_Y"],
["I/bl_E_0p000m_N_8p500m_U_0p000m", "dish0_Y", "dish1_Y"],
["I/bl_E_6p000m_N_8p500m_U_0p000m", "dish0_Y", "dish1_Y"],
]
memory_per_task = 2.0
start_chunk = 0
----------------------------------------------------------------
Expanded baseline-beam mapping:
----------------------------------------------------------------
I/bl_E_6p500m_N_0p000m_U_0p000m/dish0_X_dish1_X
I/bl_E_0p000m_N_8p500m_U_0p000m/dish0_X_dish1_X
I/bl_E_6p000m_N_8p500m_U_0p000m/dish0_X_dish1_X
I/bl_E_6p500m_N_0p000m_U_0p000m/dish0_Y_dish1_Y
I/bl_E_0p000m_N_8p500m_U_0p000m/dish0_Y_dish1_Y
I/bl_E_6p000m_N_8p500m_U_0p000m/dish0_Y_dish1_Y
----------------------------------------------------------------
Distribution Plan:
----------------------------------------------------------------
Split into 6 chunks across 1 task(s).
Chunk Statistics:
chan/chunk mem/chunk [GB] chunks/task
mean 5.00 0.07 6.00
std 0.00 0.01 NaN
min 5.00 0.06 6.00
25% 5.00 0.06 6.00
50% 5.00 0.07 6.00
75% 5.00 0.08 6.00
max 5.00 0.08 6.00
----------------------------------------------------------------
Per Baseline Statistics:
----------------------------------------------------------------
Plan generated in 0.02s.
[2026-03-15 08:28:56] [Task: 1/1] Task initialized with OMP_NUM_THREADS=unset
[2026-03-15 08:28:56] [Task: 1/1] Starting on chunk index 0 (1/6 for this task).
[2026-03-15 08:28:56] [Task: 1/1] Starting baseline group I/bl_E_0p000m_N_8p500m_U_0p000m, ifreq ∈ [0->5)
[2026-03-15 08:28:56] [Task: 1/1] Loaded integrator cache in 0.18s
[2026-03-15 08:28:57] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.38s | Visibility IO write: 0.00s
[2026-03-15 08:28:57] [Task: 1/1] Finished global chunk index 0 (local: 1/6) | MAXRSS: 0.40 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:57] [Task: 1/1] Starting on chunk index 1 (2/6 for this task).
[2026-03-15 08:28:57] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.03s | Visibility IO write: 0.00s
[2026-03-15 08:28:57] [Task: 1/1] Finished global chunk index 1 (local: 2/6) | MAXRSS: 0.40 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:57] [Task: 1/1] Starting on chunk index 2 (3/6 for this task).
[2026-03-15 08:28:57] [Task: 1/1] Starting baseline group I/bl_E_6p000m_N_8p500m_U_0p000m, ifreq ∈ [0->5)
[2026-03-15 08:28:57] [Task: 1/1] Loaded integrator cache in 0.25s
[2026-03-15 08:28:57] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.03s | Visibility IO write: 0.00s
[2026-03-15 08:28:57] [Task: 1/1] Finished global chunk index 2 (local: 3/6) | MAXRSS: 0.46 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:57] [Task: 1/1] Starting on chunk index 3 (4/6 for this task).
[2026-03-15 08:28:57] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.02s | Visibility IO write: 0.00s
[2026-03-15 08:28:57] [Task: 1/1] Finished global chunk index 3 (local: 4/6) | MAXRSS: 0.46 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:57] [Task: 1/1] Starting on chunk index 4 (5/6 for this task).
[2026-03-15 08:28:57] [Task: 1/1] Starting baseline group I/bl_E_6p500m_N_0p000m_U_0p000m, ifreq ∈ [0->5)
[2026-03-15 08:28:57] [Task: 1/1] Loaded integrator cache in 0.25s
[2026-03-15 08:28:57] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.03s | Visibility IO write: 0.01s
[2026-03-15 08:28:57] [Task: 1/1] Finished global chunk index 4 (local: 5/6) | MAXRSS: 0.47 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:57] [Task: 1/1] Starting on chunk index 5 (6/6 for this task).
[2026-03-15 08:28:57] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.03s | Visibility IO write: 0.00s
[2026-03-15 08:28:57] [Task: 1/1] Finished global chunk index 5 (local: 6/6) | MAXRSS: 0.47 GB | Est. time remaining: 00h00m.
Finished at 2026-03-15 08:28:57. Total time: 00h00m.
Visibilities from voltage beams generated successfully.
Step 5: Generate Visibilities from Power Beams (CLI)
We write a second TOML configuration file that maps each baseline / polarisation-pair
to the power beam groups in the store, then invoke
serval genvis from-power-beams.
vis_zarr = DATA_DIR / "visibilities.zarr"
genvis_toml = DATA_DIR / "genvis_config.toml"
# Build the pol/baseline-to-beam mapping for all baselines × all polarisation pairs
beam_group_names = [name for name, _, _ in beam_pairs]
pol_bl_mapping = [
[f"I/{tag}", beam_name]
for tag in bl_tags
for beam_name in beam_group_names
]
genvis_cfg = tomlkit.document()
genvis_cfg.add("store_location", str(vis_zarr))
genvis_cfg.add("integrator_cache", str(cache_zarr))
genvis_cfg.add("power_beams", str(pbeam_zarr))
genvis_cfg.add("pol_baseline_mapping", pol_bl_mapping)
genvis_cfg.add("memory_per_task", "2 GiB")
genvis_toml.write_text(tomlkit.dumps(genvis_cfg))
print(tomlkit.dumps(genvis_cfg))
store_location = "data/e2e_cli_example/visibilities.zarr"
integrator_cache = "data/e2e_cli_example/integrator_cache.zarr"
power_beams = "data/e2e_cli_example/power_beams.zarr"
pol_baseline_mapping = [["I/bl_E_6p500m_N_0p000m_U_0p000m", "dish0_X_dish1_X"], ["I/bl_E_6p500m_N_0p000m_U_0p000m", "dish0_Y_dish1_Y"], ["I/bl_E_0p000m_N_8p500m_U_0p000m", "dish0_X_dish1_X"], ["I/bl_E_0p000m_N_8p500m_U_0p000m", "dish0_Y_dish1_Y"], ["I/bl_E_6p000m_N_8p500m_U_0p000m", "dish0_X_dish1_X"], ["I/bl_E_6p000m_N_8p500m_U_0p000m", "dish0_Y_dish1_Y"]]
memory_per_task = "2 GiB"
result = subprocess.run(
[
"serval", "--dist-mode", "serial",
"genvis", "from-power-beams",
"--config", str(genvis_toml),
],
capture_output=True,
text=True,
)
print(result.stdout)
if result.returncode != 0:
print("STDERR:", result.stderr)
raise RuntimeError(f"genvis failed (return code {result.returncode})")
print("Visibilities generated successfully.")
Welcome to:
____ _____ ______ ___ _
/ ___|| ____| _ \ \ / / \ | |
\___ \| _| | |_) \ \ / / _ \ | |
___) | |___| _ < \ V / ___ \| |___
|____/|_____|_| \_\ \_/_/ \_\_____| v1.0.1b0
----------------------------------------------------------------
Distribution mode: serial (1 task).
----------------------------------------------------------------
Generating visibilities from power beams...
----------------------------------------------------------------
Parsed, converted and merged parameters:
----------------------------------------------------------------
store_location = "data/e2e_cli_example/visibilities.zarr"
integrator_cache = "data/e2e_cli_example/integrator_cache.zarr"
power_beams = "data/e2e_cli_example/power_beams.zarr"
pol_baseline_mapping = [
["I/bl_E_6p500m_N_0p000m_U_0p000m", "dish0_X_dish1_X"],
["I/bl_E_6p500m_N_0p000m_U_0p000m", "dish0_Y_dish1_Y"],
["I/bl_E_0p000m_N_8p500m_U_0p000m", "dish0_X_dish1_X"],
["I/bl_E_0p000m_N_8p500m_U_0p000m", "dish0_Y_dish1_Y"],
["I/bl_E_6p000m_N_8p500m_U_0p000m", "dish0_X_dish1_X"],
["I/bl_E_6p000m_N_8p500m_U_0p000m", "dish0_Y_dish1_Y"],
]
memory_per_task = 2.0
start_chunk = 0
----------------------------------------------------------------
Expanded baseline-beam mapping:
----------------------------------------------------------------
I/bl_E_6p500m_N_0p000m_U_0p000m/dish0_X_dish1_X
I/bl_E_6p500m_N_0p000m_U_0p000m/dish0_Y_dish1_Y
I/bl_E_0p000m_N_8p500m_U_0p000m/dish0_X_dish1_X
I/bl_E_0p000m_N_8p500m_U_0p000m/dish0_Y_dish1_Y
I/bl_E_6p000m_N_8p500m_U_0p000m/dish0_X_dish1_X
I/bl_E_6p000m_N_8p500m_U_0p000m/dish0_Y_dish1_Y
----------------------------------------------------------------
Distribution Plan:
----------------------------------------------------------------
Split into 6 chunks across 1 task(s).
Chunk Statistics:
chan/chunk mem/chunk [GB] chunks/task
mean 5.00 0.07 6.00
std 0.00 0.01 NaN
min 5.00 0.06 6.00
25% 5.00 0.06 6.00
50% 5.00 0.07 6.00
75% 5.00 0.08 6.00
max 5.00 0.08 6.00
----------------------------------------------------------------
Per Baseline Statistics:
----------------------------------------------------------------
Plan generated in 0.03s.
[2026-03-15 08:28:59] [Task: 1/1] Task initialized with OMP_NUM_THREADS=unset
[2026-03-15 08:28:59] [Task: 1/1] Starting on chunk index 0 (1/6 for this task).
[2026-03-15 08:28:59] [Task: 1/1] Starting baseline group I/bl_E_0p000m_N_8p500m_U_0p000m, ifreq ∈ [0->5)
[2026-03-15 08:28:59] [Task: 1/1] Loaded integrator cache in 0.16s
[2026-03-15 08:28:59] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.00s | Visibility IO write: 0.00s
[2026-03-15 08:28:59] [Task: 1/1] Finished global chunk index 0 (local: 1/6) | MAXRSS: 0.39 GB | Est. time remaining: 00h00m.
[2026-03-15 08:28:59] [Task: 1/1] Starting on chunk index 1 (2/6 for this task).
[2026-03-15 08:29:00] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.00s | Visibility IO write: 0.00s
[2026-03-15 08:29:00] [Task: 1/1] Finished global chunk index 1 (local: 2/6) | MAXRSS: 0.39 GB | Est. time remaining: 00h00m.
[2026-03-15 08:29:00] [Task: 1/1] Starting on chunk index 2 (3/6 for this task).
[2026-03-15 08:29:00] [Task: 1/1] Starting baseline group I/bl_E_6p000m_N_8p500m_U_0p000m, ifreq ∈ [0->5)
[2026-03-15 08:29:00] [Task: 1/1] Loaded integrator cache in 0.24s
[2026-03-15 08:29:00] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.00s | Visibility IO write: 0.00s
[2026-03-15 08:29:00] [Task: 1/1] Finished global chunk index 2 (local: 3/6) | MAXRSS: 0.41 GB | Est. time remaining: 00h00m.
[2026-03-15 08:29:00] [Task: 1/1] Starting on chunk index 3 (4/6 for this task).
[2026-03-15 08:29:00] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.00s | Visibility IO write: 0.00s
[2026-03-15 08:29:00] [Task: 1/1] Finished global chunk index 3 (local: 4/6) | MAXRSS: 0.41 GB | Est. time remaining: 00h00m.
[2026-03-15 08:29:00] [Task: 1/1] Starting on chunk index 4 (5/6 for this task).
[2026-03-15 08:29:00] [Task: 1/1] Starting baseline group I/bl_E_6p500m_N_0p000m_U_0p000m, ifreq ∈ [0->5)
[2026-03-15 08:29:00] [Task: 1/1] Loaded integrator cache in 0.22s
[2026-03-15 08:29:00] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.01s | Visibility IO write: 0.00s
[2026-03-15 08:29:00] [Task: 1/1] Finished global chunk index 4 (local: 5/6) | MAXRSS: 0.43 GB | Est. time remaining: 00h00m.
[2026-03-15 08:29:00] [Task: 1/1] Starting on chunk index 5 (6/6 for this task).
[2026-03-15 08:29:00] [Task: 1/1] Visibility generation: 0.00 s | Beam IO read: 0.00s | Visibility IO write: 0.00s
[2026-03-15 08:29:00] [Task: 1/1] Finished global chunk index 5 (local: 6/6) | MAXRSS: 0.43 GB | Est. time remaining: 00h00m.
Finished at 2026-03-15 08:29:00. Total time: 00h00m.
Visibilities generated successfully.
Step 6: Inspect Power-Beam Visibilities
The CLI stores visibilities at {pol}/{baseline_tag}/{beam_group}/vis inside the
output Zarr store. The second axis corresponds to equally-spaced east rotation angle
(ERA) samples covering \([0°, 360°)\).
# Load all Visibilities containers upfront (lazy zarr arrays, no data read yet).
vis_containers = {
(tag, beam_name): Visibilities.from_zarr_store(vis_zarr, f"I/{tag}/{beam_name}")
for tag in bl_tags
for beam_name in beam_group_names
}
n_era = vis_containers[bl_tags[0], beam_group_names[0]].vis.shape[1]
vis_cli = {} # store for later comparison
def _make_fig(title):
fig, axes = plt.subplots(
len(baselines_enu), len(beam_group_names),
figsize=(5 * len(beam_group_names), 3 * len(baselines_enu)),
squeeze=False,
)
fig.suptitle(title)
return fig, axes
def _decorate(ax, i, j, tag, beam_name, xlabel):
ax.set_title(f"{tag}\n{beam_name}")
ax.set_xlabel(xlabel)
ax.set_ylabel(r"Re[$\mathcal{V}$]")
if i == 0 and j == len(beam_group_names) - 1:
ax.legend(fontsize=7)
# ── Plot 1: vs ERA ────────────────────────────────────────────────────────────
fig1, axes1 = _make_fig("Visibilities vs Earth Rotation Angle")
for i, (tag, bl) in enumerate(zip(bl_tags, baselines_enu)):
for j, beam_name in enumerate(beam_group_names):
vis_obj = vis_containers[tag, beam_name]
vis = np.asarray(vis_obj.vis)
vis_cli[(tag, beam_name)] = vis
ax = axes1[i, j]
for k, f in enumerate(frequencies_MHz):
ax.plot(vis_obj.era_deg, vis[k].real, label=f"{f:.0f} MHz")
_decorate(ax, i, j, tag, beam_name, "ERA [deg]")
fig1.tight_layout()
plt.show()
# ── Plot 2: vs Right Ascension ────────────────────────────────────────────────
fig2, axes2 = _make_fig("Visibilities vs Right Ascension")
ra_grid = np.linspace(0, 360, n_era, endpoint=False)
for i, (tag, bl) in enumerate(zip(bl_tags, baselines_enu)):
for j, beam_name in enumerate(beam_group_names):
vis_ra = vis_containers[tag, beam_name].resample_to_ra(
ra_grid,
latitude=latitude_rad,
longitude=longitude_rad,
altitude=pointing_altitude,
azimuth=pointing_azimuth,
)
ax = axes2[i, j]
for k, f in enumerate(frequencies_MHz):
ax.plot(vis_ra.ra_deg, vis_ra.vis[k].real, label=f"{f:.0f} MHz")
_decorate(ax, i, j, tag, beam_name, "RA [deg]")
fig2.tight_layout()
plt.show()
# ── Plot 3: vs UTC time on a specific day ─────────────────────────────────────
fig3, axes3 = _make_fig("Visibilities vs UTC time on 2026-01-01")
obs_start = Time("2026-01-01T00:00:00", scale="utc")
sidereal_day_s = 86164.1 # seconds in one sidereal day
obs_times = obs_start + np.linspace(0, sidereal_day_s, n_era) * units.s
for i, (tag, bl) in enumerate(zip(bl_tags, baselines_enu)):
for j, beam_name in enumerate(beam_group_names):
vis_time = vis_containers[tag, beam_name].resample_to_times(obs_times)
time_hours = (vis_time.ut1_time - obs_start.ut1.jd) * 24
ax = axes3[i, j]
for k, f in enumerate(frequencies_MHz):
ax.plot(time_hours, vis_time.vis[k].real, label=f"{f:.0f} MHz")
_decorate(ax, i, j, tag, beam_name, "UTC [h] on 2026-01-01")
fig3.tight_layout()
plt.show()



Step 7: Compare Power-Beam and Voltage-Beam Visibilities
Both CLI paths should produce numerically identical results: from-voltage-beams calls the
same PowerBeam.from_voltage_beams code internally, so any difference is a bug.
vis_vb_root = zarr.open(str(vis_vb_zarr), mode="r")
vis_pb_root = zarr.open(str(vis_zarr), mode="r")
# Both paths use the same group-name convention, so names match directly.
beam_names = ["dish0_X_dish1_X", "dish0_Y_dish1_Y"]
for i, (tag, _) in enumerate(zip(bl_tags, baselines_enu)):
for j, b_name in enumerate(beam_names):
vis_vb = vis_vb_root[f"I/{tag}/{b_name}/vis"][:]
vis_pb = vis_pb_root[f"I/{tag}/{b_name}/vis"][:]
np.testing.assert_allclose(
vis_vb, vis_pb, rtol=1e-10,
err_msg=f"Mismatch for {tag} / {b_name}",
)
print("All visibilities agree to rtol=1e-10.")
All visibilities agree to rtol=1e-10.