Simulation Subpackage
Simulation Subpackage
This subpackage provides comprehensive simulation capabilities for CMB observations including signal, foreground, and noise components.
Modules
- cmb
CMB power spectrum calculation and map generation with cosmic birefringence effects (isotropic and time-dependent models).
- fg
Galactic foreground simulation including thermal dust and synchrotron emission with realistic spatial and spectral properties.
- mask
Mask generation for observations including Galactic plane cuts, point source masking, CO emission masking, and apodization.
- noise
Realistic noise simulation for LAT and SAT including white noise and 1/f components with proper spatial correlations.
- sky
Complete sky simulation classes that combine CMB, foregrounds, and noise for LAT and SAT observations. Includes bandpass integration and component separation via harmonic ILC.
Classes
- CMB
CMB simulation with cosmic birefringence
- synfast_pol
Fast polarization map generation
- Foreground
Foreground component simulation
- BandpassInt
Bandpass integration utilities
- HILC
Harmonic Internal Linear Combination
- Mask
Observation mask generation
- Noise
Noise realization generation
- NoiseSpectra
Noise power spectra
- LATsky
LAT observation simulation
- SATsky
SAT observation simulation
- LATskyC
LAT calibration simulation
- SATskyC
SAT calibration simulation
Overview
The simulation subpackage provides comprehensive tools for generating realistic CMB observations including signal, foregrounds, and instrumental effects.
Main Classes
CMB map generator with cosmic birefringence support. |
|
Galactic foreground generator for CMB analysis. |
|
Sky mask generator and manager for CMB analysis. |
|
Instrumental noise simulator for Simons Observatory telescopes. |
|
CMB Module
CMB Simulation Module
This module provides comprehensive functionality for simulating Cosmic Microwave Background (CMB) temperature and polarization maps with support for:
Lensed and unlensed CMB realizations
Cosmic birefringence effects (isotropic and anisotropic)
Custom power spectra from CAMB
Auxiliary scalar fields (for patchy τ models)
Beam convolution and coordinate transformations
Classes
- CMB
Main class for generating CMB maps with optional birefringence rotation.
Functions
- synfast_pol
Generate polarized CMB maps with optional auxiliary fields.
- hp_alm2map_spin
Fast spin-weighted spherical harmonic transforms.
- get_camb_cls
Compute CMB power spectra using CAMB.
Example
Generate a CMB realization with isotropic birefringence:
from cobi.simulation import CMB
cmb = CMB(
nside=512,
lmax=1500,
beta=0.35, # isotropic rotation angle in degrees
lensing=True,
libdir='./cmb_sims'
)
# Get rotated CMB map
cmb_map = cmb.get_map(idx=0, apply_beta=True)
Notes
The module uses lenspyx for efficient lensing and supports MPI parallelization for generating large simulation sets.
- synfast_pol(nside, spectra)[source]
Generate polarized CMB maps (Q, U) with optional auxiliary field A.
- class CMB(libdir, nside, model='iso', beta=None, mass=None, Acb=None, lensing=False, Acb_sim_config=None, verbose=True)[source]
Bases:
objectCMB map generator with cosmic birefringence support.
This class handles the generation of CMB temperature and polarization maps with optional cosmic birefringence rotation. Supports three birefringence models:
iso: Isotropic (constant) rotation angle β
iso_td: Time-dependent isotropic model with axion mass parameter
aniso: Anisotropic (spatially-varying) rotation from patchy reionization
- Parameters:
libdir (
str) – Directory for caching simulation products and power spectra.nside (
int) – HEALPix resolution parameter (nside = 2^n).model (
{'iso', 'iso_td', 'aniso'}, default'iso') – Birefringence model type.beta (
float, optional) – Isotropic rotation angle in degrees (required for model=’iso’).mass (
float, optional) – Axion mass parameter in units of 10^-22 eV. Must be 1, 1.5, or 6.4 (required for model=’iso_td’).Acb (
float, optional) – Amplitude of anisotropic birefringence power spectrum (required for model=’aniso’).lensing (
bool, defaultFalse) – Whether to include CMB lensing effects.Acb_sim_config (
dict, optional) – Configuration for anisotropic simulation batches with keys: - ‘alpha_vary_index’: [start, end] indices for varying α - ‘alpha_cons_index’: [start, end] indices for constant α - ‘null_alpha_index’: [start, end] indices for null test (α=0)
Examples
Generate isotropic birefringence simulation:
cmb = CMB( libdir='./cmb_sims', nside=512, model='iso', beta=0.35, lensing=True ) # Get rotated CMB map cmb_map = cmb.get_map(idx=0, apply_beta=True)
Generate anisotropic birefringence:
cmb = CMB( libdir='./cmb_sims', nside=512, model='aniso', Acb=1e-6, lensing=True ) # Get map with spatially-varying rotation cmb_map = cmb.get_map(idx=0, apply_beta=True)
Notes
Uses CAMB for power spectrum generation
Supports efficient caching of realizations
MPI-aware for parallel simulation generation
Lensing implemented via lenspyx
- __init__(libdir, nside, model='iso', beta=None, mass=None, Acb=None, lensing=False, Acb_sim_config=None, verbose=True)[source]
- __validate_Acb_sim_config__()[source]
Validate the Acb_sim_config dictionary structure.
- Return type:
None
- __get_alpha_mode__(idx)[source]
Determine which alpha mode to use for a given index.
Parameters: idx (int): The realization index.
Returns: str: One of ‘vary’, ‘constant’, or ‘null’
- __dl2cl__(arr, unit_only=False)[source]
Convert Dl to Cl.
- Parameters:
arr (numpy.ndarray)
- Return type:
- __td_eb__(dl=True)[source]
Read the E and B mode power spectra from the CMB power spectra data.
- Return type:
- get_power(dl=True)[source]
Get the CMB power spectra.
Parameters: dl (bool): If True, return the power spectra with dl factor else without dl factor.
Returns: Dict[str, np.ndarray]: A dictionary containing the CMB power spectra.
- Parameters:
dl (bool)
- Return type:
- get_unlensed_spectra(dl=True, dtype='d')[source]
Retrieve the unlensed scalar spectra from the power spectrum data.
Parameters: dl (bool, optional): If True, returns Dl (C_l * l * (l + 1) / 2π). Defaults to True. dtype (str, optional): Specifies the format of the returned spectra.
‘d’ returns a dictionary with keys ‘tt’, ‘ee’, ‘bb’, ‘te’.
‘a’ returns the full array of power spectra.
Defaults to ‘d’.
Returns: Union[Dict[str, Any], np.ndarray]:
A dictionary containing individual power spectra for ‘tt’, ‘ee’, ‘bb’, ‘te’ if dtype is ‘d’.
The full array of unlensed scalar power spectra if dtype is ‘a’.
Raises: ValueError: If dtype is not ‘d’ or ‘a’.
- get_lensed_spectra(dl=True, dtype='d')[source]
Retrieve the lensed scalar spectra from the power spectrum data.
Parameters: dl (bool, optional): If True, returns Dl (C_l * l * (l + 1) / 2π). Defaults to True. dtype (str, optional): Specifies the format of the returned spectra.
‘d’ returns a dictionary with keys ‘tt’, ‘ee’, ‘bb’, ‘te’.
‘a’ returns the full array of power spectra.
Defaults to ‘d’.
Returns: Union[Dict[str, Any], np.ndarray]:
A dictionary containing individual power spectra for ‘tt’, ‘ee’, ‘bb’, ‘te’ if dtype is ‘d’.
The full array of lensed scalar power spectra if dtype is ‘a’.
Raises: ValueError: If dtype is not ‘d’ or ‘a’.
- get_cb_unlensed_spectra(beta=0.0, dl=True, dtype='d', new=False)[source]
Calculate the cosmic birefringence (CB) unlensed spectra with a given rotation angle beta.
Parameters: beta (float, optional): The rotation angle in degrees for the cosmic birefringence effect. Defaults to 0.3 degrees. dl (bool, optional): If True, returns Dl (C_l * l * (l + 1) / 2π). Defaults to True. dtype (str, optional): Specifies the format of the returned spectra.
‘d’ returns a dictionary with keys ‘tt’, ‘ee’, ‘bb’, ‘te’, ‘eb’, ‘tb’.
‘a’ returns an array of power spectra.
Defaults to ‘d’.
- new (bool, optional): Determines the ordering of the spectra in the array if dtype is ‘a’.
If True, returns [TT, EE, BB, TE, EB, TB]. If False, returns [TT, TE, TB, EE, EB, BB]. Defaults to False.
Returns: Union[Dict[str, np.ndarray], np.ndarray]:
A dictionary containing individual CB unlensed power spectra for ‘tt’, ‘ee’, ‘bb’, ‘te’, ‘eb’, ‘tb’ if dtype is ‘d’.
An array of CB unlensed power spectra with ordering determined by the new parameter if dtype is ‘a’.
Raises: ValueError: If dtype is not ‘d’ or ‘a’.
Notes: The method applies a rotation by alpha degrees to the E and B mode spectra to account for cosmic birefringence.
- Parameters:
- Return type:
- get_cb_lensed_spectra(beta=0.0, dl=True, dtype='d', new=False)[source]
Calculate the cosmic birefringence (CB) lensed spectra with a given rotation angle beta.
Parameters: beta (float, optional): The rotation angle in degrees for the cosmic birefringence effect. Defaults to 0.3 degrees. dl (bool, optional): If True, returns Dl (C_l * l * (l + 1) / 2π). Defaults to True. dtype (str, optional): Specifies the format of the returned spectra.
‘d’ returns a dictionary with keys ‘tt’, ‘ee’, ‘bb’, ‘te’, ‘eb’, ‘tb’.
‘a’ returns an array of power spectra.
Defaults to ‘d’.
- new (bool, optional): Determines the ordering of the spectra in the array if dtype is ‘a’.
If True, returns [TT, EE, BB, TE, EB, TB]. If False, returns [TT, TE, TB, EE, EB, BB]. Defaults to False.
Returns: Union[Dict[str, np.ndarray], np.ndarray]:
A dictionary containing individual CB lensed power spectra for ‘tt’, ‘ee’, ‘bb’, ‘te’, ‘eb’, ‘tb’ if dtype is ‘d’.
An array of CB lensed power spectra with ordering determined by the new parameter if dtype is ‘a’.
Raises: ValueError: If dtype is not ‘d’ or ‘a’.
Notes: The method applies a rotation by alpha degrees to the E and B mode spectra to account for cosmic birefringence.
- Parameters:
- Return type:
- get_iso_const_cb_gaussian_lensed_QU(idx)[source]
Generate or retrieve the Q and U Stokes parameters after applying cosmic birefringence.
Parameters: idx (int): Index for the realization of the CMB map.
Returns: List[np.ndarray]: A list containing the Q and U Stokes parameter maps as NumPy arrays.
Notes: The method applies a rotation to the E and B mode spherical harmonics to simulate the effect of cosmic birefringence. If the map for the given idx exists in the specified directory, it reads the map from the file. Otherwise, it generates the Q and U maps, applies the birefringence, and saves the resulting map to a FITS file.
- Parameters:
idx (int)
- Return type:
- alpha_alm(idx)[source]
Generate the alpha alm for the anisotropic model.
Parameters: idx (int): Index for the realization of the CMB map.
Returns: np.ndarray: The alpha alm as a NumPy array.
Notes: The method generates the alpha alm for the anisotropic model. The alpha alm is generated as a random realization of the Cl_AA power spectrum. The behavior depends on Acb_sim_config: - ‘vary’: Each idx gets a unique seed (default) - ‘constant’: All indices in this range use the same fixed seed - ‘null’: Returns zeros (no rotation)
- Parameters:
idx (int)
- Return type:
- alpha_map(idx)[source]
Generate a map of the rotation angle alpha for the anisotropic model.
Parameters: idx (int): Index for the realization of the CMB map.
Returns: np.ndarray: A map of the rotation angle alpha as a NumPy array.
Notes: The method generates a map of the rotation angle alpha for the anisotropic model. The map is generated as a random realization of the Cl_AA power spectrum.
- Parameters:
idx (int)
- Return type:
Classes
- class CMB(libdir, nside, model='iso', beta=None, mass=None, Acb=None, lensing=False, Acb_sim_config=None, verbose=True)[source]
Bases:
objectCMB map generator with cosmic birefringence support.
This class handles the generation of CMB temperature and polarization maps with optional cosmic birefringence rotation. Supports three birefringence models:
iso: Isotropic (constant) rotation angle β
iso_td: Time-dependent isotropic model with axion mass parameter
aniso: Anisotropic (spatially-varying) rotation from patchy reionization
- Parameters:
libdir (
str) – Directory for caching simulation products and power spectra.nside (
int) – HEALPix resolution parameter (nside = 2^n).model (
{'iso', 'iso_td', 'aniso'}, default'iso') – Birefringence model type.beta (
float, optional) – Isotropic rotation angle in degrees (required for model=’iso’).mass (
float, optional) – Axion mass parameter in units of 10^-22 eV. Must be 1, 1.5, or 6.4 (required for model=’iso_td’).Acb (
float, optional) – Amplitude of anisotropic birefringence power spectrum (required for model=’aniso’).lensing (
bool, defaultFalse) – Whether to include CMB lensing effects.Acb_sim_config (
dict, optional) – Configuration for anisotropic simulation batches with keys: - ‘alpha_vary_index’: [start, end] indices for varying α - ‘alpha_cons_index’: [start, end] indices for constant α - ‘null_alpha_index’: [start, end] indices for null test (α=0)
Examples
Generate isotropic birefringence simulation:
cmb = CMB( libdir='./cmb_sims', nside=512, model='iso', beta=0.35, lensing=True ) # Get rotated CMB map cmb_map = cmb.get_map(idx=0, apply_beta=True)
Generate anisotropic birefringence:
cmb = CMB( libdir='./cmb_sims', nside=512, model='aniso', Acb=1e-6, lensing=True ) # Get map with spatially-varying rotation cmb_map = cmb.get_map(idx=0, apply_beta=True)
Notes
Uses CAMB for power spectrum generation
Supports efficient caching of realizations
MPI-aware for parallel simulation generation
Lensing implemented via lenspyx
- __init__(libdir, nside, model='iso', beta=None, mass=None, Acb=None, lensing=False, Acb_sim_config=None, verbose=True)[source]
- __validate_Acb_sim_config__()[source]
Validate the Acb_sim_config dictionary structure.
- Return type:
None
- __get_alpha_mode__(idx)[source]
Determine which alpha mode to use for a given index.
Parameters: idx (int): The realization index.
Returns: str: One of ‘vary’, ‘constant’, or ‘null’
- __dl2cl__(arr, unit_only=False)[source]
Convert Dl to Cl.
- Parameters:
arr (numpy.ndarray)
- Return type:
- __td_eb__(dl=True)[source]
Read the E and B mode power spectra from the CMB power spectra data.
- Return type:
- get_power(dl=True)[source]
Get the CMB power spectra.
Parameters: dl (bool): If True, return the power spectra with dl factor else without dl factor.
Returns: Dict[str, np.ndarray]: A dictionary containing the CMB power spectra.
- Parameters:
dl (bool)
- Return type:
- get_unlensed_spectra(dl=True, dtype='d')[source]
Retrieve the unlensed scalar spectra from the power spectrum data.
Parameters: dl (bool, optional): If True, returns Dl (C_l * l * (l + 1) / 2π). Defaults to True. dtype (str, optional): Specifies the format of the returned spectra.
‘d’ returns a dictionary with keys ‘tt’, ‘ee’, ‘bb’, ‘te’.
‘a’ returns the full array of power spectra.
Defaults to ‘d’.
Returns: Union[Dict[str, Any], np.ndarray]:
A dictionary containing individual power spectra for ‘tt’, ‘ee’, ‘bb’, ‘te’ if dtype is ‘d’.
The full array of unlensed scalar power spectra if dtype is ‘a’.
Raises: ValueError: If dtype is not ‘d’ or ‘a’.
- get_lensed_spectra(dl=True, dtype='d')[source]
Retrieve the lensed scalar spectra from the power spectrum data.
Parameters: dl (bool, optional): If True, returns Dl (C_l * l * (l + 1) / 2π). Defaults to True. dtype (str, optional): Specifies the format of the returned spectra.
‘d’ returns a dictionary with keys ‘tt’, ‘ee’, ‘bb’, ‘te’.
‘a’ returns the full array of power spectra.
Defaults to ‘d’.
Returns: Union[Dict[str, Any], np.ndarray]:
A dictionary containing individual power spectra for ‘tt’, ‘ee’, ‘bb’, ‘te’ if dtype is ‘d’.
The full array of lensed scalar power spectra if dtype is ‘a’.
Raises: ValueError: If dtype is not ‘d’ or ‘a’.
- get_cb_unlensed_spectra(beta=0.0, dl=True, dtype='d', new=False)[source]
Calculate the cosmic birefringence (CB) unlensed spectra with a given rotation angle beta.
Parameters: beta (float, optional): The rotation angle in degrees for the cosmic birefringence effect. Defaults to 0.3 degrees. dl (bool, optional): If True, returns Dl (C_l * l * (l + 1) / 2π). Defaults to True. dtype (str, optional): Specifies the format of the returned spectra.
‘d’ returns a dictionary with keys ‘tt’, ‘ee’, ‘bb’, ‘te’, ‘eb’, ‘tb’.
‘a’ returns an array of power spectra.
Defaults to ‘d’.
- new (bool, optional): Determines the ordering of the spectra in the array if dtype is ‘a’.
If True, returns [TT, EE, BB, TE, EB, TB]. If False, returns [TT, TE, TB, EE, EB, BB]. Defaults to False.
Returns: Union[Dict[str, np.ndarray], np.ndarray]:
A dictionary containing individual CB unlensed power spectra for ‘tt’, ‘ee’, ‘bb’, ‘te’, ‘eb’, ‘tb’ if dtype is ‘d’.
An array of CB unlensed power spectra with ordering determined by the new parameter if dtype is ‘a’.
Raises: ValueError: If dtype is not ‘d’ or ‘a’.
Notes: The method applies a rotation by alpha degrees to the E and B mode spectra to account for cosmic birefringence.
- Parameters:
- Return type:
- get_cb_lensed_spectra(beta=0.0, dl=True, dtype='d', new=False)[source]
Calculate the cosmic birefringence (CB) lensed spectra with a given rotation angle beta.
Parameters: beta (float, optional): The rotation angle in degrees for the cosmic birefringence effect. Defaults to 0.3 degrees. dl (bool, optional): If True, returns Dl (C_l * l * (l + 1) / 2π). Defaults to True. dtype (str, optional): Specifies the format of the returned spectra.
‘d’ returns a dictionary with keys ‘tt’, ‘ee’, ‘bb’, ‘te’, ‘eb’, ‘tb’.
‘a’ returns an array of power spectra.
Defaults to ‘d’.
- new (bool, optional): Determines the ordering of the spectra in the array if dtype is ‘a’.
If True, returns [TT, EE, BB, TE, EB, TB]. If False, returns [TT, TE, TB, EE, EB, BB]. Defaults to False.
Returns: Union[Dict[str, np.ndarray], np.ndarray]:
A dictionary containing individual CB lensed power spectra for ‘tt’, ‘ee’, ‘bb’, ‘te’, ‘eb’, ‘tb’ if dtype is ‘d’.
An array of CB lensed power spectra with ordering determined by the new parameter if dtype is ‘a’.
Raises: ValueError: If dtype is not ‘d’ or ‘a’.
Notes: The method applies a rotation by alpha degrees to the E and B mode spectra to account for cosmic birefringence.
- Parameters:
- Return type:
- get_iso_const_cb_gaussian_lensed_QU(idx)[source]
Generate or retrieve the Q and U Stokes parameters after applying cosmic birefringence.
Parameters: idx (int): Index for the realization of the CMB map.
Returns: List[np.ndarray]: A list containing the Q and U Stokes parameter maps as NumPy arrays.
Notes: The method applies a rotation to the E and B mode spherical harmonics to simulate the effect of cosmic birefringence. If the map for the given idx exists in the specified directory, it reads the map from the file. Otherwise, it generates the Q and U maps, applies the birefringence, and saves the resulting map to a FITS file.
- Parameters:
idx (int)
- Return type:
- alpha_alm(idx)[source]
Generate the alpha alm for the anisotropic model.
Parameters: idx (int): Index for the realization of the CMB map.
Returns: np.ndarray: The alpha alm as a NumPy array.
Notes: The method generates the alpha alm for the anisotropic model. The alpha alm is generated as a random realization of the Cl_AA power spectrum. The behavior depends on Acb_sim_config: - ‘vary’: Each idx gets a unique seed (default) - ‘constant’: All indices in this range use the same fixed seed - ‘null’: Returns zeros (no rotation)
- Parameters:
idx (int)
- Return type:
- alpha_map(idx)[source]
Generate a map of the rotation angle alpha for the anisotropic model.
Parameters: idx (int): Index for the realization of the CMB map.
Returns: np.ndarray: A map of the rotation angle alpha as a NumPy array.
Notes: The method generates a map of the rotation angle alpha for the anisotropic model. The map is generated as a random realization of the Cl_AA power spectrum.
- Parameters:
idx (int)
- Return type:
Functions
Foreground Module
Foreground Simulation Module
This module handles Galactic foreground emission modeling for CMB analysis, including dust and synchrotron radiation with frequency-dependent properties.
The module provides:
Multi-frequency foreground map generation
Bandpass integration for realistic instrument response
PySM3 sky model integration (d1, s1 models)
Spatial templates for dust and synchrotron
Coordinate transformation support
Classes
- BandpassInt
Handles bandpass profile integration for multi-frequency observations.
- Foreground
Main class for generating dust and synchrotron foreground maps.
Example
Generate foreground maps for multiple frequencies:
from cobi.simulation import Foreground
import numpy as np
fg = Foreground(
libdir='./fg_sims',
nside=512,
freqs=np.array([27, 39, 93, 145, 225, 280]), # GHz
models=['d1', 's1'], # dust and synchrotron
coord='C' # Celestial coordinates
)
# Get foreground realization
fg_maps = fg.get_map(idx=0) # Shape: (nfreqs, 3, npix)
Notes
Uses PySM3 for generating foreground templates with realistic spectral energy distributions. Supports bandpass integration for accurate modeling of detector response.
- class BandpassInt(libdir)[source]
Bases:
object- Parameters:
libdir (str)
- __init__(libdir)[source]
Initializes the BandpassInt class, loading bandpass profiles from a specified file.
- Parameters:
libdir (str)
- get_profile(band)[source]
Retrieves the frequency and bandpass profile for a specified band.
Parameters: band (str): The frequency band for which to retrieve the profile.
Returns: Tuple[np.ndarray, np.ndarray]: A tuple containing two NumPy arrays:
nu: Array of frequencies (GHz) where the bandpass is defined.
bp: Array of bandpass values corresponding to the frequencies.
- Parameters:
band (str)
- Return type:
- class Foreground(libdir, nside, dust_model, sync_model, bandpass=False, verbose=True)[source]
Bases:
objectGalactic foreground generator for CMB analysis.
This class generates multi-frequency Galactic foreground emission maps using PySM3 models for dust and synchrotron radiation. Supports bandpass integration for realistic detector response.
- Parameters:
libdir (
str) – Directory for caching foreground maps.nside (
int) – HEALPix resolution parameter.dust_model (
int) – PySM3 dust model number (e.g., 1 for d1, 6 for d6).sync_model (
int) – PySM3 synchrotron model number (e.g., 1 for s1, 3 for s3).bandpass (
bool, defaultFalse) – Whether to apply bandpass integration using detector profiles.
- bp_profile
Bandpass profile handler (if enabled).
- Type:
BandpassIntorNone
- dustQU(band)[source]
Generate dust Q/U polarization maps at given frequency.
- Parameters:
band (str)
- Return type:
- syncQU(band)[source]
Generate synchrotron Q/U polarization maps at given frequency.
- Parameters:
band (str)
- Return type:
- get_map(freqs, idx)
Get combined foreground map for multiple frequencies.
Examples
Generate foreground maps without bandpass:
fg = Foreground( libdir='./fg_sims', nside=512, dust_model=1, sync_model=1, bandpass=False ) # Get dust emission at 145 GHz dust_map = fg.dustQU('145')
With bandpass integration:
fg = Foreground( libdir='./fg_sims', nside=512, dust_model=1, sync_model=1, bandpass=True ) # Generate multi-frequency maps freqs = ['027', '039', '093', '145', '225', '280'] fg_maps = fg.get_map(freqs, idx=0)
Notes
Uses PySM3 for physically-motivated foreground models
Caches maps to disk for efficient reuse
Supports MPI parallelization
Output in μK_CMB units
- __init__(libdir, nside, dust_model, sync_model, bandpass=False, verbose=True)[source]
Initializes the Foreground class for generating and handling dust and synchrotron foreground maps.
Parameters: libdir (str): Directory where the foreground maps will be stored. nside (int): HEALPix resolution parameter. dust_model (int): Model number for the dust emission. sync_model (int): Model number for the synchrotron emission. bandpass (bool, optional): If True, bandpass integration is applied. Defaults to False.
- dustQU(band)[source]
Generates or retrieves the Q and U Stokes parameters for dust emission at a given frequency band.
Parameters: band (str): The frequency band.
Returns: np.ndarray: A NumPy array containing the Q and U maps.
- Parameters:
band (str)
- Return type:
- class HILC[source]
Bases:
objectThis class is used to perform Harmonic ILC on the foregrounds.
- harmonic_ilc_alm(alms, lbins=None)[source]
This method is used to perform Harmonic ILC on the foregrounds.
- empirical_harmonic_covariance(alms)[source]
The method empirical_harmonic_covariance is used to compute the empirical harmonic covariance.
- Parameters:
alms (np.ndarray) – The foregrounds in alm format.
- Returns:
The empirical harmonic covariance.
- Return type:
np.ndarray
Classes
- class BandpassInt(libdir)[source]
Bases:
object- Parameters:
libdir (str)
- __init__(libdir)[source]
Initializes the BandpassInt class, loading bandpass profiles from a specified file.
- Parameters:
libdir (str)
- get_profile(band)[source]
Retrieves the frequency and bandpass profile for a specified band.
Parameters: band (str): The frequency band for which to retrieve the profile.
Returns: Tuple[np.ndarray, np.ndarray]: A tuple containing two NumPy arrays:
nu: Array of frequencies (GHz) where the bandpass is defined.
bp: Array of bandpass values corresponding to the frequencies.
- Parameters:
band (str)
- Return type:
- class Foreground(libdir, nside, dust_model, sync_model, bandpass=False, verbose=True)[source]
Bases:
objectGalactic foreground generator for CMB analysis.
This class generates multi-frequency Galactic foreground emission maps using PySM3 models for dust and synchrotron radiation. Supports bandpass integration for realistic detector response.
- Parameters:
libdir (
str) – Directory for caching foreground maps.nside (
int) – HEALPix resolution parameter.dust_model (
int) – PySM3 dust model number (e.g., 1 for d1, 6 for d6).sync_model (
int) – PySM3 synchrotron model number (e.g., 1 for s1, 3 for s3).bandpass (
bool, defaultFalse) – Whether to apply bandpass integration using detector profiles.
- bp_profile
Bandpass profile handler (if enabled).
- Type:
BandpassIntorNone
- dustQU(band)[source]
Generate dust Q/U polarization maps at given frequency.
- Parameters:
band (str)
- Return type:
- syncQU(band)[source]
Generate synchrotron Q/U polarization maps at given frequency.
- Parameters:
band (str)
- Return type:
- get_map(freqs, idx)
Get combined foreground map for multiple frequencies.
Examples
Generate foreground maps without bandpass:
fg = Foreground( libdir='./fg_sims', nside=512, dust_model=1, sync_model=1, bandpass=False ) # Get dust emission at 145 GHz dust_map = fg.dustQU('145')
With bandpass integration:
fg = Foreground( libdir='./fg_sims', nside=512, dust_model=1, sync_model=1, bandpass=True ) # Generate multi-frequency maps freqs = ['027', '039', '093', '145', '225', '280'] fg_maps = fg.get_map(freqs, idx=0)
Notes
Uses PySM3 for physically-motivated foreground models
Caches maps to disk for efficient reuse
Supports MPI parallelization
Output in μK_CMB units
- __init__(libdir, nside, dust_model, sync_model, bandpass=False, verbose=True)[source]
Initializes the Foreground class for generating and handling dust and synchrotron foreground maps.
Parameters: libdir (str): Directory where the foreground maps will be stored. nside (int): HEALPix resolution parameter. dust_model (int): Model number for the dust emission. sync_model (int): Model number for the synchrotron emission. bandpass (bool, optional): If True, bandpass integration is applied. Defaults to False.
- dustQU(band)[source]
Generates or retrieves the Q and U Stokes parameters for dust emission at a given frequency band.
Parameters: band (str): The frequency band.
Returns: np.ndarray: A NumPy array containing the Q and U maps.
- Parameters:
band (str)
- Return type:
- class HILC[source]
Bases:
objectThis class is used to perform Harmonic ILC on the foregrounds.
- harmonic_ilc_alm(alms, lbins=None)[source]
This method is used to perform Harmonic ILC on the foregrounds.
- empirical_harmonic_covariance(alms)[source]
The method empirical_harmonic_covariance is used to compute the empirical harmonic covariance.
- Parameters:
alms (np.ndarray) – The foregrounds in alm format.
- Returns:
The empirical harmonic covariance.
- Return type:
np.ndarray
Mask Module
Mask Generation Module
This module provides tools for creating and managing sky masks used in CMB analysis, including galactic masks, point source masks, and survey footprints.
Features:
Multiple mask types (LAT, SAT, CO, point sources, galactic)
Apodization with multiple methods (C1, C2, Gaussian)
Galactic cut options
Combination of multiple mask types
Automatic mask caching for efficiency
Classes
- Mask
Main class for creating and managing sky masks with apodization.
Example
Create a LAT survey mask with galactic cut:
from cobi.simulation import Mask
mask = Mask(
libdir='./masks',
nside=512,
select='lat',
apo_scale=1.0, # degrees
apo_method='C2',
gal_cut=20 # degrees from galactic plane
)
# Get the mask
mask_map = mask.mask
Combine multiple masks:
mask = Mask(
libdir='./masks',
nside=512,
select='lat+gal+ps', # LAT + galactic + point sources
apo_scale=1.0
)
Notes
Masks are cached to disk for reuse. The module supports NaMaster apodization methods for optimal mode-coupling correction in power spectrum estimation.
- class Mask(libdir, nside, select, apo_scale=0.0, apo_method='C2', gal_cut=0, verbose=True)[source]
Bases:
objectSky mask generator and manager for CMB analysis.
This class creates and handles various types of sky masks including survey footprints (LAT/SAT), galactic plane cuts, point source masks, and CO line emission masks. Supports apodization for optimal power spectrum estimation.
- Parameters:
libdir (
str) – Directory for caching mask files.nside (
int) – HEALPix resolution parameter.select (
str) – Mask type selector. Can combine multiple masks with ‘+’: - ‘lat’: LAT survey footprint - ‘sat’: SAT survey footprint - ‘co’: CO line emission mask - ‘ps’: Point source mask - ‘gal’: Galactic plane mask Example: ‘lat+gal+ps’ combines LAT, galactic, and point source masks.apo_scale (
float, default0.0) – Apodization scale in degrees. If 0, no apodization applied.apo_method (
{'C1', 'C2', 'Gaussian'}, default'C2') – Apodization method compatible with NaMaster.gal_cut (
float,int, orstr, default0) – Galactic cut specification: - float < 1: f_sky fraction (e.g., 0.4 for 40% sky) - int > 1: percentage (e.g., 40 for 40% sky) - str: direct percentage ‘40’, ‘60’, ‘70’, ‘80’, ‘90’ Only used when ‘gal’ or ‘GAL’ in select.
- mask
The combined mask array (values 0-1).
- Type:
ndarray
- get_mask()
Returns the mask array.
Examples
Create LAT mask with galactic cut:
mask = Mask( libdir='./masks', nside=512, select='lat+gal', gal_cut=20, # 20% sky apo_scale=1.0, apo_method='C2' ) mask_array = mask.mask print(f"Sky fraction: {mask.fsky:.3f}")
Simple point source mask:
ps_mask = Mask( libdir='./masks', nside=512, select='ps', apo_scale=0.5 )
Combined mask for full analysis:
full_mask = Mask( libdir='./masks', nside=512, select='lat+gal+ps+co', gal_cut=40, apo_scale=1.0 )
Notes
Masks are cached to disk for efficient reuse
Apodization reduces mode-coupling in power spectra
Multiple masks are combined via multiplication
Compatible with NaMaster (pymaster) workflows
- __init__(libdir, nside, select, apo_scale=0.0, apo_method='C2', gal_cut=0, verbose=True)[source]
Initializes the Mask class for handling and generating sky masks.
Parameters: nside (int): HEALPix resolution parameter. libdir (Optional[str], optional): Directory where the mask may be saved or loaded from. Defaults to None.
- __load_mask_healper__()[source]
Loads a mask from a file.
Returns: np.ndarray: The mask array.
- Return type:
Classes
- class Mask(libdir, nside, select, apo_scale=0.0, apo_method='C2', gal_cut=0, verbose=True)[source]
Bases:
objectSky mask generator and manager for CMB analysis.
This class creates and handles various types of sky masks including survey footprints (LAT/SAT), galactic plane cuts, point source masks, and CO line emission masks. Supports apodization for optimal power spectrum estimation.
- Parameters:
libdir (
str) – Directory for caching mask files.nside (
int) – HEALPix resolution parameter.select (
str) – Mask type selector. Can combine multiple masks with ‘+’: - ‘lat’: LAT survey footprint - ‘sat’: SAT survey footprint - ‘co’: CO line emission mask - ‘ps’: Point source mask - ‘gal’: Galactic plane mask Example: ‘lat+gal+ps’ combines LAT, galactic, and point source masks.apo_scale (
float, default0.0) – Apodization scale in degrees. If 0, no apodization applied.apo_method (
{'C1', 'C2', 'Gaussian'}, default'C2') – Apodization method compatible with NaMaster.gal_cut (
float,int, orstr, default0) – Galactic cut specification: - float < 1: f_sky fraction (e.g., 0.4 for 40% sky) - int > 1: percentage (e.g., 40 for 40% sky) - str: direct percentage ‘40’, ‘60’, ‘70’, ‘80’, ‘90’ Only used when ‘gal’ or ‘GAL’ in select.
- mask
The combined mask array (values 0-1).
- Type:
ndarray
- get_mask()
Returns the mask array.
Examples
Create LAT mask with galactic cut:
mask = Mask( libdir='./masks', nside=512, select='lat+gal', gal_cut=20, # 20% sky apo_scale=1.0, apo_method='C2' ) mask_array = mask.mask print(f"Sky fraction: {mask.fsky:.3f}")
Simple point source mask:
ps_mask = Mask( libdir='./masks', nside=512, select='ps', apo_scale=0.5 )
Combined mask for full analysis:
full_mask = Mask( libdir='./masks', nside=512, select='lat+gal+ps+co', gal_cut=40, apo_scale=1.0 )
Notes
Masks are cached to disk for efficient reuse
Apodization reduces mode-coupling in power spectra
Multiple masks are combined via multiplication
Compatible with NaMaster (pymaster) workflows
- __init__(libdir, nside, select, apo_scale=0.0, apo_method='C2', gal_cut=0, verbose=True)[source]
Initializes the Mask class for handling and generating sky masks.
Parameters: nside (int): HEALPix resolution parameter. libdir (Optional[str], optional): Directory where the mask may be saved or loaded from. Defaults to None.
- __load_mask_healper__()[source]
Loads a mask from a file.
Returns: np.ndarray: The mask array.
- Return type:
Noise Module
Noise Simulation Module
This module simulates realistic instrumental noise for Simons Observatory observations, including white noise, atmospheric noise, and correlated noise between detector pairs.
Features:
SO LAT and SAT noise models (v3.1.1)
Sensitivity modes (baseline, goal)
Atmospheric and 1/f noise components
Cross-frequency correlations
Map-based and harmonic-space realizations
Anisotropic noise from hit maps
Classes
- Noise
Main class for generating noise realizations and spectra.
Functions
- NoiseSpectra
Compute noise power spectra for SO telescopes.
Example
Generate noise maps for LAT baseline sensitivity:
from cobi.simulation import Noise
import numpy as np
noise = Noise(
libdir='./noise_sims',
nside=512,
lmax=3000,
freqs=np.array([27, 39, 93, 145, 225, 280]),
fwhm=np.array([7.4, 5.1, 2.2, 1.4, 1.0, 0.9]),
telescope='LAT',
sensitivity_mode='baseline',
fsky=0.4,
atm_noise=True
)
# Get noise realization
noise_maps = noise.get_map(idx=0) # Shape: (nfreqs, 3, npix)
Anisotropic noise from hit maps:
noise_aniso = Noise(
libdir='./noise_sims',
nside=512,
nside_aniso=4096, # high-res for hit maps
isotropic=False,
hitmap_file='hitmaps.h5'
)
Notes
Uses SO noise models v3.1.1 for realistic noise properties. Supports both isotropic and anisotropic noise based on actual survey hit maps.
- class Noise(nside, fsky, telescope, sim='NC', atm_noise=False, nsplits=2, aso=False, verbose=True)[source]
Bases:
objectInstrumental noise simulator for Simons Observatory telescopes.
This class generates realistic noise realizations for SO LAT and SAT observations based on the v3.1.1 noise models. Supports white noise, atmospheric/1f noise, and correlated noise between detector pairs within the same optical tube.
- Parameters:
nside (
int) – HEALPix resolution parameter.fsky (
float) – Sky fraction for noise power calculation (0 < fsky ≤ 1).telescope (
{'LAT', 'SAT'}) – Telescope identifier (Large Aperture or Small Aperture).sim (
{'NC', 'TOD'}, default'NC') – Simulation type: - ‘NC’: Noise curves from SO models (analytic) - ‘TOD’: Time-ordered data based simulations (uses SO products)atm_noise (
bool, defaultFalse) – Include atmospheric and 1/f noise components.nsplits (
int, default2) – Number of data splits (for split-based null tests).aso (bool)
- get_map(idx, freqs, fwhm, tube)
Generate noise map realization for given configuration.
- get_spectra(idx, freqs, fwhm, tube)
Get noise power spectra for given configuration.
Examples
LAT baseline noise with atmosphere:
noise = Noise( nside=512, fsky=0.4, telescope='LAT', sim='NC', atm_noise=True, nsplits=2 ) # Generate noise for 6 LAT frequencies freqs = np.array([27, 39, 93, 145, 225, 280]) fwhm = np.array([7.4, 5.1, 2.2, 1.4, 1.0, 0.9]) tube = np.array([0, 0, 1, 1, 2, 2]) noise_maps = noise.get_map(idx=0, freqs=freqs, fwhm=fwhm, tube=tube)
SAT white noise only:
noise_sat = Noise( nside=512, fsky=0.1, telescope='SAT', sim='NC', atm_noise=False ) freqs = np.array([27, 39, 93, 145, 225, 280]) noise_maps = noise_sat.get_map(idx=0, freqs=freqs, fwhm=None, tube=None)
Notes
Uses SO v3.1.1 noise models (baseline sensitivity)
Correlated noise between frequency pairs in same tube
Supports independent realizations via idx parameter
MPI-aware for parallel generation
Output in μK_CMB units
- __init__(nside, fsky, telescope, sim='NC', atm_noise=False, nsplits=2, aso=False, verbose=True)[source]
Initializes the Noise class for generating noise maps with or without atmospheric noise.
Parameters: nside (int): HEALPix resolution parameter. atm_noise (bool, optional): If True, includes atmospheric noise. Defaults to False. nhits (bool, optional): If True, includes hit count map. Defaults to False. nsplits (int, optional): Number of data splits to consider. Defaults to 2.
- property rand_alm: numpy.ndarray
Generates random spherical harmonic coefficients (alm) with a specified power spectrum.
Returns: np.ndarray: A complex array of spherical harmonic coefficients.
- property cholesky_matrix_elements: Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray]
Computes the Cholesky matrix elements for the noise model.
Returns: Tuple of nine np.ndarray elements representing the Cholesky matrix components.
- white_noise_maps_freq(freq)[source]
Generates white noise maps for a specific frequency band.
Returns: np.ndarray: An array of white noise maps for the specified frequency band.
- Parameters:
freq (str)
- Return type:
- atm_noise_maps_freq(idx, freq)[source]
Generates atmospheric noise maps using Cholesky decomposition.
Returns: np.ndarray: An array of atmospheric noise maps for different frequency bands.
- Parameters:
- Return type:
- atm_noise_maps(split, idx)[source]
Generates atmospheric noise maps using Cholesky decomposition.
Returns: np.ndarray: An array of atmospheric noise maps for different frequency bands.
- Return type:
- noiseQU_NC(idx)[source]
Generates Q and U polarization noise maps based on the noise model.
Returns: np.ndarray: An array of Q and U noise maps.
- Return type:
- noiseQU_NC_freq(idx, freq)[source]
Generates Q and U polarization noise maps for a specific frequency band.
Returns: np.ndarray: An array of Q and U noise maps for the specified frequency band.
- Parameters:
- Return type:
- noiseQU(idx=None)[source]
Generates Q and U polarization noise maps based on the noise model.
Returns: np.ndarray: An array of Q and U noise maps.
- Parameters:
idx (int | None)
- Return type:
Classes
- class Noise(nside, fsky, telescope, sim='NC', atm_noise=False, nsplits=2, aso=False, verbose=True)[source]
Bases:
objectInstrumental noise simulator for Simons Observatory telescopes.
This class generates realistic noise realizations for SO LAT and SAT observations based on the v3.1.1 noise models. Supports white noise, atmospheric/1f noise, and correlated noise between detector pairs within the same optical tube.
- Parameters:
nside (
int) – HEALPix resolution parameter.fsky (
float) – Sky fraction for noise power calculation (0 < fsky ≤ 1).telescope (
{'LAT', 'SAT'}) – Telescope identifier (Large Aperture or Small Aperture).sim (
{'NC', 'TOD'}, default'NC') – Simulation type: - ‘NC’: Noise curves from SO models (analytic) - ‘TOD’: Time-ordered data based simulations (uses SO products)atm_noise (
bool, defaultFalse) – Include atmospheric and 1/f noise components.nsplits (
int, default2) – Number of data splits (for split-based null tests).aso (bool)
- get_map(idx, freqs, fwhm, tube)
Generate noise map realization for given configuration.
- get_spectra(idx, freqs, fwhm, tube)
Get noise power spectra for given configuration.
Examples
LAT baseline noise with atmosphere:
noise = Noise( nside=512, fsky=0.4, telescope='LAT', sim='NC', atm_noise=True, nsplits=2 ) # Generate noise for 6 LAT frequencies freqs = np.array([27, 39, 93, 145, 225, 280]) fwhm = np.array([7.4, 5.1, 2.2, 1.4, 1.0, 0.9]) tube = np.array([0, 0, 1, 1, 2, 2]) noise_maps = noise.get_map(idx=0, freqs=freqs, fwhm=fwhm, tube=tube)
SAT white noise only:
noise_sat = Noise( nside=512, fsky=0.1, telescope='SAT', sim='NC', atm_noise=False ) freqs = np.array([27, 39, 93, 145, 225, 280]) noise_maps = noise_sat.get_map(idx=0, freqs=freqs, fwhm=None, tube=None)
Notes
Uses SO v3.1.1 noise models (baseline sensitivity)
Correlated noise between frequency pairs in same tube
Supports independent realizations via idx parameter
MPI-aware for parallel generation
Output in μK_CMB units
- __init__(nside, fsky, telescope, sim='NC', atm_noise=False, nsplits=2, aso=False, verbose=True)[source]
Initializes the Noise class for generating noise maps with or without atmospheric noise.
Parameters: nside (int): HEALPix resolution parameter. atm_noise (bool, optional): If True, includes atmospheric noise. Defaults to False. nhits (bool, optional): If True, includes hit count map. Defaults to False. nsplits (int, optional): Number of data splits to consider. Defaults to 2.
- property rand_alm: numpy.ndarray
Generates random spherical harmonic coefficients (alm) with a specified power spectrum.
Returns: np.ndarray: A complex array of spherical harmonic coefficients.
- property cholesky_matrix_elements: Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray]
Computes the Cholesky matrix elements for the noise model.
Returns: Tuple of nine np.ndarray elements representing the Cholesky matrix components.
- white_noise_maps_freq(freq)[source]
Generates white noise maps for a specific frequency band.
Returns: np.ndarray: An array of white noise maps for the specified frequency band.
- Parameters:
freq (str)
- Return type:
- atm_noise_maps_freq(idx, freq)[source]
Generates atmospheric noise maps using Cholesky decomposition.
Returns: np.ndarray: An array of atmospheric noise maps for different frequency bands.
- Parameters:
- Return type:
- atm_noise_maps(split, idx)[source]
Generates atmospheric noise maps using Cholesky decomposition.
Returns: np.ndarray: An array of atmospheric noise maps for different frequency bands.
- Return type:
- noiseQU_NC(idx)[source]
Generates Q and U polarization noise maps based on the noise model.
Returns: np.ndarray: An array of Q and U noise maps.
- Return type:
- noiseQU_NC_freq(idx, freq)[source]
Generates Q and U polarization noise maps for a specific frequency band.
Returns: np.ndarray: An array of Q and U noise maps for the specified frequency band.
- Parameters:
- Return type:
- noiseQU(idx=None)[source]
Generates Q and U polarization noise maps based on the noise model.
Returns: np.ndarray: An array of Q and U noise maps.
- Parameters:
idx (int | None)
- Return type:
Functions
Sky Module
Sky Simulation Module
This module provides a unified interface for generating complete multi-frequency sky simulations combining CMB, foregrounds, and noise with realistic instrumental effects and data processing.
The SkySimulation class orchestrates:
CMB signal with cosmic birefringence
Galactic foregrounds (dust, synchrotron)
Instrumental noise (LAT/SAT models)
Beam convolution
Component separation (ILC)
Survey masks
Coordinate transformations
Classes
- SkySimulation
Main class for end-to-end sky simulation and processing.
- HILC
Harmonic-space Internal Linear Combination for component separation.
Example
Generate a complete LAT simulation:
from cobi.simulation import SkySimulation
import numpy as np
sky = SkySimulation(
libdir='./sky_sims',
nside=512,
freqs=np.array([27, 39, 93, 145, 225, 280]),
fwhm=np.array([7.4, 5.1, 2.2, 1.4, 1.0, 0.9]),
tube=np.array([0, 0, 1, 1, 2, 2]), # tube assignments
cb_model='iso',
beta=0.35,
lensing=True,
sensitivity_mode='baseline',
telescope='LAT',
fsky=0.4
)
# Get full simulation (all components)
full_sky = sky.get_sim(idx=0, component='all')
# Get CMB-only
cmb_only = sky.get_sim(idx=0, component='cmb')
# Get ILC-cleaned map
cleaned = sky.get_ilc(idx=0)
Multi-component analysis:
# Individual components
cmb = sky.get_sim(idx=0, component='cmb')
fg = sky.get_sim(idx=0, component='fg')
noise = sky.get_sim(idx=0, component='noise')
# Verify: all = cmb + fg + noise
all_sky = sky.get_sim(idx=0, component='all')
Notes
This module provides the highest-level interface for simulation generation. It handles caching, parallelization, and ensures consistency across components. Supports both isotropic and anisotropic cosmic birefringence models.
- class SkySimulation(libdir, nside, freqs, fwhm, tube, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, Acb_sim_config=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Bases:
objectUnified multi-frequency sky simulation framework.
This is the highest-level class that orchestrates complete end-to-end simulations combining CMB, foregrounds, noise, and instrumental effects for SO LAT/SAT observations. Includes component separation via harmonic ILC.
- Parameters:
libdir (
str) – Base directory for all simulation products and caching.nside (
int) – HEALPix resolution parameter.freqs (
ndarray) – Array of observation frequencies in GHz (e.g., [27, 39, 93, 145, 225, 280]).fwhm (
ndarray) – Beam FWHM in arcminutes for each frequency.tube (
ndarray) – Optical tube assignment for each frequency (for correlated noise).cb_model (
{'iso', 'iso_td', 'aniso'}, default'iso') – Cosmic birefringence model type.beta (
float, default0.35) – Isotropic birefringence angle in degrees (for cb_model=’iso’).mass (
float, default1.5) – Axion mass in 10^-22 eV units (for cb_model=’iso_td’). Must be 1, 1.5, or 6.4.Acb (
float, default1e-6) – Anisotropic birefringence amplitude (for cb_model=’aniso’).dust_model (
int, default10) – PySM3 dust model number.sync_model (
int, default5) – PySM3 synchrotron model number.bandpass (
bool, defaultTrue) – Apply bandpass integration for foregrounds.alpha (
floatorarray, default0.0) – Miscalibration polarization angle(s) in degrees. - float: same angle for all frequencies - array: per-frequency angles (must match len(freqs))alpha_err (
float, default0.0) – Gaussian uncertainty on miscalibration angle (degrees).miscal_samples (
int, optional) – Number of miscalibration angle realizations. Default: max index from Acb_sim_config or 300.Acb_sim_config (
dict, optional) – Anisotropic birefringence simulation batch configuration. Keys: ‘alpha_vary_index’, ‘alpha_cons_index’, ‘null_alpha_index’ Values: [start, end] index ranges.noise_model (
{'NC', 'TOD'}, default'NC') – Noise simulation type (curves vs. time-ordered data).atm_noise (
bool, defaultTrue) – Include atmospheric/1f noise.nsplits (
int, default2) – Number of data splits for null tests.gal_cut (
int, default0) – Galactic mask cut (0, 40, 60, 70, 80, 90).hilc_bins (
int, default10) – Number of multipole bins for harmonic ILC.deconv_maps (
bool, defaultFalse) – Deconvolve beam from output maps.fldname_suffix (
str, default'') – Additional suffix for output directory naming.sht_backend (
{'healpy', 'ducc0'}, default'healpy') – Spherical harmonic transform backend.aso (bool)
- freqs
Observation frequencies.
- Type:
ndarray
- cmb
CMB simulation object.
- Type:
CMB
- foreground
Foreground simulation object.
- Type:
Foreground
- noise
Noise simulation object.
- Type:
Noise
- mask
Combined survey mask.
- Type:
ndarray
- get_sim(idx, component='all')
Get simulation realization (cmb/fg/noise/all).
- get_ilc(idx, return_noise=False)
Get ILC-cleaned CMB map.
- get_obs(idx, freqs_idx=None)
Get observed maps at specific frequencies.
Examples
Complete LAT simulation pipeline:
from cobi.simulation import SkySimulation import numpy as np sky = SkySimulation( libdir='./sims', nside=512, freqs=np.array([27, 39, 93, 145, 225, 280]), fwhm=np.array([7.4, 5.1, 2.2, 1.4, 1.0, 0.9]), tube=np.array([0, 0, 1, 1, 2, 2]), cb_model='iso', beta=0.35, lensing=True, dust_model=10, sync_model=5, alpha=0.0, # no miscalibration atm_noise=True, gal_cut=40 ) # Get full simulation full_map = sky.get_sim(idx=0, component='all') # Get CMB-only cmb_map = sky.get_sim(idx=0, component='cmb') # Get ILC-cleaned map cleaned_map = sky.get_ilc(idx=0)
Multi-frequency miscalibration:
sky = SkySimulation( libdir='./sims', nside=512, freqs=np.array([93, 145, 225]), fwhm=np.array([2.2, 1.4, 1.0]), tube=np.array([1, 1, 2]), alpha=[0.2, -0.1, 0.15], # per-frequency alpha_err=0.05 # with uncertainty )
Notes
Automatically manages caching and directory structure
MPI-aware for parallel processing
Supports batch processing with Acb_sim_config
ILC uses harmonic-space cleaning
All outputs in μK_CMB units
- __init__(libdir, nside, freqs, fwhm, tube, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, Acb_sim_config=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
freqs (numpy.ndarray)
fwhm (numpy.ndarray)
tube (numpy.ndarray)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha_err (float)
miscal_samples (int | None)
Acb_sim_config (dict | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- class LATsky(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, Acb_sim_config=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Bases:
SkySimulation- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
Acb_sim_config (dict | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- freqs = ['27', '39', '93', '145', '225', '280']
- fwhm = [7.4, 5.1, 2.2, 1.4, 1.0, 0.9]
- tube = ['LF', 'LF', 'MF', 'MF', 'HF', 'HF']
- __init__(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, Acb_sim_config=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
Acb_sim_config (dict | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- class SATsky(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Bases:
SkySimulation- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- freqs = ['27', '39', '93', '145', '225', '280']
- fwhm = [91, 63, 30, 17, 11, 9]
- tube = ['S1', 'S1', 'S2', 'S2', 'S3', 'S3']
- __init__(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- class LATskyC(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', verbose=True)[source]
Bases:
SkySimulation- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
verbose (bool)
- freqs = ['93', '145']
- fwhm = [2.2, 1.4]
- tube = ['MF', 'MF']
- __init__(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
verbose (bool)
- class SATskyC(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', verbose=True)[source]
Bases:
SkySimulation- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
verbose (bool)
- freqs = ['93', '145']
- fwhm = [30, 17]
- tube = ['S2', 'S2']
- __init__(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
verbose (bool)
Classes
- class SkySimulation(libdir, nside, freqs, fwhm, tube, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, Acb_sim_config=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Bases:
objectUnified multi-frequency sky simulation framework.
This is the highest-level class that orchestrates complete end-to-end simulations combining CMB, foregrounds, noise, and instrumental effects for SO LAT/SAT observations. Includes component separation via harmonic ILC.
- Parameters:
libdir (
str) – Base directory for all simulation products and caching.nside (
int) – HEALPix resolution parameter.freqs (
ndarray) – Array of observation frequencies in GHz (e.g., [27, 39, 93, 145, 225, 280]).fwhm (
ndarray) – Beam FWHM in arcminutes for each frequency.tube (
ndarray) – Optical tube assignment for each frequency (for correlated noise).cb_model (
{'iso', 'iso_td', 'aniso'}, default'iso') – Cosmic birefringence model type.beta (
float, default0.35) – Isotropic birefringence angle in degrees (for cb_model=’iso’).mass (
float, default1.5) – Axion mass in 10^-22 eV units (for cb_model=’iso_td’). Must be 1, 1.5, or 6.4.Acb (
float, default1e-6) – Anisotropic birefringence amplitude (for cb_model=’aniso’).dust_model (
int, default10) – PySM3 dust model number.sync_model (
int, default5) – PySM3 synchrotron model number.bandpass (
bool, defaultTrue) – Apply bandpass integration for foregrounds.alpha (
floatorarray, default0.0) – Miscalibration polarization angle(s) in degrees. - float: same angle for all frequencies - array: per-frequency angles (must match len(freqs))alpha_err (
float, default0.0) – Gaussian uncertainty on miscalibration angle (degrees).miscal_samples (
int, optional) – Number of miscalibration angle realizations. Default: max index from Acb_sim_config or 300.Acb_sim_config (
dict, optional) – Anisotropic birefringence simulation batch configuration. Keys: ‘alpha_vary_index’, ‘alpha_cons_index’, ‘null_alpha_index’ Values: [start, end] index ranges.noise_model (
{'NC', 'TOD'}, default'NC') – Noise simulation type (curves vs. time-ordered data).atm_noise (
bool, defaultTrue) – Include atmospheric/1f noise.nsplits (
int, default2) – Number of data splits for null tests.gal_cut (
int, default0) – Galactic mask cut (0, 40, 60, 70, 80, 90).hilc_bins (
int, default10) – Number of multipole bins for harmonic ILC.deconv_maps (
bool, defaultFalse) – Deconvolve beam from output maps.fldname_suffix (
str, default'') – Additional suffix for output directory naming.sht_backend (
{'healpy', 'ducc0'}, default'healpy') – Spherical harmonic transform backend.aso (bool)
- freqs
Observation frequencies.
- Type:
ndarray
- cmb
CMB simulation object.
- Type:
CMB
- foreground
Foreground simulation object.
- Type:
Foreground
- noise
Noise simulation object.
- Type:
Noise
- mask
Combined survey mask.
- Type:
ndarray
- get_sim(idx, component='all')
Get simulation realization (cmb/fg/noise/all).
- get_ilc(idx, return_noise=False)
Get ILC-cleaned CMB map.
- get_obs(idx, freqs_idx=None)
Get observed maps at specific frequencies.
Examples
Complete LAT simulation pipeline:
from cobi.simulation import SkySimulation import numpy as np sky = SkySimulation( libdir='./sims', nside=512, freqs=np.array([27, 39, 93, 145, 225, 280]), fwhm=np.array([7.4, 5.1, 2.2, 1.4, 1.0, 0.9]), tube=np.array([0, 0, 1, 1, 2, 2]), cb_model='iso', beta=0.35, lensing=True, dust_model=10, sync_model=5, alpha=0.0, # no miscalibration atm_noise=True, gal_cut=40 ) # Get full simulation full_map = sky.get_sim(idx=0, component='all') # Get CMB-only cmb_map = sky.get_sim(idx=0, component='cmb') # Get ILC-cleaned map cleaned_map = sky.get_ilc(idx=0)
Multi-frequency miscalibration:
sky = SkySimulation( libdir='./sims', nside=512, freqs=np.array([93, 145, 225]), fwhm=np.array([2.2, 1.4, 1.0]), tube=np.array([1, 1, 2]), alpha=[0.2, -0.1, 0.15], # per-frequency alpha_err=0.05 # with uncertainty )
Notes
Automatically manages caching and directory structure
MPI-aware for parallel processing
Supports batch processing with Acb_sim_config
ILC uses harmonic-space cleaning
All outputs in μK_CMB units
- __init__(libdir, nside, freqs, fwhm, tube, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, Acb_sim_config=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
freqs (numpy.ndarray)
fwhm (numpy.ndarray)
tube (numpy.ndarray)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha_err (float)
miscal_samples (int | None)
Acb_sim_config (dict | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- class LATsky(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, Acb_sim_config=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Bases:
SkySimulation- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
Acb_sim_config (dict | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- freqs = ['27', '39', '93', '145', '225', '280']
- fwhm = [7.4, 5.1, 2.2, 1.4, 1.0, 0.9]
- tube = ['LF', 'LF', 'MF', 'MF', 'HF', 'HF']
- __init__(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, Acb_sim_config=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
Acb_sim_config (dict | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- class SATsky(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Bases:
SkySimulation- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- freqs = ['27', '39', '93', '145', '225', '280']
- fwhm = [91, 63, 30, 17, 11, 9]
- tube = ['S1', 'S1', 'S2', 'S2', 'S3', 'S3']
- __init__(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', sht_backend='healpy', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
sht_backend (str)
verbose (bool)
- class LATskyC(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', verbose=True)[source]
Bases:
SkySimulation- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
verbose (bool)
- freqs = ['93', '145']
- fwhm = [2.2, 1.4]
- tube = ['MF', 'MF']
- __init__(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', aso=False, atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
aso (bool)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
verbose (bool)
- class SATskyC(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', verbose=True)[source]
Bases:
SkySimulation- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
verbose (bool)
- freqs = ['93', '145']
- fwhm = [30, 17]
- tube = ['S2', 'S2']
- __init__(libdir, nside, cb_model='iso', beta=0.35, mass=1.5, Acb=1e-06, lensing=True, dust_model=10, sync_model=5, bandpass=True, alpha=0.0, alpha_err=0.0, miscal_samples=None, noise_model='NC', atm_noise=True, nsplits=2, gal_cut=0, hilc_bins=10, deconv_maps=False, fldname_suffix='', verbose=True)[source]
Initialize the SkySimulation.
See class docstring for detailed parameter descriptions.
- Parameters:
libdir (str)
nside (int)
cb_model (str)
beta (float)
mass (float)
Acb (float)
lensing (bool)
dust_model (int)
sync_model (int)
bandpass (bool)
alpha (float)
alpha_err (float)
miscal_samples (int | None)
noise_model (str)
atm_noise (bool)
nsplits (int)
gal_cut (int)
hilc_bins (int)
deconv_maps (bool)
fldname_suffix (str)
verbose (bool)