cobi.simulation.CMB
- class CMB(libdir, nside, model='iso', beta=None, mass=None, Acb=None, lensing=False, Acb_sim_config=None, verbose=True)[source]
CMB 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]
Methods
__init__(libdir, nside[, model, beta, mass, ...])alpha_alm(idx)Generate the alpha alm for the anisotropic model.
alpha_map(idx)Generate a map of the rotation angle alpha for the anisotropic model.
cl_aa()Compute the Cl_AA power spectrum for the anisotropic model.
cl_pp()compute the CMB power spectra using CAMB.
get_cb_lensed_QU(idx)get_cb_lensed_mass_spectra([dl, dtype, new])get_cb_lensed_spectra([beta, dl, dtype, new])Calculate the cosmic birefringence (CB) lensed spectra with a given rotation angle beta.
get_cb_unlensed_mass_spectra([dl, dtype, new])get_cb_unlensed_spectra([beta, dl, dtype, new])Calculate the cosmic birefringence (CB) unlensed spectra with a given rotation angle beta.
Generate or retrieve the Q and U Stokes parameters after applying cosmic birefringence.
get_lensed_spectra([dl, dtype])Retrieve the lensed scalar spectra from the power spectrum data.
get_power([dl])Get the CMB power spectra.
get_unlensed_spectra([dl, dtype])Retrieve the unlensed scalar spectra from the power spectrum data.
grad_phi_alm(idx)phi_alm(idx)scale_invariant(Acb)- __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: