InterpolatingComponent

class pysm.InterpolatingComponent(path, input_units, nside, interpolation_kind='linear', has_polarization=True, map_dist=None, verbose=False)[source] [edit on github]

Bases: pysm.Model

PySM component interpolating between precomputed maps

In order to save memory, maps are converted to float32, if this is not acceptable, please open an issue on the PySM repository. When you create the model, PySM checks the folder of the templates and stores a list of available frequencies. Once you call get_emission, maps are read, ud_graded to the target nside and stored for future use. This is useful if you are running many channels with a similar bandpass. If not, you can call cached_maps.clear() to remove the cached maps.

Parameters:
path : str

Path should contain maps named as the frequency in GHz e.g. 20.fits or 20.5.fits or 00100.fits

input_units : str

Any unit available in PySM3 e.g. “uK_RJ”, “uK_CMB”

nside : int

HEALPix NSIDE of the output maps

interpolation_kind : string

Currently only linear is implemented

has_polarization : bool

whether or not to simulate also polarization maps

map_dist : pysm.MapDistribution

Required for partial sky or MPI, see the PySM docs

verbose : bool

Control amount of output

Methods Summary

get_emission(self, freqs[, weights]) This function evaluates the component model at a either a single frequency, an array of frequencies, or over a bandpass.
get_filenames(self, path)
read_map_by_frequency(self, freq)
read_map_file(self, freq, filename)

Methods Documentation

get_emission(self, freqs: Unit("GHz"), weights=None) -> Unit("uK_RJ")[source] [edit on github]

This function evaluates the component model at a either a single frequency, an array of frequencies, or over a bandpass.

Parameters:
freqs: scalar or array astropy.units.Quantity

Frequency at which the model should be evaluated, in a frequency which can be converted to GHz using astropy.units. If an array of frequencies is provided, integrate using trapz with a equal weighting, i.e. simulate a top-hat bandpass.

weights: np.array, optional

Array of weights describing the frequency response of the instrument, i.e. the bandpass. Weights are normalized and applied in Jy/sr.

Returns:
output : astropy.units.Quantity

Simulated map at the given frequency or integrated over the given bandpass. The shape of the output is (3,npix) for polarized components, (1,npix) for temperature-only components. Output is in uK_RJ.

get_filenames(self, path)[source] [edit on github]
read_map_by_frequency(self, freq)[source] [edit on github]
read_map_file(self, freq, filename)[source] [edit on github]