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