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