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