Sky model

The Model class is the main interface to cosmoglobe. The __call__ method of the SkyModel object returns the simulated sky emission.

class cosmoglobe.sky.model.SkyModel(nside, components, version=None)[source]

Sky model representing an initialized instance of the Cosmoglobe Sky Model.

This class acts as a container for the various components making up the Cosmoglobe Sky Model and provides methods to simulate the sky. The primary use case of this class is to call its __call__ method, which simulates the sky at a single frequency \(\nu\), or integrated over a bandpass \(\tau\).

Examples

>>> import cosmoglobe
>>> chain = cosmoglobe.get_test_chain()
>>> model = cosmoglobe.SkyModel.from_chain(chain, nside=256)
>>> print(model)
SkyModel(
    version: BeyondPlanck
    nside: 256
    components(
        (ame): SpinningDust(freq_peak)
        (cmb): CMB()
        (dust): ModifiedBlackbody(beta, T)
        (ff): LinearOpticallyThin(T_e)
        (radio): AGNPowerLaw(alpha)
        (synch): PowerLaw(beta)
    )
)

Simulating the full sky emission at some frequency, given a beam FWHM:

>>> import astropy.units as u
>>> model(50*u.GHz, fwhm=30*u.arcmin)
[[ 2.25809786e+03  2.24380103e+03  2.25659060e+03 ... -2.34783682e+03
  -2.30464421e+03 -2.30387946e+03]
 [-1.64627550e+00  2.93583564e-01 -1.06788937e+00 ... -1.64354585e+01
   1.60621841e+01 -1.05506092e+01]
 [-4.15682825e+00  3.08881971e-01 -1.02012415e+00 ...  5.44745701e+00
  -4.71776995e+00  4.39850830e+00]] uK_RJ
__call__(freqs, weights=None, *, components=None, fwhm=<Quantity 0. arcmin>, output_unit='uK_RJ')[source]

Simulates and returns the full sky model emission.

Parameters
  • freqs (Quantity) – A frequency, or a list of frequencies corresponding to the bandpass weights.

  • weights (Optional[Quantity]) – Bandpass weights corresponding to the frequencies in freqs. The units of this quantity must match the units of the related detector, i.e, they must be compatible with K_RJ, K_CMB, or Jy/sr. Even if the bandpass weights are already normalized and now have units of 1/Hz, redefine the quantities unites to be that of the detector. If bandpass is None and freqs is a single frequency, a delta peak is assumed. Defaults to None.

  • components (Optional[List[str]]) – List of component labels. If None, all components in the sky model is included. Defaults to None.

  • fwhm (Quantity) – The full width half max parameter of the Gaussian. Defaults to 0.0, which indicates no smoothing of output maps.

  • output_unit (str) – The output units of the emission. For instance ‘uK_RJ’, ‘uK_CMB’, or ‘MJy/sr’. If Kelvin, Rayleigh-Jeans or CMB must be specified.

Returns

The simulated emission given the Cosmoglobe Sky Model.

Return type

emission

__init__(nside, components, version=None)[source]

Initializes an instance of the Cosmoglobe Sky Model.

Parameters
  • nside (int) – Healpix resolution of the maps in sky model.

  • components (Dict[str, SkyComponent]) – A list of pre-initialized `SkyComponent`to include in the model.

  • version (Optional[str]) – The version of the Cosmoglobe Model used in the sky model.

classmethod from_chain(chain, nside, components=None, model='BeyondPlanck', samples=- 1, burn_in=None)[source]

Initializes the SkyModel from a Cosmoglobe chain.

Parameters
  • chain (Union[str, Chain]) – Path to a Cosmoglobe chainfile or a Chain object.

  • nside (int) – Model HEALPIX map resolution parameter.

  • components (Optional[List[str]]) – List of components to include in the model.

  • model (str) – String representing which sky model to use. Defaults to BeyondPlanck.

  • samples (Union[range, int, Literal[‘all’]]) – The sample number for which to extract the model. If the input is ‘all’, then the model will an average of all samples in the chain. Defaults to the last sample in the chain.

  • burn_in (Optional[int]) – Burn in sample for which all previous samples are disregarded.

Returns

Initialized sky model.

Return type

sky_model

classmethod from_hub()[source]

Initializes the SkyModel from a public cosmoglobe data release.

Return type

SkyModel

property nside: int

HEALPIX map resolution of the maps in the Sky Model.

Return type

int

remove_dipole(gal_cut=<Quantity 10. deg>, return_dipole=False)[source]

Removes the Solar dipole (and monopole), from the CMB sky component.

Parameters
  • gal_cut (Quantity) – Ignores pixels in range [-gal_cut;+gal_cut] when fitting and subtracting the dipole.

  • return_dipole (bool) – If True, returns the dipole map.

Returns

The dipole of the CMB sky component.

Return type

dipole