CMBDipole¶
- class pysm3.CMBDipole(nside: int, amp, T_cmb, dip_lon, dip_lat, map_dist=None, quadrupole_correction: bool = False)[source] [edit on github]¶
Bases:
objectSimulate the CMB dipole anisotropy as a full-sky HEALPix map.
The dipole is modeled as a temperature fluctuation due to the observer’s motion with respect to the CMB rest frame, following the relativistic Doppler effect.
- Parameters:
- nsideint
HEALPix NSIDE parameter for the output map (dimensionless).
- ampfloat
Amplitude of the dipole in micro-Kelvin (uK_CMB).
- T_cmbfloat
CMB monopole temperature in Kelvin (K_CMB).
- dip_lonfloat
Galactic longitude of the dipole direction, in degrees (deg).
- dip_latfloat
Galactic latitude of the dipole direction, in degrees (deg).
- Returns:
- dipole_mapastropy.units.Quantity
Full-sky HEALPix map (array) of the dipole temperature anisotropy in Kelvin.
Notes
The dipole amplitude is calculated as:
ΔT/T = (v/c) * cos(θ)
where θ is the angle between the dipole direction and the line of sight.
References
For the best-fit parameters, see Table 1 of “Planck intermediate results. LVII. Joint Planck LFI and HFI data processing” https://arxiv.org/pdf/2007.04997.pdf
Methods Summary
get_emission(freqs[, weights])Return the dipole emission map, integrating over the bandpass if needed.
Methods Documentation
- get_emission(freqs: Annotated[Quantity, Unit('GHz')], weights=None) Annotated[Quantity, Unit('uK_RJ')][source] [edit on github]¶
Return the dipole emission map, integrating over the bandpass if needed.
- Parameters:
- freqsQuantity
Frequency or array of frequencies (for bandpass integration).
- weightsarray-like, optional
Integration weights for the bandpass.
- Returns:
- dipole_mapastropy.units.Quantity
Full-sky HEALPix map (array) of the dipole temperature anisotropy in uK_RJ.