Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1"""Footprints: Some relevant LSST footprints, including utilities to build them. 

2 

3The goal here is to make it easy to build typical target maps and then their associated combined 

4survey inputs (maps in each filter, including scaling between filters; the associated cloud and 

5sky brightness maps that would have limits for WFD, etc.). 

6 

7For generic use for defining footprints from scratch, there is also a utility that simply generates 

8the healpix points across the sky, along with their corresponding RA/Dec/Galactic l,b/Ecliptic l,b values. 

9""" 

10 

11import os 

12import numpy as np 

13import healpy as hp 

14from astropy.coordinates import SkyCoord 

15from astropy import units as u 

16from .utils import set_default_nside, int_rounded 

17from lsst.sims.utils import _hpid2RaDec, _angularSeparation 

18from lsst.sims.utils import Site 

19from lsst.utils import getPackageDir 

20 

21__all__ = ['ra_dec_hp_map', 'generate_all_sky', 'get_dustmap', 

22 'WFD_healpixels', 'WFD_no_gp_healpixels', 'WFD_bigsky_healpixels', 'WFD_no_dust_healpixels', 

23 'SCP_healpixels', 'NES_healpixels', 

24 'galactic_plane_healpixels', #'low_lat_plane_healpixels', 'bulge_healpixels', 

25 'magellanic_clouds_healpixels', 

26 'generate_goal_map', 'standard_goals', 

27 'calc_norm_factor', 'filter_count_ratios', 'Step_line', 'Footprints', 'Footprint', 

28 'Step_slopes'] 

29 

30 

31class Base_pixel_evolution(object): 

32 """Helper class that can be used to describe the time evolution of a HEALpix in a footprint 

33 """ 

34 

35 def __init__(self, period=365.25, rise=1., t_start=0.): 

36 self.period = period 

37 self.rise = rise 

38 self.t_start = t_start 

39 

40 def __call__(self, mjd_in, phase): 

41 pass 

42 

43 

44class Step_line(Base_pixel_evolution): 

45 """ 

46 Parameters 

47 ---------- 

48 period : float (365.25) 

49 The period to use 

50 rise : float (1.) 

51 How much the curve should rise every period 

52 """ 

53 def __call__(self, mjd_in, phase): 

54 t = mjd_in+phase - self.t_start 

55 n_periods = np.floor(t/(self.period)) 

56 result = n_periods*self.rise 

57 tphased = t % self.period 

58 step_area = np.where(tphased > self.period/2.)[0] 

59 result[step_area] += (tphased[step_area] - self.period/2)*self.rise/(0.5*self.period) 

60 result[np.where(t < 0)] = 0 

61 return result 

62 

63 

64class Step_slopes(Base_pixel_evolution): 

65 """ 

66 Parameters 

67 ---------- 

68 period : float (365.25) 

69 The period to use 

70 rise : np.array-like 

71 How much the curve should rise each period. 

72 """ 

73 def __call__(self, mjd_in, phase): 

74 steps = np.array(self.rise) 

75 t = mjd_in+phase - self.t_start 

76 season = np.floor(t/(self.period)) 

77 season = season.astype(int) 

78 plateus = np.cumsum(steps) 

79 result = plateus[season] 

80 tphased = t % self.period 

81 step_area = np.where(tphased > self.period/2.)[0] 

82 result[step_area] += (tphased[step_area] - self.period/2)*steps[season+1][step_area]/(0.5*self.period) 

83 result[np.where(t < 0)] = 0 

84 

85 return result 

86 

87 

88class Footprint(object): 

89 """An object to compute the desired survey footprint at a given time 

90 

91 Parameters 

92 ---------- 

93 mjd_start : float 

94 The MJD the survey starts on 

95 sun_RA_start : float 

96 The RA of the sun at the start of the survey (radians) 

97 

98 """ 

99 def __init__(self, mjd_start, sun_RA_start=0, nside=32, 

100 filters={'u': 0, 'g': 1, 'r': 2, 'i': 3, 'z': 4, 'y': 5}, 

101 period=365.25, step_func=None): 

102 self.period = period 

103 self.nside = nside 

104 if step_func is None: 

105 step_func = Step_line() 

106 self.step_func = step_func 

107 self.mjd_start = mjd_start 

108 self.sun_RA_start = sun_RA_start 

109 self.npix = hp.nside2npix(nside) 

110 self.filters = filters 

111 self.ra, self.dec = _hpid2RaDec(self.nside, np.arange(self.npix)) 

112 # Set the phase of each healpixel. If RA to sun is zero, we are at phase np.pi/2. 

113 self.phase = (-self.ra + self.sun_RA_start + np.pi/2) % (2.*np.pi) 

114 self.phase = self.phase * (self.period/2./np.pi) 

115 # Empty footprints to start 

116 self.out_dtype = list(zip(filters, [float]*len(filters))) 

117 self.footprints = np.zeros((len(filters), self.npix), dtype=float) 

118 self.estimate = np.zeros((len(filters), self.npix), dtype=float) 

119 self.current_footprints = np.zeros((len(filters), self.npix), dtype=float) 

120 self.zero = self.step_func(0., self.phase) 

121 self.mjd_current = None 

122 

123 def set_footprint(self, filtername, values): 

124 self.footprints[self.filters[filtername], :] = values 

125 

126 def get_footprint(self, filtername): 

127 return self.footprints[self.filters[filtername], :] 

128 

129 def _update_mjd(self, mjd, norm=True): 

130 if mjd != self.mjd_current: 

131 self.mjd_current = mjd 

132 t_elapsed = mjd - self.mjd_start 

133 

134 norm_coverage = self.step_func(t_elapsed, self.phase) 

135 norm_coverage -= self.zero 

136 max_coverage = np.max(norm_coverage) 

137 if max_coverage != 0: 

138 norm_coverage = norm_coverage/max_coverage 

139 self.current_footprints = self.footprints * norm_coverage 

140 c_sum = np.sum(self.current_footprints) 

141 if norm: 

142 if c_sum != 0: 

143 self.current_footprints = self.current_footprints/c_sum 

144 

145 def arr2struc(self, inarr): 

146 """take an array and convert it to labled struc array 

147 """ 

148 result = np.empty(self.npix, dtype=self.out_dtype) 

149 for key in self.filters: 

150 result[key] = inarr[self.filters[key]] 

151 # Argle bargel, why doesn't this view work? 

152 # struc = inarr.view(dtype=self.out_dtype).squeeze() 

153 return result 

154 

155 def estimate_counts(self, mjd, nvisits=2.2e6, fov_area=9.6): 

156 """Estimate the counts we'll get after some time and visits 

157 """ 

158 pix_area = hp.nside2pixarea(self.nside, degrees=True) 

159 pix_per_visit = fov_area/pix_area 

160 self._update_mjd(mjd, norm=True) 

161 self.estimate = self.current_footprints * pix_per_visit * nvisits 

162 return self.arr2struc(self.estimate) 

163 

164 def __call__(self, mjd, array=False): 

165 """ 

166 Returns 

167 ------- 

168 a numpy array with the normalized number of observations that should be at each HEALpix. 

169 Multiply by the number of HEALpix observations (all filters), to convert to the number of observations 

170 desired. 

171 """ 

172 self._update_mjd(mjd) 

173 #if array: 

174 # return self.current_footprints 

175 #else: 

176 return self.arr2struc(self.current_footprints) 

177 

178 

179class Footprints(Footprint): 

180 """An object to combine multiple Footprint objects. 

181 """ 

182 def __init__(self, footprint_list): 

183 self.footprint_list = footprint_list 

184 self.mjd_current = None 

185 self.current_footprints = 0 

186 # Should probably run a check that all the footprints are compatible (same nside, etc) 

187 self.npix = footprint_list[0].npix 

188 self.out_dtype = footprint_list[0].out_dtype 

189 self.filters = footprint_list[0].filters 

190 self.nside = footprint_list[0].nside 

191 

192 self.footprints = np.zeros((len(self.filters), self.npix), dtype=float) 

193 for Fp in self.footprint_list: 

194 self.footprints += Fp.footprints 

195 

196 def set_footprint(self, filtername, values): 

197 pass 

198 

199 def _update_mjd(self, mjd, norm=True): 

200 if mjd != self.mjd_current: 

201 self.mjd_current = mjd 

202 self.current_footprints = 0. 

203 for fp in self.footprint_list: 

204 fp._update_mjd(mjd, norm=False) 

205 self.current_footprints += fp.current_footprints 

206 c_sum = np.sum(self.current_footprints) 

207 if norm: 

208 if c_sum != 0: 

209 self.current_footprints = self.current_footprints/c_sum 

210 

211 

212def ra_dec_hp_map(nside=None): 

213 """ 

214 Return all the RA,dec points for the centers of a healpix map, in radians. 

215 """ 

216 if nside is None: 

217 nside = set_default_nside() 

218 ra, dec = _hpid2RaDec(nside, np.arange(hp.nside2npix(nside))) 

219 return ra, dec 

220 

221 

222def get_dustmap(nside=None): 

223 if nside is None: 

224 nside = set_default_nside() 

225 ebvDataDir = getPackageDir('sims_maps') 

226 filename = 'DustMaps/dust_nside_%i.npz' % nside 

227 dustmap = np.load(os.path.join(ebvDataDir, filename))['ebvMap'] 

228 return dustmap 

229 

230 

231def generate_all_sky(nside=None, elevation_limit=20, mask=hp.UNSEEN): 

232 """Set up a healpix map over the entire sky. 

233 Calculate RA & Dec, Galactic l & b, Ecliptic l & b, for all healpixels. 

234 Calculate max altitude, to set to areas which LSST cannot reach (set these to hp.unseen). 

235 

236 This is intended to be a useful tool to use to set up target maps, beyond the standard maps 

237 provided in these utilities. Masking based on RA, Dec, Galactic or Ecliptic lat and lon is easier. 

238 

239 Parameters 

240 ---------- 

241 nside : int, opt 

242 Resolution for the healpix maps. 

243 Default None uses lsst.sims.featureScheduler.utils.set_default_nside to set default (often 32). 

244 elevation_limit : float, opt 

245 Elevation limit for map. 

246 Parts of the sky which do not reach this elevation limit will be set to mask. 

247 mask : float, opt 

248 Mask value for 'unreachable' parts of the sky, defined as elevation < 20. 

249 Note that the actual limits will be set elsewhere, using the observatory model. 

250 This limit is for use when understanding what the maps could look like. 

251 

252 Returns 

253 ------- 

254 dict of np.ndarray 

255 Returns map, RA/Dec, Gal l/b, Ecl l/b (each an np.ndarray IN RADIANS) in a dictionary. 

256 """ 

257 if nside is None: 

258 nside = set_default_nside() 

259 

260 # Calculate coordinates of everything. 

261 skymap = np.zeros(hp.nside2npix(nside), float) 

262 ra, dec = ra_dec_hp_map(nside=nside) 

263 coord = SkyCoord(ra=ra * u.rad, dec=dec * u.rad, frame='icrs') 

264 eclip_lat = coord.barycentrictrueecliptic.lat.deg 

265 eclip_lon = coord.barycentrictrueecliptic.lon.deg 

266 gal_lon = coord.galactic.l.deg 

267 gal_lat = coord.galactic.b.deg 

268 

269 # Calculate max altitude (when on meridian). 

270 lsst_site = Site('LSST') 

271 elev_max = np.pi / 2. - np.abs(dec - lsst_site.latitude_rad) 

272 skymap = np.where(int_rounded(elev_max) >= int_rounded(np.radians(elevation_limit), skymap, mask)) 

273 

274 return {'map': skymap, 'ra': np.degrees(ra), 'dec': np.degrees(dec), 

275 'eclip_lat': eclip_lat, 'eclip_lon': eclip_lon, 

276 'gal_lat': gal_lat, 'gal_lon': gal_lon} 

277 

278 

279def WFD_healpixels(nside=None, dec_min=-62.5, dec_max=3.6): 

280 """ 

281 Define a region based on declination limits only. 

282 

283 Parameters 

284 ---------- 

285 nside : int, opt 

286 Resolution for the healpix maps. 

287 Default None uses lsst.sims.featureScheduler.utils.set_default_nside to set default (often 32). 

288 dec_min : float, opt 

289 Minimum declination of the region (deg). Default -62.5. 

290 dec_max : float, opt 

291 Maximum declination of the region (deg). Default 3.6. 

292 

293 Returns 

294 ------- 

295 np.ndarray 

296 Healpix map with regions in declination-limited 'wfd' region as 1. 

297 """ 

298 if nside is None: 

299 nside = set_default_nside() 

300 

301 ra, dec = ra_dec_hp_map(nside=nside) 

302 result = np.zeros(ra.size, float) 

303 dec = int_rounded(dec) 

304 good = np.where((dec >= int_rounded(np.radians(dec_min))) & 

305 (dec <= int_rounded(np.radians(dec_max)))) 

306 result[good] = 1 

307 return result 

308 

309 

310def WFD_no_gp_healpixels(nside, dec_min=-62.5, dec_max=3.6, 

311 center_width=10., end_width=4., gal_long1=290., gal_long2=70.): 

312 """ 

313 Define a wide fast deep region with a galactic plane limit. 

314 

315 Parameters 

316 ---------- 

317 nside : int, opt 

318 Resolution for the healpix maps. 

319 Default None uses lsst.sims.featureScheduler.utils.set_default_nside to set default (often 32). 

320 dec_min : float, opt 

321 Minimum declination of the region (deg). 

322 dec_max : float, opt 

323 Maximum declination of the region (deg). 

324 center_width : float, opt 

325 Width across the central part of the galactic plane region. 

326 end_width : float, opt 

327 Width across the remainder of the galactic plane region. 

328 gal_long1 : float, opt 

329 Longitude at which to start tapering from center_width to end_width. 

330 gal_long2 : float, opt 

331 Longitude at which to stop tapering from center_width to end_width. 

332 

333 Returns 

334 ------- 

335 np.ndarray 

336 Healpix map with regions in declination-limited 'wfd' region as 1. 

337 """ 

338 wfd_dec = WFD_healpixels(nside, dec_min=dec_min, dec_max=dec_max) 

339 gp = galactic_plane_healpixels(nside=nside, center_width=center_width, end_width=end_width, 

340 gal_long1=gal_long1, gal_long2=gal_long2) 

341 sky = np.where(wfd_dec - gp > 0, wfd_dec - gp, 0) 

342 return sky 

343 

344 

345def WFD_bigsky_healpixels(nside): 

346 sky = WFD_no_gp_healpixels(nside, dec_min=-72.25, dec_max=12.4, center_width=14.9, 

347 gal_long1=0, gal_long2=360) 

348 return sky 

349 

350 

351def WFD_no_dust_healpixels(nside, dec_min=-72.25, dec_max=12.4, dust_limit=0.19): 

352 """Define a WFD region with a dust extinction limit. 

353 

354 Parameters 

355 ---------- 

356 nside : int, opt 

357 Resolution for the healpix maps. 

358 Default None uses lsst.sims.featureScheduler.utils.set_default_nside to set default (often 32). 

359 dec_min : float, opt 

360 Minimum dec of the region (deg). Default -72.5 deg. 

361 dec_max : float, opt. 

362 Maximum dec of the region (deg). Default 12.5 deg. 

363 1.75 is the FOV radius in deg. 

364 dust_limit : float, None 

365 Remove pixels with E(B-V) values greater than dust_limit from the footprint. 

366 

367 Returns 

368 ------- 

369 result : numpy array 

370 """ 

371 if nside is None: 

372 nside = set_default_nside() 

373 

374 ra, dec = ra_dec_hp_map(nside=nside) 

375 dustmap = get_dustmap(nside) 

376 

377 result = np.zeros(ra.size, float) 

378 # First set based on dec range. 

379 dec = int_rounded(dec) 

380 good = np.where((dec >= int_rounded(np.radians(dec_min))) & 

381 (dec <= int_rounded(np.radians(dec_max)))) 

382 result[good] = 1 

383 # Now remove areas with dust extinction beyond the limit. 

384 result = np.where(dustmap >= dust_limit, 0, result) 

385 return result 

386 

387 

388def SCP_healpixels(nside=None, dec_max=-60.): 

389 """ 

390 Define the South Celestial Pole region. Return a healpix map with SCP pixels as 1. 

391 """ 

392 if nside is None: 

393 nside = set_default_nside() 

394 

395 ra, dec = ra_dec_hp_map(nside=nside) 

396 result = np.zeros(ra.size, float) 

397 good = np.where(int_rounded(dec) < int_rounded(np.radians(dec_max))) 

398 result[good] = 1 

399 return result 

400 

401 

402def NES_healpixels(nside=None, min_EB=-30.0, max_EB = 10.0, dec_min=2.8): 

403 """ 

404 Define the North Ecliptic Spur region. Return a healpix map with NES pixels as 1. 

405 

406 Parameters 

407 ---------- 

408 nside : int 

409 A valid healpix nside 

410 min_EB : float (-30.) 

411 Minimum barycentric true ecliptic latitude (deg) 

412 max_EB : float (10.) 

413 Maximum barycentric true ecliptic latitude (deg) 

414 dec_min : float (2.8) 

415 Minimum dec in degrees 

416 

417 Returns 

418 ------- 

419 result : numpy array 

420 """ 

421 if nside is None: 

422 nside = set_default_nside() 

423 

424 ra, dec = ra_dec_hp_map(nside=nside) 

425 result = np.zeros(ra.size, float) 

426 coord = SkyCoord(ra=ra*u.rad, dec=dec*u.rad) 

427 eclip_lat = coord.barycentrictrueecliptic.lat.radian 

428 eclip_lat = int_rounded(eclip_lat) 

429 dec = int_rounded(dec) 

430 good = np.where((eclip_lat > int_rounded(np.radians(min_EB))) & 

431 (eclip_lat < int_rounded(np.radians(max_EB))) & 

432 (dec > int_rounded(np.radians(dec_min)))) 

433 result[good] = 1 

434 

435 return result 

436 

437 

438def galactic_plane_healpixels(nside=None, center_width=10., end_width=4., 

439 gal_long1=290., gal_long2=70.): 

440 """ 

441 Define a Galactic Plane region. 

442 

443 Parameters 

444 ---------- 

445 nside : int, opt 

446 Resolution for the healpix maps. 

447 Default None uses lsst.sims.featureScheduler.utils.set_default_nside to set default (often 32). 

448 center_width : float, opt 

449 Width at the center of the galactic plane region. 

450 end_width : float, opt 

451 Width at the remainder of the galactic plane region. 

452 gal_long1 : float, opt 

453 Longitude at which to start the GP region. 

454 gal_long2 : float, opt 

455 Longitude at which to stop the GP region. 

456 Order matters for gal_long1 / gal_long2! 

457 

458 Returns 

459 ------- 

460 np.ndarray 

461 Healpix map with galactic plane regions set to 1. 

462 """ 

463 if nside is None: 

464 nside = set_default_nside() 

465 ra, dec = ra_dec_hp_map(nside=nside) 

466 

467 coord = SkyCoord(ra=ra*u.rad, dec=dec*u.rad) 

468 gal_lon, gal_lat = coord.galactic.l.deg, coord.galactic.b.deg 

469 # Reject anything beyond the central width. 

470 sky = np.where(np.abs(gal_lat) < center_width, 1, 0) 

471 # Apply the galactic longitude cuts, so that plane goes between gal_long1 to gal_long2. 

472 # This is NOT the shortest distance between the angles. 

473 gp_length = (gal_long2 - gal_long1) % 360 

474 # If the length is greater than 0 then we can add additional cuts. 

475 if gp_length > 0: 

476 # First, remove anything outside the gal_long1/gal_long2 region. 

477 sky = np.where(int_rounded((gal_lon - gal_long1) % 360) < int_rounded(gp_length), sky, 0) 

478 # Add the tapers. 

479 # These slope from the center (gp_center @ center_width) 

480 # to the edge (gp_center + gp_length/2 @ end_width). 

481 half_width = gp_length / 2. 

482 slope = (center_width - end_width) / half_width 

483 gp_center = (gal_long1 + half_width) % 360 

484 gp_dist = gal_lon - gp_center 

485 gp_dist = np.abs(np.where(int_rounded(gp_dist) > int_rounded(180), (180 - gp_dist) % 180, gp_dist)) 

486 lat_limit = np.abs(center_width - slope * gp_dist) 

487 sky = np.where(int_rounded(np.abs(gal_lat)) < int_rounded(lat_limit), sky, 0) 

488 return sky 

489 

490 

491def magellanic_clouds_healpixels(nside=None, lmc_radius=10, smc_radius=5): 

492 """ 

493 Define the Galactic Plane region. Return a healpix map with GP pixels as 1. 

494 """ 

495 if nside is None: 

496 nside = set_default_nside() 

497 ra, dec = ra_dec_hp_map(nside=nside) 

498 result = np.zeros(hp.nside2npix(nside)) 

499 

500 lmc_ra = np.radians(80.893860) 

501 lmc_dec = np.radians(-69.756126) 

502 lmc_radius = np.radians(lmc_radius) 

503 

504 smc_ra = np.radians(13.186588) 

505 smc_dec = np.radians(-72.828599) 

506 smc_radius = np.radians(smc_radius) 

507 

508 dist_to_lmc = _angularSeparation(lmc_ra, lmc_dec, ra, dec) 

509 lmc_pix = np.where(int_rounded(dist_to_lmc) < int_rounded(lmc_radius)) 

510 result[lmc_pix] = 1 

511 

512 dist_to_smc = _angularSeparation(smc_ra, smc_dec, ra, dec) 

513 smc_pix = np.where(int_rounded(dist_to_smc) < int_rounded(smc_radius)) 

514 result[smc_pix] = 1 

515 return result 

516 

517 

518def generate_goal_map(nside=None, NES_fraction = .3, WFD_fraction = 1., 

519 SCP_fraction=0.4, GP_fraction = 0.2, 

520 NES_min_EB = -30., NES_max_EB = 10, NES_dec_min = 3.6, 

521 SCP_dec_max=-62.5, gp_center_width=10., 

522 gp_end_width=4., gp_long1=290., gp_long2=70., 

523 wfd_dec_min=-62.5, wfd_dec_max=3.6, 

524 generate_id_map=False): 

525 """ 

526 Handy function that will put together a target map in the proper order. 

527 """ 

528 if nside is None: 

529 nside = set_default_nside() 

530 

531 # Note, some regions overlap, thus order regions are added is important. 

532 result = np.zeros(hp.nside2npix(nside), dtype=float) 

533 id_map = np.zeros(hp.nside2npix(nside), dtype=int) 

534 pid = 1 

535 prop_name_dict = dict() 

536 

537 if NES_fraction > 0.: 

538 nes = NES_healpixels(nside=nside, min_EB = NES_min_EB, max_EB = NES_max_EB, 

539 dec_min=NES_dec_min) 

540 result[np.where(nes != 0)] = 0 

541 result += NES_fraction*nes 

542 id_map[np.where(nes != 0)] = 1 

543 pid += 1 

544 prop_name_dict[1] = 'NorthEclipticSpur' 

545 

546 if WFD_fraction > 0.: 

547 wfd = WFD_healpixels(nside=nside, dec_min=wfd_dec_min, dec_max=wfd_dec_max) 

548 result[np.where(wfd != 0)] = 0 

549 result += WFD_fraction*wfd 

550 id_map[np.where(wfd != 0)] = 3 

551 pid += 1 

552 prop_name_dict[3] = 'WideFastDeep' 

553 

554 if SCP_fraction > 0.: 

555 scp = SCP_healpixels(nside=nside, dec_max=SCP_dec_max) 

556 result[np.where(scp != 0)] = 0 

557 result += SCP_fraction*scp 

558 id_map[np.where(scp != 0)] = 2 

559 pid += 1 

560 prop_name_dict[2] = 'SouthCelestialPole' 

561 

562 if GP_fraction > 0.: 

563 gp = galactic_plane_healpixels(nside=nside, center_width=gp_center_width, 

564 end_width=gp_end_width, gal_long1=gp_long1, 

565 gal_long2=gp_long2) 

566 result[np.where(gp != 0)] = 0 

567 result += GP_fraction*gp 

568 id_map[np.where(gp != 0)] = 4 

569 pid += 1 

570 prop_name_dict[4] = 'GalacticPlane' 

571 

572 if generate_id_map: 

573 return result, id_map, prop_name_dict 

574 else: 

575 return result 

576 

577 

578def standard_goals(nside=None): 

579 """ 

580 A quick function to generate the "standard" goal maps. This is the traditional WFD/mini survey footprint. 

581 """ 

582 if nside is None: 

583 nside = set_default_nside() 

584 

585 result = {} 

586 result['u'] = generate_goal_map(nside=nside, NES_fraction=0., 

587 WFD_fraction=0.31, SCP_fraction=0.15, 

588 GP_fraction=0.15, 

589 wfd_dec_min=-62.5, wfd_dec_max=3.6) 

590 result['g'] = generate_goal_map(nside=nside, NES_fraction=0.2, 

591 WFD_fraction=0.44, SCP_fraction=0.15, 

592 GP_fraction=0.15, 

593 wfd_dec_min=-62.5, wfd_dec_max=3.6) 

594 result['r'] = generate_goal_map(nside=nside, NES_fraction=0.46, 

595 WFD_fraction=1.0, SCP_fraction=0.15, 

596 GP_fraction=0.15, 

597 wfd_dec_min=-62.5, wfd_dec_max=3.6) 

598 result['i'] = generate_goal_map(nside=nside, NES_fraction=0.46, 

599 WFD_fraction=1.0, SCP_fraction=0.15, 

600 GP_fraction=0.15, 

601 wfd_dec_min=-62.5, wfd_dec_max=3.6) 

602 result['z'] = generate_goal_map(nside=nside, NES_fraction=0.4, 

603 WFD_fraction=0.9, SCP_fraction=0.15, 

604 GP_fraction=0.15, 

605 wfd_dec_min=-62.5, wfd_dec_max=3.6) 

606 result['y'] = generate_goal_map(nside=nside, NES_fraction=0., 

607 WFD_fraction=0.9, SCP_fraction=0.15, 

608 GP_fraction=0.15, 

609 wfd_dec_min=-62.5, wfd_dec_max=3.6) 

610 return result 

611 

612 

613def calc_norm_factor(goal_dict, radius=1.75): 

614 """Calculate how to normalize a Target_map_basis_function. 

615 This is basically: 

616 the area of the fov / area of a healpixel / 

617 the sum of all of the weighted-healpix values in the footprint. 

618 

619 Parameters 

620 ----------- 

621 goal_dict : dict of healpy maps 

622 The target goal map(s) being used 

623 radius : float (1.75) 

624 Radius of the FoV (degrees) 

625 

626 Returns 

627 ------- 

628 Value to use as Target_map_basis_function norm_factor kwarg 

629 """ 

630 all_maps_sum = 0 

631 for key in goal_dict: 

632 good = np.where(goal_dict[key] > 0) 

633 all_maps_sum += goal_dict[key][good].sum() 

634 nside = hp.npix2nside(goal_dict[key].size) 

635 hp_area = hp.nside2pixarea(nside, degrees=True) 

636 norm_val = radius**2*np.pi/hp_area/all_maps_sum 

637 return norm_val 

638 

639 

640def filter_count_ratios(target_maps): 

641 """Given the goal maps, compute the ratio of observations we want in each filter. 

642 This is basically: 

643 per filter, sum the number of pixels in each map and return this per filter value, normalized 

644 so that the total sum across all filters is 1. 

645 """ 

646 results = {} 

647 all_norm = 0. 

648 for key in target_maps: 

649 good = target_maps[key] > 0 

650 results[key] = np.sum(target_maps[key][good]) 

651 all_norm += results[key] 

652 for key in results: 

653 results[key] /= all_norm 

654 return results