lsst.meas.astrom  13.0-14-g9415442+26
 All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups Pages
matchPessimisticB.py
Go to the documentation of this file.
1 from __future__ import absolute_import, division, print_function
2 
3 import numpy as np
4 from scipy.spatial import cKDTree
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 
20 class MatchTolerancePessimistic(MatchTolerance):
21  """Stores match tolerances for use in AstrometryTask and later
22  iterations of the matcher.
23 
24  Attributes
25  ----------
26  maxMatchDist : lsst.afw.geom.Angle
27  autoMaxMatchDist : lsst.afw.geom.Angle
28  maxShift : lsst.afw.geom.Angle
29  lastMatchedPattern : int
30  failedPatternList : list of ints
31  """
32 
33  def __init__(self, maxMatchDist=None, autoMaxMatchDist=None,
34  maxShift=None, lastMatchedPattern=None,
35  failedPatternList=None):
36  """Construct a MatchPessimisticTolerance
37 
38  MatchPessimisticBTask relies on several state variables to be
39  preserved over different iterations in the
40  AstrometryTask.matchAndFitWcs loop of AstrometryTask.
41 
42  Parameters
43  ----------
44  maxMatchDist : afw.geom.Angle
45  Current 2 sigma scatter from the previous matched wcs (if it
46  exists. It is None if this is the first iteration.)
47  autoMatxMatchDist : afw.geom.Angle
48  Result of the automated match tolerance generation.
49  maxShift : afw.geom.Angle
50  None for the first iteration or is the magnitude of the previous
51  iteration's wcs shift.
52  lastMatchedPattern : int
53  Reference to the position in themagnitude sorted source array where
54  a successful pattern match was found.
55  failedPatternList : list of ints
56  List of ints specifying indicies in the magnitude sourced source
57  array to skip. These are skipped are previous iterations that are
58  likely false positives due to the code having to soften after a
59  pattern is matched.
60  """
61  self.maxMatchDist = maxMatchDist
62  self.autoMaxMatchDist = autoMaxMatchDist
63  self.maxShift = maxShift
64  self.lastMatchedPattern = lastMatchedPattern
65  if failedPatternList is None:
67 
68 
69 class MatchPessimisticBConfig(pexConfig.Config):
70  """Configuration for MatchPessimisticBTask
71  """
72  numBrightStars = pexConfig.RangeField(
73  doc="Number of bright stars to use. Sets the max number of patterns "
74  "that can be tested.",
75  dtype=int,
76  default=200,
77  min=2,
78  )
79  minMatchedPairs = pexConfig.RangeField(
80  doc="Minimum number of matched pairs; see also minFracMatchedPairs.",
81  dtype=int,
82  default=30,
83  min=2,
84  )
85  minFracMatchedPairs = pexConfig.RangeField(
86  doc="Minimum number of matched pairs as a fraction of the smaller of "
87  "the number of reference stars or the number of good sources; "
88  "the actual minimum is the smaller of this value or "
89  "minMatchedPairs.",
90  dtype=float,
91  default=0.3,
92  min=0,
93  max=1,
94  )
95  matcherIterations = pexConfig.RangeField(
96  doc="Number of softening iterations in matcher.",
97  dtype=int,
98  default=5,
99  min=1,
100  )
101  maxOffsetPix = pexConfig.RangeField(
102  doc="Maximum allowed shift of WCS, due to matching (pixel).",
103  dtype=int,
104  default=300,
105  max=4000,
106  )
107  maxRotationDeg = pexConfig.RangeField(
108  doc="Rotation angle allowed between sources and position reference "
109  "objects (degrees).",
110  dtype=float,
111  default=1.0,
112  max=6.0,
113  )
114  numPointsForShape = pexConfig.Field(
115  doc="Number of points to define a shape for matching.",
116  dtype=int,
117  default=6,
118  )
119  numPointsForShapeAttempt = pexConfig.Field(
120  doc="Number of points to try for creating a shape. This value should "
121  "be greater than or equal to numPointsForShape.",
122  dtype=int,
123  default=7,
124  )
125  numPatternConsensus = pexConfig.Field(
126  doc="Number of implied shift/rotations from patterns that must agree "
127  "before it a given shift/rotation is accepted. This is only used "
128  "after the first softening iteration fails and if both the "
129  "number of reference and source objects is greater than "
130  "numBrightStars.",
131  dtype=int,
132  default=3,
133  )
134  sourceSelector = sourceSelectorRegistry.makeField(
135  doc="How to select sources for cross-matching. The default "
136  "matcherSourceSelector removes objects with low S/N, bad "
137  "saturated objects, edge objects, and interpolated objects.",
138  default="matcherPessimistic"
139  )
140 
141  def setDefaults(self):
142  sourceSelector = self.sourceSelector["matcherPessimistic"]
143  sourceSelector.setDefaults()
144 
145  def validate(self):
146  pexConfig.Config.validate(self)
148  raise ValueError("numPointsForShapeAttempt must be greater than "
149  "or equal to numPointsForShape.")
150 
151 
152 # The following block adds links to this task from the Task Documentation page.
153 # \addtogroup LSST_task_documentation
154 # \{
155 # \page measAstrom_MatchPessimisticBTask
156 # \ref MatchPessimisticBTask "MatchPessimisticBTask"
157 # Match sources to reference objects
158 # \}
159 
160 
161 class MatchPessimisticBTask(pipeBase.Task):
162  """!Match sources to reference objects
163 
164  @anchor MatchPessimisticBTask_
165 
166  @section meas_astrom_MatchPessimisticB_Contents Contents
167 
168  - @ref meas_astrom_MatchPessimisticB_Purpose
169  - @ref meas_astrom_MatchPessimisticB_Initialize
170  - @ref meas_astrom_MatchPessimisticB_IO
171  - @ref meas_astrom_MatchPessimisticB_Config
172  - @ref meas_astrom_MatchPessimisticB_Example
173  - @ref meas_astrom_MatchPessimisticB_Debug
174 
175  @section meas_astrom_MatchPessimisticB_Purpose Description
176 
177  Match sources to reference objects. This is often done as a preliminary
178  step to fitting an astrometric or photometric solution. For details about
179  the matching algorithm see pessimistic_pattern_matcher_b_3D.py
180 
181  @section meas_astrom_MatchPessimisticB_Initialize Task initialisation
182 
183  @copydoc \_\_init\_\_
184 
185  @section meas_astrom_MatchPessimisticB_IO Invoking the Task
186 
187  @copydoc matchObjectsToSources
188 
189  @section meas_astrom_MatchPessimisticB_Config Configuration
190  parameters
191 
192  See @ref MatchPessimisticBConfig
193 
194  To modify the tests for good sources for matching, create a new
195  sourceSelector class in meas_algorithms and use it in the config.
196 
197  @section meas_astrom_MatchPessimisticB_Example A complete example of
198  using MatchPessimisticBTask
199 
200  MatchPessimisticBTask is a subtask of AstrometryTask, which is called by
201  PhotoCalTask. See \ref meas_photocal_photocal_Example.
202 
203  @section meas_astrom_MatchPessimisticB_Debug Debug variables
204 
205  The @link lsst.pipe.base.cmdLineTask.CmdLineTask command line task@endlink
206  interface supports a flag @c -d to import @b debug.py from your
207  @c PYTHONPATH; see @ref baseDebug for more about @b debug.py files.
208 
209  The available variables in MatchPessimisticBTask are:
210  <DL>
211  <DT> @c verbose (bool)
212  <DD> If True then the matcher prints debug messages to stdout
213  </DL>
214 
215  To investigate the @ref meas_astrom_MatchPessimisticB_Debug, put something
216  like
217  @code{.py}
218  import lsstDebug
219  def DebugInfo(name):
220  # N.b. lsstDebug.Info(name) would call us recursively
221  debug = lsstDebug.getInfo(name)
222  if name == "lsst.pipe.tasks.astrometry":
223  debug.verbose = True
224 
225  return debug
226 
227  lsstDebug.Info = DebugInfo
228  @endcode
229  into your debug.py file and run this task with the @c --debug flag.
230  """
231 
232  ConfigClass = MatchPessimisticBConfig
233  _DefaultName = "matchObjectsToSources"
234 
235  def __init__(self, **kwargs):
236  pipeBase.Task.__init__(self, **kwargs)
237  self.makeSubtask("sourceSelector")
238 
239  @pipeBase.timeMethod
240  def matchObjectsToSources(self, refCat, sourceCat, wcs, refFluxField,
241  match_tolerance=None):
242  """!Match sources to position reference stars
243 
244  @param[in] refCat catalog of reference objects that overlap the
245  exposure; reads fields for:
246  - coord
247  - the specified flux field
248  @param[in] sourceCat catalog of sources found on an exposure;
249  Please check the required fields of your specified source selector
250  that the correct flags are present.
251  @param[in] wcs estimated WCS
252  @param[in] refFluxField field of refCat to use for flux
253  @param[in] match_tolerance is a MatchTolerance class object or None.
254  This this class is used to comunicate state between AstrometryTask
255  and MatcherTask. AstrometryTask will also set the MatchTolerance
256  class variable maxMatchDist based on the scatter AstrometryTask has
257  found after fitting for the wcs.
258  @return an lsst.pipe.base.Struct with fields:
259  - matches a list of matches, each instance of
260  lsst.afw.table.ReferenceMatch
261  - usableSourcCat a catalog of sources potentially usable for
262  matching.
263  - match_tolerance a MatchTolerance object containing the resulting
264  state variables from the match.
265  """
266  import lsstDebug
267  debug = lsstDebug.Info(__name__)
268 
269  # If we get an empty tolerance struct create the variables we need for
270  # this matcher.
271  if match_tolerance is None:
272  match_tolerance = MatchTolerancePessimistic()
273 
274  # usableSourceCat: sources that are good but may be saturated
275  numSources = len(sourceCat)
276  selectedSources = self.sourceSelector.selectSources(sourceCat)
277  goodSourceCat = selectedSources.sourceCat
278  numUsableSources = len(goodSourceCat)
279  self.log.info("Purged %d sources, leaving %d good sources" %
280  (numSources - numUsableSources, numUsableSources))
281 
282  if len(goodSourceCat) == 0:
283  raise pipeBase.TaskError("No sources are good")
284 
285  # avoid accidentally using sourceCat; use goodSourceCat from now on
286  del sourceCat
287 
288  minMatchedPairs = min(self.config.minMatchedPairs,
289  int(self.config.minFracMatchedPairs *
290  min([len(refCat), len(goodSourceCat)])))
291 
292  doMatchReturn = self._doMatch(
293  refCat=refCat,
294  sourceCat=goodSourceCat,
295  wcs=wcs,
296  refFluxField=refFluxField,
297  numUsableSources=numUsableSources,
298  minMatchedPairs=minMatchedPairs,
299  match_tolerance=match_tolerance,
300  sourceFluxField=self.sourceSelector.fluxField,
301  verbose=debug.verbose,
302  )
303  matches = doMatchReturn.matches
304  match_tolerance = doMatchReturn.match_tolerance
305 
306  if len(matches) == 0:
307  raise RuntimeError("Unable to match sources")
308 
309  self.log.info("Matched %d sources" % len(matches))
310  if len(matches) < minMatchedPairs:
311  self.log.warn("Number of matches is smaller than request")
312 
313  return pipeBase.Struct(
314  matches=matches,
315  usableSourceCat=goodSourceCat,
316  match_tolerance=match_tolerance,
317  )
318 
319  @pipeBase.timeMethod
320  def _doMatch(self, refCat, sourceCat, wcs, refFluxField, numUsableSources,
321  minMatchedPairs, match_tolerance, sourceFluxField, verbose):
322  """!Implementation of matching sources to position reference stars
323 
324  Unlike matchObjectsToSources, this method does not check if the sources
325  are suitable.
326 
327  @param[in] refCat catalog of position reference stars that overlap an
328  exposure
329  @param[in] sourceCat catalog of sources found on the exposure
330  @param[in] wcs estimated WCS of exposure
331  @param[in] refFluxField field of refCat to use for flux
332  @param[in] numUsableSources number of usable sources (sources with
333  known centroid that are not near the edge, but may be saturated)
334  @param[in] minMatchedPairs minimum number of matches
335  @param[in] match_tolerance a MatchTolerance object containing
336  variables specifying matcher tolerances and state from possible
337  previous runs.
338  @param[in] sourceInfo SourceInfo for the sourceCat
339  @param[in] verbose true to print diagnostic information to std::cout
340 
341  @return a list of matches, an instance of
342  lsst.afw.table.ReferenceMatch, a MatchTolerance object
343  """
344 
345  # Load the source and reference catalog as spherical points
346  # in numpy array. We do this rather than relying on internal
347  # lsst C objects for simplicity and because we require
348  # objects contiguous in memory. We need to do these slightly
349  # differently for the reference and source cats as they are
350  # different catalog objects with different fields.
351  ref_array = np.empty((len(refCat), 4), dtype=np.float64)
352  for ref_idx, refObj in enumerate(refCat):
353  theta = np.pi / 2 - refObj.getDec().asRadians()
354  phi = refObj.getRa().asRadians()
355  flux = refObj[refFluxField]
356  ref_array[ref_idx, :] = \
357  self._latlong_flux_to_xyz_mag(theta, phi, flux)
358 
359  src_array = np.empty((len(sourceCat), 4), dtype=np.float64)
360  for src_idx, srcObj in enumerate(sourceCat):
361  coord = wcs.pixelToSky(srcObj.getCentroid())
362  theta = np.pi / 2 - coord.getLatitude().asRadians()
363  phi = coord.getLongitude().asRadians()
364  flux = srcObj.getPsfFlux()
365  src_array[src_idx, :] = \
366  self._latlong_flux_to_xyz_mag(theta, phi, flux)
367 
368  # Set configurable defaults when we encounter None type or set
369  # state based on previous run of AstrometryTask._matchAndFitWcs.
370  if match_tolerance.maxShift is None:
371  maxShiftArcseconds = (self.config.maxOffsetPix *
372  wcs.pixelScale().asArcseconds())
373  else:
374  # We don't want to clamp down too hard on the allowed shift so
375  # we test that the smallest we ever allow is the pixel scale.
376  maxShiftArcseconds = np.max(
377  (match_tolerance.maxShift.asArcseconds(),
378  wcs.pixelScale().asArcseconds()))
379 
380  # If our tolerances are not set from a previous run, estiamte a
381  # starting tolerance guess from statistics of patterns we can create
382  # on both the source and reference catalog. We use the smaller of the
383  # two.
384  if match_tolerance.maxMatchDist is None:
385  self.log.debug("Computing source statistics...")
386  maxMatchDistArcSecSrc = self._get_pair_pattern_statistics(
387  src_array)
388  self.log.debug("Computing reference statistics...")
389  maxMatchDistArcSecRef = self._get_pair_pattern_statistics(
390  ref_array)
391  maxMatchDistArcSec = np.min((maxMatchDistArcSecSrc,
392  maxMatchDistArcSecRef))
393  match_tolerance.autoMaxDist = afwgeom.Angle(maxMatchDistArcSec,
394  afwgeom.arcseconds)
395  else:
396  maxMatchDistArcSec = np.max(
397  (wcs.pixelScale().asArcseconds(),
398  np.min((match_tolerance.maxMatchDist.asArcseconds(),
399  match_tolerance.autoMaxDist.asArcseconds()))))
400 
401  # Make sure the data we are considering is dense enough to require
402  # the consensus mode of the matcher. If not default to Optimistic
403  # pattern matcher behavior.
404  numConsensus = self.config.numPatternConsensus
405  minObjectsForConsensus = \
406  self.config.numBrightStars + self.config.numPointsForShapeAttempt
407  if ref_array.shape[0] < minObjectsForConsensus or \
408  src_array.shape[0] < minObjectsForConsensus:
409  numConsensus = 1
410 
411  self.log.debug("Current tol maxDist: %.4f arcsec" %
412  maxMatchDistArcSec)
413  self.log.debug("Current shift: %.4f arcsec" %
414  maxShiftArcseconds)
415 
416  # Create our matcher object.
417  pyPPMb = PessimisticPatternMatcherB(ref_array[:, :3], self.log)
418 
419  # Start the ineration over our tolerances.
420  for try_idx in range(self.config.matcherIterations):
421  if try_idx == 0 and \
422  match_tolerance.lastMatchedPattern is not None:
423  # If we are on the first, most stringent tolerance,
424  # and have already found a match, the matcher should behave
425  # like an optimistic pattern matcher. Exiting at the first
426  # match.
427  matcher_struct = pyPPMb.match(
428  source_array=src_array,
429  n_check=self.config.numPointsForShapeAttempt,
430  n_match=self.config.numPointsForShape,
431  n_agree=1,
432  max_n_patterns=self.config.numBrightStars,
433  max_shift=maxShiftArcseconds,
434  max_rotation=self.config.maxRotationDeg,
435  max_dist=maxMatchDistArcSec,
436  min_matches=minMatchedPairs,
437  pattern_skip_array=np.array(
438  match_tolerance.failedPatternList)
439  )
440  else:
441  # Once we fail and start softening tolerences we switch to
442  # pessimistic or consenus mode where 3 patterns must agree
443  # on a rotation before exiting. We double the match dist
444  # tolerance each round and add 1 to the pattern complexiy and
445  # two to the number of candidate spokes to check.
446  matcher_struct = pyPPMb.match(
447  source_array=src_array,
448  n_check=self.config.numPointsForShapeAttempt + 2 * try_idx,
449  n_match=self.config.numPointsForShape + try_idx,
450  n_agree=numConsensus,
451  max_n_patterns=self.config.numBrightStars,
452  max_shift=(self.config.maxOffsetPix *
453  wcs.pixelScale().asArcseconds()),
454  max_rotation=self.config.maxRotationDeg,
455  max_dist=maxMatchDistArcSec * 2. ** try_idx,
456  min_matches=minMatchedPairs,
457  pattern_skip_array=np.array(
458  match_tolerance.failedPatternList)
459  )
460  if try_idx == 0 and \
461  len(matcher_struct.match_ids) == 0 and \
462  match_tolerance.lastMatchedPattern is not None:
463  # If we found a pattern on a previous match-fit iteration and
464  # can't find an optimistic match on our first try with the
465  # tolerances as found in the previous match-fit,
466  # the match we found in the last iteration was likely bad. We
467  # append the bad match's index to the a list of
468  # patterns/matches to skip on subsequent iterations.
469  match_tolerance.failedPatternList.append(
470  match_tolerance.lastMatchedPattern)
471  match_tolerance.lastMatchedPattern = None
472  elif len(matcher_struct.match_ids) > 0:
473  # Match found, save a bit a state regarding this pattern
474  # in the struct and exit.
475  match_tolerance.maxShift = afwgeom.Angle(matcher_struct.shift,
476  afwgeom.arcseconds)
477  match_tolerance.lastMatchedPattern = \
478  matcher_struct.pattern_idx
479  break
480 
481  # The matcher returns all the nearest neighbors that agree between
482  # the reference and source catalog. For the currrent astrometry solver
483  # we need to remove as many false positives as possible before sending
484  # the matches off to the astronomery solver.
485  distances_arcsec = np.degrees(matcher_struct.distances_rad) * 3600
486  # Rather than have a hard cut on the max distance we allow for
487  # objects we preform a iterative sigma clipping. We input a
488  # starting point for the clipping based on the largest of the
489  # distances we have. The hope is this distance will not be so
490  # large as to have the sigma_clip return a large unphysical value
491  # and not so small that we don't allow for some distortion in the
492  # WCS between source and reference.
493  clipped_dist_arcsec = self._sigma_clip(
494  values=distances_arcsec, n_sigma=2,
495  input_cut=np.max((10 * wcs.pixelScale().asArcseconds(),
496  maxMatchDistArcSec,
497  maxShiftArcseconds)))
498  # We pick the largest of the sigma clipping, the unsofted
499  # maxMatchDistArcSec or 2 pixels. This pevents
500  # the AstrometryTask._matchAndFitWCS from overfitting to a small
501  # number of objects and also allows the WCS fitter to bring in more
502  # matches as the WCS fit improves.
503  dist_cut_arcsec = np.max(
504  (2 * wcs.pixelScale().asArcseconds(),
505  clipped_dist_arcsec, maxMatchDistArcSec))
506 
507  # A match has been found, return our list of matches and
508  # return.
509  matches = []
510  for match_id_pair, dist_arcsec in zip(matcher_struct.match_ids,
511  distances_arcsec):
512  if dist_arcsec < dist_cut_arcsec:
513  match = afwTable.ReferenceMatch()
514  match.first = refCat[match_id_pair[1]]
515  match.second = sourceCat[match_id_pair[0]]
516  # We compute the true distance along and sphere instead
517  # and store it in units of arcseconds. The previous
518  # distances we used were aproximate.
519  match.distance = match.first.getCoord().angularSeparation(
520  match.second.getCoord()).asArcseconds()
521  matches.append(match)
522 
523  return pipeBase.Struct(
524  matches=matches,
525  match_tolerance=match_tolerance,
526  )
527 
528  def _sigma_clip(self, values, n_sigma, input_cut=None):
529  """
530  Find the distance cut by recursively cliping by n_sigma.
531 
532  The algorithm converges once the length of the array on the
533  nth iteration is equal to the the length of the array on the
534  n+1th interation.
535 
536  Parameters
537  ----------
538  values : float array
539  float array of values to caculate the cut on.
540  n_sigma : float
541  Number of sigma to use for the recersive clipping.
542  input_cut : float
543  Optional starting cut for the clipping.
544 
545  Returns
546  -------
547  float
548  Value at which the clipping algorithm converged or the
549  value after a maximum of 100 iterations.
550  """
551  # Initialize the variables that will persist over the different
552  # iteration steps.
553  hold_total = 0
554  if input_cut is None:
555  hold_cut = np.max(values)
556  else:
557  hold_cut = input_cut
558  # Start clipping
559  for clip_idx in range(100):
560  # Find the mask and current total objects.
561  mask = (values < hold_cut)
562  total = len(values[mask])
563  # Test to see if the current total is equal to the previous.
564  # If it is we exit returning the previous cut.
565  if total == hold_total:
566  return hold_cut
567  # If we have not converged store the current total.
568  hold_tot = total
569  # Recompute the cut using the current mask.
570  hold_cut = (np.mean(values[mask]) +
571  n_sigma * np.std(values[mask], ddof=1))
572  return hold_cut
573 
574  def _latlong_flux_to_xyz_mag(self, theta, phi, flux):
575  r"""Convert angles theta and phi plux a flux into unit sphere
576  x, y, z, and a relative magnitude.
577 
578  Takes in a afw catalog objet and converts the catalog object RA, DECs
579  to points on the unit sphere. Also convets the flux into a simple,
580  non-zeropointed magnitude for relative sorting.
581 
582  Parameters
583  ----------
584  theta : float
585  Angle from the north pole (z axis) of the sphere
586  phi : float
587  Rotation around the sphere
588 
589  Return
590  ------
591  float array
592  Spherical unit vector x, y, z on the unitsphere.
593  """
594  output_array = np.empty(4, dtype=np.float64)
595  output_array[0] = np.sin(theta)*np.cos(phi)
596  output_array[1] = np.sin(theta)*np.sin(phi)
597  output_array[2] = np.cos(theta)
598  if flux > 0:
599  output_array[3] = -2.5 * np.log10(flux)
600  else:
601  # Set flux to a very faint mag if its for some reason it
602  # does not exist
603  output_array[3] = 99.
604 
605  return output_array
606 
607  def _get_pair_pattern_statistics(self, cat_array):
608  """ Compute the tolerances for the matcher automatically by comparing
609  pinwheel patterns as we would in the matcher.
610 
611  We test how simplar the patterns we can create from a given set of
612  objects by computing the spoke lengths for each pattern and sorting
613  them from smallest to largest. The match tolerance is the average
614  distance per spoke between the closest two patterns in the sorted
615  spoke length space.
616 
617  Parameters
618  ----------
619  cat_array : float array
620  array of 3 vectors representing the x, y, z position of catalog
621  objedts on the unit sphere.
622 
623  Returns
624  -------
625  float
626  Suggested max match tolerance distance calculated from comparisons
627  between pinwheel patterns used in optimistic/pessimsitic pattern
628  matcher.
629  """
630 
631  self.log.debug("Starting automated tolerance calculation...")
632 
633  # Create an empty array of all the patterns we possibly make
634  # sorting from brightest to faintest.
635  pattern_array = np.empty(
636  (cat_array.shape[0] - self.config.numPointsForShape,
637  self.config.numPointsForShape - 1))
638  flux_args_array = np.argsort(cat_array[:, -1])
639 
640  # Sort our input array.
641  tmp_sort_array = cat_array[flux_args_array]
642 
643  # Start making patterns.
644  for start_idx in range(cat_array.shape[0] -
645  self.config.numPointsForShape):
646  pattern_points = tmp_sort_array[start_idx:start_idx +
647  self.config.numPointsForShape, :-1]
648  pattern_delta = pattern_points[1:, :] - pattern_points[0, :]
649  pattern_array[start_idx, :] = np.sqrt(
650  pattern_delta[:, 0] ** 2 +
651  pattern_delta[:, 1] ** 2 +
652  pattern_delta[:, 2] ** 2)
653 
654  # When we store the length of each spoke in our pattern we
655  # sort from shortest to longest so we have a defined space
656  # to compare them in.
657  pattern_array[start_idx, :] = pattern_array[
658  start_idx, np.argsort(pattern_array[start_idx, :])]
659 
660  # Create a searchable tree object of the patterns and find
661  # for any given pattern the closest pattern in the sorted
662  # spoke length space.
663  dist_tree = cKDTree(
664  pattern_array[:, :(self.config.numPointsForShape - 1)])
665  dist_nearest_array, ids = dist_tree.query(
666  pattern_array[:, :(self.config.numPointsForShape - 1)], k=2)
667  dist_nearest_array = dist_nearest_array[:, 1]
668  dist_nearest_array.sort()
669 
670  # We use the two closest patterns to set our tolerance.
671  dist_idx = 0
672  dist_tol = (np.degrees(dist_nearest_array[dist_idx]) * 3600. /
673  (self.config.numPointsForShape - 1.))
674 
675  self.log.debug("Automated tolerance")
676  self.log.debug("\tdistance/match tol: %.4f [arcsec]" % dist_tol)
677 
678  return dist_tol
def matchObjectsToSources
Match sources to position reference stars.
def _doMatch
Implementation of matching sources to position reference stars.