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