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