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

1import numpy as np 

2from lsst.sims.featureScheduler import features 

3from lsst.sims.featureScheduler import utils 

4import healpy as hp 

5import matplotlib.pylab as plt 

6import warnings 

7from lsst.sims.featureScheduler.basis_functions import Base_basis_function 

8from lsst.sims.utils import _hpid2RaDec 

9 

10 

11 

12__all__ = ["Target_map_modulo_basis_function", "Footprint_basis_function", "Footprint_rolling_basis_function"] 

13 

14 

15 

16class Footprint_basis_function(Base_basis_function): 

17 """Basis function that tries to maintain a uniformly covered footprint 

18 

19 Parameters 

20 ---------- 

21 filtername : str ('r') 

22 The filter for this footprint 

23 footprint : lsst.sims.featureScheduler.utils.Footprint object 

24 The desired footprint. 

25 all_footprints_sum : float (None) 

26 If using multiple filters, the sum of all the footprints. Needed to make sure basis functions are 

27 normalized properly across all fitlers. 

28 out_of_bounds_val : float (-10) 

29 The value to set the basis function for regions that are not in the footprint (default -10, np.nan is 

30 another good value to use) 

31 

32 """ 

33 def __init__(self, filtername='r', nside=None, footprint=None, 

34 out_of_bounds_val=-10., window_size=6.): 

35 

36 super(Footprint_basis_function, self).__init__(nside=nside, filtername=filtername) 

37 self.footprint = footprint 

38 

39 self.survey_features = {} 

40 # All the observations in all filters 

41 self.survey_features['N_obs_all'] = features.N_observations(nside=nside, filtername=None) 

42 self.survey_features['N_obs'] = features.N_observations(nside=nside, filtername=filtername) 

43 

44 # should probably actually loop over all the target maps? 

45 self.out_of_bounds_area = np.where(footprint.get_footprint(self.filtername) <= 0)[0] 

46 self.out_of_bounds_val = out_of_bounds_val 

47 

48 def _calc_value(self, conditions, indx=None): 

49 

50 # Find out what the footprint object thinks we should have been observed 

51 desired_footprint_normed = self.footprint(conditions.mjd)[self.filtername] 

52 

53 # Compute how many observations we should have on the sky 

54 desired = desired_footprint_normed * np.sum(self.survey_features['N_obs_all'].feature) 

55 result = desired - self.survey_features['N_obs'].feature 

56 result[self.out_of_bounds_area] = self.out_of_bounds_val 

57 return result 

58 

59 

60class Footprint_rolling_basis_function(Base_basis_function): 

61 """Let's get the rolling really right. 

62 

63 Parameters 

64 ---------- 

65 footprints : list of np.array 

66 List of HEALpix arrays. The footprints should all have matching sums and have the same nside. 

67 all_footprints_sum : float 

68 The sum of footprints over all filters. 

69 all_rolling_sum : float 

70 The sum (over all filters) of the region of the maps that changs. 

71 season_modulo : int (2) 

72 The modulo to pass to utils.season_calc. 

73 season_length : float (365.25) 

74 How long a season should be (days). 

75 max_season : int (None) 

76 If set, the season calc will return -1 for values greater than max_season 

77 day_offset : np.array (None) 

78 Offset to pass to utils.season_calc (days). 

79 

80 """ 

81 

82 def __init__(self, filtername='r', nside=None, footprints=None, all_footprints_sum=None, all_rolling_sum=None, out_of_bounds_val=-10, 

83 season_modulo=2, season_length=365.25, max_season=None, day_offset=None, window_size=6.): 

84 super(Footprint_rolling_basis_function, self).__init__(nside=nside, filtername=filtername) 

85 

86 # OK, going to find the parts of the map that are the same everywhere, and compute the 

87 # basis function the same as usual for those. 

88 same_footprint = np.ones(footprints[0].size, dtype=bool) 

89 for footprint in footprints[0:-1]: 

90 same_footprint = same_footprint & (footprint == footprints[-1]) 

91 

92 sum_footprints = footprints[0]*0 

93 for footprint in footprints: 

94 sum_footprints += footprint 

95 self.constant_footprint_indx = np.where((same_footprint == True) & (sum_footprints > 0))[0] 

96 self.rolling_footprint_indx = np.where((same_footprint == False) & (sum_footprints > 0))[0] 

97 

98 self.season_modulo = season_modulo 

99 self.season_length = season_length 

100 self.max_season = max_season 

101 self.day_offset = day_offset 

102 self.footprints = footprints 

103 

104 self.max_day_offset = np.max(self.day_offset) 

105 

106 self.all_footprints_sum = all_footprints_sum 

107 self.all_rolling_sum = all_rolling_sum 

108 

109 self.survey_features = {} 

110 # Set a season for -1 (for before rolling or after max_season) 

111 # All the observations in the given filter 

112 self.survey_features['N_obs_%i' % -1] = features.N_observations(nside=nside, filtername=filtername) 

113 # All the observations in all filters 

114 self.survey_features['N_obs_all_%i' % -1] = features.N_observations(nside=nside, filtername=None) 

115 

116 for i, temp in enumerate(footprints[0:-1]): 

117 # Observation in a season, in filtername 

118 self.survey_features['N_obs_%i' % i] = features.N_observations_season(i, filtername=filtername, 

119 nside=self.nside, 

120 modulo=season_modulo, 

121 offset=day_offset, 

122 max_season=max_season, 

123 season_length=season_length) 

124 # Count of all the observations taken in a season 

125 self.survey_features['N_obs_all_%i' % i] = features.N_observations_season(i, filtername=None, 

126 modulo=season_modulo, 

127 offset=day_offset, 

128 nside=self.nside, 

129 max_season=max_season, 

130 season_length=season_length) 

131 

132 # Now I need to track the observations taken in each season. 

133 self.out_of_bounds_area = np.where(footprint <= 0)[0] 

134 self.out_of_bounds_val = out_of_bounds_val 

135 

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

137 

138 def _calc_value(self, conditions, indx=None): 

139 result = self.result.copy() 

140 

141 # Compute what season it is at each pixel 

142 seasons = conditions.season(modulo=self.season_modulo, 

143 max_season=self.max_season, season_length=self.season_length) 

144 

145 # Update the coverage of the sky so far 

146 # If RA to sun is zero, we are at phase np.pi/2. 

147 coverage_map_phase = (conditions.ra - conditions.sun_RA_start+np.pi/2.) % (2.*np.pi) 

148 zero = utils.step_line(0., 1., 365.25, phase=coverage_map_phase*365.25/2/np.pi) 

149 t_elapsed = conditions.mjd - conditions.mjd_start 

150 norm_coverage_raw = utils.step_line(t_elapsed, 1., 365.25, phase=coverage_map_phase*365.25/2/np.pi) 

151 norm_coverage_raw -= zero 

152 

153 norm_coverage = norm_coverage_raw/np.max(norm_coverage_raw) 

154 norm_footprint = self.footprints[-1] * norm_coverage 

155 

156 # Compute the constant parts of the footprint like before 

157 desired = norm_footprint / self.all_footprints_sum * np.sum(self.survey_features['N_obs_all_-1'].feature) 

158 result[self.constant_footprint_indx] = desired[self.constant_footprint_indx] - self.survey_features['N_obs_-1'].feature[self.constant_footprint_indx] 

159 

160 # Now for the rolling sections 

161 norm_coverage = 0 

162 ns = np.floor((t_elapsed+self.max_day_offset)/365.25) 

163 if np.max(ns) > -1: 

164 ns_on = ns/self.season_modulo 

165 ns_off = ns - ns_on 

166 norm_coverage = norm_coverage_raw - ns_off 

167 norm_coverage = norm_coverage/np.max(norm_coverage) 

168 

169 for season in np.unique(seasons[self.rolling_footprint_indx]): 

170 if season == -1: 

171 nf = norm_footprint 

172 else: 

173 nf = self.footprints[season] * norm_coverage 

174 season_indx = np.where(seasons[self.rolling_footprint_indx] == season)[0] 

175 desired = nf[self.rolling_footprint_indx[season_indx]] / self.all_rolling_sum * np.sum(self.survey_features['N_obs_all_%i' % season].feature[self.rolling_footprint_indx]) 

176 result[self.rolling_footprint_indx[season_indx]] = desired - self.survey_features['N_obs_%i' % season].feature[self.rolling_footprint_indx][season_indx] 

177 

178 result[self.out_of_bounds_area] = self.out_of_bounds_val 

179 return result 

180 

181 

182class Target_map_modulo_basis_function(Base_basis_function): 

183 """Basis function that tracks number of observations and tries to match a specified spatial distribution 

184 can enter multiple maps that will be used at different times in the survey 

185 

186 Parameters 

187 ---------- 

188 day_offset : np.array 

189 Healpix map that has the offset to be applied to each pixel when computing what season it is on. 

190 filtername : (string 'r') 

191 The name of the filter for this target map. 

192 nside: int (default_nside) 

193 The healpix resolution. 

194 target_maps : list of numpy array (None) 

195 healpix maps showing the ratio of observations desired for all points on the sky. Last map will be used 

196 for season -1. Probably shouldn't support going to season less than -1. 

197 norm_factor : float (0.00010519) 

198 for converting target map to number of observations. Should be the area of the camera 

199 divided by the area of a healpixel divided by the sum of all your goal maps. Default 

200 value assumes LSST foV has 1.75 degree radius and the standard goal maps. If using 

201 mulitple filters, see lsst.sims.featureScheduler.utils.calc_norm_factor for a utility 

202 that computes norm_factor. 

203 out_of_bounds_val : float (-10.) 

204 Reward value to give regions where there are no observations requested (unitless). 

205 season_modulo : int (2) 

206 The value to modulate the season by (years). 

207 max_season : int (None) 

208 For seasons higher than this value (pre-modulo), the final target map is used. 

209 

210 """ 

211 def __init__(self, day_offset=None, filtername='r', nside=None, target_maps=None, 

212 norm_factor=None, out_of_bounds_val=-10., season_modulo=2, max_season=None, 

213 season_length=365.25): 

214 

215 super(Target_map_modulo_basis_function, self).__init__(nside=nside, filtername=filtername) 

216 

217 if norm_factor is None: 

218 warnings.warn('No norm_factor set, use utils.calc_norm_factor if using multiple filters.') 

219 self.norm_factor = 0.00010519 

220 else: 

221 self.norm_factor = norm_factor 

222 

223 self.survey_features = {} 

224 # Map of the number of observations in filter 

225 

226 for i, temp in enumerate(target_maps[0:-1]): 

227 self.survey_features['N_obs_%i' % i] = features.N_observations_season(i, filtername=filtername, 

228 nside=self.nside, 

229 modulo=season_modulo, 

230 offset=day_offset, 

231 max_season=max_season, 

232 season_length=season_length) 

233 # Count of all the observations taken in a season 

234 self.survey_features['N_obs_count_all_%i' % i] = features.N_obs_count_season(i, filtername=None, 

235 season_modulo=season_modulo, 

236 offset=day_offset, 

237 nside=self.nside, 

238 max_season=max_season, 

239 season_length=season_length) 

240 # Set the final one to be -1 

241 self.survey_features['N_obs_%i' % -1] = features.N_observations_season(-1, filtername=filtername, 

242 nside=self.nside, 

243 modulo=season_modulo, 

244 offset=day_offset, 

245 max_season=max_season, 

246 season_length=season_length) 

247 self.survey_features['N_obs_count_all_%i' % -1] = features.N_obs_count_season(-1, filtername=None, 

248 season_modulo=season_modulo, 

249 offset=day_offset, 

250 nside=self.nside, 

251 max_season=max_season, 

252 season_length=season_length) 

253 if target_maps is None: 

254 self.target_maps = utils.generate_goal_map(filtername=filtername, nside=self.nside) 

255 else: 

256 self.target_maps = target_maps 

257 # should probably actually loop over all the target maps? 

258 self.out_of_bounds_area = np.where(self.target_maps[0] == 0)[0] 

259 self.out_of_bounds_val = out_of_bounds_val 

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

261 self.all_indx = np.arange(self.result.size) 

262 

263 # For computing what day each healpix is on 

264 if day_offset is None: 

265 self.day_offset = np.zeros(hp.nside2npix(self.nside), dtype=float) 

266 else: 

267 self.day_offset = day_offset 

268 

269 self.season_modulo = season_modulo 

270 self.max_season = max_season 

271 self.season_length = season_length 

272 

273 def _calc_value(self, conditions, indx=None): 

274 """ 

275 Parameters 

276 ---------- 

277 indx : list (None) 

278 Index values to compute, if None, full map is computed 

279 Returns 

280 ------- 

281 Healpix reward map 

282 """ 

283 

284 result = self.result.copy() 

285 if indx is None: 

286 indx = self.all_indx 

287 

288 # Compute what season it is at each pixel 

289 seasons = utils.season_calc(conditions.night, offset=self.day_offset, 

290 modulo=self.season_modulo, max_season=self.max_season, 

291 season_length=self.season_length) 

292 

293 composite_target = self.result.copy()[indx] 

294 composite_nobs = self.result.copy()[indx] 

295 

296 composite_goal_N = self.result.copy()[indx] 

297 

298 for season in np.unique(seasons): 

299 season_indx = np.where(seasons == season)[0] 

300 composite_target[season_indx] = self.target_maps[season][season_indx] 

301 composite_nobs[season_indx] = self.survey_features['N_obs_%i' % season].feature[season_indx] 

302 composite_goal_N[season_indx] = composite_target[season_indx] * self.survey_features['N_obs_count_all_%i' % season].feature * self.norm_factor 

303 

304 result[indx] = composite_goal_N - composite_nobs[indx] 

305 result[self.out_of_bounds_area] = self.out_of_bounds_val 

306 

307 return result