SONoiseSimulator

class mapsims.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, homogenous=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')[source] [edit on github]

Bases: object

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
nsideint

nside of HEALPix map. If None, uses rectangular pixel geometry specified through shape and wcs.

shapetuple of ints

shape of ndmap array (see pixell.pixell.enmap). Must also specify wcs.

wcsastropy.wcs.wcs.WCS instance

World Coordinate System for geometry of map (see pixell.pixell.enmap). Must also specify shape.

ell_maxint

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_CMBbool

True, output is in microK_CMB, False output is in microK_RJ

sensitivity_modestr

Value should be threshold, baseline or goal to use predefined sensitivities

apply_beam_correctionbool

Include the effect of the beam in the noise angular power spectrum

apply_kludge_correctionbool

If True, reduce the hitcount by a factor of 0.85 to account for not-uniformity in the scanning

homogenousbool

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_ellint

The input spectra have significant power at low \(\ell\), we can zero that power specifying an integer \(\ell\) value here. The power spectra at \(\ell < \ell_0\) are set to zero.

rolloff_ellint

Low ell power damping, see the docstring of so_noise_models.so_models_v3.SO_Noise_Calculator_Public_v3_1_1.rolloff

survey_efficiencyfloat

Fraction of calendar time that may be used to compute map depth.

full_covariancebool

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_yearsint

Total number of years for the Large Aperture telescopes survey

LA_noise_modelstr

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

elevationfloat

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_yearsint

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_fractionoptional,float

If homogenous is True, this sky_fraction is used for the noise curves.

cache_hitmapsbool

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_parametersPath or str

See the help of MapSims

Methods Summary

get_beam_fwhm(tube[, band])

Get beam FWHMs in arcminutes corresponding to the tueb.

get_hitmaps([tube, hitmap])

Get and process hitmaps and sky fractions for the provided tube or provided an external one.

get_inverse_variance(tube[, output_units, ...])

Get the inverse noise variance in each pixel for the requested tube.

get_noise_indices(tube[, band])

Gets indices in the so_noise_model package of a channel or the 2 channels of a tube

get_noise_spectra(tube[, ...])

Get the noise power spectra corresponding to the requested tube from the SO noise model code.

get_white_noise_power(tube, sky_fraction[, ...])

Get white noise power in uK^2-sr (units='sr') or uK^2-arcmin^2 (units='arcmin2') corresponding to the tube name tube.

simulate(tube[, output_units, seed, ...])

Create a random realization of the noise power spectrum

Methods Documentation

get_beam_fwhm(tube, band=None)[source] [edit on github]

Get beam FWHMs in arcminutes corresponding to the tueb. This is useful if non-beam-deconvolved sims are requested and you want to know what beam to apply to your signal simulation.

Parameters
tubestr

Specify a tube.

bandstr,optional

Optionally specify the band name within the tube to get just its white noise.

Returns
beamtuple of floats

The beam FWHM in arcminutes either as a tuple for the pair of bands in the tube, or just for the specific band requested.

get_hitmaps(tube=None, hitmap=None)[source] [edit on github]

Get and process hitmaps and sky fractions for the provided tube or provided an external one.

Parameters
tubestr

Specify a tube.

hitmapstring 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.

Returns
hitmapsndarray or ndmap

Processed hitmaps.

sky_fractionsfloat

The sky fraction covered by the survey determined from the hitmaps.

get_inverse_variance(tube, output_units='uK_CMB', hitmap=None, white_noise_rms=None)[source] [edit on github]

Get the inverse noise variance in each pixel for the requested tube.

Parameters
tubestr

Specify a specific tube.

output_unitsstr

Output unit supported by PySM.units, e.g. uK_CMB or K_RJ

hitmapstring 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_rmsfloat 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
ivar_mapndarray or ndmap

Numpy array with the HEALPix or CAR map of the inverse variance in each pixel. The default units are uK^(-2). This is an extensive quantity that depends on the size of pixels.

get_noise_indices(tube, band=None)[source] [edit on github]

Gets indices in the so_noise_model package of a channel or the 2 channels of a tube

get_noise_spectra(tube, ncurve_sky_fraction=1, return_corr=False)[source] [edit on github]

Get the noise power spectra corresponding to the requested tube from the SO noise model code.

Parameters
tubestr

Specify a tube.

ncurve_sky_fractionfloat,optional

The sky fraction to report to the noise simulator code. In the current implementation, the default is to pass a sky fraction of 1, and scale the result by the corresponding sky fraction determined from each band’s hitmap.

return_corrbool

If True, returns cross-correlation N_XY / sqrt(N_XX * N_YY) coeffient instead of cross-correlation power N_XY in the third row of the returned arrays. This is more convenient sometimes, e.g. when you need to scale the auto-correlation power by some factor.

Returns
ell(nells,) ndarray

Array of nells multipoles starting at ell=0 and spaced by delta_ell=1 corresponding to the noise power spectra nells_T and nells_P

nells_T(3,nells) ndarray

The first two rows contain the temperature auto-correlation of the noise power spectra of each band in the tube. The third row contains the correlation power between the two bands by default, but you can get the cross-correlation coefficient instead by setting return_corr=True.

nells_P(3,nells) ndarray

Same as for nells_T but for polarization.

get_white_noise_power(tube, sky_fraction, band=None, units='sr')[source] [edit on github]

Get white noise power in uK^2-sr (units=’sr’) or uK^2-arcmin^2 (units=’arcmin2’) corresponding to the tube name tube. This is useful if you want to generate your own simulations that do not have the atmospheric component.

Parameters
tubestr

Specify a tube.

sky_fractionfloat

The sky fraction covered by the survey.

bandstr,optional

Optionally specify the band name within the tube to get just its white noise.

units: str

‘sr’ for white noise power in uK^2-steradian and ‘arcmin2’ for the same in uK^2-arcmin^2 units.

Returns
wnoisetuple of floats

The white noise variance in the requested units either as a tuple for the pair of bands in the tube, or just for the specific band requested.

simulate(tube, output_units='uK_CMB', seed=None, nsplits=1, mask_value=None, atmosphere=True, hitmap=None, white_noise_rms=None)[source] [edit on github]

Create a random realization of the noise power spectrum

Parameters
tubestr

Specify a specific tube.

output_unitsstr

Output unit supported by PySM.units, e.g. uK_CMB or K_RJ

seedinteger 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.

nsplitsinteger, 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_valuefloat, 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.

atmospherebool, 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.

hitmapstring 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_rmsfloat 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_mapndarray 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. The second dimension corresponds the three polarization Stokes components I,Q,U and the third dimension corresponds to independent split realizations of the noise.