lsst.meas.astrom  16.0-18-g51a54b3+1
matchPessimisticB.py
Go to the documentation of this file.
1 
2 import numpy as np
3 from scipy.spatial import cKDTree
4 from scipy.stats import sigmaclip
5 
6 import lsst.pex.config as pexConfig
7 import lsst.pipe.base as pipeBase
8 import lsst.afw.geom as afwgeom
9 import lsst.afw.table as afwTable
10 from lsst.meas.algorithms.sourceSelector import sourceSelectorRegistry
11 
12 from .matchOptimisticB import MatchTolerance
13 
14 from .pessimistic_pattern_matcher_b_3D import PessimisticPatternMatcherB
15 
16 __all__ = ["MatchPessimisticBTask", "MatchPessimisticBConfig",
17  "MatchTolerancePessimistic"]
18 
19 
21  """Stores match tolerances for use in AstrometryTask and later
22  iterations of the matcher.
23 
24  MatchPessimisticBTask relies on several state variables to be
25  preserved over different iterations in the
26  AstrometryTask.matchAndFitWcs loop of AstrometryTask.
27 
28  Parameters
29  ----------
30  maxMatchDist : `lsst.geom.Angle`
31  Maximum distance to consider a match from the previous match/fit
32  iteration.
33  autoMaxMatchDist : `lsst.geom.Angle`
34  Automated estimation of the maxMatchDist from the sky statistics of the
35  source and reference catalogs.
36  maxShift : `lsst.geom.Angle`
37  Maximum shift found in the previous match/fit cycle.
38  lastMatchedPattern : `int`
39  Index of the last source pattern that was matched into the reference
40  data.
41  failedPatternList : `list` of `int`
42  Previous matches were found to be false positives.
43  """
44 
45  def __init__(self, maxMatchDist=None, autoMaxMatchDist=None,
46  maxShift=None, lastMatchedPattern=None,
47  failedPatternList=None):
48  self.maxMatchDist = maxMatchDist
49  self.autoMaxMatchDist = autoMaxMatchDist
50  self.maxShift = maxShift
51  self.lastMatchedPattern = lastMatchedPattern
52  if failedPatternList is None:
54 
55 
56 class MatchPessimisticBConfig(pexConfig.Config):
57  """Configuration for MatchPessimisticBTask
58  """
59  numBrightStars = pexConfig.RangeField(
60  doc="Number of bright stars to use. Sets the max number of patterns "
61  "that can be tested.",
62  dtype=int,
63  default=200,
64  min=2,
65  )
66  minMatchedPairs = pexConfig.RangeField(
67  doc="Minimum number of matched pairs; see also minFracMatchedPairs.",
68  dtype=int,
69  default=30,
70  min=2,
71  )
72  minFracMatchedPairs = pexConfig.RangeField(
73  doc="Minimum number of matched pairs as a fraction of the smaller of "
74  "the number of reference stars or the number of good sources; "
75  "the actual minimum is the smaller of this value or "
76  "minMatchedPairs.",
77  dtype=float,
78  default=0.3,
79  min=0,
80  max=1,
81  )
82  matcherIterations = pexConfig.RangeField(
83  doc="Number of softening iterations in matcher.",
84  dtype=int,
85  default=5,
86  min=1,
87  )
88  maxOffsetPix = pexConfig.RangeField(
89  doc="Maximum allowed shift of WCS, due to matching (pixel). "
90  "When changing this value, the "
91  "LoadReferenceObjectsConfig.pixelMargin should also be updated.",
92  dtype=int,
93  default=300,
94  max=4000,
95  )
96  maxRotationDeg = pexConfig.RangeField(
97  doc="Rotation angle allowed between sources and position reference "
98  "objects (degrees).",
99  dtype=float,
100  default=1.0,
101  max=6.0,
102  )
103  numPointsForShape = pexConfig.Field(
104  doc="Number of points to define a shape for matching.",
105  dtype=int,
106  default=6,
107  )
108  numPointsForShapeAttempt = pexConfig.Field(
109  doc="Number of points to try for creating a shape. This value should "
110  "be greater than or equal to numPointsForShape. Besides "
111  "loosening the signal to noise cut in the matcherSourceSelector, "
112  "increasing this number will solve CCDs where no match was found.",
113  dtype=int,
114  default=6,
115  )
116  minMatchDistPixels = pexConfig.RangeField(
117  doc="Distance in units of pixels to always consider a source-"
118  "reference pair a match. This prevents the astrometric fitter "
119  "from over-fitting and removing stars that should be matched and "
120  "allows for inclusion of new matches as the wcs improves.",
121  dtype=float,
122  default=1.0,
123  min=0.0,
124  max=6.0,
125  )
126  numPatternConsensus = pexConfig.Field(
127  doc="Number of implied shift/rotations from patterns that must agree "
128  "before it a given shift/rotation is accepted. This is only used "
129  "after the first softening iteration fails and if both the "
130  "number of reference and source objects is greater than "
131  "numBrightStars.",
132  dtype=int,
133  default=3,
134  )
135  sourceSelector = sourceSelectorRegistry.makeField(
136  doc="How to select sources for cross-matching. The default "
137  "matcherSourceSelector removes objects with low S/N, bad "
138  "saturated objects, edge objects, and interpolated objects.",
139  default="matcherPessimistic"
140  )
141 
142  def setDefaults(self):
143  sourceSelector = self.sourceSelector["matcherPessimistic"]
144  sourceSelector.setDefaults()
145 
146  def validate(self):
147  pexConfig.Config.validate(self)
149  raise ValueError("numPointsForShapeAttempt must be greater than "
150  "or equal to numPointsForShape.")
151 
152 
153 # The following block adds links to this task from the Task Documentation page.
154 # \addtogroup LSST_task_documentation
155 # \{
156 # \page measAstrom_MatchPessimisticBTask
157 # \ref MatchPessimisticBTask "MatchPessimisticBTask"
158 # Match sources to reference objects
159 # \}
160 
161 
162 class MatchPessimisticBTask(pipeBase.Task):
163  """Match sources to reference objects.
164  """
165 
166  ConfigClass = MatchPessimisticBConfig
167  _DefaultName = "matchObjectsToSources"
168 
169  def __init__(self, **kwargs):
170  pipeBase.Task.__init__(self, **kwargs)
171  self.makeSubtask("sourceSelector")
172 
173  @pipeBase.timeMethod
174  def matchObjectsToSources(self, refCat, sourceCat, wcs, refFluxField,
175  match_tolerance=None):
176  """Match sources to position reference stars
177 
178  refCat : `lsst.afw.table.SimpleCatalog`
179  catalog of reference objects that overlap the exposure; reads
180  fields for:
181 
182  - coord
183  - the specified flux field
184 
185  sourceCat : `lsst.afw.table.SourceCatalog`
186  catalog of sources found on an exposure; Please check the required
187  fields of your specified source selector that the correct flags are present.
188  wcs : `lsst.afw.geom.SkyWcs`
189  estimated WCS
190  refFluxField : `str`
191  field of refCat to use for flux
192  match_tolerance : `lsst.meas.astrom.MatchTolerancePessimistic`
193  is a MatchTolerance class object or `None`. This this class is used
194  to communicate state between AstrometryTask and MatcherTask.
195  AstrometryTask will also set the MatchTolerance class variable
196  maxMatchDist based on the scatter AstrometryTask has found after
197  fitting for the wcs.
198 
199  Returns
200  -------
201  result : `lsst.pipe.base.Struct`
202  Result struct with components:
203 
204  - ``matches`` : source to reference matches found (`list` of
205  `lsst.afw.table.ReferenceMatch`)
206  - ``usableSourcCat`` : a catalog of sources potentially usable for
207  matching and WCS fitting (`lsst.afw.table.SourceCatalog`).
208  - ``match_tolerance`` : a MatchTolerance object containing the
209  resulting state variables from the match
210  (`lsst.meas.astrom.MatchTolerancePessimistic`).
211  """
212  import lsstDebug
213  debug = lsstDebug.Info(__name__)
214 
215  # If we get an empty tolerance struct create the variables we need for
216  # this matcher.
217  if match_tolerance is None:
218  match_tolerance = MatchTolerancePessimistic()
219 
220  # usableSourceCat: sources that are good but may be saturated
221  numSources = len(sourceCat)
222  selectedSources = self.sourceSelector.run(sourceCat)
223  goodSourceCat = selectedSources.sourceCat
224  numUsableSources = len(goodSourceCat)
225  self.log.info("Purged %d sources, leaving %d good sources" %
226  (numSources - numUsableSources, numUsableSources))
227 
228  if len(goodSourceCat) == 0:
229  raise pipeBase.TaskError("No sources are good")
230 
231  # avoid accidentally using sourceCat; use goodSourceCat from now on
232  del sourceCat
233 
234  minMatchedPairs = min(self.config.minMatchedPairs,
235  int(self.config.minFracMatchedPairs *
236  min([len(refCat), len(goodSourceCat)])))
237 
238  doMatchReturn = self._doMatch(
239  refCat=refCat,
240  sourceCat=goodSourceCat,
241  wcs=wcs,
242  refFluxField=refFluxField,
243  numUsableSources=numUsableSources,
244  minMatchedPairs=minMatchedPairs,
245  match_tolerance=match_tolerance,
246  sourceFluxField=self.sourceSelector.fluxField,
247  verbose=debug.verbose,
248  )
249  matches = doMatchReturn.matches
250  match_tolerance = doMatchReturn.match_tolerance
251 
252  if len(matches) == 0:
253  raise RuntimeError("Unable to match sources")
254 
255  self.log.info("Matched %d sources" % len(matches))
256  if len(matches) < minMatchedPairs:
257  self.log.warn("Number of matches is smaller than request")
258 
259  return pipeBase.Struct(
260  matches=matches,
261  usableSourceCat=goodSourceCat,
262  match_tolerance=match_tolerance,
263  )
264 
265  @pipeBase.timeMethod
266  def _doMatch(self, refCat, sourceCat, wcs, refFluxField, numUsableSources,
267  minMatchedPairs, match_tolerance, sourceFluxField, verbose):
268  """Implementation of matching sources to position reference stars
269 
270  Unlike matchObjectsToSources, this method does not check if the sources
271  are suitable.
272 
273  Parameters
274  ----------
275  refCat : `lsst.afw.table.SimpleCatalog`
276  catalog of position reference stars that overlap an exposure
277  sourceCat : `lsst.afw.table.SourceCatalog`
278  catalog of sources found on the exposure
279  wcs : `lsst.afw.geom.SkyWcs`
280  estimated WCS of exposure
281  refFluxField : `str`
282  field of refCat to use for flux
283  numUsableSources : `int`
284  number of usable sources (sources with known centroid that are not
285  near the edge, but may be saturated)
286  minMatchedPairs : `int`
287  minimum number of matches
288  match_tolerance : `lsst.meas.astrom.MatchTolerancePessimistic`
289  a MatchTolerance object containing variables specifying matcher
290  tolerances and state from possible previous runs.
291  sourceFluxField : `str`
292  Name of the flux field in the source catalog.
293  verbose : `bool`
294  Set true to print diagnostic information to std::cout
295 
296  Returns
297  -------
298  result :
299  Results struct with components:
300 
301  - ``matches`` : a list the matches found
302  (`list` of `lsst.afw.table.ReferenceMatch`).
303  - ``match_tolerance`` : MatchTolerance containing updated values from
304  this fit iteration (`lsst.meas.astrom.MatchTolerancePessimistic`)
305  """
306 
307  # Load the source and reference catalog as spherical points
308  # in numpy array. We do this rather than relying on internal
309  # lsst C objects for simplicity and because we require
310  # objects contiguous in memory. We need to do these slightly
311  # differently for the reference and source cats as they are
312  # different catalog objects with different fields.
313  ref_array = np.empty((len(refCat), 4), dtype=np.float64)
314  for ref_idx, refObj in enumerate(refCat):
315  theta = np.pi / 2 - refObj.getDec().asRadians()
316  phi = refObj.getRa().asRadians()
317  flux = refObj[refFluxField]
318  ref_array[ref_idx, :] = \
319  self._latlong_flux_to_xyz_mag(theta, phi, flux)
320 
321  src_array = np.empty((len(sourceCat), 4), dtype=np.float64)
322  for src_idx, srcObj in enumerate(sourceCat):
323  coord = wcs.pixelToSky(srcObj.getCentroid())
324  theta = np.pi / 2 - coord.getLatitude().asRadians()
325  phi = coord.getLongitude().asRadians()
326  flux = srcObj.getPsfInstFlux()
327  src_array[src_idx, :] = \
328  self._latlong_flux_to_xyz_mag(theta, phi, flux)
329 
330  # Set configurable defaults when we encounter None type or set
331  # state based on previous run of AstrometryTask._matchAndFitWcs.
332  if match_tolerance.maxShift is None:
333  maxShiftArcseconds = (self.config.maxOffsetPix *
334  wcs.getPixelScale().asArcseconds())
335  else:
336  # We don't want to clamp down too hard on the allowed shift so
337  # we test that the smallest we ever allow is the pixel scale.
338  maxShiftArcseconds = np.max(
339  (match_tolerance.maxShift.asArcseconds(),
340  self.config.minMatchDistPixels *
341  wcs.getPixelScale().asArcseconds()))
342 
343  # If our tolerances are not set from a previous run, estimate a
344  # starting tolerance guess from the statistics of patterns we can
345  # create on both the source and reference catalog. We use the smaller
346  # of the two.
347  if match_tolerance.maxMatchDist is None:
348  self.log.debug("Computing source statistics...")
349  maxMatchDistArcSecSrc = self._get_pair_pattern_statistics(
350  src_array)
351  self.log.debug("Computing reference statistics...")
352  maxMatchDistArcSecRef = self._get_pair_pattern_statistics(
353  ref_array)
354  maxMatchDistArcSec = np.min((maxMatchDistArcSecSrc,
355  maxMatchDistArcSecRef))
356  match_tolerance.autoMaxDist = afwgeom.Angle(maxMatchDistArcSec,
357  afwgeom.arcseconds)
358  else:
359  maxMatchDistArcSec = np.max(
360  (self.config.minMatchDistPixels *
361  wcs.getPixelScale().asArcseconds(),
362  np.min((match_tolerance.maxMatchDist.asArcseconds(),
363  match_tolerance.autoMaxDist.asArcseconds()))))
364 
365  # Make sure the data we are considering is dense enough to require
366  # the consensus mode of the matcher. If not default to Optimistic
367  # pattern matcher behavior.
368  numConsensus = self.config.numPatternConsensus
369  minObjectsForConsensus = \
370  self.config.numBrightStars + self.config.numPointsForShapeAttempt
371  if ref_array.shape[0] < minObjectsForConsensus or \
372  src_array.shape[0] < minObjectsForConsensus:
373  numConsensus = 1
374 
375  self.log.debug("Current tol maxDist: %.4f arcsec" %
376  maxMatchDistArcSec)
377  self.log.debug("Current shift: %.4f arcsec" %
378  maxShiftArcseconds)
379 
380  # Create our matcher object.
381  pyPPMb = PessimisticPatternMatcherB(ref_array[:, :3], self.log)
382 
383  match_found = False
384  # Start the iteration over our tolerances.
385  for soften_pattern in range(self.config.matcherIterations):
386  for soften_dist in range(self.config.matcherIterations):
387  if soften_pattern == 0 and soften_dist == 0 and \
388  match_tolerance.lastMatchedPattern is not None:
389  # If we are on the first, most stringent tolerance,
390  # and have already found a match, the matcher should behave
391  # like an optimistic pattern matcher. Exiting at the first
392  # match.
393  run_n_consent = 1
394  else:
395  # If we fail or this is the first match attempt, set the
396  # pattern consensus to the specified config value.
397  run_n_consent = numConsensus
398  # We double the match dist tolerance each round and add 1 to the
399  # to the number of candidate spokes to check.
400  matcher_struct = pyPPMb.match(
401  source_array=src_array,
402  n_check=self.config.numPointsForShapeAttempt + soften_pattern,
403  n_match=self.config.numPointsForShape,
404  n_agree=run_n_consent,
405  max_n_patterns=self.config.numBrightStars,
406  max_shift=maxShiftArcseconds,
407  max_rotation=self.config.maxRotationDeg,
408  max_dist=maxMatchDistArcSec * 2. ** soften_dist,
409  min_matches=minMatchedPairs,
410  pattern_skip_array=np.array(
411  match_tolerance.failedPatternList)
412  )
413 
414  if soften_pattern == 0 and soften_dist == 0 and \
415  len(matcher_struct.match_ids) == 0 and \
416  match_tolerance.lastMatchedPattern is not None:
417  # If we found a pattern on a previous match-fit iteration and
418  # can't find an optimistic match on our first try with the
419  # tolerances as found in the previous match-fit,
420  # the match we found in the last iteration was likely bad. We
421  # append the bad match's index to the a list of
422  # patterns/matches to skip on subsequent iterations.
423  match_tolerance.failedPatternList.append(
424  match_tolerance.lastMatchedPattern)
425  match_tolerance.lastMatchedPattern = None
426  maxShiftArcseconds = \
427  self.config.maxOffsetPix * wcs.getPixelScale().asArcseconds()
428  elif len(matcher_struct.match_ids) > 0:
429  # Match found, save a bit a state regarding this pattern
430  # in the match tolerance class object and exit.
431  match_tolerance.maxShift = \
432  matcher_struct.shift * afwgeom.arcseconds
433  match_tolerance.lastMatchedPattern = \
434  matcher_struct.pattern_idx
435  match_found = True
436  break
437  if match_found:
438  break
439 
440  # If we didn't find a match, exit early.
441  if not match_found:
442  return pipeBase.Struct(
443  matches=[],
444  match_tolerance=match_tolerance,
445  )
446 
447  # The matcher returns all the nearest neighbors that agree between
448  # the reference and source catalog. For the current astrometric solver
449  # we need to remove as many false positives as possible before sending
450  # the matches off to the solver. The low value of 100 and high value of
451  # 2 are the low number of sigma and high respectively. The exact values
452  # were found after testing on data of various reference/source
453  # densities and astrometric distortion quality, specifically the
454  # visits: HSC (3358), DECam (406285, 410827),
455  # CFHT (793169, 896070, 980526).
456  distances_arcsec = np.degrees(matcher_struct.distances_rad) * 3600
457  clip_max_dist = np.max(
458  (sigmaclip(distances_arcsec, low=100, high=2)[-1],
459  self.config.minMatchDistPixels * wcs.getPixelScale().asArcseconds())
460  )
461  # Assuming the number of matched objects surviving the clip_max_dist
462  # cut if greater the requested min number of pairs, we select the
463  # smaller of the current maxMatchDist or the sigma clipped distance.
464  if not np.isfinite(clip_max_dist):
465  clip_max_dist = maxMatchDistArcSec
466 
467  if clip_max_dist < maxMatchDistArcSec and \
468  len(distances_arcsec[distances_arcsec < clip_max_dist]) < \
469  minMatchedPairs:
470  dist_cut_arcsec = maxMatchDistArcSec
471  else:
472  dist_cut_arcsec = np.min((clip_max_dist, maxMatchDistArcSec))
473 
474  # A match has been found, return our list of matches and
475  # return.
476  matches = []
477  for match_id_pair, dist_arcsec in zip(matcher_struct.match_ids,
478  distances_arcsec):
479  if dist_arcsec < dist_cut_arcsec:
480  match = afwTable.ReferenceMatch()
481  match.first = refCat[match_id_pair[1]]
482  match.second = sourceCat[match_id_pair[0]]
483  # We compute the true distance along and sphere instead
484  # and store it in units of arcseconds. The previous
485  # distances we used were approximate.
486  match.distance = match.first.getCoord().separation(
487  match.second.getCoord()).asArcseconds()
488  matches.append(match)
489 
490  return pipeBase.Struct(
491  matches=matches,
492  match_tolerance=match_tolerance,
493  )
494 
495  def _latlong_flux_to_xyz_mag(self, theta, phi, flux):
496  """Convert angles theta and phi and a flux into unit sphere
497  x, y, z, and a relative magnitude.
498 
499  Takes in a afw catalog object and converts the catalog object RA, DECs
500  to points on the unit sphere. Also converts the flux into a simple,
501  non-zero-pointed magnitude for relative sorting.
502 
503  Parameters
504  ----------
505  theta : `float`
506  Angle from the north pole (z axis) of the sphere
507  phi : `float`
508  Rotation around the sphere
509 
510  Return
511  ------
512  output_array : `numpy.ndarray`, (N, 4)
513  Spherical unit vector x, y, z with flux.
514  """
515  output_array = np.empty(4, dtype=np.float64)
516  output_array[0] = np.sin(theta)*np.cos(phi)
517  output_array[1] = np.sin(theta)*np.sin(phi)
518  output_array[2] = np.cos(theta)
519  if flux > 0:
520  output_array[3] = -2.5 * np.log10(flux)
521  else:
522  # Set flux to a very faint mag if its for some reason it
523  # does not exist
524  output_array[3] = 99.
525 
526  return output_array
527 
528  def _get_pair_pattern_statistics(self, cat_array):
529  """ Compute the tolerances for the matcher automatically by comparing
530  pinwheel patterns as we would in the matcher.
531 
532  We test how similar the patterns we can create from a given set of
533  objects by computing the spoke lengths for each pattern and sorting
534  them from smallest to largest. The match tolerance is the average
535  distance per spoke between the closest two patterns in the sorted
536  spoke length space.
537 
538  Parameters
539  ----------
540  cat_array : `numpy.ndarray`, (N, 3)
541  array of 3 vectors representing the x, y, z position of catalog
542  objects on the unit sphere.
543 
544  Returns
545  -------
546  dist_tol : `float`
547  Suggested max match tolerance distance calculated from comparisons
548  between pinwheel patterns used in optimistic/pessimistic pattern
549  matcher.
550  """
551 
552  self.log.debug("Starting automated tolerance calculation...")
553 
554  # Create an empty array of all the patterns we possibly make
555  # sorting from brightest to faintest.
556  pattern_array = np.empty(
557  (cat_array.shape[0] - self.config.numPointsForShape,
558  self.config.numPointsForShape - 1))
559  flux_args_array = np.argsort(cat_array[:, -1])
560 
561  # Sort our input array.
562  tmp_sort_array = cat_array[flux_args_array]
563 
564  # Start making patterns.
565  for start_idx in range(cat_array.shape[0] -
566  self.config.numPointsForShape):
567  pattern_points = tmp_sort_array[start_idx:start_idx +
568  self.config.numPointsForShape, :-1]
569  pattern_delta = pattern_points[1:, :] - pattern_points[0, :]
570  pattern_array[start_idx, :] = np.sqrt(
571  pattern_delta[:, 0] ** 2 +
572  pattern_delta[:, 1] ** 2 +
573  pattern_delta[:, 2] ** 2)
574 
575  # When we store the length of each spoke in our pattern we
576  # sort from shortest to longest so we have a defined space
577  # to compare them in.
578  pattern_array[start_idx, :] = pattern_array[
579  start_idx, np.argsort(pattern_array[start_idx, :])]
580 
581  # Create a searchable tree object of the patterns and find
582  # for any given pattern the closest pattern in the sorted
583  # spoke length space.
584  dist_tree = cKDTree(
585  pattern_array[:, :(self.config.numPointsForShape - 1)])
586  dist_nearest_array, ids = dist_tree.query(
587  pattern_array[:, :(self.config.numPointsForShape - 1)], k=2)
588  dist_nearest_array = dist_nearest_array[:, 1]
589  dist_nearest_array.sort()
590 
591  # We use the two closest patterns to set our tolerance.
592  dist_idx = 0
593  dist_tol = (np.degrees(dist_nearest_array[dist_idx]) * 3600. /
594  (self.config.numPointsForShape - 1.))
595 
596  self.log.debug("Automated tolerance")
597  self.log.debug("\tdistance/match tol: %.4f [arcsec]" % dist_tol)
598 
599  return dist_tol
def _doMatch(self, refCat, sourceCat, wcs, refFluxField, numUsableSources, minMatchedPairs, match_tolerance, sourceFluxField, verbose)
def __init__(self, maxMatchDist=None, autoMaxMatchDist=None, maxShift=None, lastMatchedPattern=None, failedPatternList=None)
def matchObjectsToSources(self, refCat, sourceCat, wcs, refFluxField, match_tolerance=None)