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