lsst.meas.astrom  13.0-16-g30f6da5+16
pessimistic_pattern_matcher_b_3D.py
Go to the documentation of this file.
1 
2 from __future__ import division, print_function, absolute_import
3 
4 import numpy as np
5 from scipy.optimize import leastsq
6 from scipy.spatial import cKDTree
7 
8 import lsst.pipe.base as pipeBase
9 
10 
11 def _rotation_matrix_chi_sq(flattened_rot_matrix, src_pattern, ref_pattern):
12  """Compute the squared differences for least squares fitting.
13 
14  Given a flattened rotation matrix, a N point pattern from the source
15  catalog and the reference pattern the sources match to, compute
16  the squared differences between the points in the two patterns after
17  rotation.
18 
19  Parameters
20  ----------
21  flattened_rot_matrix : float array
22  A flattened array representing a 3x3 rotation matrix. The array is
23  flattened to comply with the API of scipy.optimize.leastsq. Flattened
24  elements are [[0, 0], [0, 1], [0, 2], [1, 0]...]
25  src_patterm : float array of 3 vectors
26  A array containing N, 3 vectors
27  ref_pattern : float array of 3 vectors
28  A array containing N, 3 vectors
29 
30  Returns
31  -------
32  float array
33  Array of differences between the vectors representing of the source
34  pattern rotated into the reference frame and the converse. This is
35  used to minimize in a least squares fitter.
36  """
37  # Unflatten the rotation matrix
38  rot_matrix = flattened_rot_matrix.reshape((3, 3))
39  # Compare the rotated source pattern to the references.
40  diff_2_ref = np.dot(rot_matrix, src_pattern.transpose()).transpose() -\
41  ref_pattern
42  # Compare the rotated reference pattern to the sources. This second
43  # test forces the fitted rotation matrix to be aproximately unitary and
44  # reversable.
45  # That is np.dot(rot_matrix.transpose(), rot_matrix) ~= idenity
46  diff_2_src = np.dot(rot_matrix.transpose(),
47  ref_pattern.transpose()).transpose() -\
48  src_pattern
49  # Return the flattened differences for use in a least squares fitter.
50  return np.concatenate((diff_2_src.flatten(), diff_2_ref.flatten()))
51 
52 
54  """ Class implementing a pessimistic version of Optimsitic Pattern Matcher
55  B (OPMb) from Tabur 2007. The class loads and stores the reference object
56  in a convienent data structure for matching any set of source obejcts that
57  are assumed to contain each other. The pessimistic nature of the algorithm
58  comes from requiring that it discovers at least two patterns that agree on
59  the correct shift and rotation for matching before exiting. The original
60  behavior of OPMb can be recovered simply. Patterns matched between the
61  input datasets are n-spoked pinwheels created from n+1 points. Refer to
62  DMTN #031 for more details. http://github.com/lsst-dm/dmtn-031
63  --------------------------------------------------------------------------
64  Attributes:
65  reference_array : float array
66  spherical points x, y, z of to use as reference objects for
67  pattern matching.
68  log : an lsst.log instance
69  pair_id_array : int array
70  Internal lookup table. Given an id in the reference array, return
71  an array of the id pair that contains this object's id sorted on
72  the distance to the pairs.
73  pair_delta_array : float array
74  Internal lookup table. Given an id in the reference array, return
75  an array of the 3 vector deltas of all other pairs sorted on their
76  distance.
77  pair_dist_array : float array
78  Internal lookup table. Given an id in the reference return an array
79  of pair distances of all other pairs sorted on distance.
80  dist_array : float array
81  Array of all pairs of objects in the reference array sorted on
82  distance.
83  id_array : int array
84  Array of id pairs that lookup into the reference array sorted
85  on pair distance.
86  delta_array : float array
87  Array of 3 vector deltas for each pair in the reference array
88  sorted on pair distance.
89  """
90 
91  def __init__(self, reference_array, log):
92  """
93  Arguments
94  ---------
95  reference_array : float array
96  Array of spherical points x, y, z to use as reference objects.
97  log : lsst.log
98  logger object for reporting warnings and failures.
99  """
100  self._reference_array = reference_array
101  self._n_reference = len(self._reference_array)
102  self.log = log
103 
105 
106  def _build_distances_and_angles(self):
107  """ Create the data structures we will use to search for our pattern
108  match in.
109 
110  Throughout this function and the rest of the
111 
112  class we use id to reference the position in the input reference
113  catalog and index to 'index' into the arrays sorted on distance.
114 
115  Parameters
116  ----------
117  self
118 
119  Returns
120  -------
121  None
122  """
123 
124  # Initialize the arrays we will need for quick look up of pairs once
125  # have a candidate spoke center.
126  self._pair_id_array = np.empty(
127  (self._n_reference, self._n_reference - 1),
128  dtype=np.uint32)
129  self._pair_delta_array = np.empty(
130  (self._n_reference, self._n_reference - 1, 3),
131  dtype=np.float64)
132  self._pair_dist_array = np.empty(
133  (self._n_reference, self._n_reference - 1),
134  dtype=np.float64)
135 
136  # Create empty lists to temporarially store our pair information per
137  # reference object. These will be concatentated into our final arrays.
138  sub_id_array_list = []
139  sub_delta_array_list = []
140  sub_dist_array_list = []
141 
142  # Loop over reference objects and store pair distances, ids, and
143  # 3 vector deltas.
144  for ref_id, ref_obj in enumerate(self._reference_array):
145 
146  # Reserve and fill the ids of each reference object pair.
147  sub_id_array = np.zeros((self._n_reference - 1 - ref_id, 2),
148  dtype=np.uint32)
149  sub_id_array[:, 0] = ref_id
150  sub_id_array[:, 1] = np.arange(ref_id + 1, self._n_reference,
151  dtype=np.uint32)
152 
153  # Compute the vector deltas for each pair of reference objects
154  # and compute and store the distances.
155  sub_delta_array = (self._reference_array[ref_id + 1:, :] -
156  ref_obj)
157  sub_dist_array = np.sqrt(sub_delta_array[:, 0] ** 2 +
158  sub_delta_array[:, 1] ** 2 +
159  sub_delta_array[:, 2] ** 2)
160 
161  # Append to our arrays to the output lists for later
162  # concatenation.
163  sub_id_array_list.append(sub_id_array)
164  sub_delta_array_list.append(sub_delta_array)
165  sub_dist_array_list.append(sub_dist_array)
166 
167  # Fill the pair look up arrays row wise and then column wise.
168  self._pair_id_array[ref_id, ref_id:] = sub_id_array[:, 1]
169  self._pair_delta_array[ref_id, ref_id:, :] = sub_delta_array
170  self._pair_dist_array[ref_id, ref_id:] = sub_dist_array
171 
172  # Don't fill the array column wise if we are on the last as
173  # these pairs have already been filled by previous
174  # iterations.
175  if ref_id < self._n_reference - 1:
176  self._pair_id_array[ref_id + 1:, ref_id] = ref_id
177  self._pair_delta_array[ref_id + 1:, ref_id, :] = \
178  sub_delta_array
179  self._pair_dist_array[ref_id + 1:, ref_id] = sub_dist_array
180 
181  # Sort each row on distance for fast look up of pairs given
182  # the id of one of the objects in the pair.
183  sorted_pair_dist_args = self._pair_dist_array[ref_id, :].argsort()
184  self._pair_dist_array[ref_id, :] = self._pair_dist_array[
185  ref_id, sorted_pair_dist_args]
186  self._pair_id_array[ref_id, :] = self._pair_id_array[
187  ref_id, sorted_pair_dist_args]
188  self._pair_delta_array[ref_id, :, :] = self._pair_delta_array[
189  ref_id, sorted_pair_dist_args, :]
190 
191  # Concantent our arrays together.
192  unsorted_id_array = np.concatenate(sub_id_array_list)
193  unsorted_delta_array = np.concatenate(sub_delta_array_list)
194  unsorted_dist_array = np.concatenate(sub_dist_array_list)
195 
196  # Sort each array on the pair distances for the inital
197  # optimistic pattern matcher lookup.
198  sorted_dist_args = unsorted_dist_array.argsort()
199  self._dist_array = unsorted_dist_array[sorted_dist_args]
200  self._id_array = unsorted_id_array[sorted_dist_args]
201  self._delta_array = unsorted_delta_array[sorted_dist_args]
202 
203  return None
204 
205  def match(self, source_array, n_check, n_match, n_agree,
206  max_n_patterns, max_shift, max_rotation, max_dist,
207  min_matches, pattern_skip_array=None):
208  """Match a given source catalog into the loaded reference catalog.
209 
210  Given array of points on the unit sphere and tolerances, we
211  attempt to match a pinwheel like pattern between these input sources
212  and the reference objects this class was created with. This pattern
213  informs of the shift and rotation needed to align the input source
214  objects into the frame of the refernces.
215 
216  Parameters
217  ----------
218  source_array: float array
219  An array of spherical x,y,z coordinates and a magnitude in units
220  of objects having a lower value for sorting. The array should be
221  of shape (N, 4).
222  n_check : int value
223  Number of sources to create a pattern from. Not all objects may be
224  checked if n_match criteria is before looping through all n_check
225  objects.
226  n_match : int value
227  Number of objects to use in constructing a pattern to match.
228  n_agree: int value
229  Number of found patterns that must agree on their shift and
230  rotation before exiting. Set this value to 1 to recover the
231  expected behavior of Optimistic Pattern Matcher B.
232  max_n_patters : int value
233  Number of patterns to create from the input source objects to
234  attempt to match into the reference objects.
235  max_shift: float value
236  Maximum allowed shift to match patterns in arcseconds.
237  max_rotation: float value
238  Maximum allowed rotation between patterns in degrees.
239  max_dist: float value
240  Maximum distance in arcseconds allowed between candidate spokes in
241  the source and reference objects. Also sets that maximum distance
242  in the intermediate verify, pattern shift/rotation agreement, and
243  final verify steps.
244  pattern_skip_array: int array
245  Patterns we would like to skip. This could be due to the pattern
246  being matched on a previous iteration that we now consider invalid.
247  This assumes the ordering of the source objects is the same
248  between different runs of the matcher which, assuming no object
249  has been inserted or the magnitudes have changed, it should be.
250 
251  Returns
252  -------
253  output_struct : pipe.base.struct
254  A pipebase struct containing the following outputs.
255 
256  matches : int array
257  (N, 2) array of matched ids for pairs. Empty list if no
258  match found.
259  distances_rad : float array
260  Radian distances between the matched objects. Empty list
261  if no match found.
262  pattern_idx : int
263  Index of matched pattern. None if no match found.
264  shift : float
265  Magnitude for the shift between the source and
266  reference objects in arcseconds. None if no match found.
267  """
268 
269  # Given our input source_array we sort on magnitude.
270  sorted_source_array = source_array[source_array[:, -1].argsort(), :3]
271  n_source = len(sorted_source_array)
272 
273  # Initialize output struct.
274  output_match_struct = pipeBase.Struct(
275  match_ids=[],
276  distances_rad=[],
277  pattern_idx=None,
278  shift=None,)
279 
280  if n_source <= 0:
281  self.log.warn("Source object array is empty. Unable to match. "
282  "Exiting matcher.")
283  return None
284 
285  # To test if the shifts and rotations we find agree with each other,
286  # we first create two test points situated at the top and bottom of
287  # where the z axis on the sphere bisects the source catalog.
288  test_vect_list = self._compute_test_vectors(source_array[:, :3])
289 
290  # We now create an empty list of our resultant rotated vectors to
291  # compare the different rotations we find.
292  rot_vect_list = []
293 
294  # Convert the tolerances to values we will use in the code.
295  max_cos_shift = np.cos(np.radians(max_shift / 3600.))
296  max_cos_rot_sq = np.cos(np.radians(max_rotation)) ** 2
297  max_dist_rad = np.radians(max_dist / 3600.)
298 
299  # Loop through the sources from brightest to faintest, grabbing a
300  # chunk of n_check each time.
301  for pattern_idx in range(np.min((max_n_patterns,
302  n_source - n_match))):
303 
304  # If this pattern is one that we matched on the past but we
305  # now want to skip, we do so here.
306  if pattern_skip_array is not None and \
307  np.any(pattern_skip_array == pattern_idx):
308  self.log.debug(
309  "Skipping previously matched bad pattern %i..." %
310  pattern_idx)
311  continue
312  # Grab the sources to attempt to create this pattern.
313  pattern = sorted_source_array[
314  pattern_idx: np.min((pattern_idx + n_check, n_source)), :3]
315 
316  # Construct a pattern given the number of points defining the
317  # pattern complexity. This is the start of the primary tests to
318  # match our source pattern into the reference objects.
319  construct_return_struct = \
321  pattern, n_match, max_cos_shift, max_cos_rot_sq,
322  max_dist_rad)
323 
324  # Our struct is None if we could not match the pattern.
325  if construct_return_struct.ref_candidates is None or \
326  construct_return_struct.shift_rot_matrix is None or\
327  construct_return_struct.cos_shift is None or \
328  construct_return_struct.sin_rot is None:
329  continue
330 
331  # Grab the output data from the Struct object.
332  ref_candidates = construct_return_struct.ref_candidates
333  shift_rot_matrix = construct_return_struct.shift_rot_matrix
334  cos_shift = construct_return_struct.cos_shift
335  sin_rot = construct_return_struct.sin_rot
336 
337  # If we didn't match enough candidates we continue to the next
338  # pattern.
339  if len(ref_candidates) < n_match:
340  continue
341 
342  # Now that we know our pattern and shift/rotation are valid we
343  # store the the rotated versions of our test points for later
344  # use.
345  tmp_rot_vect_list = []
346  for test_vect in test_vect_list:
347  tmp_rot_vect_list.append(np.dot(shift_rot_matrix, test_vect))
348  tmp_rot_vect_list.append(pattern_idx)
349  rot_vect_list.append(tmp_rot_vect_list)
350 
351  # Test if we have enough rotations, which agree, or if we
352  # are in optimistic mode.
353  if self._test_rotation_agreement(rot_vect_list, max_dist_rad) < \
354  n_agree - 1:
355  continue
356 
357  # Perform final verify.
358  match_sources_struct = self._match_sources(
359  source_array[:, :3], shift_rot_matrix)
360 
361  # Check that we have enough matches.
362  if len(match_sources_struct.match_ids) >= min_matches:
363  # Convert the observed shift to arcseconds
364  shift = np.degrees(np.arccos(cos_shift)) * 3600.
365  # Print information to the logger.
366  self.log.debug("Succeeded after %i patterns." % pattern_idx)
367  self.log.debug("\tShift %.4f arcsec" % shift)
368  self.log.debug("\tRotation: %.4f deg" %
369  np.degrees(np.arcsin(sin_rot)))
370 
371  # Fill the struct and return.
372  output_match_struct.match_ids = \
373  match_sources_struct.match_ids
374  output_match_struct.distances_rad = \
375  match_sources_struct.distances_rad
376  output_match_struct.pattern_idx = pattern_idx
377  output_match_struct.shift = shift
378  return output_match_struct
379 
380  self.log.warn("Failed after %i patterns." % (pattern_idx + 1))
381  return output_match_struct
382 
383  def _compute_test_vectors(self, source_array):
384  """Compute spherical 3 vectors at the edges of the x, y, z extent
385  of the input source catalog.
386 
387  Parameters
388  ----------
389  source_array : float array (N, 3)
390  array of 3 vectors representing possitions on the unit
391  sphere.
392 
393  Returns
394  -------
395  list of 3 vectors
396  list of vectors representing the maxinum extents in x, y, z
397  of the input source array. These are used with the rotations
398  the code finds to test for agreement from different patterns
399  when the code is running in pessimistic mode.
400  """
401 
402  # Get the center of source_array.
403  if np.any(np.logical_not(np.isfinite(source_array))):
404  self.log.warn("Input source objects contain non-finite values. "
405  "This could end badly.")
406  center_vect = np.nanmean(source_array, axis=0)
407 
408  # So that our rotation test works over the full sky we compute
409  # the max extent in each cartesian direction x,y,z.
410  xbtm_vect = np.array([np.min(source_array[:, 0]), center_vect[1],
411  center_vect[2]], dtype=np.float64)
412  xtop_vect = np.array([np.max(source_array[:, 0]), center_vect[1],
413  center_vect[2]], dtype=np.float64)
414  xbtm_vect /= np.sqrt(np.dot(xbtm_vect, xbtm_vect))
415  xtop_vect /= np.sqrt(np.dot(xtop_vect, xtop_vect))
416 
417  ybtm_vect = np.array([center_vect[0], np.min(source_array[:, 1]),
418  center_vect[2]], dtype=np.float64)
419  ytop_vect = np.array([center_vect[0], np.max(source_array[:, 1]),
420  center_vect[2]], dtype=np.float64)
421  ybtm_vect /= np.sqrt(np.dot(ybtm_vect, ybtm_vect))
422  ytop_vect /= np.sqrt(np.dot(ytop_vect, ytop_vect))
423 
424  zbtm_vect = np.array([center_vect[0], center_vect[1],
425  np.min(source_array[:, 2])], dtype=np.float64)
426  ztop_vect = np.array([center_vect[0], center_vect[1],
427  np.max(source_array[:, 2])], dtype=np.float64)
428  zbtm_vect /= np.sqrt(np.dot(zbtm_vect, zbtm_vect))
429  ztop_vect /= np.sqrt(np.dot(ztop_vect, ztop_vect))
430 
431  # Return our list of vectors for later rotation testing.
432  return [xbtm_vect, xtop_vect, ybtm_vect, ytop_vect,
433  zbtm_vect, ztop_vect]
434 
435  def _construct_pattern_and_shift_rot_matrix(self, src_pattern_array,
436  n_match, max_cos_theta_shift,
437  max_cos_rot_sq, max_dist_rad):
438  """Test an input source pattern against the reference catalog.
439 
440  Returns the candidate matched patterns and their
441  implied rotation matrices or None.
442 
443  Parameters
444  ----------
445  src_pattern_array : float array
446  Sub selection of source 3 vectors to create a pattern from
447  n_match : int
448  Number of points to attempt to create a pattern from. Must be
449  >= len(src_pattern_array)
450  max_cos_theta_shift : float
451  Maximum shift allowed between two patterns' centers.
452  max_cos_rot_sq : float
453  Maximum rotation beteween two patterns that have been shifted
454  to have their centers on top of each other.
455  max_dist_rad : float
456  Maximum delta distance allowed between the source and reference
457  pair distances to consider the reference pair a candidate for
458  the source pair. Also sets the tolerance between the opening
459  angles of the spokes when compared to the reference.
460 
461  Return
462  -------
463  lsst.pipe.base.Struct
464  Return a Struct containing the following data:
465 
466  ref_candidates : list of ints
467  ids of the matched pattern in the internal reference_array
468  object.
469  src_candidates : list of ints
470  Pattern ids of the sources matched.
471  shift_rot_matrix : float array
472  3x3 matrix specifying the full shift
473  and rotation between the reference and source objects. Rotates
474  source into reference frame. None if match is not found.
475  cos_shift : float
476  Magnitude of the shift found between the two patten
477  centers. None if match is not found.
478  sin_rot : float value of the rotation to align the already shifted
479  source pattern to the reference pattern. None if no match
480  found.
481  """
482 
483  # Create our place holder variables for the matched sources and
484  # referces. The source list starts with the 0th and first indexed
485  # objects as we are garunteed to use those and these define both
486  # the shift and rotation of the final pattern.
487  output_matched_pattern = pipeBase.Struct(
488  ref_candidates=[],
489  src_candidates=[],
490  shift_rot_matrix=None,
491  cos_shift=None,
492  sin_rot=None)
493 
494  # Create the delta vectors and distances we will need to assemble the
495  # spokes of the pattern.
496  src_delta_array = np.empty((len(src_pattern_array) - 1, 3))
497  src_delta_array[:, 0] = (src_pattern_array[1:, 0] -
498  src_pattern_array[0, 0])
499  src_delta_array[:, 1] = (src_pattern_array[1:, 1] -
500  src_pattern_array[0, 1])
501  src_delta_array[:, 2] = (src_pattern_array[1:, 2] -
502  src_pattern_array[0, 2])
503  src_dist_array = np.sqrt(src_delta_array[:, 0]**2 +
504  src_delta_array[:, 1]**2 +
505  src_delta_array[:, 2]**2)
506 
507  # Our first test. We search the reference dataset for pairs
508  # that have the same length as our first source pairs to with
509  # plus/minus the max_dist tolerance.
510  ref_dist_index_array = self._find_candidate_reference_pairs(
511  src_dist_array[0], self._dist_array, max_dist_rad)
512 
513  # Start our loop over the candidate reference objects.
514  for ref_dist_idx in ref_dist_index_array:
515  # We have two candidates for which reference object corresponds
516  # with the source at the center of our pattern. As such we loop
517  # over and test both possiblities.
518  tmp_ref_pair_list = self._id_array[ref_dist_idx]
519  for pair_idx, ref_id in enumerate(tmp_ref_pair_list):
520  src_candidates = [0, 1]
521  ref_candidates = []
522  shift_rot_matrix = None
523  cos_shift = None
524  sin_rot = None
525  # Test the angle between our candidate ref center and the
526  # source center of our pattern. This angular distance also
527  # defines the shift we will later use.
528  ref_center = self._reference_array[ref_id]
529  cos_shift = np.dot(src_pattern_array[0], ref_center)
530  if cos_shift < max_cos_theta_shift:
531  continue
532 
533  # We can now append this one as a candaite.
534  ref_candidates.append(ref_id)
535  ref_delta = self._delta_array[ref_dist_idx]
536  # If the candidate reference center we found is second in
537  # this pair we need to reverse the direction of the
538  # corresponding pair's delta vector.
539  if pair_idx == 0:
540  ref_candidates.append(
541  tmp_ref_pair_list[1])
542  else:
543  ref_candidates.append(
544  tmp_ref_pair_list[0])
545  ref_delta *= -1
546 
547  # For dense fields it will be faster to compute the absolute
548  # rotation this pair suggests first rather than saving it
549  # after all the spokes are found. We then compute the cos^2
550  # of the rotation and first part of the rotation matrix from
551  # source to reference frame.
552  test_rot_struct = self._test_rotation(
553  src_pattern_array[0], ref_center, src_delta_array[0],
554  ref_delta, cos_shift, max_cos_rot_sq)
555  if test_rot_struct.cos_rot_sq is None or \
556  test_rot_struct.shift_matrix is None:
557  continue
558 
559  # Get the data from the return struct.
560  cos_rot_sq = test_rot_struct.cos_rot_sq
561  shift_matrix = test_rot_struct.shift_matrix
562 
563  # Now that we have a candidate first spoke and reference
564  # pattern center, we mask our future search to only those
565  # pairs that contain our candidate reference center.
566  ref_dist = self._dist_array[ref_dist_idx]
567  tmp_ref_delta_array = self._pair_delta_array[ref_id]
568  tmp_ref_dist_arary = self._pair_dist_array[ref_id]
569  tmp_ref_id_array = self._pair_id_array[ref_id]
570 
571  # Now we feed this sub data to match the spokes of
572  # our pattern.
573  pattern_spoke_struct = self._create_pattern_spokes(
574  src_pattern_array[0], src_delta_array, src_dist_array,
575  self._reference_array[ref_id], ref_id, ref_delta,
576  ref_dist, tmp_ref_delta_array, tmp_ref_dist_arary,
577  tmp_ref_id_array, max_dist_rad,
578  n_match)
579 
580  # If we don't find enough candidates we can continue to the
581  # next reference center pair.
582  if len(pattern_spoke_struct.ref_spoke_list) < n_match - 2 or \
583  len(pattern_spoke_struct.src_spoke_list) < n_match - 2:
584  continue
585 
586  # If we have the right number of matched ids we store these.
587  ref_candidates.extend(pattern_spoke_struct.ref_spoke_list)
588  src_candidates.extend(pattern_spoke_struct.src_spoke_list)
589 
590  # We can now create our full rotation matrix for both the
591  # shift and rotation. Reminder shift, aligns the pattern
592  # centers, rotation rotates the spokes on top of each other.
593  shift_rot_struct = self._create_shift_rot_matrix(
594  cos_rot_sq, shift_matrix, src_delta_array[0],
595  self._reference_array[ref_id], ref_delta)
596  # If we fail to create the rotation matrix, continue to the
597  # next objects.
598  if shift_rot_struct.sin_rot is None or \
599  shift_rot_struct.shift_rot_matrix is None:
600  continue
601 
602  # Get the data from the return struct.
603  sin_rot = shift_rot_struct.sin_rot
604  shift_rot_matrix = shift_rot_struct.shift_rot_matrix
605 
606  # Now that we have enough candidates we test to see if it
607  # passes intermediate verify. This shifts and rotates the
608  # source pattern into the reference frame and tests that each
609  # source/reference object pair is within max_dist.
610  fit_shift_rot_matrix = self._intermediate_verify(
611  src_pattern_array[src_candidates],
612  self._reference_array[ref_candidates],
613  shift_rot_matrix, max_dist_rad)
614  if fit_shift_rot_matrix is not None:
615  # Fill the struct and return.
616  output_matched_pattern.ref_candidates = ref_candidates
617  output_matched_pattern.src_candidates = src_candidates
618  output_matched_pattern.shift_rot_matrix = \
619  fit_shift_rot_matrix
620  output_matched_pattern.cos_shift = cos_shift
621  output_matched_pattern.sin_rot = sin_rot
622  return output_matched_pattern
623 
624  return output_matched_pattern
625 
626  def _find_candidate_reference_pairs(self, src_dist, ref_dist_array,
627  max_dist_rad):
628  """Wrap numpy.searchsorted to find the range of reference spokes
629  within a spoke distance tolerance of our source spoke.
630 
631  Returns an array sorted from the smallest absolute delta distance
632  between source and reference spoke length. This sorting increases the
633  speed for the pattern search greatly.
634 
635  Parameters
636  ----------
637  src_dist : float radians
638  float value of the distance we would like to search for in
639  the reference array in radians.
640  ref_dist_array : float array
641  sorted array of distances in radians.
642  max_dist_rad : float
643  maximum plus/minus search to find in the reference array in
644  radians.
645 
646  Return
647  ------
648  int array
649  indices lookup into the input ref_dist_array sorted by the
650  difference in value to the src_dist from absolute value
651  smallest to largest.
652  """
653  # Find the index of the minimum and maximum values that satisfy
654  # the tolerance.
655  start_idx = np.searchsorted(ref_dist_array, src_dist - max_dist_rad)
656  end_idx = np.searchsorted(ref_dist_array, src_dist + max_dist_rad,
657  side='right')
658 
659  # If these are equal there are no candidates and we exit.
660  if start_idx == end_idx:
661  return []
662 
663  # Make sure the endpoints of the input array are respected.
664  if start_idx < 0:
665  start_idx = 0
666  if end_idx > ref_dist_array.shape[0]:
667  end_idx = ref_dist_array.shape[0]
668 
669  # Now we sort the indices from smallest absolute delta dist difference
670  # to the largest and return the vector. This step greatly increases
671  # the speed of the algorithm.
672  tmp_diff_array = np.fabs(ref_dist_array[start_idx:end_idx] - src_dist)
673  return tmp_diff_array.argsort() + start_idx
674 
675  def _test_rotation(self, src_center, ref_center, src_delta, ref_delta,
676  cos_shift, max_cos_rot_sq):
677  """ Test if the rotation implied between the source
678  pattern and reference pattern is within tolerance. To test this
679  we need to create the first part of our spherical rotation matrix
680  which we also return for use later.
681 
682  Parameters
683  ----------
684  src_center : float array3
685  pattern.
686  ref_center : float array
687  3 vector defining the center of the candidate refence pinwheel
688  pattern.
689  src_delta : float array
690  3 vector delta between the source pattern center and the end of
691  the pinwheel spoke.
692  ref_delta : float array
693  3 vector delta of the candidate matched reference pair
694  cos_shift : float
695  Cosine of the angle between the source and reference candidate
696  centers.
697  max_cos_rot_sq : float
698  candidate reference pair after shifting the centers on top of each
699  other. The function will return None if the rotation implied is
700  greater than max_cos_rot_sq.
701 
702  Returns
703  -------
704  lsst.pipe.base.Struct
705  Return a pipe.base.Struct containing the following data.
706 
707  cos_rot_sq : float
708  magnitude of the rotation needed to align the two patterns
709  after their centers are shifted on top of each other.
710  None if rotation test fails.
711  shift_matrix : float array
712  3x3 rotation matrix describing the shift needed to align
713  the source and candidate reference center.
714  None if rotation test fails.
715  """
716 
717  # Make sure the sine is a real number.
718  if cos_shift > 1.0:
719  cos_shift = 1.
720  elif cos_shift < -1.0:
721  cos_shift = -1.
722  sin_shift = np.sqrt(1 - cos_shift ** 2)
723 
724  # If the sine of our shift is zero we only need to use the idenity
725  # matrix for the shift. Else we construct the rotation matrix for
726  # shift.
727  if sin_shift > 0:
728  rot_axis = np.cross(src_center, ref_center)
729  rot_axis /= sin_shift
730  shift_matrix = self._create_spherical_rotation_matrix(
731  rot_axis, cos_shift, sin_shift)
732  else:
733  shift_matrix = np.identity(3)
734 
735  # Now that we have our shift we apply it to the src delta vector
736  # and check the rotation.
737  rot_src_delta = np.dot(shift_matrix, src_delta)
738  cos_rot_sq = (np.dot(rot_src_delta, ref_delta)**2 /
739  (np.dot(rot_src_delta, rot_src_delta) *
740  np.dot(ref_delta, ref_delta)))
741  # If the rotation isn't in tolerance return None.
742  if cos_rot_sq < max_cos_rot_sq:
743  return pipeBase.Struct(
744  cos_rot_sq=None,
745  shift_matrix=None,)
746  return pipeBase.Struct(
747  cos_rot_sq=cos_rot_sq,
748  shift_matrix=shift_matrix,)
749 
750  def _create_spherical_rotation_matrix(self, rot_axis, cos_rotation,
751  sin_rotion):
752  """Construct a generalized 3D rotation matrix about a given
753  axis.
754 
755  Parameters
756  ----------
757  rot_axis : float array
758  3 vector defining the axis to rotate about.
759  cos_rotation : float
760  cosine of the rotation angle.
761  sin_rotation : float
762  sine of the rotation angle.
763 
764  Return
765  ------
766  float array
767  3x3 spherical, rotation matrix.
768  """
769 
770  rot_cross_matrix = np.array(
771  [[0., -rot_axis[2], rot_axis[1]],
772  [rot_axis[2], 0., -rot_axis[0]],
773  [-rot_axis[1], rot_axis[0], 0.]], dtype=np.float64)
774  shift_matrix = (cos_rotation*np.identity(3) +
775  sin_rotion*rot_cross_matrix +
776  (1. - cos_rotation)*np.outer(rot_axis, rot_axis))
777 
778  return shift_matrix
779 
780  def _create_pattern_spokes(self, src_ctr, src_delta_array, src_dist_array,
781  ref_ctr, ref_ctr_id, ref_delta, ref_dist,
782  ref_delta_array, ref_dist_array,
783  ref_id_array, max_dist_rad, n_match):
784  """ Create the individual spokes that make up the pattern now that the
785  shift and rotation are within tolerance.
786 
787  If we can't create a valid pattern we exit early.
788 
789  Parameters
790  ----------
791  src_ctr : float array
792  3 vector of the source pinwheel center
793  src_delta_array : float array
794  Array of 3 vector deltas between the source center and the pairs
795  that make up the remaining spokes of the pinwheel
796  src_dist_array : float array
797  Array of the distances of each src_delta in the pinwheel
798  ref_ctr : float array
799  3 vector of the candidate refenerce center
800  ref_ctr_id : int
801  id of the ref_ctr in the master reference array
802  ref_delta : float array
803  3 vector of the first candidate pair of the pinwheel. This is
804  the candidate pair that was matched in the
805  main _construct_pattern_and_shift_rot_matrix loop
806  ref_dist : float
807  Radian distance of the first candidate reference pair
808  ref_delta_array : float array
809  Array of 3 vector deltas that are have the current candidate
810  reference center as part of the pair
811  ref_dist_array : float array
812  Array of vector distances for each of the reference pairs
813  ref_id_array : int array
814  Array of id lookups into the master reference array that our
815  center id object is paired with.
816  max_dist_rad : float
817  Maximum search distance
818  n_match : int
819  Number of source deltas that must be matched into the reference
820  deltas in order to consider this a successful pattern match.
821 
822  Returns
823  -------
824  lsst.pipe.base.Struct
825  The Struct contains the following data:
826 
827  ref_spoke_list : list of ints specifying ids into the master
828  reference array
829  src_spoke_list : list of ints specifying indices into the current
830  source pattern that is being tested.
831  """
832  # Struct where we will be putting our results.
833  output_spokes = pipeBase.Struct(
834  ref_spoke_list=[],
835  src_spoke_list=[],)
836 
837  # Counter for number of spokes we failed to find a reference
838  # candidate for. We break the loop if we haven't found enough.
839  n_fail = 0
840  ref_spoke_list = []
841  src_spoke_list = []
842  # Loop over the source pairs.
843  for src_idx in range(1, len(src_dist_array)):
844  if n_fail > len(src_dist_array) - (n_match - 1):
845  break
846 
847  # Given our length tolerance we can use it to compute a tolernace
848  # on the angle between our spoke.
849  src_sin_tol = (max_dist_rad /
850  (src_dist_array[src_idx] + max_dist_rad))
851  # Test if the small angle approximation will still hold. This is
852  # defined as when sin(theta) ~= theta to within 0.1% of each
853  # other. This also implicitly sets a minimum spoke length that we
854  # can use.
855  if src_sin_tol > 0.0447:
856  n_fail += 1
857  continue
858 
859  # Find the reference pairs that include our candidate pattern
860  # center and sort them in increasing delta
861  ref_dist_idx_array = self._find_candidate_reference_pairs(
862  src_dist_array[src_idx], ref_dist_array, max_dist_rad)
863 
864  # Test the spokes and return the id of the reference object.
865  # Return None if no match is found.
866  ref_id = self._test_spoke(
867  src_ctr, src_delta_array[src_idx], src_dist_array[src_idx],
868  src_delta_array[0], src_dist_array[0], ref_ctr, ref_ctr_id,
869  ref_delta, ref_dist, ref_dist_idx_array, ref_delta_array,
870  ref_dist_array,
871  ref_id_array, src_sin_tol)
872  if ref_id is None:
873  n_fail += 1
874  continue
875 
876  # Append the successful indices to our list. The src_idx needs
877  # an extra iteration to skip the first and second source objects.
878  ref_spoke_list.append(ref_id)
879  src_spoke_list.append(src_idx + 1)
880  # If we found enough reference objects we can return early. This is
881  # n_match - 2 as we already have 2 source objects matched into the
882  # reference data.
883  if len(ref_spoke_list) >= n_match - 2:
884  # Set the struct data and return the struct.
885  output_spokes.ref_spoke_list = ref_spoke_list
886  output_spokes.src_spoke_list = src_spoke_list
887  return output_spokes
888 
889  return output_spokes
890 
891  def _test_spoke(self, src_ctr, src_delta, src_dist, src_ctr_delta,
892  src_ctr_dist, ref_ctr, ref_ctr_id, ref_delta, ref_dist,
893  ref_dist_idx_array, ref_delta_array, ref_dist_array,
894  ref_id_array, src_sin_tol):
895  """Test the opening angle between the first spoke of our pattern
896  for the soruce object against the reference object.
897 
898  This method makes heavy use of the small angle approximation to perform
899  the comparison.
900 
901  Parameters
902  ----------
903  src_ctr : float array
904  3 vector of the source pinwheel center
905  src_delta : float array
906  3 vector delta from the source center and the source object that
907  makes up the current spoke of the pinwheel we are testing.
908  src_dist : float array
909  Distance of the current spoke we are testing
910  src_ctr_delta : float array
911  3 vector delta between the center of the pattern and the first
912  spoke of the pattern. Used to test compute the opening angle
913  between the current spoke and the first spoke.
914  src_ctr_dist : float
915  Distance between the pairs that make up src_ctr_delta
916  ref_ctr : float array
917  3 vector of the candidate reference center
918  ref_ctr_id : int
919  id lookup of the ref_ctr into the master reference array
920  ref_delta : float array
921  3 vector of the first candidate pair of the pinwheel. That is
922  the candidate pair that was matched in the
923  main _construct_pattern_and_shift_rot_matrix loop
924  ref_dist : float
925  Radian distance of the first candidate reference pair
926  ref_dist_idx_array : int array
927  Indices sorted by the delta distance between the source
928  spoke we are trying to test and the candidate reference
929  spokes.
930  ref_delta_array : float array
931  Array of 3 vector deltas that are have the current candidate
932  reference center as part of the pair
933  ref_dist_array : float array
934  Array of vector distances for each of the reference pairs
935  ref_id_array : int array
936  Array of id lookups into the master reference array that our
937  center id object is paired with.
938  src_sin_tol : float
939  Sine of tolerance allowed between source and reference spoke
940  opening angles.
941 
942  Returns
943  -------
944  None or int
945  If we can not find a candidate spoke we return None else we
946  return an int id into the master reference array.
947  """
948 
949  # Precompute all of the source only cross and dot products so we don't
950  # have to do it for each interation in the reference loop.
951  cos_theta_src = (np.dot(src_delta, src_ctr_delta) /
952  (src_dist * src_ctr_dist))
953  cross_src = (np.cross(src_delta, src_ctr_delta) /
954  (src_dist * src_ctr_dist))
955  dot_cross_src = np.dot(cross_src, src_ctr)
956 
957  # Loop over our candidate reference objects.
958  for ref_dist_idx in ref_dist_idx_array:
959  # Check the direction of the delta vector.
960  ref_sign = 1
961  if ref_id_array[ref_dist_idx] < ref_ctr_id:
962  ref_sign = -1
963 
964  # Compute the cos between our "center" reference vector and the
965  # current reference candidate.
966  cos_theta_ref = ref_sign * (
967  np.dot(ref_delta_array[ref_dist_idx], ref_delta) /
968  (ref_dist_array[ref_dist_idx] * ref_dist))
969 
970  # Make sure we can safely make the comparison in case
971  # our "center" and candidate vectors are mostly aligned.
972  if cos_theta_ref ** 2 < (1 - src_sin_tol ** 2):
973  cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2 /
974  (1 - cos_theta_ref ** 2))
975  else:
976  cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2 /
977  src_sin_tol ** 2)
978  # Test the difference of the cosine of the reference angle against
979  # the source angle. Assumes that the delta between the two is
980  # small.
981  if cos_sq_comparison > src_sin_tol ** 2:
982  continue
983 
984  # The cosine tests the manitude of the angle but not
985  # its direction. To do that we need to know the sine as well.
986  # This cross product calculation does that.
987  cross_ref = ref_sign * (
988  np.cross(ref_delta_array[ref_dist_idx], ref_delta) /
989  (ref_dist_array[ref_dist_idx] * ref_dist))
990  dot_cross_ref = np.dot(cross_ref, ref_ctr)
991 
992  # Check the value of the cos again to make sure that it is not
993  # near zero.
994  if abs(cos_theta_src) < src_sin_tol:
995  sin_comparison = (dot_cross_src - dot_cross_ref) / src_sin_tol
996  else:
997  sin_comparison = \
998  (dot_cross_src - dot_cross_ref) / cos_theta_ref
999 
1000  if abs(sin_comparison) > src_sin_tol:
1001  continue
1002 
1003  # Return the correct id of the candidate we found.
1004  return ref_id_array[ref_dist_idx]
1005 
1006  return None
1007 
1008  def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta,
1009  ref_ctr, ref_delta):
1010  """ Create the final part of our spherical rotation matrix.
1011 
1012  Parameters
1013  ----------
1014  cos_rot_sq : float
1015  cosine of the rotation needed to align our source and reference
1016  candidate patterns.
1017  shift_matrix : float array
1018  3x3 rotation matrix for shifting the source pattern center on top
1019  of the candidate reference pattern center.
1020  src_delta : float array
1021  3 vector delta of representing the first spoke of the source
1022  pattern
1023  ref_ctr : float array
1024  3 vector on the unitsphere representing the center of our
1025  reference pattern.
1026  ref_delta : float array
1027  3 vector delta made by the first pair of the reference pattern.
1028 
1029  Returns
1030  -------
1031  lsst.pipe.base.Struct
1032  Struct object containing the following data:
1033 
1034  sin_rot : float sine of the amount of rotation between the
1035  source and reference pattern. We use sine here as it is
1036  signed and tells us the chirality of the rotation.
1037  shift_rot_matrix : float array representing the 3x3 rotation
1038  matrix that takes the source patern and shifts and rotates
1039  it to align with the reference pattern.
1040  """
1041  cos_rot = np.sqrt(cos_rot_sq)
1042  rot_src_delta = np.dot(shift_matrix, src_delta)
1043  delta_dot_cross = np.dot(np.cross(rot_src_delta, ref_delta), ref_ctr)
1044 
1045  sin_rot = np.sign(delta_dot_cross) * np.sqrt(1 - cos_rot_sq)
1046  rot_matrix = self._create_spherical_rotation_matrix(
1047  ref_ctr, cos_rot, sin_rot)
1048 
1049  shift_rot_matrix = np.dot(rot_matrix, shift_matrix)
1050 
1051  return pipeBase.Struct(
1052  sin_rot=sin_rot,
1053  shift_rot_matrix=shift_rot_matrix,)
1054 
1055  def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix,
1056  max_dist_rad):
1057  """ Perform an intermediate verify step. Rotate the matches references
1058  into the source frame and test their distances against tolerance. Only
1059  return true if all points are within tolerance.
1060 
1061  Parameters
1062  ----------
1063  src_pattern : float array
1064  Array of 3 vectors representing the points that make up our source
1065  pinwheel pattern.
1066  ref_pattern : float array
1067  Array of 3 vectors representing our candidate reference pinwheel
1068  pattern.
1069  shift_rot_matrix : float array
1070  3x3 rotation matrix that takes the source objects and rotates them
1071  onto the frame of the reference objects
1072  max_dist_rad : float
1073  Maximum distance allowed to consider two objects the same.
1074 
1075  Returns
1076  -------
1077  float array or None
1078  Return the fitted shift/rotation matrix if all of the points in our
1079  source pattern are within max_dist_rad of their matched reference
1080  objects. Returns None if this criteria is not satisfied.
1081  """
1082  if len(src_pattern) != len(ref_pattern):
1083  raise ValueError(
1084  "Source pattern length does not match ref pattern.\n"
1085  "\t source pattern len=%i, reference pattern len=%i" %
1086  (len(src_pattern), len(ref_pattern)))
1087 
1089  src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad):
1090  # Now that we know our initial shift and rot matrix is valid we
1091  # want to fit the implied shift/rotation using all points from
1092  # our pattern. This is a more robust rotation matrix as our
1093  # intial matrix only used the first 2 points from the source
1094  # pattern to estimate the shift and rotation. The matrix we
1095  # fit is forced to be unitary so no determainant tolerance is
1096  # required.
1097  fit_shift_rot_matrix = leastsq(
1098  _rotation_matrix_chi_sq, x0=shift_rot_matrix.flatten(),
1099  args=(src_pattern, ref_pattern,))[0].reshape((3, 3))
1100  # Do another verify just in case the fit wondered off.
1102  src_pattern, ref_pattern, fit_shift_rot_matrix,
1103  max_dist_rad):
1104  return fit_shift_rot_matrix
1105 
1106  return None
1107 
1108  def _intermediate_verify_comparison(self, src_pattern, ref_pattern,
1109  shift_rot_matrix, max_dist_rad):
1110  """Test the input rotation matrix against the input source and
1111  ref patterns.
1112 
1113  If every point in the source patern after rotation is within a
1114  distance of max_dist_rad to its candidate reference point, we
1115  return True.
1116 
1117  Parameters
1118  ----------
1119  src_pattern : float array
1120  Array of 3 vectors representing the points that make up our source
1121  pinwheel pattern.
1122  ref_pattern : float array
1123  Array of 3 vectors representing our candidate reference pinwheel
1124  pattern.
1125  shift_rot_matrix : float array
1126  3x3 rotation matrix that takes the source objects and rotates them
1127  onto the frame of the reference objects
1128  max_dist_rad : float
1129  Maximum distance allowed to consider two objects the same.
1130 
1131 
1132  Returns
1133  -------
1134  bool
1135  True if all rotated source points are within max_dist_rad of
1136  the candidate references matches.
1137  """
1138  shifted_ref_pattern = np.dot(shift_rot_matrix.transpose(),
1139  ref_pattern.transpose()).transpose()
1140  tmp_delta_array = src_pattern - shifted_ref_pattern
1141  tmp_dist_array = (tmp_delta_array[:, 0] ** 2 +
1142  tmp_delta_array[:, 1] ** 2 +
1143  tmp_delta_array[:, 2] ** 2)
1144  return np.all(tmp_dist_array < max_dist_rad ** 2)
1145 
1146  def _test_rotation_agreement(self, rot_vect_list, max_dist_rad):
1147  """ Test this rotation against the previous N found and return
1148  the number that a agree within tolerance to where our test
1149  points are.
1150 
1151  Parameters
1152  ----------
1153  rot_vect_list : list of lists of float arrays
1154  list of lists of rotated 3 vectors representing the maximum x, y,
1155  z extent on the unit sphere of the input source objects roated by
1156  the candidate rotations into the reference frame.
1157  max_dist_rad : float
1158  maximum distance in radians to consider two points "agreeing" on
1159  a rotation
1160 
1161  Returns
1162  -------
1163  int
1164  Number of candidate rotations that agree for all of the rotated
1165  test 3 vectors.
1166  """
1167 
1168  self.log.debug("Comparing pattern %i to previous %i rotations..." %
1169  (rot_vect_list[-1][-1], len(rot_vect_list) - 1))
1170 
1171  tot_consent = 0
1172  for rot_idx in range(max((len(rot_vect_list) - 1), 0)):
1173  tmp_dist_list = []
1174  for vect_idx in range(len(rot_vect_list[rot_idx]) - 1):
1175  tmp_delta_vect = (rot_vect_list[rot_idx][vect_idx] -
1176  rot_vect_list[-1][vect_idx])
1177  tmp_dist_list.append(
1178  np.dot(tmp_delta_vect, tmp_delta_vect))
1179  if np.all(np.array(tmp_dist_list) < max_dist_rad ** 2):
1180  tot_consent += 1
1181  return tot_consent
1182 
1183  def _match_sources(self, source_array, shift_rot_matrix):
1184  """ Shift both the reference and source catalog to the the respective
1185  frames and find their nearest neighbor using a kdTree. Removes all
1186  matches who do not agree when either the refernce or source catalog is
1187  rotated and removes all matches greated than the requested distance.
1188 
1189  Parameters
1190  ----------
1191  source_array : float array
1192  array of 3 vectors representing the source objects we are trying
1193  to match into the source catalog.
1194  shift_rot_matrix : float array
1195  3x3 rotation matrix that performs the spherical rotation from the
1196  source frame into the reference frame.
1197 
1198  Returns
1199  -------
1200  lsst.pipe.base.Struct
1201  A Struct object containing the following data
1202  matches : a (N, 2) array of integer ids into the source and
1203  reference arrays. Matches are only returned for those that
1204  satisfy the distance and handshake criteria.
1205  distances : float array of the distance between each match in
1206  radians after the shift and rotation is applied.
1207  """
1208  shifted_references = np.dot(
1209  shift_rot_matrix.transpose(),
1210  self._reference_array.transpose()).transpose()
1211  shifted_sources = np.dot(
1212  shift_rot_matrix,
1213  source_array.transpose()).transpose()
1214 
1215  ref_matches = np.empty((len(shifted_references), 2),
1216  dtype=np.uint32)
1217  src_matches = np.empty((len(shifted_sources), 2),
1218  dtype=np.uint32)
1219 
1220  ref_matches[:, 1] = np.arange(len(shifted_references),
1221  dtype=np.uint32)
1222  src_matches[:, 0] = np.arange(len(shifted_sources),
1223  dtype=np.uint32)
1224 
1225  ref_kdtree = cKDTree(self._reference_array)
1226  src_kdtree = cKDTree(source_array)
1227 
1228  tmp_src_dist, tmp_src_idx = src_kdtree.query(shifted_references)
1229  tmp_ref_dist, tmp_ref_idx = ref_kdtree.query(shifted_sources)
1230 
1231  ref_matches[:, 0] = tmp_src_idx
1232  src_matches[:, 1] = tmp_ref_idx
1233 
1234  handshake_mask = self._handshake_match(ref_matches, src_matches)
1235  return pipeBase.Struct(
1236  match_ids=ref_matches[handshake_mask],
1237  distances_rad=tmp_src_dist[handshake_mask],)
1238 
1239  def _handshake_match(self, matches_ref, matches_src):
1240  """Return only those matches where both the source
1241  and reference objects agree they they are each others'
1242  nearest neighbor.
1243 
1244  Parameters
1245  ----------
1246  matches_ref : int array
1247  (N, 2) int array of nearest neighbor matches between shifted and
1248  rotated source objects matched into the references.
1249  matches_src : int array
1250  (M, 2) int array of nearest neighbor matches between shifted and
1251  rotated reference objects matched into the sources.
1252 
1253  Return
1254  ------
1255  bool array
1256  Return the array positions where the two match catalogs agree.
1257  """
1258  handshake_mask_array = np.zeros(len(matches_ref), dtype=np.bool)
1259 
1260  for ref_match_idx, match in enumerate(matches_ref):
1261  src_match_idx = np.searchsorted(matches_src[:, 0], match[0])
1262  if match[1] == matches_src[src_match_idx, 1]:
1263  handshake_mask_array[ref_match_idx] = True
1264  return handshake_mask_array
def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta, ref_ctr, ref_delta)
def _test_spoke(self, src_ctr, src_delta, src_dist, src_ctr_delta, src_ctr_dist, ref_ctr, ref_ctr_id, ref_delta, ref_dist, ref_dist_idx_array, ref_delta_array, ref_dist_array, ref_id_array, src_sin_tol)
def _intermediate_verify_comparison(self, src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad)
def _construct_pattern_and_shift_rot_matrix(self, src_pattern_array, n_match, max_cos_theta_shift, max_cos_rot_sq, max_dist_rad)
def _test_rotation(self, src_center, ref_center, src_delta, ref_delta, cos_shift, max_cos_rot_sq)
def match(self, source_array, n_check, n_match, n_agree, max_n_patterns, max_shift, max_rotation, max_dist, min_matches, pattern_skip_array=None)
def _create_pattern_spokes(self, src_ctr, src_delta_array, src_dist_array, ref_ctr, ref_ctr_id, ref_delta, ref_dist, ref_delta_array, ref_dist_array, ref_id_array, max_dist_rad, n_match)
def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad)