On-the-fly simulations in CAR pixelization

We'll now switch to working with maps that are in the Plate-Caree (CAR) projection stored in 2d grids of rectangular pixels. The 2d rectangular aspect allows us to work with maps as numpy arrays, which is especially efficient at storing partial skies. The CAR projection aspect means that spherical harmonic transforms are possible through the libsharp library.

We use the pixell library to work with these maps. Please watch the tutorial on using pixell here: https://www.simonsfoundation.org/event/so-theory-tutorial/

In [1]:
# We load the mapsims library as before
import mapsims
# pixell.enmap is for general map manipulation
# pixell.curvedsky is for SHT-related operations
# pixell.utils is for general utilities; we use it for units
# pixell.enplot lets us make high-resolution visualizations of the map using the Python Image Library
from pixell import enmap, curvedsky as cs, utils, enplot
import numpy as np
# We'll use healpy to do an alm -> cl operation
import healpy as hp
import matplotlib.pyplot as plt

Noise class recap

The core class for noise simulations is SONoiseSimulator. It wraps around the SO noise forecast code, which provides noise power spectra (including atmospheric 1/f and white noise) assuming uniform coverage in a specified sky fraction. In addition, it knows about the "hit-count" information for various SO detector tubes (produced by the scan strategy / time domain working groups) through interfaces with remote locations for maps that encode these. It caches these hitmaps locally on disk on first use, and in memory during repeated calls for realizations of the noise.

The hit-counts are used to modulate the SO noise forecast power spectra to produce Gaussian realization maps that are inhomogenous.

Let's examine the class to see how to initialize it:

In [2]:
mapsims.noise.SONoiseSimulator?
Init signature:
mapsims.noise.SONoiseSimulator(
    nside=None,
    shape=None,
    wcs=None,
    ell_max=None,
    return_uK_CMB=True,
    sensitivity_mode='baseline',
    apply_beam_correction=False,
    apply_kludge_correction=True,
    homogeneous=False,
    no_power_below_ell=None,
    rolloff_ell=50,
    survey_efficiency=0.2,
    full_covariance=True,
    LA_years=5,
    LA_noise_model='SOLatV3point1',
    elevation=50,
    SA_years=5,
    SA_one_over_f_mode='pessimistic',
    sky_fraction=None,
    cache_hitmaps=True,
    boolean_sky_fraction=False,
    instrument_parameters='simonsobs_instrument_parameters_2020.06',
)
Docstring:      <no docstring>
Init docstring:
Simulate noise maps for Simons Observatory

Simulate the noise power spectrum in spherical harmonics domain and then generate a map
in microK_CMB or microK_RJ (based on return_uK_CMB)

In the constructor, this object calls the published 20180822 noise simulator and generates
the expected noise power spectra for all channels.
Then you need to call the `simulate` method with a channel identifier to create a simulated map.

Parameters
----------

nside : int
    nside of HEALPix map. If None, uses
    rectangular pixel geometry specified through shape and wcs.
shape : tuple of ints
    shape of ndmap array (see pixell.pixell.enmap). Must also specify wcs.
wcs : astropy.wcs.wcs.WCS instance
    World Coordinate System for geometry of map (see pixell.pixell.enmap). Must
    also specify shape.
ell_max : int
    Maximum ell for the angular power spectrum, if not provided set to 3 * nside when using healpix
    or 10000 * (1.0 / pixel_height_arcmin) when using CAR, corresponding roughly to the Nyquist
    frequency.
return_uK_CMB : bool
    True, output is in microK_CMB, False output is in microK_RJ
sensitivity_mode : str
    Value should be threshold, baseline or goal to use predefined sensitivities
apply_beam_correction : bool
    Include the effect of the beam in the noise angular power spectrum
apply_kludge_correction : bool
    If True, reduce the hitcount by a factor of 0.85 to account for not-uniformity in the scanning
homogeneous : bool
    Set to True to generate full-sky maps with no hit-count variation, with noise curves
    corresponding to a survey that covers a sky fraction of sky_fraction (defaults to 1).
no_power_below_ell : int
    The input spectra have significant power at low :math:`\ell`,
    we can zero that power specifying an integer :math:`\ell` value here.
    The power spectra at :math:`\ell < \ell_0` are set to zero.
rolloff_ell : int
    Low ell power damping, see the docstring of
    `so_noise_models.so_models_v3.SO_Noise_Calculator_Public_v3_1_1.rolloff`
survey_efficiency : float
    Fraction of calendar time that may be used to compute map depth.
full_covariance : bool
    Whether or not to include the intra-tube covariance between bands.
    If white noise (atmosphere=False) sims are requested, no
    covariance is included regardless of the value of full_covariance.
LA_years : int
    Total number of years for the Large Aperture telescopes survey
LA_noise_model : str
    Noise model among the ones available in `so_noise_model`, "SOLatV3point1" is default, "SOLatV3" is
    the model released in 2018 which had a bug in the atmosphere contribution
elevation : float
    Elevation of the scans in degrees, the V3.1.1 noise model includes elevation
    dependence for the LAT. This should reproduced original V3 results at the
    reference elevation of 50 degrees.
SA_years : int
    Total number of years for the Small Aperture telescopes survey
SA_one_over_f_mode : {"pessimistic", "optimistic", "none"}
    Correlated noise performance of the detectors on the Small Aperture telescopes
sky_fraction : optional,float
    If homogeneous is True, this sky_fraction is used for the noise curves.
cache_hitmaps : bool
    If True, caches hitmaps.
boolean_sky_fraction: bool
    If True, determines sky fraction based on fraction of hitmap that is zero. If False,
    determines sky_fraction from <Nhits>.
instrument_parameters : Path or str
    See the help of MapSims
File:           /global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/noise.py
Type:           type
Subclasses:     

We see that we need to specify either HEALPIX nside or the shape,wcs pair for a CAR map. Let's request simulated maps with 16 arcmin wide pixels since for this demo we'll work with SAT maps. We ask for our output maps to be projected to a full-sky template with that resolution.

CAR geometry noise maps

In [3]:
shape,wcs = enmap.fullsky_geometry(res=16.0 * utils.arcmin)
noise_sim = mapsims.noise.SONoiseSimulator(shape=shape,wcs=wcs)

We've instantiated the noise simulator object. The main member function we need is simulate, we generates a Gaussian random field with the default inhomogenity and noise power spectra. Let's examine its signature:

In [4]:
noise_sim.simulate?
Signature:
noise_sim.simulate(
    tube,
    output_units='uK_CMB',
    seed=None,
    nsplits=1,
    mask_value=None,
    atmosphere=True,
    hitmap=None,
    white_noise_rms=None,
)
Docstring:
Create a random realization of the noise power spectrum

Parameters
----------

tube : str
    Specify a tube (for SO: ST0-ST3, LT0-LT6) see the `tubes` attribute
output_units : str
    Output unit supported by PySM.units, e.g. uK_CMB or K_RJ
seed : integer or tuple of integers, optional
    Specify a seed. The seed is converted to a tuple if not already
    one and appended to (0,0,6,tube_id) to avoid collisions between
    tubes, with the signal sims and with ACT noise sims, where
    tube_id is the integer ID of the tube.
nsplits : integer, optional
    Number of splits to generate. The splits will have independent noise
    realizations, with noise power scaled by a factor of nsplits, i.e. atmospheric
    noise is assumed to average down with observing time the same way
    the white noise does. By default, only one split (the coadd) is generated.
mask_value : float, optional
    The value to set in masked (unobserved) regions. By default, it uses
    the value in default_mask_value, which for healpix is healpy.UNSEEN
    and for CAR is numpy.nan.
atmosphere : bool, optional
    Whether to include the correlated 1/f from the noise model. This is
    True by default. If it is set to False, then a pure white noise map
    is generated from the white noise power in the noise model, and
    the covariance between arrays is ignored.
hitmap : string or map, optional
    Provide the path to a hitmap to override the default used for
    the tube. You could also provide the hitmap as an array
    directly.
white_noise_rms : float or tuple of floats, optional
    Optionally scale the simulation so that the small-scale limit white noise
    level is white_noise_rms in uK-arcmin (either a single number or
    a pair for the dichroic array).

Returns
-------

output_map : ndarray or ndmap
    Numpy array with the HEALPix or CAR map realization of noise.
    The shape of the returned array is (2,3,nsplits,)+oshape, where
    oshape is (npix,) for HEALPix and (Ny,Nx) for CAR.
    The first dimension of size 2 corresponds to the two different
    bands within a dichroic tube.
    See the `band_id` attribute of the Channel class
    to identify which is the index of a Channel in the array.

    The second dimension corresponds to independent split realizations
    of the noise, e.g. it is 1 for full mission.

    The third dimension corresponds to the three polarization
    Stokes components I,Q,U

    The last dimension is the number of pixels
File:      /global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/noise.py
Type:      method
In [5]:
print(noise_sim.tubes)
print(noise_sim.tubes['ST0'][0].center_frequency)
defaultdict(<class 'list'>, {'LT0': [Channel LT0_UHF1, Channel LT0_UHF2], 'LT1': [Channel LT1_UHF1, Channel LT1_UHF2], 'LT2': [Channel LT2_MFF1, Channel LT2_MFF2], 'LT3': [Channel LT3_MFF1, Channel LT3_MFF2], 'LT4': [Channel LT4_MFS1, Channel LT4_MFS2], 'LT5': [Channel LT5_MFS1, Channel LT5_MFS2], 'LT6': [Channel LT6_LF1, Channel LT6_LF2], 'ST0': [Channel ST0_UHF1, Channel ST0_UHF2], 'ST1': [Channel ST1_MFF1, Channel ST1_MFF2], 'ST2': [Channel ST2_MFS1, Channel ST2_MFS2], 'ST3': [Channel ST3_LF1, Channel ST3_LF2]})
225.7 GHz

We need to specify the name of an SO tube. Let's go with the SAT's ST0 tube. This is a dichroic array containing 150 GHz and 90 GHz channels. We specify a mask value in the unobserved region of zero instead of the default numpy.nan:

In [6]:
tube = 'ST0'
omap = noise_sim.simulate(tube,mask_value=0)
/global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/utils.py:71: UserWarning: Access data from /global/cfs/cdirs/sobs/www/so_mapsims_data/v0.2/car/ST0_UHF1_01_of_20.nominal_telescope_all_time_all_hmap_CAR_12.00_arcmin.fits.gz
  warnings.warn(f"Access data from {full_path}")
/global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/utils.py:71: UserWarning: Access data from /global/cfs/cdirs/sobs/www/so_mapsims_data/v0.2/car/ST0_UHF2_01_of_20.nominal_telescope_all_time_all_hmap_CAR_12.00_arcmin.fits.gz
  warnings.warn(f"Access data from {full_path}")
/global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/noise.py:406: UserWarning: WCS of hitmap with nearest pixel-size is not compatible, so interpolating hitmap
  "WCS of hitmap with nearest pixel-size is not compatible, so interpolating hitmap"
/global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/noise.py:406: UserWarning: WCS of hitmap with nearest pixel-size is not compatible, so interpolating hitmap
  "WCS of hitmap with nearest pixel-size is not compatible, so interpolating hitmap"

Let's examine the shape and wcs of the returned map:

In [7]:
print(omap.shape)
print(omap.wcs)
(2, 1, 3, 676, 1350)
car:{cdelt:[-0.2667,0.2667],crval:[0.1333,0],crpix:[675.5,339]}

Why is the map 5-dimensional? The dimension convention in this noise generator is:

(narray,nsplit,ncomp,Ny,Nx)

Ny and Nx are the usual vertical (declination) and horizontal (RA) pixel lengths of the underlying 2d maps. ncomp=3 since there are three maps in each array for each of the I,Q,U Stokes components. This code allows us to generate indepndent splits of the map, where each split has nsplits larger noise power than specified, with the default being no splitting (nsplits=1). Finally, the leading dimensions of the map always have length narray=2 corresponding to the two correlated dichroic arrays within a tube, in this case the 90 GHz and 150 GHz channels of the LAT's LT3 tube.

Let's focus on the intensity I of the 90 GHz map, and plot it:

In [8]:
tmap1 = omap[0,0,0]
enplot.show(enplot.plot(tmap1,downgrade=1,grid=True,ticks=20,mask=0,colorbar=True,range=3))

That looks like noise. Noice. It's fairly inhomogeneous, with the noise blowing up near the edges. We don't see any stripiness -- because the SO map-based noise model is too simplistic to capture that! We also find that compared to previous versions of the SO map-based noise sims, we're not just seeing noise modes at the largest scales dominate, since the noise power has been rolled off on large scales to simulate what happens in a realistic map-maker run that doesn't converge fully on the largest scales (and the atmospheric noise doesn't keep arbitrarily increasing anyway).

We also see that we waste quite a lot of space in memory storing all those zero pixels. Let's utilize the cut-sky benefits of working with CAR pixels and define a smaller template limited to the -80 to 30 degree declination region, and re-initialize our noise simulator with that.

In [9]:
shape,wcs = enmap.band_geometry(dec_cut=np.asarray((-80,30))*utils.degree,res=16.0 * utils.arcmin)
noise_sim = mapsims.noise.SONoiseSimulator(shape=shape,wcs=wcs)
In [10]:
omap = noise_sim.simulate(tube,mask_value=0)
/global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/utils.py:71: UserWarning: Access data from /global/cfs/cdirs/sobs/www/so_mapsims_data/v0.2/car/ST0_UHF1_01_of_20.nominal_telescope_all_time_all_hmap_CAR_12.00_arcmin.fits.gz
  warnings.warn(f"Access data from {full_path}")
/global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/utils.py:71: UserWarning: Access data from /global/cfs/cdirs/sobs/www/so_mapsims_data/v0.2/car/ST0_UHF2_01_of_20.nominal_telescope_all_time_all_hmap_CAR_12.00_arcmin.fits.gz
  warnings.warn(f"Access data from {full_path}")
/global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/noise.py:406: UserWarning: WCS of hitmap with nearest pixel-size is not compatible, so interpolating hitmap
  "WCS of hitmap with nearest pixel-size is not compatible, so interpolating hitmap"
/global/common/software/sobs/mbs/lib/python3.7/site-packages/mapsims/noise.py:406: UserWarning: WCS of hitmap with nearest pixel-size is not compatible, so interpolating hitmap
  "WCS of hitmap with nearest pixel-size is not compatible, so interpolating hitmap"
In [11]:
tmap1 = omap[0,0,0]
plot = lambda x,**kwargs: enplot.show(enplot.plot(x,downgrade=1,grid=True,ticks=20,mask=0,**kwargs))
plot(tmap1,range=3)