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