HI Intensity Mapping generator with LHS sampling & multi-CPU batching
Project description
HIcrafterHIcrafter generates realistic HI intensity maps across redshift ranges. Features: Gaussian beam/noise/masks/zebras, single maps or Latin Hypercube parameter suites for ML/SBI training, multi-CPU parallel batching. Seeds logged 1→N. Perfect for cosmological simulations. |
|
|
HIcrafter is designed for rapid experimentation and reproducibility, enabling large-scale simulation campaigns while preserving precise control over astrophysical and instrumental systematics. |
This README explains the package structure, configuration options, and provides guidance for running simulations and using the included tutorials.
Table of Contents
Overview
HIcrafter generates full-sky neutral hydrogen (HI) intensity maps by combining realizations of the large-scale matter density field with an effective model for the HI brightness temperature, followed by the application of instrumental and observational effects.
HI intensity map construction
The HI brightness temperature fluctuation along a direction n̂ is modeled as
T_HI(n̂) = ∫_{z_min}^{z_max} T̄_HI(z) · b_HI(z) · δ_m(n̂, z) dz
where δ_m(n̂, z) is the matter overdensity field, b_HI(z) is the HI bias, and T̄_HI(z) is the mean
HI brightness temperature.
In practice, the redshift integral is discretized into radial shells, yielding
T_HI(n̂) ≈ Σ_i ⟨T̄_HI⟩_i · b_HI(z_i) · δ_i(n̂)
where δ_i(n̂) are lognormal realizations of the matter field on shell i.
The mean HI brightness temperature is parametrized as
T̄_HI(z) = 189 × 4 × 10^-4 · (1 + z)^2.6 · h / H(z)
and the HI bias is modeled as
b_HI(z) = 0.6 + 0.3(1 + z)
Beam smoothing
Finite angular resolution is modeled by convolution with a Gaussian beam. In spherical harmonic space, the observed coefficients are
a_obs(ℓ, m) = b_ℓ · a(ℓ, m)
with beam transfer function
b_ℓ = exp[ -½ · ℓ(ℓ + 1) · σ_b^2 ]
σ_b = θ_FWHM / sqrt(8 ln 2)
Thermal noise
Instrumental noise is modeled as additive, uncorrelated Gaussian noise in pixel space,
n(n̂) ~ N(0, σ_n^2)
such that the noisy map is
T'_HI(n̂) = T_HI(n̂) + n(n̂)
Zebra-striping systematics
Scan-synchronous striping systematics are modeled as a deterministic modulation added to the map,
S(θ, φ) = A · sin[ (2π / λ) · (sinθ cosφ) ]
where (θ, φ) are spherical coordinates, λ sets the stripe wavelength, and A is the stripe amplitude.
The final observed map is therefore
T_obs(n̂) = (T_HI * B)(n̂) + n(n̂) + S(n̂)
with optional masking applied multiplicatively in pixel space.
Features
- Lognormal HI intensity maps — full-sky HI brightness temperature maps generated from lognormal matter-field realizations via GLASS.
- Gaussian beam smoothing — convolve the map with a Gaussian beam of arbitrary FWHM (degrees).
- Thermal noise — add pixel-level white Gaussian noise with a configurable noise level.
- Zebra-striping systematics — inject scan-synchronous striping artefacts.
- Sky masking — apply an external mask (HEALPix FITS or
.npyfile). - Latin Hypercube suites — generate large batches of maps spanning a cosmological parameter space for ML / SBI training.
- Reproducible seeds — every map is generated from a deterministic seed; seeds are logged 1 → N across a batch.
- Multi-CPU batching — produce N maps efficiently with
generate_batch(). - Flexible cosmology — five configurable ΛCDM parameters (h, As, Ωc, Ωb, ns) via CAMB.
Installation
We recommend installing HIcrafter in a new conda environment.
Install HIcrafter using pip:
pip install hicrafter
To access tutorials and examples, clone the repository:
git clone https://github.com/PaulineGorbatchev/HIcrafter.git
Dependencies are listed in requirements.txt:
numpy>=1.21.0
healpy>=1.15.0
camb>=0.5.0
glass>=0.5.0
scipy>=1.8.0
Usage
Import HIGenerator from the hicrafter package and call generate_map() to produce a
single HEALPix map, or generate_batch() to produce a numbered set of maps saved to disk.
from hicrafter import HIGenerator
gen = HIGenerator(nside=64, z_min=0.40, z_max=0.45, seed=42)
hi_map = gen.generate_map() # returns numpy array of shape (hp.nside2npix(nside),)
To generate a batch of N maps (each with an independent seed offset), use:
maps = gen.generate_batch(n_maps=100, output_dir="output/maps")
See the basic pipeline tutorial and Latin Hypercube tutorial for full worked examples.
Configuration
All parameters are passed to the HIGenerator constructor.
| Parameter | Type | Default | Description |
|---|---|---|---|
h | float | 0.7 | Dimensionless Hubble constant. |
As | float | 2×10−9 | Scalar amplitude of primordial power spectrum. |
Oc | float | 0.25 | Cold dark matter density Ωc. |
Ob | float | 0.05 | Baryon density Ωb. |
ns | float | 0.965 | Scalar spectral index. |
nside | int | 32 | HEALPix resolution parameter (capped at 512 for safety). Pixel count = 12 × nside². |
z_min | float | 0.4 | Minimum redshift of the survey volume. |
z_max | float | 0.45 | Maximum redshift of the survey volume. |
nbins | int | 1 | Number of tomographic redshift bins. |
sigmaz0 | float | 1×10−4 | Photometric redshift uncertainty σz,0. |
beam_deg | float or None | None | Gaussian beam FWHM in degrees. Set to None to skip smoothing. |
noise | bool | True | Whether to add thermal noise. |
noise_level | float | 1.0 | Standard deviation of the Gaussian noise (same units as THI). |
mask_file | str or None | None | Path to a sky mask (.npy or HEALPix FITS). Masked pixels are set to zero. |
zebras | bool | False | Add scan-synchronous zebra-striping systematics. |
seed | int | 1 | Base random seed for reproducibility. |
Examples
Baseline HI map
Generate a baseline HI intensity map without instrumental or observational effects.
gen = HIGenerator(
nside=32,
z_min=0.40,
z_max=0.45,
nbins=1,
sigmaz0=1e-4,
beam_deg=None,
noise=False,
zebras=False,
seed=1,
)
hi_map = gen.generate_map()
Beam smoothing
Apply Gaussian beam smoothing to the HI intensity map.
gen_beam = HIGenerator(
nside=32,
z_min=0.40,
z_max=0.45,
nbins=1,
sigmaz0=1e-4,
beam_deg=1.5,
noise=False,
zebras=False,
seed=1,
)
hi_beam = gen_beam.generate_map()
Thermal noise
Add thermal noise on top of the beam-smoothed HI map.
gen_noise = HIGenerator(
nside=32,
z_min=0.40,
z_max=0.45,
nbins=1,
sigmaz0=1e-4,
beam_deg=1.5,
noise=True,
noise_level=1.0,
zebras=False,
seed=1,
)
hi_noise = gen_noise.generate_map()
Zebra-striping systematics
Add zebra-striping systematics in addition to beam smoothing and noise.
gen_zebra = HIGenerator(
nside=32,
z_min=0.40,
z_max=0.45,
nbins=1,
sigmaz0=1e-4,
beam_deg=1.5,
noise=True,
noise_level=1.0,
zebras=True,
seed=1,
)
hi_zebra = gen_zebra.generate_map()
Citation
If you use HIcrafter in your research, please cite the repository:
@software{hicrafter,
author = {Gorbatchev, Pauline},
title = {HIcrafter: HI Intensity Map Generator},
url = {https://github.com/PaulineGorbatchev/HIcrafter},
version = {0.3.0}
}
Project details
Release history Release notifications | RSS feed
Download files
Download the file for your platform. If you're not sure which to choose, learn more about installing packages.
Source Distribution
Built Distribution
Filter files by name, interpreter, ABI, and platform.
If you're not sure about the file name format, learn more about wheel file names.
Copy a direct link to the current filters
File details
Details for the file hicrafter-0.3.0.tar.gz.
File metadata
- Download URL: hicrafter-0.3.0.tar.gz
- Upload date:
- Size: 9.6 kB
- Tags: Source
- Uploaded using Trusted Publishing? No
- Uploaded via: twine/6.2.0 CPython/3.11.8
File hashes
| Algorithm | Hash digest | |
|---|---|---|
| SHA256 |
997e8eece779971baa6745935f06391d2113d34b8aef6393215a09b56f408d54
|
|
| MD5 |
dbda82077412f5c7ef50a729f09f7c11
|
|
| BLAKE2b-256 |
6a047effb13eff0295752478f2c23d2a9024fd78da208a177f1366c870983021
|
File details
Details for the file hicrafter-0.3.0-py3-none-any.whl.
File metadata
- Download URL: hicrafter-0.3.0-py3-none-any.whl
- Upload date:
- Size: 8.9 kB
- Tags: Python 3
- Uploaded using Trusted Publishing? No
- Uploaded via: twine/6.2.0 CPython/3.11.8
File hashes
| Algorithm | Hash digest | |
|---|---|---|
| SHA256 |
6ccb8aa668d9290e4a33e2aaf9cf056b3ba696b4bd485ba83f22825b4e7176bb
|
|
| MD5 |
025906aa338f45b547cfd8e535afaa6d
|
|
| BLAKE2b-256 |
32c25c69a54d9b3a5416af75ea99a0fe2f095b744f5be3f3dfe35f42372811e8
|