|
| 1 | +import numpy as np |
| 2 | +import healpy as hp |
| 3 | +from numba import njit |
| 4 | +from .. import utils |
| 5 | + |
| 6 | + |
| 7 | +# from astropy import constants as const |
| 8 | +# |
| 9 | +from .. import units as u |
| 10 | + |
| 11 | +from .template import Model |
| 12 | + |
| 13 | +import h5py |
| 14 | + |
| 15 | + |
| 16 | +@njit |
| 17 | +def fwhm2sigma(fwhm): |
| 18 | + """Converts the Full Width Half Maximum of a Gaussian beam to its standard deviation""" |
| 19 | + return fwhm / (2.0 * np.sqrt(2.0 * np.log(2.0))) |
| 20 | + |
| 21 | + |
| 22 | +@njit |
| 23 | +def flux2amp(flux, fwhm): |
| 24 | + """Converts the total flux of a radio source to the peak amplitude of its Gaussian |
| 25 | + beam representation, taking into account the width of the beam as specified |
| 26 | + by its FWHM |
| 27 | +
|
| 28 | + Parameters |
| 29 | + ---------- |
| 30 | + flux: float |
| 31 | + Total flux of the radio source |
| 32 | + fwhm: float |
| 33 | + Full Width Half Maximum of the beam in radians |
| 34 | +
|
| 35 | + Returns |
| 36 | + ------- |
| 37 | + amp: float |
| 38 | + Peak amplitude of the Gaussian beam representation of the radio source""" |
| 39 | + sigma = fwhm2sigma(fwhm) |
| 40 | + return flux / (2 * np.pi * sigma**2) |
| 41 | + |
| 42 | + |
| 43 | +@njit |
| 44 | +def evaluate_poly(p, x): |
| 45 | + """Low level polynomial evaluation, both input are 1D |
| 46 | + same interface of np.polyval. |
| 47 | + Having this implemented in numba should allow numba |
| 48 | + to provide better optimization. If not, just use |
| 49 | + np.polyval directly.""" |
| 50 | + |
| 51 | + out = 0 |
| 52 | + N = len(p) |
| 53 | + for i in range(N): |
| 54 | + out += p[i] * x ** (N - 1 - i) |
| 55 | + return out |
| 56 | + |
| 57 | + |
| 58 | +@njit |
| 59 | +def evaluate_model(freqs, weights, coeff): |
| 60 | + """Integrate log polynomial model across the bandpass for |
| 61 | + each source in the catalog |
| 62 | +
|
| 63 | + Parameters |
| 64 | + ---------- |
| 65 | + freqs: np.array |
| 66 | + Array of frequencies in GHz |
| 67 | + weights: np.array |
| 68 | + Array of relative bandpass weights already normalized |
| 69 | + Same length of freqs |
| 70 | + coeff: 2D np.array (n_sources, n_coeff) |
| 71 | + Array of log polynomial coefficients for each source |
| 72 | +
|
| 73 | + Returns |
| 74 | + ------- |
| 75 | + flux: np.array |
| 76 | + Array of the flux of each source integrated over the band |
| 77 | + """ |
| 78 | + n_sources = coeff.shape[0] |
| 79 | + logfreqs = np.log(freqs) |
| 80 | + out = np.zeros(n_sources, dtype=np.float64) |
| 81 | + assert len(freqs) == len(weights) |
| 82 | + if len(freqs) == 1: |
| 83 | + for i_source in range(n_sources): |
| 84 | + out[i_source] = evaluate_poly(coeff[i_source, :], logfreqs[0]) |
| 85 | + else: |
| 86 | + flux = np.zeros(len(freqs), dtype=np.float64) |
| 87 | + for i_source in range(n_sources): |
| 88 | + for i_freq in range(len(freqs)): |
| 89 | + flux[i_freq] = evaluate_poly(coeff[i_source, :], logfreqs[i_freq]) |
| 90 | + out[i_source] = np.trapz(flux * weights, x=freqs) |
| 91 | + return out |
| 92 | + |
| 93 | + |
| 94 | +class PointSourceCatalog(Model): |
| 95 | + """Model for a Catalog of point sources defined with their coordinates and |
| 96 | + a model of their emission based on a logpolynomial of frequency. |
| 97 | + The beam convolution is performed in map domain with `pixell`. |
| 98 | +
|
| 99 | + The catalog should be in HDF5 format, with the fields: |
| 100 | + theta: colatitude in radians |
| 101 | + phi: longitude in radians |
| 102 | + logpolycoefflux and logpolycoefpolflux: polynomial coefficients in natural |
| 103 | + logaritm (`np.log`) of the frequency, typically 4th order, but accepts |
| 104 | + any order. (source_index, pol_order). Unit needs to be Jy |
| 105 | + each field should have an attribute units which is checked when loading |
| 106 | + a model. No conversion is performed. |
| 107 | + See the documentation and the unit tests for examples on how to create a |
| 108 | + catalog file with `xarray`. |
| 109 | +
|
| 110 | + Parameters |
| 111 | + ---------- |
| 112 | + catalog_filename: str or Path |
| 113 | + Path to the catalog HDF5 file |
| 114 | + """ |
| 115 | + |
| 116 | + def __init__( |
| 117 | + self, |
| 118 | + catalog_filename, |
| 119 | + nside=None, |
| 120 | + target_wcs=None, |
| 121 | + map_dist=None, |
| 122 | + ): |
| 123 | + self.catalog_filename = catalog_filename |
| 124 | + self.nside = nside |
| 125 | + self.shape = (3, hp.nside2npix(nside)) |
| 126 | + self.wcs = target_wcs |
| 127 | + |
| 128 | + with h5py.File(self.catalog_filename) as f: |
| 129 | + assert f["theta"].attrs["units"].decode("UTF-8") == "rad" |
| 130 | + assert f["phi"].attrs["units"].decode("UTF-8") == "rad" |
| 131 | + assert f["logpolycoefflux"].attrs["units"].decode("UTF-8") == "Jy" |
| 132 | + assert f["logpolycoefpolflux"].attrs["units"].decode("UTF-8") == "Jy" |
| 133 | + |
| 134 | + assert map_dist is None, "Distributed execution not supported" |
| 135 | + |
| 136 | + def get_fluxes(self, freqs: u.GHz, coeff="logpolycoefflux", weights=None): |
| 137 | + """Get catalog fluxes in Jy integrated over a bandpass""" |
| 138 | + weights /= np.trapz(weights, x=freqs.to_value(u.GHz)) |
| 139 | + with h5py.File(self.catalog_filename) as f: |
| 140 | + flux = evaluate_model(freqs.to_value(u.GHz), weights, np.array(f[coeff])) |
| 141 | + return flux * u.Jy |
| 142 | + |
| 143 | + @u.quantity_input |
| 144 | + def get_emission( |
| 145 | + self, |
| 146 | + freqs: u.GHz, |
| 147 | + fwhm: [u.arcmin, None] = None, |
| 148 | + weights=None, |
| 149 | + output_units=u.uK_RJ, |
| 150 | + car_map_resolution: [u.arcmin, None] = None, |
| 151 | + return_car=False, |
| 152 | + ): |
| 153 | + """Generate a HEALPix or CAR map of the catalog emission integrated on the bandpass |
| 154 | + and convolved with the beam |
| 155 | +
|
| 156 | + Parameters |
| 157 | + ---------- |
| 158 | + freqs: np.array |
| 159 | + Array of frequencies in GHz |
| 160 | + fwhm: float or None |
| 161 | + Full Width Half Maximum of the beam in arcminutes, if None, each source is assigned |
| 162 | + to a single pixel |
| 163 | + weights: np.array |
| 164 | + Array of relative bandpass weights already normalized |
| 165 | + Same length of freqs, if None, uniform weights are assumed |
| 166 | + output_units: astropy.units |
| 167 | + Output units of the map |
| 168 | + car_map_resolution: float |
| 169 | + Resolution of the CAR map used by pixell to generate the map, if None, |
| 170 | + it is set to half of the resolution of the HEALPix map given by `self.nside` |
| 171 | + return_car: bool |
| 172 | + If True return a CAR map, if False return a HEALPix map |
| 173 | +
|
| 174 | + Returns |
| 175 | + ------- |
| 176 | + output_map: np.array |
| 177 | + Output HEALPix or CAR map""" |
| 178 | + with h5py.File(self.catalog_filename) as f: |
| 179 | + pix = hp.ang2pix(self.nside, f["theta"], f["phi"]) |
| 180 | + scaling_factor = utils.bandpass_unit_conversion( |
| 181 | + freqs, weights, output_unit=output_units, input_unit=u.Jy / u.sr |
| 182 | + ) |
| 183 | + pix_size = hp.nside2pixarea(self.nside) * u.sr |
| 184 | + if car_map_resolution is None: |
| 185 | + car_map_resolution = (hp.nside2resol(self.nside) * u.rad) / 2 |
| 186 | + |
| 187 | + # Make sure the resolution evenly divides the map vertically |
| 188 | + if (car_map_resolution.to_value(u.rad) % np.pi) > 1e-8: |
| 189 | + car_map_resolution = ( |
| 190 | + np.pi / np.round(np.pi / car_map_resolution.to_value(u.rad)) |
| 191 | + ) * u.rad |
| 192 | + fluxes_I = self.get_fluxes(freqs, weights=weights, coeff="logpolycoefflux") |
| 193 | + |
| 194 | + if fwhm is None: |
| 195 | + output_map = np.zeros(self.shape, dtype=np.float32) * output_units |
| 196 | + # sum, what if we have 2 sources on the same pixel? |
| 197 | + output_map[0, pix] += fluxes_I / pix_size * scaling_factor |
| 198 | + else: |
| 199 | + |
| 200 | + from pixell import ( |
| 201 | + enmap, |
| 202 | + pointsrcs, |
| 203 | + ) |
| 204 | + |
| 205 | + shape, wcs = enmap.fullsky_geometry( |
| 206 | + car_map_resolution.to_value(u.radian), |
| 207 | + dims=(3,), |
| 208 | + variant="fejer1", |
| 209 | + ) |
| 210 | + output_map = enmap.enmap(np.zeros(shape, dtype=np.float32), wcs) |
| 211 | + r, p = pointsrcs.expand_beam(fwhm2sigma(fwhm.to_value(u.rad))) |
| 212 | + with h5py.File(self.catalog_filename) as f: |
| 213 | + pointing = np.column_stack( |
| 214 | + (np.pi / 2 - np.array(f["theta"]), np.array(f["phi"])) |
| 215 | + ) |
| 216 | + output_map[0] = pointsrcs.sim_objects( |
| 217 | + shape, |
| 218 | + wcs, |
| 219 | + pointing, |
| 220 | + flux2amp( |
| 221 | + fluxes_I.to_value(u.Jy) * scaling_factor.value, |
| 222 | + fwhm.to_value(u.rad), |
| 223 | + ), # to peak amplitude and to output units |
| 224 | + ((r, p)), |
| 225 | + ) |
| 226 | + |
| 227 | + del fluxes_I |
| 228 | + fluxes_P = self.get_fluxes(freqs, weights=weights, coeff="logpolycoefpolflux") |
| 229 | + # set seed so that the polarization angle is always the same for each run |
| 230 | + # could expose to the interface if useful |
| 231 | + np.random.seed(56567) |
| 232 | + psirand = np.random.uniform( |
| 233 | + low=-np.pi / 2.0, high=np.pi / 2.0, size=len(fluxes_P) |
| 234 | + ) |
| 235 | + if fwhm is None: |
| 236 | + output_map[1, pix] += ( |
| 237 | + fluxes_P / pix_size * scaling_factor * np.cos(2 * psirand) |
| 238 | + ) |
| 239 | + output_map[2, pix] += ( |
| 240 | + fluxes_P / pix_size * scaling_factor * np.sin(2 * psirand) |
| 241 | + ) |
| 242 | + else: |
| 243 | + pols = [(1, np.cos)] |
| 244 | + pols.append((2, np.sin)) |
| 245 | + for i_pol, sincos in pols: |
| 246 | + output_map[i_pol] = pointsrcs.sim_objects( |
| 247 | + shape, |
| 248 | + wcs, |
| 249 | + pointing, |
| 250 | + flux2amp( |
| 251 | + fluxes_P.to_value(u.Jy) |
| 252 | + * scaling_factor.value |
| 253 | + * sincos(2 * psirand), |
| 254 | + fwhm.to_value(u.rad), |
| 255 | + ), |
| 256 | + ((r, p)), |
| 257 | + ) |
| 258 | + if return_car: |
| 259 | + return output_map |
| 260 | + else: |
| 261 | + from pixell import reproject |
| 262 | + |
| 263 | + return reproject.map2healpix(output_map, self.nside) |
0 commit comments