LSSTApplications  17.0+103,17.0+11,17.0+61,18.0.0+13,18.0.0+25,18.0.0+5,18.0.0+52,18.0.0-4-g68ffd23,18.1.0-1-g0001055+8,18.1.0-1-g03d53ef+1,18.1.0-1-g1349e88+28,18.1.0-1-g2505f39+22,18.1.0-1-g380d4d4+27,18.1.0-1-g5315e5e+1,18.1.0-1-g5e4b7ea+10,18.1.0-1-g7e8fceb+1,18.1.0-1-g85f8cd4+23,18.1.0-1-g9a6769a+13,18.1.0-1-ga1a4c1a+22,18.1.0-1-gd55f500+17,18.1.0-12-g42eabe8e+10,18.1.0-14-gd04256d+15,18.1.0-16-g430f6a53+1,18.1.0-17-gd2166b6e4,18.1.0-18-gb5d19ff+1,18.1.0-2-gfbf3545+7,18.1.0-2-gfefb8b5+16,18.1.0-3-g52aa583+13,18.1.0-3-g62b5e86+14,18.1.0-3-g8f4a2b1+17,18.1.0-3-g9bc06b8+7,18.1.0-3-gb69f684+9,18.1.0-4-g1ee41a7+1,18.1.0-5-g6dbcb01+13,18.1.0-5-gc286bb7+3,18.1.0-6-g48bdcd3+2,18.1.0-6-gd05e160+9,18.1.0-7-gc4d902b+2,18.1.0-7-gebc0338+8,18.1.0-9-gae7190a+10,w.2019.38
LSSTDataManagementBasePackage
anetBasicAstrometry.py
Go to the documentation of this file.
1 #
2 # LSST Data Management System
3 #
4 # Copyright 2008-2015 AURA/LSST.
5 #
6 # This product includes software developed by the
7 # LSST Project (http://www.lsst.org/).
8 #
9 # This program is free software: you can redistribute it and/or modify
10 # it under the terms of the GNU General Public License as published by
11 # the Free Software Foundation, either version 3 of the License, or
12 # (at your option) any later version.
13 #
14 # This program is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 # GNU General Public License for more details.
18 #
19 # You should have received a copy of the LSST License Statement and
20 # the GNU General Public License along with this program. If not,
21 # see <https://www.lsstcorp.org/LegalNotices/>.
22 #
23 __all__ = ["InitialAstrometry", "ANetBasicAstrometryConfig", "ANetBasicAstrometryTask"]
24 
25 import math
26 import sys
27 
28 import numpy as np
29 
30 import lsst.daf.base as dafBase
31 from lsst.pex.config import Field, RangeField, ListField
32 import lsst.pex.exceptions as pexExceptions
33 import lsst.pipe.base as pipeBase
34 import lsst.afw.image as afwImage
35 import lsst.afw.math as afwMath
36 import lsst.afw.table as afwTable
37 import lsst.geom as geom
38 import lsst.meas.algorithms.utils as maUtils
39 from .loadAstrometryNetObjects import LoadAstrometryNetObjectsTask, LoadMultiIndexes
40 from lsst.meas.astrom import displayAstrometry, makeMatchStatisticsInRadians
41 import lsst.meas.astrom.sip as astromSip
42 from . import cleanBadPoints
43 
44 
46  """
47  Object returned by Astrometry.determineWcs
48 
49  getWcs(): sipWcs or tanWcs
50  getMatches(): sipMatches or tanMatches
51 
52  Other fields are:
53  solveQa (PropertyList)
54  tanWcs (Wcs)
55  tanMatches (MatchList)
56  sipWcs (Wcs)
57  sipMatches (MatchList)
58  refCat (lsst::afw::table::SimpleCatalog)
59  matchMeta (PropertyList)
60  """
61 
62  def __init__(self):
63  self.tanWcs = None
64  self.tanMatches = None
65  self.sipWcs = None
66  self.sipMatches = None
67  self.refCat = None
69  self.solveQa = None
70 
71  def getMatches(self):
72  """
73  Get "sipMatches" -- MatchList using the SIP WCS solution, if it
74  exists, or "tanMatches" -- MatchList using the TAN WCS solution
75  otherwise.
76  """
77  return self.sipMatches or self.tanMatches
78 
79  def getWcs(self):
80  """
81  Returns the SIP WCS, if one was found, or a TAN WCS
82  """
83  return self.sipWcs or self.tanWcs
84 
85  matches = property(getMatches)
86  wcs = property(getWcs)
87 
88  # "Not very pythonic!" complains Paul.
89  # Consider these methods deprecated; if you want these elements, just
90  # .grab them.
91  def getSipWcs(self):
92  return self.sipWcs
93 
94  def getTanWcs(self):
95  return self.tanWcs
96 
97  def getSipMatches(self):
98  return self.sipMatches
99 
100  def getTanMatches(self):
101  return self.tanMatches
102 
103  def getMatchMetadata(self):
104  return self.matchMeta
105 
107  return self.solveQa
108 
109 
110 class ANetBasicAstrometryConfig(LoadAstrometryNetObjectsTask.ConfigClass):
111 
112  maxCpuTime = RangeField(
113  doc="Maximum CPU time to spend solving, in seconds",
114  dtype=float,
115  default=0.,
116  min=0.,
117  )
118  matchThreshold = RangeField(
119  doc="Matching threshold for Astrometry.net solver (log-odds)",
120  dtype=float,
121  default=math.log(1e12),
122  min=math.log(1e6),
123  )
124  maxStars = RangeField(
125  doc="Maximum number of stars to use in Astrometry.net solving",
126  dtype=int,
127  default=50,
128  min=10,
129  )
130  useWcsPixelScale = Field(
131  doc="Use the pixel scale from the input exposure\'s WCS headers?",
132  dtype=bool,
133  default=True,
134  )
135  useWcsRaDecCenter = Field(
136  doc="Use the RA,Dec center information from the input exposure\'s WCS headers?",
137  dtype=bool,
138  default=True,
139  )
140  useWcsParity = Field(
141  doc="Use the parity (flip / handedness) of the image from the input exposure\'s WCS headers?",
142  dtype=bool,
143  default=True,
144  )
145  raDecSearchRadius = RangeField(
146  doc="When useWcsRaDecCenter=True, this is the radius, in degrees, around the RA,Dec center "
147  r"specified in the input exposure\'s WCS to search for a solution.",
148  dtype=float,
149  default=1.0,
150  min=0.0,
151  )
152  pixelScaleUncertainty = RangeField(
153  doc="Range of pixel scales, around the value in the WCS header, to search. "
154  "If the value of this field is X and the nominal scale is S, "
155  "the range searched will be S/X to S*X",
156  dtype=float,
157  default=1.1,
158  min=1.001,
159  )
160  catalogMatchDist = RangeField(
161  doc="Matching radius (arcsec) for matching sources to reference objects",
162  dtype=float,
163  default=1.0,
164  min=0.0,
165  )
166  cleaningParameter = RangeField(
167  doc="Sigma-clipping parameter in cleanBadPoints.py",
168  dtype=float,
169  default=3.0,
170  min=0.0,
171  )
172  calculateSip = Field(
173  doc="Compute polynomial SIP distortion terms?",
174  dtype=bool,
175  default=True,
176  )
177  sipOrder = RangeField(
178  doc="Polynomial order of SIP distortion terms",
179  dtype=int,
180  default=4,
181  min=2,
182  )
183  badFlags = ListField(
184  doc="List of flags which cause a source to be rejected as bad",
185  dtype=str,
186  default=[
187  "slot_Centroid_flag", # bad centroids
188  "base_PixelFlags_flag_edge",
189  "base_PixelFlags_flag_saturated",
190  "base_PixelFlags_flag_crCenter", # cosmic rays
191  ],
192  )
193  allFluxes = Field(
194  doc="Retrieve all available fluxes (and errors) from catalog?",
195  dtype=bool,
196  default=True,
197  )
198  maxIter = RangeField(
199  doc="maximum number of iterations of match sources and fit WCS"
200  "ignored if not fitting a WCS",
201  dtype=int,
202  default=5,
203  min=1,
204  )
205  matchDistanceSigma = RangeField(
206  doc="The match and fit loop stops when maxMatchDist minimized: "
207  " maxMatchDist = meanMatchDist + matchDistanceSigma*stdDevMatchDistance "
208  " (where the mean and std dev are computed using outlier rejection);"
209  " ignored if not fitting a WCS",
210  dtype=float,
211  default=2,
212  min=0,
213  )
214 
215 
216 class ANetBasicAstrometryTask(pipeBase.Task):
217  """!Basic implemeentation of the astrometry.net astrometrical fitter
218 
219  A higher-level class ANetAstrometryTask takes care of dealing with the fact
220  that the initial WCS is probably only a pure TAN SIP, yet we may have
221  significant distortion and a good estimate for that distortion.
222 
223 
224  About Astrometry.net index files (astrometry_net_data):
225 
226  There are three components of an index file: a list of stars
227  (stored as a star kd-tree), a list of quadrangles of stars ("quad
228  file") and a list of the shapes ("codes") of those quadrangles,
229  stored as a code kd-tree.
230 
231  Each index covers a region of the sky, defined by healpix nside
232  and number, and a range of angular scales. In LSST, we share the
233  list of stars in a part of the sky between multiple indexes. That
234  is, the star kd-tree is shared between multiple indices (quads and
235  code kd-trees). In the astrometry.net code, this is called a
236  "multiindex".
237 
238  It is possible to "unload" and "reload" multiindex (and index)
239  objects. When "unloaded", they consume no FILE or mmap resources.
240 
241  The multiindex object holds the star kd-tree and gives each index
242  object it holds a pointer to it, so it is necessary to
243  multiindex_reload_starkd() before reloading the indices it holds.
244  The multiindex_unload() method, on the other hand, unloads its
245  starkd and unloads each index it holds.
246  """
247  ConfigClass = ANetBasicAstrometryConfig
248  _DefaultName = "aNetBasicAstrometry"
249 
250  def __init__(self,
251  config,
252  andConfig=None,
253  **kwargs):
254  r"""!Construct an ANetBasicAstrometryTask
255 
256  @param[in] config configuration (an instance of self.ConfigClass)
257  @param[in] andConfig astrometry.net data config (an instance of AstromNetDataConfig, or None);
258  if None then use andConfig.py in the astrometry_net_data product (which must be setup)
259  @param[in] kwargs additional keyword arguments for pipe_base Task.\_\_init\_\_
260 
261  @throw RuntimeError if andConfig is None and the configuration cannot be found,
262  either because astrometry_net_data is not setup in eups
263  or because the setup version does not include the file "andConfig.py"
264  """
265  pipeBase.Task.__init__(self, config=config, **kwargs)
266  self.config = config
267  # this is not a subtask because it cannot safely be retargeted
269  config=self.config,
270  andConfig=andConfig,
271  log=self.log,
272  name="loadAN",
273  )
274  self.refObjLoader._readIndexFiles()
275 
276  def memusage(self, prefix=''):
277  # Not logging at DEBUG: do nothing
278  if self.log.getLevel() > self.log.DEBUG:
279  return
280  from astrometry.util.ttime import get_memusage
281  mu = get_memusage()
282  ss = []
283  for key in ['VmPeak', 'VmSize', 'VmRSS', 'VmData']:
284  if key in mu:
285  ss.append(key + ': ' + ' '.join(mu[key]))
286  if 'mmaps' in mu:
287  ss.append('Mmaps: %i' % len(mu['mmaps']))
288  if 'mmaps_total' in mu:
289  ss.append('Mmaps: %i kB' % (mu['mmaps_total'] / 1024))
290  self.log.debug(prefix + 'Memory: ' + ', '.join(ss))
291 
292  def _getImageParams(self, exposure=None, bbox=None, wcs=None, filterName=None, wcsRequired=True):
293  """Get image parameters
294 
295  @param[in] exposure exposure (an afwImage.Exposure) or None
296  @param[in] bbox bounding box (an geom.Box2I) or None; if None then bbox must be specified
297  @param[in] wcs WCS (an afwImage.Wcs) or None; if None then exposure must be specified
298  @param[in] filterName filter name, a string, or None; if None exposure must be specified
299  @param[in] wcsRequired if True then either wcs must be specified or exposure must contain a wcs;
300  if False then the returned wcs may be None
301  @return these items:
302  - bbox bounding box; guaranteed to be set
303  - wcs WCS if known, else None
304  - filterName filter name if known, else None
305  @throw RuntimeError if bbox cannot be determined, or wcs cannot be determined and wcsRequired True
306  """
307  if exposure is not None:
308  if bbox is None:
309  bbox = exposure.getBBox()
310  self.log.debug("Setting bbox = %s from exposure metadata", bbox)
311  if wcs is None:
312  self.log.debug("Setting wcs from exposure metadata")
313  wcs = exposure.getWcs()
314  if filterName is None:
315  filterName = exposure.getFilter().getName()
316  self.log.debug("Setting filterName = %r from exposure metadata", filterName)
317  if bbox is None:
318  raise RuntimeError("bbox or exposure must be specified")
319  if wcs is None and wcsRequired:
320  raise RuntimeError("wcs or exposure (with a WCS) must be specified")
321  return bbox, wcs, filterName
322 
323  def useKnownWcs(self, sourceCat, wcs=None, exposure=None, filterName=None, bbox=None, calculateSip=None):
324  """!Return an InitialAstrometry object, just like determineWcs,
325  but assuming the given input WCS is correct.
326 
327  This involves searching for reference sources within the WCS
328  area, and matching them to the given 'sourceCat'. If
329  'calculateSip' is set, we will try to compute a TAN-SIP
330  distortion correction.
331 
332  @param[in] sourceCat list of detected sources in this image.
333  @param[in] wcs your known WCS, or None to get from exposure
334  @param[in] exposure the exposure holding metadata for this image;
335  if None then you must specify wcs, filterName and bbox
336  @param[in] filterName string, filter name, eg "i", or None to get from exposure`
337  @param[in] bbox bounding box of image, or None to get from exposure
338  @param[in] calculateSip calculate SIP distortion terms for the WCS? If None
339  then use self.config.calculateSip. To disable WCS fitting set calculateSip=False
340 
341  @note this function is also called by 'determineWcs' (via 'determineWcs2'), since the steps are all
342  the same.
343  """
344  # return value:
345  astrom = InitialAstrometry()
346 
347  if calculateSip is None:
348  calculateSip = self.config.calculateSip
349 
350  bbox, wcs, filterName = self._getImageParams(
351  exposure=exposure,
352  bbox=bbox,
353  wcs=wcs,
354  filterName=filterName,
355  wcsRequired=True,
356  )
357  refCat = self.refObjLoader.loadPixelBox(
358  bbox=bbox,
359  wcs=wcs,
360  filterName=filterName,
361  photoCalib=None,
362  ).refCat
363  astrom.refCat = refCat
364  catids = [src.getId() for src in refCat]
365  uids = set(catids)
366  self.log.debug('%i reference sources; %i unique IDs', len(catids), len(uids))
367  matches = self._getMatchList(sourceCat, refCat, wcs)
368  uniq = set([sm.second.getId() for sm in matches])
369  if len(matches) != len(uniq):
370  self.log.warn('The list of matched stars contains duplicate reference source IDs '
371  '(%i sources, %i unique ids)', len(matches), len(uniq))
372  if len(matches) == 0:
373  self.log.warn('No matches found between input sources and reference catalogue.')
374  return astrom
375 
376  self.log.debug('%i reference objects match input sources using input WCS', len(matches))
377  astrom.tanMatches = matches
378  astrom.tanWcs = wcs
379 
380  srcids = [s.getId() for s in sourceCat]
381  for m in matches:
382  assert(m.second.getId() in srcids)
383  assert(m.second in sourceCat)
384 
385  if calculateSip:
386  sipwcs, matches = self._calculateSipTerms(wcs, refCat, sourceCat, matches, bbox=bbox)
387  if sipwcs == wcs:
388  self.log.debug('Failed to find a SIP WCS better than the initial one.')
389  else:
390  self.log.debug('%i reference objects match input sources using SIP WCS',
391  len(matches))
392  astrom.sipWcs = sipwcs
393  astrom.sipMatches = matches
394 
395  wcs = astrom.getWcs()
396  # _getMatchList() modifies the source list RA,Dec coordinates.
397  # Here, we make them consistent with the WCS we are returning.
398  for src in sourceCat:
399  src.updateCoord(wcs)
400  astrom.matchMeta = _createMetadata(bbox, wcs, filterName)
401  return astrom
402 
403  def determineWcs(self, sourceCat, exposure, **kwargs):
404  """Find a WCS solution for the given 'sourceCat' in the given
405  'exposure', getting other parameters from config.
406 
407  Valid kwargs include:
408 
409  'radecCenter', an afw.geom.SpherePoint giving the ICRS RA,Dec position
410  of the center of the field. This is used to limit the
411  search done by Astrometry.net (to make it faster and load
412  fewer index files, thereby using less memory). Defaults to
413  the RA,Dec center from the exposure's WCS; turn that off
414  with the boolean kwarg 'useRaDecCenter' or config option
415  'useWcsRaDecCenter'
416 
417  'useRaDecCenter', a boolean. Don't use the RA,Dec center from
418  the exposure's initial WCS.
419 
420  'searchRadius', in degrees, to search for a solution around
421  the given 'radecCenter'; default from config option
422  'raDecSearchRadius'.
423 
424  'useParity': parity is the 'flip' of the image. Knowing it
425  reduces the search space (hence time) for Astrometry.net.
426  The parity can be computed from the exposure's WCS (the
427  sign of the determinant of the CD matrix); this option
428  controls whether we do that or force Astrometry.net to
429  search both parities. Default from config.useWcsParity.
430 
431  'pixelScale': geom.Angle, estimate of the angle-per-pixel
432  (ie, arcseconds per pixel). Defaults to a value derived
433  from the exposure's WCS. If enabled, this value, plus or
434  minus config.pixelScaleUncertainty, will be used to limit
435  Astrometry.net's search.
436 
437  'usePixelScale': boolean. Use the pixel scale to limit
438  Astrometry.net's search? Defaults to config.useWcsPixelScale.
439 
440  'filterName', a string, the filter name of this image. Will
441  be mapped through the 'filterMap' config dictionary to a
442  column name in the astrometry_net_data index FITS files.
443  Defaults to the exposure.getFilter() value.
444 
445  'bbox', bounding box of exposure; defaults to exposure.getBBox()
446 
447  """
448  assert(exposure is not None)
449 
450  margs = kwargs.copy()
451  if 'searchRadius' not in margs:
452  margs.update(searchRadius=self.config.raDecSearchRadius * geom.degrees)
453  if 'usePixelScale' not in margs:
454  margs.update(usePixelScale=self.config.useWcsPixelScale)
455  if 'useRaDecCenter' not in margs:
456  margs.update(useRaDecCenter=self.config.useWcsRaDecCenter)
457  if 'useParity' not in margs:
458  margs.update(useParity=self.config.useWcsParity)
459  margs.update(exposure=exposure)
460  return self.determineWcs2(sourceCat=sourceCat, **margs)
461 
462  def determineWcs2(self, sourceCat, **kwargs):
463  """Get a blind astrometric solution for the given catalog of sources.
464 
465  We need:
466  -the image size;
467  -the filter
468 
469  And if available, we can use:
470  -an initial Wcs estimate;
471  --> RA,Dec center
472  --> pixel scale
473  --> "parity"
474 
475  (all of which are metadata of Exposure).
476 
477  filterName: string
478  imageSize: (W,H) integer tuple/iterable
479  pixelScale: lsst.geom.Angle per pixel.
480  radecCenter: lsst.afw.coord.Coord
481  """
482  wcs, qa = self.getBlindWcsSolution(sourceCat, **kwargs)
483  kw = {}
484  # Keys passed to useKnownWcs
485  for key in ['exposure', 'bbox', 'filterName']:
486  if key in kwargs:
487  kw[key] = kwargs[key]
488  astrom = self.useKnownWcs(sourceCat, wcs=wcs, **kw)
489  astrom.solveQa = qa
490  astrom.tanWcs = wcs
491  return astrom
492 
493  def getBlindWcsSolution(self, sourceCat,
494  exposure=None,
495  wcs=None,
496  bbox=None,
497  radecCenter=None,
498  searchRadius=None,
499  pixelScale=None,
500  filterName=None,
501  doTrim=False,
502  usePixelScale=True,
503  useRaDecCenter=True,
504  useParity=True,
505  searchRadiusScale=2.):
506  if not useRaDecCenter and radecCenter is not None:
507  raise RuntimeError('radecCenter is set, but useRaDecCenter is False. Make up your mind!')
508  if not usePixelScale and pixelScale is not None:
509  raise RuntimeError('pixelScale is set, but usePixelScale is False. Make up your mind!')
510 
511  bbox, wcs, filterName = self._getImageParams(
512  exposure=exposure,
513  bbox=bbox,
514  wcs=wcs,
515  filterName=filterName,
516  wcsRequired=False,
517  )
518 
519  bboxD = geom.Box2D(bbox)
520  xc, yc = bboxD.getCenter()
521  parity = None
522 
523  if wcs is not None:
524  if pixelScale is None:
525  if usePixelScale:
526  pixelScale = wcs.getPixelScale()
527  self.log.debug('Setting pixel scale estimate = %.3f from given WCS estimate',
528  pixelScale.asArcseconds())
529 
530  if radecCenter is None:
531  if useRaDecCenter:
532  radecCenter = wcs.pixelToSky(xc, yc)
533  self.log.debug('Setting RA,Dec center estimate = (%.3f, %.3f) from given WCS '
534  'estimate, using pixel center = (%.1f, %.1f)',
535  radecCenter.getLongitude().asDegrees(),
536  radecCenter.getLatitude().asDegrees(), xc, yc)
537 
538  if searchRadius is None:
539  if useRaDecCenter:
540  assert(pixelScale is not None)
541  pixRadius = math.hypot(*bboxD.getDimensions()) / 2
542  searchRadius = (pixelScale * pixRadius * searchRadiusScale)
543  self.log.debug('Using RA,Dec search radius = %.3f deg, from pixel scale, '
544  'image size, and searchRadiusScale = %g',
545  searchRadius, searchRadiusScale)
546  if useParity:
547  parity = wcs.isFlipped
548  self.log.debug('Using parity = %s' % (parity and 'True' or 'False'))
549 
550  if doTrim:
551  n = len(sourceCat)
552  if exposure is not None:
553  exposureBBoxD = geom.Box2D(exposure.getMaskedImage().getBBox())
554  else:
555  exposureBBoxD = bboxD
556  sourceCat = self._trimBadPoints(sourceCat, exposureBBoxD)
557  self.log.debug("Trimming: kept %i of %i sources", n, len(sourceCat))
558 
559  wcs, qa = self._solve(
560  sourceCat=sourceCat,
561  wcs=wcs,
562  bbox=bbox,
563  pixelScale=pixelScale,
564  radecCenter=radecCenter,
565  searchRadius=searchRadius,
566  parity=parity,
567  filterName=filterName,
568  )
569  if wcs is None:
570  raise RuntimeError("Unable to match sources with catalog.")
571  self.log.info('Got astrometric solution from Astrometry.net')
572 
573  rdc = wcs.pixelToSky(xc, yc)
574  self.log.debug('New WCS says image center pixel (%.1f, %.1f) -> RA,Dec (%.3f, %.3f)',
575  xc, yc, rdc.getLongitude().asDegrees(), rdc.getLatitude().asDegrees())
576  return wcs, qa
577 
578  def getSipWcsFromWcs(self, wcs, bbox, ngrid=20, linearizeAtCenter=True):
579  """!Get a TAN-SIP WCS, starting from an existing WCS.
580 
581  It uses your WCS to compute a fake grid of corresponding "stars" in pixel and sky coords,
582  and feeds that to the regular SIP code.
583 
584  @param[in] wcs initial WCS
585  @param[in] bbox bounding box of image
586  @param[in] ngrid number of grid points along x and y for fitting (fit at ngrid^2 points)
587  @param[in] linearizeAtCenter if True, get a linear approximation of the input
588  WCS at the image center and use that as the TAN initialization for
589  the TAN-SIP solution. You probably want this if your WCS has its
590  CRPIX outside the image bounding box.
591  """
592  # Ugh, build src and ref tables
593  srcSchema = afwTable.SourceTable.makeMinimalSchema()
594  key = srcSchema.addField("centroid", type="PointD")
595  srcTable = afwTable.SourceTable.make(srcSchema)
596  srcTable.defineCentroid("centroid")
597  srcs = srcTable
598  refs = afwTable.SimpleTable.make(afwTable.SimpleTable.makeMinimalSchema())
599  cref = []
600  csrc = []
601  (W, H) = bbox.getDimensions()
602  x0, y0 = bbox.getMin()
603  for xx in np.linspace(0., W, ngrid):
604  for yy in np.linspace(0, H, ngrid):
605  src = srcs.makeRecord()
606  src.set(key.getX(), x0 + xx)
607  src.set(key.getY(), y0 + yy)
608  csrc.append(src)
609  rd = wcs.pixelToSky(xx + x0, yy + y0)
610  ref = refs.makeRecord()
611  ref.setCoord(rd)
612  cref.append(ref)
613 
614  if linearizeAtCenter:
615  # Linearize the original WCS around the image center to create a
616  # TAN WCS.
617  # Reference pixel in LSST coords
618  crpix = geom.Box2D(bbox).getCenter()
619  crval = wcs.pixelToSky(crpix)
620  crval = crval.getPosition(geom.degrees)
621  # Linearize *AT* crval to get effective CD at crval.
622  # (we use the default skyUnit of degrees as per WCS standard)
623  aff = wcs.linearizePixelToSky(crval)
624  cd = aff.getLinear().getMatrix()
625  wcs = afwImage.Wcs(crval, crpix, cd)
626 
627  return self.getSipWcsFromCorrespondences(wcs, cref, csrc, (W, H), x0=x0, y0=y0)
628 
629  def getSipWcsFromCorrespondences(self, origWcs, refCat, sourceCat, bbox):
630  """Produce a SIP solution given a list of known correspondences.
631 
632  Unlike _calculateSipTerms, this does not iterate the solution;
633  it assumes you have given it a good sets of corresponding stars.
634 
635  NOTE that "refCat" and "sourceCat" are assumed to be the same length;
636  entries "refCat[i]" and "sourceCat[i]" are assumed to be correspondences.
637 
638  @param[in] origWcs the WCS to linearize in order to get the TAN part of the TAN-SIP WCS.
639  @param[in] refCat reference source catalog
640  @param[in] sourceCat source catalog
641  @param[in] bbox bounding box of image
642  """
643  sipOrder = self.config.sipOrder
644  matches = []
645  for ci, si in zip(refCat, sourceCat):
646  matches.append(afwTable.ReferenceMatch(ci, si, 0.))
647 
648  sipObject = astromSip.makeCreateWcsWithSip(matches, origWcs, sipOrder, bbox)
649  return sipObject.getNewWcs()
650 
651  def _calculateSipTerms(self, origWcs, refCat, sourceCat, matches, bbox):
652  """!Iteratively calculate SIP distortions and regenerate matches based on improved WCS.
653 
654  @param[in] origWcs original WCS object, probably (but not necessarily) a TAN WCS;
655  this is used to set the baseline when determining whether a SIP
656  solution is any better; it will be returned if no better SIP solution
657  can be found.
658  @param[in] refCat reference source catalog
659  @param[in] sourceCat sources in the image to be solved
660  @param[in] matches list of supposedly matched sources, using the "origWcs".
661  @param[in] bbox bounding box of image, which is used when finding reverse SIP coefficients.
662  """
663  sipOrder = self.config.sipOrder
664  wcs = origWcs
665 
666  lastMatchSize = len(matches)
667  lastMatchStats = self._computeMatchStatsOnSky(wcs=wcs, matchList=matches)
668  for i in range(self.config.maxIter):
669  # fit SIP terms
670  try:
671  sipObject = astromSip.makeCreateWcsWithSip(matches, wcs, sipOrder, bbox)
672  proposedWcs = sipObject.getNewWcs()
673  self.plotSolution(matches, proposedWcs, bbox.getDimensions())
674  except pexExceptions.Exception as e:
675  self.log.warn('Failed to calculate distortion terms. Error: ', str(e))
676  break
677 
678  # update the source catalog
679  for source in sourceCat:
680  skyPos = proposedWcs.pixelToSky(source.getCentroid())
681  source.setCoord(skyPos)
682 
683  # use new WCS to get new matchlist.
684  proposedMatchlist = self._getMatchList(sourceCat, refCat, proposedWcs)
685  proposedMatchSize = len(proposedMatchlist)
686  proposedMatchStats = self._computeMatchStatsOnSky(wcs=proposedWcs, matchList=proposedMatchlist)
687 
688  self.log.debug(
689  "SIP iteration %i: %i objects match, previous = %i;" %
690  (i, proposedMatchSize, lastMatchSize) +
691  " clipped mean scatter = %s arcsec, previous = %s; " %
692  (proposedMatchStats.distMean.asArcseconds(), lastMatchStats.distMean.asArcseconds()) +
693  " max match dist = %s arcsec, previous = %s" %
694  (proposedMatchStats.maxMatchDist.asArcseconds(),
695  lastMatchStats.maxMatchDist.asArcseconds())
696  )
697 
698  if lastMatchStats.maxMatchDist <= proposedMatchStats.maxMatchDist:
699  self.log.debug(
700  "Fit WCS: use iter %s because max match distance no better in next iter: " % (i-1,) +
701  " %g < %g arcsec" % (lastMatchStats.maxMatchDist.asArcseconds(),
702  proposedMatchStats.maxMatchDist.asArcseconds()))
703  break
704 
705  wcs = proposedWcs
706  matches = proposedMatchlist
707  lastMatchSize = proposedMatchSize
708  lastMatchStats = proposedMatchStats
709 
710  return wcs, matches
711 
712  def plotSolution(self, matches, wcs, imageSize):
713  """Plot the solution, when debugging is turned on.
714 
715  @param matches The list of matches
716  @param wcs The Wcs
717  @param imageSize 2-tuple with the image size (W,H)
718  """
719  import lsstDebug
720  display = lsstDebug.Info(__name__).display
721  if not display:
722  return
723 
724  try:
725  import matplotlib.pyplot as plt
726  except ImportError as e:
727  self.log.warn("Unable to import matplotlib: %s", e)
728  return
729 
730  fig = plt.figure(1)
731  fig.clf()
732  try:
733  fig.canvas._tkcanvas._root().lift() # == Tk's raise, but raise is a python reserved word
734  except Exception: # protect against API changes
735  pass
736 
737  num = len(matches)
738  x = np.zeros(num)
739  y = np.zeros(num)
740  dx = np.zeros(num)
741  dy = np.zeros(num)
742  for i, m in enumerate(matches):
743  x[i] = m.second.getX()
744  y[i] = m.second.getY()
745  pixel = wcs.skyToPixel(m.first.getCoord())
746  dx[i] = x[i] - pixel.getX()
747  dy[i] = y[i] - pixel.getY()
748 
749  subplots = maUtils.makeSubplots(fig, 2, 2, xgutter=0.1, ygutter=0.1, pygutter=0.04)
750 
751  def plotNext(x, y, xLabel, yLabel, xMax):
752  ax = next(subplots)
753  ax.set_autoscalex_on(False)
754  ax.set_xbound(lower=0, upper=xMax)
755  ax.scatter(x, y)
756  ax.set_xlabel(xLabel)
757  ax.set_ylabel(yLabel)
758  ax.axhline(0.0)
759 
760  plotNext(x, dx, "x", "dx", imageSize[0])
761  plotNext(x, dy, "x", "dy", imageSize[0])
762  plotNext(y, dx, "y", "dx", imageSize[1])
763  plotNext(y, dy, "y", "dy", imageSize[1])
764 
765  fig.show()
766 
767  while True:
768  try:
769  reply = input("Pausing for inspection, enter to continue... [hpQ] ").strip()
770  except EOFError:
771  reply = "n"
772 
773  reply = reply.split()
774  if reply:
775  reply = reply[0]
776  else:
777  reply = ""
778 
779  if reply in ("", "h", "p", "Q"):
780  if reply == "h":
781  print("h[elp] p[db] Q[uit]")
782  continue
783  elif reply == "p":
784  import pdb
785  pdb.set_trace()
786  elif reply == "Q":
787  sys.exit(1)
788  break
789 
790  def _computeMatchStatsOnSky(self, wcs, matchList):
791  """Compute on-sky radial distance statistics for a match list
792 
793  @param[in] wcs WCS for match list; an lsst.afw.image.Wcs
794  @param[in] matchList list of matches between reference object and sources;
795  a list of lsst.afw.table.ReferenceMatch;
796  the source centroid and reference object coord are read
797 
798  @return a pipe_base Struct containing these fields:
799  - distMean clipped mean of on-sky radial separation
800  - distStdDev clipped standard deviation of on-sky radial separation
801  - maxMatchDist distMean + self.config.matchDistanceSigma*distStdDev
802  """
803  distStatsInRadians = makeMatchStatisticsInRadians(wcs, matchList,
804  afwMath.MEANCLIP | afwMath.STDEVCLIP)
805  distMean = distStatsInRadians.getValue(afwMath.MEANCLIP)*geom.radians
806  distStdDev = distStatsInRadians.getValue(afwMath.STDEVCLIP)*geom.radians
807  return pipeBase.Struct(
808  distMean=distMean,
809  distStdDev=distStdDev,
810  maxMatchDist=distMean + self.config.matchDistanceSigma*distStdDev,
811  )
812 
813  def _getMatchList(self, sourceCat, refCat, wcs):
814  dist = self.config.catalogMatchDist * geom.arcseconds
815  clean = self.config.cleaningParameter
816  matcher = astromSip.MatchSrcToCatalogue(refCat, sourceCat, wcs, dist)
817  matches = matcher.getMatches()
818  if matches is None:
819  # Produce debugging stats...
820  X = [src.getX() for src in sourceCat]
821  Y = [src.getY() for src in sourceCat]
822  R1 = [src.getRa().asDegrees() for src in sourceCat]
823  D1 = [src.getDec().asDegrees() for src in sourceCat]
824  R2 = [src.getRa().asDegrees() for src in refCat]
825  D2 = [src.getDec().asDegrees() for src in refCat]
826  # for src in sourceCat:
827  # self.log.debug("source: x,y (%.1f, %.1f), RA,Dec (%.3f, %.3f)" %
828  # (src.getX(), src.getY(), src.getRa().asDegrees(), src.getDec().asDegrees()))
829  # for src in refCat:
830  # self.log.debug("ref: RA,Dec (%.3f, %.3f)" %
831  # (src.getRa().asDegrees(), src.getDec().asDegrees()))
832  self.loginfo('_getMatchList: %i sources, %i reference sources' % (len(sourceCat), len(refCat)))
833  if len(sourceCat):
834  self.loginfo(
835  'Source range: x [%.1f, %.1f], y [%.1f, %.1f], RA [%.3f, %.3f], Dec [%.3f, %.3f]' %
836  (min(X), max(X), min(Y), max(Y), min(R1), max(R1), min(D1), max(D1)))
837  if len(refCat):
838  self.loginfo('Reference range: RA [%.3f, %.3f], Dec [%.3f, %.3f]' %
839  (min(R2), max(R2), min(D2), max(D2)))
840  raise RuntimeError('No matches found between image and catalogue')
841  matches = cleanBadPoints.clean(matches, wcs, nsigma=clean)
842  return matches
843 
844  def getColumnName(self, filterName, columnMap, default=None):
845  """
846  Returns the column name in the astrometry_net_data index file that will be used
847  for the given filter name.
848 
849  @param filterName Name of filter used in exposure
850  @param columnMap Dict that maps filter names to column names
851  @param default Default column name
852  """
853  filterName = self.config.filterMap.get(filterName, filterName) # Exposure filter --> desired filter
854  try:
855  return columnMap[filterName] # Desired filter --> a_n_d column name
856  except KeyError:
857  self.log.warn("No column in configuration for filter '%s'; using default '%s'" %
858  (filterName, default))
859  return default
860 
861  def _solve(self, sourceCat, wcs, bbox, pixelScale, radecCenter, searchRadius, parity, filterName=None):
862  """
863  @param[in] parity True for flipped parity, False for normal parity, None to leave parity unchanged
864  """
865  solver = self.refObjLoader._getSolver()
866 
867  imageSize = bbox.getDimensions()
868  x0, y0 = bbox.getMin()
869 
870  # select sources with valid x, y, flux
871  xybb = geom.Box2D()
872  goodsources = afwTable.SourceCatalog(sourceCat.table)
873  badkeys = [goodsources.getSchema().find(name).key for name in self.config.badFlags]
874 
875  for s in sourceCat:
876  if np.isfinite(s.getX()) and np.isfinite(s.getY()) and np.isfinite(s.getPsfInstFlux()) \
877  and self._isGoodSource(s, badkeys):
878  goodsources.append(s)
879  xybb.include(geom.Point2D(s.getX() - x0, s.getY() - y0))
880  self.log.info("Number of selected sources for astrometry : %d" % (len(goodsources)))
881  if len(goodsources) < len(sourceCat):
882  self.log.debug('Keeping %i of %i sources with finite X,Y positions and PSF flux',
883  len(goodsources), len(sourceCat))
884  self.log.debug('Feeding sources in range x=[%.1f, %.1f], y=[%.1f, %.1f] ' +
885  '(after subtracting x0,y0 = %.1f,%.1f) to Astrometry.net',
886  xybb.getMinX(), xybb.getMaxX(), xybb.getMinY(), xybb.getMaxY(), x0, y0)
887  # setStars sorts them by PSF flux.
888  solver.setStars(goodsources, x0, y0)
889  solver.setMaxStars(self.config.maxStars)
890  solver.setImageSize(*imageSize)
891  solver.setMatchThreshold(self.config.matchThreshold)
892  raDecRadius = None
893  if radecCenter is not None:
894  raDecRadius = (radecCenter.getLongitude().asDegrees(), radecCenter.getLatitude().asDegrees(),
895  searchRadius.asDegrees())
896  solver.setRaDecRadius(*raDecRadius)
897  self.log.debug('Searching for match around RA,Dec = (%g, %g) with radius %g deg' %
898  raDecRadius)
899 
900  if pixelScale is not None:
901  dscale = self.config.pixelScaleUncertainty
902  scale = pixelScale.asArcseconds()
903  lo = scale / dscale
904  hi = scale * dscale
905  solver.setPixelScaleRange(lo, hi)
906  self.log.debug(
907  'Searching for matches with pixel scale = %g +- %g %% -> range [%g, %g] arcsec/pix',
908  scale, 100.*(dscale-1.), lo, hi)
909 
910  if parity is not None:
911  solver.setParity(parity)
912  self.log.debug('Searching for match with parity = %s', str(parity))
913 
914  # Find and load index files within RA,Dec range and scale range.
915  if radecCenter is not None:
916  multiInds = self.refObjLoader._getMIndexesWithinRange(radecCenter, searchRadius)
917  else:
918  multiInds = self.refObjLoader.multiInds
919  qlo, qhi = solver.getQuadSizeRangeArcsec()
920 
921  toload_multiInds = set()
922  toload_inds = []
923  for mi in multiInds:
924  for i in range(len(mi)):
925  ind = mi[i]
926  if not ind.overlapsScaleRange(qlo, qhi):
927  continue
928  toload_multiInds.add(mi)
929  toload_inds.append(ind)
930 
931  import lsstDebug
932  if lsstDebug.Info(__name__).display:
933  # Use separate context for display, since astrometry.net can segfault if we don't...
934  with LoadMultiIndexes(toload_multiInds):
935  displayAstrometry(refCat=self.refObjLoader.loadPixelBox(bbox, wcs, filterName).refCat,
936  frame=lsstDebug.Info(__name__).frame, pause=lsstDebug.Info(__name__).pause)
937 
938  with LoadMultiIndexes(toload_multiInds):
939  solver.addIndices(toload_inds)
940  self.memusage('Index files loaded: ')
941 
942  cpulimit = self.config.maxCpuTime
943  solver.run(cpulimit)
944 
945  self.memusage('Solving finished: ')
946 
947  self.memusage('Index files unloaded: ')
948 
949  if solver.didSolve():
950  self.log.debug('Solved!')
951  wcs = solver.getWcs()
952 
953  if x0 != 0 or y0 != 0:
954  wcs = wcs.copyAtShiftedPixelOrigin(geom.Extent2D(x0, y0))
955 
956  else:
957  self.log.warn('Did not get an astrometric solution from Astrometry.net')
958  wcs = None
959  # Gather debugging info...
960 
961  # -are there any reference stars in the proposed search area?
962  # log the number found and discard the results
963  if radecCenter is not None:
964  self.refObjLoader.loadSkyCircle(radecCenter, searchRadius, filterName)
965 
966  qa = solver.getSolveStats()
967  self.log.debug('qa: %s', qa.toString())
968  return wcs, qa
969 
970  def _isGoodSource(self, candsource, keys):
971  for k in keys:
972  if candsource.get(k):
973  return False
974  return True
975 
976  @staticmethod
977  def _trimBadPoints(sourceCat, bbox, wcs=None):
978  """Remove elements from catalog whose xy positions are not within the given bbox.
979 
980  sourceCat: a Catalog of SimpleRecord or SourceRecord objects
981  bbox: an afwImage.Box2D
982  wcs: if not None, will be used to compute the xy positions on-the-fly;
983  this is required when sources actually contains SimpleRecords.
984 
985  Returns:
986  a list of Source objects with xAstrom, yAstrom within the bbox.
987  """
988  keep = type(sourceCat)(sourceCat.table)
989  for s in sourceCat:
990  point = s.getCentroid() if wcs is None else wcs.skyToPixel(s.getCoord())
991  if bbox.contains(point):
992  keep.append(s)
993  return keep
994 
995 
996 def _createMetadata(bbox, wcs, filterName):
997  """
998  Create match metadata entries required for regenerating the catalog
999 
1000  @param bbox bounding box of image (pixels)
1001  @param filterName Name of filter, used for magnitudes
1002  @return Metadata
1003  """
1004  meta = dafBase.PropertyList()
1005 
1006  bboxD = geom.Box2D(bbox)
1007  cx, cy = bboxD.getCenter()
1008  radec = wcs.pixelToSky(cx, cy)
1009  meta.add('RA', radec.getRa().asDegrees(), 'field center in degrees')
1010  meta.add('DEC', radec.getDec().asDegrees(), 'field center in degrees')
1011  pixelRadius = math.hypot(*bboxD.getDimensions())/2.0
1012  skyRadius = wcs.getPixelScale() * pixelRadius
1013  meta.add('RADIUS', skyRadius.asDegrees(),
1014  'field radius in degrees, approximate')
1015  meta.add('SMATCHV', 1, 'SourceMatchVector version number')
1016  if filterName is not None:
1017  meta.add('FILTER', str(filterName), 'LSST filter name for tagalong data')
1018  return meta
def getSipWcsFromWcs(self, wcs, bbox, ngrid=20, linearizeAtCenter=True)
Get a TAN-SIP WCS, starting from an existing WCS.
A floating-point coordinate rectangle geometry.
Definition: Box.h:305
Class for storing ordered metadata with comments.
Definition: PropertyList.h:68
def useKnownWcs(self, sourceCat, wcs=None, exposure=None, filterName=None, bbox=None, calculateSip=None)
Return an InitialAstrometry object, just like determineWcs, but assuming the given input WCS is corre...
Provides consistent interface for LSST exceptions.
Definition: Exception.h:107
daf::base::PropertySet * set
Definition: fits.cc:884
int min
def getBlindWcsSolution(self, sourceCat, exposure=None, wcs=None, bbox=None, radecCenter=None, searchRadius=None, pixelScale=None, filterName=None, doTrim=False, usePixelScale=True, useRaDecCenter=True, useParity=True, searchRadiusScale=2.)
def _solve(self, sourceCat, wcs, bbox, pixelScale, radecCenter, searchRadius, parity, filterName=None)
Basic implemeentation of the astrometry.net astrometrical fitter.
def __init__(self, config, andConfig=None, kwargs)
Construct an ANetBasicAstrometryTask.
table::Key< int > type
Definition: Detector.cc:167
Lightweight representation of a geometric match between two records.
Definition: fwd.h:101
def _getImageParams(self, exposure=None, bbox=None, wcs=None, filterName=None, wcsRequired=True)
def getLevel(loggername)
int max
def displayAstrometry(refCat=None, sourceCat=None, distortedCentroidKey=None, bbox=None, exposure=None, matches=None, frame=1, title="", pause=True)
Definition: display.py:34
afw::math::Statistics makeMatchStatisticsInRadians(afw::geom::SkyWcs const &wcs, std::vector< MatchT > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl=afw::math::StatisticsControl())
Compute statistics of on-sky radial separation for a match list, in radians.
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects...
bool strip
Definition: fits.cc:883
def _calculateSipTerms(self, origWcs, refCat, sourceCat, matches, bbox)
Iteratively calculate SIP distortions and regenerate matches based on improved WCS.