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# (c) Stefan Countryman, 2016-2018 

2 

3""" 

4FileHandlers associated with HEALPix skymap loading and manipulation. Provides 

5a subclass of FileHandler that specializes in loading and saving HEALPix FITS 

6file skymaps and putting arbitrary skymap pixelizations and functions on the 

7sphere into identical HEALPix pixelization schemes to simplify mathematical 

8analyses. 

9 

10Some documentation: 

11http://healpix.sourceforge.net/html/intronode4.htm 

12https://healpy.readthedocs.io/en/latest/healpy_pix.html 

13""" 

14 

15import sys 

16import types 

17import logging 

18import operator 

19import textwrap 

20from functools import reduce 

21from llama.filehandler import JSONFile, FileHandler 

22from llama.utils import archive_figs 

23from llama.files.healpix.skymap import ( 

24 HEALPixSkyMap, 

25 HEALPixRepresentableFileHandler, 

26 HEALPixPSF, 

27) 

28 

29LOGGER = logging.getLogger(__name__) 

30 

31# smallest possible feature size in radians 

32MIN_SCALE_RAD = 4e-9 

33 

34 

35def pixel_radii_rad(): 

36 """Get the maximum pixel radii in radians for each value of ``k``, where 

37 ``nside = 2**k``. 

38 """ 

39 import healpy as hp 

40 return hp.pixelfunc.max_pixrad([2**k for k in range(30)]) 

41 

42 

43def scale2nside(scale, degrees=True): 

44 """Find the correct HEALPix scale NSIDE for a given feature resolution. For 

45 example, to capture features of size 0.2 degrees across, you want a HEALPix 

46 skymap with NSIDE large enough so that the max inter-pixel distance is less 

47 than 0.1 degrees. Note that HEALPix skymaps cannot have NSIDE > 2**30; if 

48 higher resolution is required, a ``ValueError`` will be raised. 

49 

50 Parameters 

51 ---------- 

52 scale : int 

53 The angular size of the smallest feature you need to capture (in the 

54 specified units). 

55 degrees : bool, optional 

56 If true, the feature size is assumed to be in degrees. Otherwise, it is 

57 taken to be in radians. 

58 

59 Returns 

60 ------- 

61 nside : int 

62 The minimum NSIDE parameter required to capture features with scale 

63 larger than ``scale``. Will always be a power of 2 less than or equal 

64 to ``2**30``. 

65 

66 Raises 

67 ------ 

68 ValueError 

69 If the feature scale is too small to represent with a HEALPix map. 

70 """ 

71 import numpy as np 

72 if degrees: 

73 scale = np.radians(scale) # pylint: disable=assignment-from-no-return 

74 if scale < MIN_SCALE_RAD: 

75 raise ValueError('scale is too small to represent in HEALPix.') 

76 return 2**(np.argwhere(pixel_radii_rad() < scale / 2)[0][0] + 1) 

77 

78 

79def subpixels(indices, inputnside, nside, nest=True): 

80 """Take a list of indices into a HEALPix skymap with ``NSIDE=inputnside`` 

81 and return a list of indices into a HEALPix skymap with ``NSIDE=nside`` 

82 where the returned list of indices correspond to subpixels of the pixels 

83 indexed from the original skymap size. Useful for getting higher-resolution 

84 subpixels for a list of pixels of interest. 

85 

86 Returns 

87 ------- 

88 subindices : array 

89 Indices for a HEALPix skymap with ``NSIDE=nside`` that correspond to 

90 subpixels from the input pixel indices, sorted in ascending order. 

91 

92 Parameters 

93 ---------- 

94 indices : array 

95 Indices for a HEALPix skymap with ``NSIDE=inputnside``. Repeated 

96 indices are ignored. 

97 inputnside : int 

98 HEALPix NSIDE (i.e. resolution) for the input pixels. 

99 nside : int 

100 NSIDE for the output pixels (must be greater than ``inputnside``) 

101 nest : bool, optional 

102 If true, input and output skymaps are interpreted as havin NEST 

103 ordering. If false, they are interpreted as having RING ordering. 

104 

105 Raises 

106 ------ 

107 ValueError 

108 If ``nside`` or ``inputnside`` are not valid NSIDE HEALPix parameters 

109 (i.e. not equal to 2**k for some integer k between 0 and 29), if 

110 ``inputnside`` is greater than ``nside``, or if the ``indices`` are 

111 invalid. 

112 """ 

113 import numpy as np 

114 import healpy as hp 

115 if not (hp.isnsideok(nside) and hp.isnsideok(inputnside)): 

116 raise ValueError("Invalid nside ({}) or inputnside ({}).".format( 

117 nside, inputnside 

118 )) 

119 if inputnside > nside: 

120 raise ValueError("inputside ({}) greater than nside ({}).".format( 

121 inputnside, nside 

122 )) 

123 maxind = hp.nside2npix(inputnside) - 1 

124 inds = np.array(indices) 

125 if np.any(np.logical_or(inds < 0, inds > maxind)): 

126 raise ValueError( 

127 "Input indices out of bounds for inputnside={}:\n{}".format( 

128 inputnside, inds 

129 ) 

130 ) 

131 if not nest: 

132 inds = hp.ring2nest(inputnside, inds) 

133 subdiv = int(nside) // int(inputnside) 

134 subindices = np.array( 

135 sorted( 

136 sum( 

137 [range(subdiv**2 * i, subdiv**2 * (i+1)) for i in set(inds)], 

138 [] 

139 ) 

140 ) 

141 ) 

142 if not nest: 

143 subindices = hp.nest2ring(nside, subindices) 

144 subindices.sort() 

145 return subindices 

146 

147 

148# pylint: disable=W0223 

149class HEALPixSkyMapFileHandler(HEALPixRepresentableFileHandler): 

150 """ 

151 A FileHandler that corresponds to a file object containing HEALPix 

152 skymaps. Provides methods for reading HEALPix settings (i.e. Nsides, 

153 nesting method, and coordinate system) as well as reading and writing 

154 HEALPix files from ``numpy`` arrays. Methods do NOT assume a single skymap, 

155 and in general HEALPix-specific methods return an object for each skymap 

156 contained in the file described by this file handler in the form of a 

157 (possibly empty!) list. 

158 """ 

159 

160 @staticmethod 

161 def read_healpix_from_file(infile, nest=True): 

162 """Read the HEALPix skymap stored in ``infile`` and return a list 

163 containing a single ``HEALPixSkyMap`` instance containing the skymap 

164 and metadata (returns a list so that other implementations of this 

165 class can optionally return many skymaps). By default, return in 

166 HEALPix NEST ordering; optionally, return in RING ordering by 

167 specifying ``nest=False``. ``infile`` can be either a FITS image file 

168 (inferred from file extension ".fits" or ".fits.gz") or a 

169 specially-formatted HDF5 file (inferred from file extension ".hdf5") of 

170 the sort produced by ``write_healpix_to_file``. Note that if loading 

171 from HDF5 files, ``nest`` is ignored. 

172 

173 >>> import tempfile, os 

174 >>> f = tempfile.NamedTemporaryFile(suffix='.hdf5', delete=False) 

175 >>> fname = f.name 

176 >>> f.close() 

177 >>> skymap0 = HEALPixSkyMap(range(12)) 

178 >>> skymap1 = HEALPixSkyMap([1.]*12) 

179 >>> HEALPixSkyMapFileHandler.write_healpix_to_file( 

180 ... [skymap0, skymap1], 

181 ... fname 

182 ... ) 

183 >>> skymaps = HEALPixSkyMapFileHandler.read_healpix_from_file( 

184 ... fname 

185 ... ) 

186 >>> skymaps[0] == skymap0 

187 True 

188 >>> skymaps[1] == skymap1 

189 True 

190 """ 

191 infilelower = infile.lower() 

192 if infilelower.endswith('.fits.gz') or infilelower.endswith('.fits'): 

193 import healpy as hp 

194 return [HEALPixSkyMap(hp.read_map(infile, nest=nest, 

195 verbose=False), nest=nest)] 

196 elif infilelower.endswith('.hdf5'): 

197 import h5py 

198 with h5py.File(infile, 'r') as h5file: 

199 pixels = h5file['skymaps/pixels'][()] 

200 return [HEALPixSkyMap(pixels[i], nest=nst) 

201 for i, nst in enumerate(h5file['skymaps/nest'][()])] 

202 

203 def get_healpix(self, nest=True, template_skymap=None): 

204 """Get a list of HEALPix skymaps stored in the file corresponding to 

205 this ``FileHandler``. 

206 

207 Parameters 

208 ---------- 

209 nest : bool, optional 

210 If ``True``, return skymaps in HEALPix NEST ordering; if ``False``, 

211 return in RING ordering. 

212 template_skymap : HEALPixSkyMap, optional 

213 If provided, force the loaded skymaps to have the same HEALPix 

214 parameters as ``template_skymap``. Overrides the ``nest`` 

215 parameter. 

216 

217 Returns 

218 ------- 

219 skymaps : list 

220 ``HEALPixSkyMap`` instances with the data contained in this file. 

221 """ 

222 skymap = self.read_healpix_from_file(self.fullpath, nest=nest) 

223 if template_skymap is not None: 

224 skymap = template_skymap.conform_skymaps(skymap)[0] 

225 return skymap 

226 

227 def get_healpix_descriptions(self, skymap_indices=None): 

228 """Get a list of descriptive strings for each HEALPix skymap contained 

229 in this file. 

230 

231 Parameters 

232 ---------- 

233 skymap_indices : list 

234 A list of tuples describing the index of each skymap used to 

235 produce this output skymap, with each entry in the tuple 

236 corresponding to a detector. In this case, the descriptive strings 

237 will be generated based on this ``FileHandler`` instance's 

238 ``DETECTORS`` property and the indices of the skymap_indices in 

239 ``skymap_indices``. The ordering of the returned descriptive 

240 strings will match the ordering of skymap_indices. If 

241 ``skymap_indices`` is not specified, the HEALPix skymaps and their 

242 descriptions will be read from file. Note that for single-detector 

243 skymaps, the elements of the ``skymap_indices`` list should just be 

244 single-element tuples corresponding to the index of each skymap. 

245 """ 

246 import h5py 

247 if skymap_indices is None: 

248 with h5py.File(self.fullpath, 'r') as h5file: 

249 return list(h5file['skymaps/descriptions'][()]) 

250 return super(HEALPixSkyMapFileHandler, 

251 self).get_healpix_descriptions(skymap_indices) 

252 

253 @staticmethod 

254 def write_healpix_to_file(skymaps, outfile, descriptions=None): 

255 """Write ``skymaps`` to ``outfile``. 

256 

257 Parameters 

258 ---------- 

259 skymaps : HEALPixSkyMap or list 

260 ``HEALPixSkyMap`` instance or list thereof to be written to file. 

261 outfile : str 

262 Path to write skymap to. File extension must conform to either a 

263 FITS image file (inferred from file extension ".fits" or 

264 ".fits.gz") or a specially-formatted HDF5 file (inferred from file 

265 extension ".hdf5"). 

266 descriptions : list, optional 

267 A list of descriptive strings, each corresponding to the skymap in 

268 ``skymaps`` of the same index. 

269 """ 

270 import numpy as np 

271 outfilelower = outfile.lower() 

272 if outfilelower.endswith('.fits.gz') or outfilelower.endswith('.fits'): 

273 import healpy as hp 

274 if isinstance(skymaps, list): 

275 if len(skymaps) != 1: 

276 raise ValueError('Can only save 1 skymap to .FITS file') 

277 skymaps = skymaps[0] 

278 hp.write_map(outfile, skymaps.pixels, nest=skymaps.nest, 

279 dtype=np.float) 

280 elif outfilelower.endswith('.hdf5'): 

281 import h5py 

282 if not isinstance(skymaps, list): 

283 skymaps = [skymaps] 

284 with h5py.File(outfile, 'w') as h5file: 

285 pixels = np.array([m.pixels for m in skymaps]) 

286 nest = np.array([m.nest for m in skymaps]) 

287 maps = h5file.create_group('skymaps') 

288 maps.create_dataset('pixels', data=pixels, compression='gzip') 

289 maps.create_dataset('nest', data=nest) 

290 if descriptions is not None: 

291 maps.create_dataset( 

292 'descriptions', 

293 data=[d.encode('utf-8') for d in descriptions], 

294 compression='gzip' 

295 ) 

296 

297 def write_healpix(self, skymaps, descriptions=None): 

298 """Write ``skymap``, which should be a ``HEALPixSkyMap`` instance or a 

299 list containing a *single* ``HEALPixSkyMap`` instance, to the full path 

300 specified by this ``FileHandler`` instance. ``outfile`` can be either a 

301 FITS image file (inferred from file extension ".fits" or ".fits.gz") or 

302 a specially-formatted HDF5 file (inferred from file extension 

303 ".hdf5").""" 

304 if descriptions is None: 

305 descriptions = self.get_healpix_descriptions( 

306 skymap_indices=[(i,) for i in range(len(skymaps))] 

307 ) 

308 self.write_healpix_to_file(skymaps, self.fullpath, 

309 descriptions=descriptions) 

310 

311 @property 

312 def num_triggers(self): 

313 """The number of triggers described by this file. Useful mostly for 

314 quickly determining if this trigger list is empty.""" 

315 fnamelower = self.FILENAME.lower() 

316 if fnamelower.endswith('.fits.gz') or fnamelower.endswith('.fits'): 

317 return 1 

318 elif fnamelower.endswith('.hdf5'): 

319 return len(self.get_healpix_descriptions()) 

320 

321 

322class LvcHEALPixSkyMapFileHandler(HEALPixSkyMapFileHandler): 

323 """ 

324 A ``HEALPixSkyMapFileHandler`` using ``astropy.table.Table`` to read 

325 the input file (in order to accomodate the NUNIQ pixel 

326 ordering/multi-resolution skymaps). 

327 """ 

328 

329 @staticmethod 

330 def read_healpix_from_file(infile, nest=True): 

331 """Read the HEALPix skymap stored in ``infile`` and return a list 

332 containing a single ``HEALPixSkyMap`` instance containing the skymap 

333 and metadata (returns a list so that other implementations of this 

334 class can optionally return many skymaps). By default, return in 

335 HEALPix NEST ordering; optionally, return in RING ordering by 

336 specifying ``nest=False``. ``infile`` can be either a FITS image file 

337 (inferred from file extension ".fits" or ".fits.gz") or a 

338 specially-formatted HDF5 file (inferred from file extension ".hdf5") of 

339 the sort produced by ``write_healpix_to_file``. Note that if loading 

340 from HDF5 files, ``nest`` is ignored. 

341 

342 >>> import tempfile, os 

343 >>> f = tempfile.NamedTemporaryFile(suffix='.hdf5', delete=False) 

344 >>> fname = f.name 

345 >>> f.close() 

346 >>> skymap0 = HEALPixSkyMap(range(12)) 

347 >>> skymap1 = HEALPixSkyMap([1.]*12) 

348 >>> HEALPixSkyMapFileHandler.write_healpix_to_file( 

349 ... [skymap0, skymap1], 

350 ... fname 

351 ... ) 

352 >>> skymaps = HEALPixSkyMapFileHandler.read_healpix_from_file( 

353 ... fname 

354 ... ) 

355 >>> skymaps[0] == skymap0 

356 True 

357 >>> skymaps[1] == skymap1 

358 True 

359 """ 

360 import h5py 

361 from ligo.skymap.io import fits 

362 infilelower = infile.lower() 

363 if infilelower.endswith('.fits.gz') or infilelower.endswith('.fits'): 

364 skymap, metadata = fits.read_sky_map(infile, nest=None) 

365 hpskymap = HEALPixSkyMap(skymap, nest=metadata['nest']) 

366 if metadata['nest'] and (not nest): 

367 return [hpskymap.to_ring()] 

368 elif (not metadata['nest']) and nest: 

369 return [hpskymap.to_nest()] 

370 return [hpskymap] 

371 elif infilelower.endswith('.hdf5'): 

372 with h5py.File(infile, 'r') as h5file: 

373 pixels = h5file['skymaps/pixels'][()] 

374 return [HEALPixSkyMap(pixels[i], nest=nst) 

375 for i, nst in enumerate(h5file['skymaps/nest'][()])] 

376 

377 

378class HEALPixSkyMapAnalysis(HEALPixSkyMapFileHandler): 

379 """ 

380 Create output skymaps that are outputs of functions on input HEALPix 

381 skymaps. 

382 """ 

383 

384 FILENAME_PREFIX = None 

385 

386 _REQUIRED = ("FILENAME_PREFIX",) 

387 

388 @staticmethod 

389 def analysis_kernel(*skymaps): 

390 """A function that operates on an arbitrary number of skymaps and 

391 returns a single skymap. This is where the actual analysis algorithm is 

392 defined. Defaults to a simple elementwise product of skymaps, but can 

393 be overridden to do more complex analyses.""" 

394 return reduce(operator.mul, skymaps) 

395 

396 def analyze_skymaps(self, *skymap_lists): 

397 """A function that takes an arbitrary number N of lists of input 

398 skymaps as its arguments and returns ``(results, indices)``, a tuple of 

399 lists where ``results`` is a list of results of combinations of skymaps 

400 from each skymap list (as output by ``analysis_kernel``) and 

401 ``indices`` is a list of tuples of indices indicating *which* skymap 

402 was used from each skymap list.""" 

403 if 0 in [len(l) for l in skymap_lists]: 

404 return ([], []) 

405 skymap_indices = [tuple()] 

406 results = list() 

407 for skymap_list in skymap_lists: 

408 new_indices = list() 

409 for skymap_index in skymap_indices: 

410 new_indices += [skymap_index + (i,) 

411 for i in range(len(skymap_list))] 

412 skymap_indices = new_indices 

413 for ituple in skymap_indices: 

414 skymaps = list() 

415 for list_number, index in enumerate(ituple): 

416 skymaps.append(skymap_lists[list_number][index]) 

417 results.append(self.analysis_kernel(*skymaps)) 

418 return (results, skymap_indices) 

419 

420 @classmethod 

421 def set_class_attributes(cls, subclass): 

422 """See ``FileHandler.set_class_attributes``; this method additionally 

423 sets the ``FILENAME`` and ``DETECTORS`` attributes based on 

424 ``subclass.DEPENDENCIES`` and ``subclass.FILENAME_PREFIX``.""" 

425 subclass.DETECTORS = tuple(sorted({d for fh in subclass.DEPENDENCIES 

426 for d in fh.DETECTORS})) 

427 detstr = '-'.join(d.abbrev.lower() for d in subclass.DETECTORS) 

428 subclass.FILENAME = '{}_{}.hdf5'.format(subclass.FILENAME_PREFIX, 

429 detstr) 

430 super().set_class_attributes(subclass) 

431 return subclass 

432 

433 def write_healpix(self, skymaps, descriptions=None): 

434 """Write ``skymap``, which should be a ``HEALPixSkyMap`` instance or a 

435 list containing a *single* ``HEALPixSkyMap`` instance, to the full path 

436 specified by this ``FileHandler`` instance. ``outfile`` must be an HDF5 

437 file (inferred from file extension ".hdf5").""" 

438 self.write_healpix_to_file(skymaps, self.fullpath, 

439 descriptions=descriptions) 

440 

441 def _generate(self): # pylint: disable=W0221 

442 """Run the ``analysis_kernel`` on all permutations of input skymap 

443 lists and save the resulting skymaps to a compressed HDF5 file.""" 

444 assert all(len(dep(self).DETECTORS) == 1 for dep in self.DEPENDENCIES) 

445 deps = sorted(self.DEPENDENCIES, 

446 key=lambda f: f(self).DETECTORS[0].name) 

447 if 0 in [d(self).num_triggers for d in deps]: 

448 results = [] 

449 descriptions = [] 

450 else: 

451 skymap_lists = [fh(self).get_healpix() for fh in deps] 

452 results, inds = self.analyze_skymaps(*skymap_lists) 

453 descriptions = self.get_healpix_descriptions(skymap_indices=inds) 

454 self.write_healpix(results, descriptions=descriptions) 

455 

456 

457class HEALPixSkyMapStats(JSONFile): 

458 """ 

459 Load a skymap from a ``HEALPixRepresentableFileHandler`` instance and 

460 generate stats on it, e.g. the likelihood integral for a 

461 ``HEALPixSkyMapAnalysis`` instance, and store them in JSON format. By 

462 default, only has one skymap file as a dependency and stores the integral 

463 of that skymap as ``integrated_likelihood``. 

464 """ 

465 

466 def _generate(self): # pylint: disable=W0221 

467 skymap_fh = self.DEPENDENCIES[0](self) 

468 skymaps = skymap_fh.get_healpix() 

469 descriptions = skymap_fh.get_healpix_descriptions() 

470 assert len(descriptions) == len(skymaps) 

471 stats = dict() 

472 for i, desc in enumerate(descriptions): 

473 stats[desc] = {'integrated_likelihood': skymaps[i].integrate()} 

474 self._write_json(stats) 

475 

476 @classmethod 

477 def set_class_attributes(cls, subclass): 

478 """See ``FileHandler.set_class_attributes``; this method sets 

479 programatically-generated filenames for this class.""" 

480 subclass.FILENAME = subclass.DEPENDENCIES[0].FILENAME+'-stats.json' 

481 super().set_class_attributes(subclass) 

482 

483 

484class HEALPixPlots(FileHandler): 

485 """ 

486 Load HEALPix skymaps and save Mollweide-projection plots in PDF and PNG 

487 format to a compressed tarfile. 

488 """ 

489 

490 def _generate(self): # pylint: disable=W0221 

491 skymap_fh = self.DEPENDENCIES[0](self) 

492 skymaps = skymap_fh.get_healpix() 

493 descriptions = skymap_fh.get_healpix_descriptions() 

494 plots = [m.mollview(title=descriptions[i]) 

495 for i, m in enumerate(skymaps)] 

496 archive_figs(plots, self.fullpath, fname_list=descriptions) 

497 

498 @classmethod 

499 def set_class_attributes(cls, subclass): 

500 """See ``FileHandler.set_class_attributes``; this method sets 

501 programatically-generated filenames for this class.""" 

502 subclass.FILENAME = subclass.DEPENDENCIES[0].FILENAME+'-plots.tar.gz' 

503 super().set_class_attributes(subclass) 

504 

505 

506# pylint: disable=R0913 

507def healpix_skymap_analysis_factory(dependencies, 

508 bases=(HEALPixSkyMapAnalysis,), 

509 make_plotters=False, 

510 make_stats=False, 

511 noclobber=False): 

512 """Programmatically create analysis classes based on input skymaps and 

513 analysis algorithm. **New classes are saved to the calling namespace** 

514 (consider this a shortcut macro for making ``FileHandler`` subclasses 

515 programatically). 

516 

517 Parameters 

518 ---------- 

519 dependencies : list 

520 ``FileHandler`` subclasses that implement the ``get_healpix`` method 

521 and can thus be used to run an analysis using some 

522 ``HEALPixSkyMapAnalysis`` implementation's ``analysis_kernel``. 

523 bases : tuple 

524 Base classes for the class being generated in the form of a tuple, e.g. 

525 ``(base1, base2, ...)`` 

526 make_plotters : bool, optional 

527 A boolean specifying whether ``HEALPixPlots`` class implementations 

528 should be created to accompany each ``HEALPixSkyMapAnalysis`` class 

529 generated. 

530 make_stats : bool, optional 

531 A boolean specifying whether ``HEALPixSkyMapStats`` implementations 

532 should be created to accompany each ``HEALPixSkyMapAnalysis`` class 

533 generated. 

534 noclobber : bool, optional 

535 If true, don't overwrite existing variable names; instead, throw 

536 a ``ValueError``. 

537 

538 Returns 

539 ------- 

540 new_classes : dict 

541 The new classes generated, where each key is the name of the new class 

542 and each value is the new class itself. 

543 """ 

544 namespace = sys._getframe(1).f_globals 

545 module = namespace['__name__'] 

546 for deps in dependencies: 

547 detstr = ''.join(sorted( 

548 { 

549 d.abbrev.capitalize() for fh in deps 

550 for d in fh.DETECTORS 

551 } 

552 )) 

553 name = 'CoincSkymap{}HDF5'.format(detstr) 

554 docs = _healpix_skymap_analysis_factory_docstring(deps, bases) 

555 newclassdict = { 

556 'FILENAME_PREFIX': 'likelihood', 

557 'DEPENDENCIES': tuple(deps), 

558 '__module__': module, 

559 '__doc__': docs[0], 

560 } 

561 # info on making new classes dynamically: 

562 # https://stackoverflow.com/questions/15247075 

563 newclass = type(name, bases, newclassdict) 

564 setters = [b for b in bases if hasattr(b, 'set_class_attributes')] 

565 if setters: 

566 setters[0].set_class_attributes(newclass) 

567 if noclobber and name in namespace: 

568 raise ValueError('{} already defined, noclobber on'.format(name)) 

569 namespace[name] = newclass 

570 if make_plotters: 

571 # TODO stop archiving figures as TarGz files, instead, use subplots 

572 # and ``hold=True`` with mollview to put them all on single figures 

573 name = 'CoincLikelihoodPlts{}'.format(detstr) 

574 newclassdict = { 

575 'DEPENDENCIES': (newclass,), 

576 '__module__': module, 

577 '__doc__': docs[1] 

578 } 

579 if noclobber and name in namespace: 

580 raise ValueError(f'{name} already defined, noclobber on') 

581 newplots = type(name, (HEALPixPlots,), newclassdict) 

582 HEALPixPlots.set_class_attributes(newplots) 

583 namespace[name] = newplots 

584 if make_stats: 

585 name = 'CoincLikelihoodStats{}'.format(detstr) 

586 newclassdict = { 

587 'DEPENDENCIES': (newclass,), 

588 '__module__': module, 

589 '__doc__': docs[2] 

590 } 

591 if noclobber and name in namespace: 

592 raise ValueError(f'{name} already defined, noclobber on') 

593 newstats = type(name, (HEALPixSkyMapStats,), newclassdict) 

594 HEALPixSkyMapStats.set_class_attributes(newstats) 

595 namespace[name] = newstats 

596 return namespace 

597 

598 

599def _healpix_skymap_analysis_factory_docstring(deps, bases): 

600 """Get the docstrings for classes generated from the given args. 

601 

602 Parameters 

603 ---------- 

604 deps : list 

605 Dependencies for the generated ``FileHandler`` subclasses. Each must 

606 implement ``get_healpix``. 

607 bases : tuple 

608 Base classes for the class being generated in the form of a tuple, e.g. 

609 ``(base1, base2, ...)`` 

610 

611 Returns 

612 ------- 

613 analysis_doc : str 

614 Docstring for a ``HEALPixSkyMapAnalysis`` instance that would be 

615 generated from the given args. 

616 plot_doc : str 

617 Same, but for a ``HEALPixPlots`` instance. 

618 stats_doc : str 

619 Same, but for a ``HEALPixSkyMapStats`` instance. 

620 """ 

621 docstring_fmt = """{} analysis products of {} skymaps generated using the 

622 ``{}.analysis_kernel`` algorithm.""" 

623 if len(deps) == 1: 

624 fmt = '``{}``' 

625 elif len(deps) == 2: 

626 fmt = '``{}`` and ``{}``' 

627 else: 

628 fmt = ('``{}``, '*(len(deps)-1) + 'and ``{}``') 

629 input_names = fmt.format(*[d.__name__ for d in deps]) 

630 for base in type('', bases, {}).mro()[1:]: 

631 if hasattr(base, 'analysis_kernel'): 

632 if isinstance(base.analysis_kernel, types.FunctionType): 

633 kernel_class = base 

634 break 

635 analysis_doc = textwrap.fill( 

636 docstring_fmt.format('Skymap', input_names, kernel_class.__name__), 

637 width=76 

638 ) 

639 plot_doc = textwrap.fill( 

640 docstring_fmt.format( 

641 'Skymap plots of', 

642 input_names, 

643 kernel_class.__name__ 

644 ), 

645 width=76 

646 ) 

647 stats_doc = textwrap.fill( 

648 docstring_fmt.format( 

649 'Integrated likelihoods calculated from', 

650 input_names, 

651 kernel_class.__name__ 

652 ), 

653 width=76 

654 ) 

655 return (analysis_doc, plot_doc, stats_doc)