3 from scipy.optimize
import least_squares
4 from scipy.spatial
import cKDTree
5 from scipy.stats
import sigmaclip
7 from .pessimisticPatternMatcherUtils
import check_spoke
11 def _rotation_matrix_chi_sq(flattened_rot_matrix,
15 """Compute the squared differences for least squares fitting.
17 Given a flattened rotation matrix, one N point pattern and another N point
18 pattern to transform into to, compute the squared differences between the
19 points in the two patterns after the rotation.
23 flattened_rot_matrix : `numpy.ndarray`, (9, )
24 A flattened array representing a 3x3 rotation matrix. The array is
25 flattened to comply with the API of scipy.optimize.least_squares.
26 Flattened elements are [[0, 0], [0, 1], [0, 2], [1, 0]...]
27 pattern_a : `numpy.ndarray`, (N, 3)
28 A array containing N, 3 vectors representing the objects we would like
29 to transform into the frame of pattern_b.
30 pattern_b : `numpy.ndarray`, (N, 3)
31 A array containing N, 3 vectors representing the reference frame we
32 would like to transform pattern_a into.
33 max_dist_rad : `float`
34 The maximum distance allowed from the pattern matching. This value is
35 used as the standard error for the resultant chi values.
39 noralized_diff : `numpy.ndarray`, (9,)
40 Array of differences between the vectors representing of the source
41 pattern rotated into the reference frame and the converse. This is
42 used to minimize in a least squares fitter.
45 rot_matrix = flattened_rot_matrix.reshape((3, 3))
47 rot_pattern_a = np.dot(rot_matrix, pattern_a.transpose()).transpose()
48 diff_pattern_a_to_b = rot_pattern_a - pattern_b
51 return diff_pattern_a_to_b.flatten() / max_dist_rad
55 """Class implementing a pessimistic version of Optimistic Pattern Matcher
56 B (OPMb) from Tabur 2007. See `DMTN-031 <http://ls.st/DMTN-031`_
60 reference_array : `numpy.ndarray`, (N, 3)
61 spherical points x, y, z of to use as reference objects for
63 log : `lsst.log.Log` or `logging.Logger`
64 Logger for outputting debug info.
68 The class loads and stores the reference object
69 in a convenient data structure for matching any set of source objects that
70 are assumed to contain each other. The pessimistic nature of the algorithm
71 comes from requiring that it discovers at least two patterns that agree on
72 the correct shift and rotation for matching before exiting. The original
73 behavior of OPMb can be recovered simply. Patterns matched between the
74 input datasets are n-spoked pinwheels created from n+1 points. Refer to
75 DMTN #031 for more details. http://github.com/lsst-dm/dmtn-031
86 raise ValueError(
"No reference objects supplied")
88 def _build_distances_and_angles(self):
89 """Create the data structures we will use to search for our pattern
92 Throughout this function and the rest of the class we use id to
93 reference the position in the input reference catalog and index to
94 'index' into the arrays sorted on distance.
110 endIdx = startIdx + self.
_n_reference_n_reference - 1 - ref_id
115 self.
_id_array_id_array[startIdx:endIdx, 0] = ref_id
116 self.
_id_array_id_array[startIdx:endIdx, 1] = np.arange(ref_id + 1,
122 self.
_dist_array_dist_array[startIdx:endIdx] = np.sqrt(
124 - ref_obj) ** 2).sum(axis=1))
130 sorted_dist_args = self.
_dist_array_dist_array.argsort()
134 def match(self, source_array, n_check, n_match, n_agree,
135 max_n_patterns, max_shift, max_rotation, max_dist,
136 min_matches, pattern_skip_array=None):
137 """Match a given source catalog into the loaded reference catalog.
139 Given array of points on the unit sphere and tolerances, we
140 attempt to match a pinwheel like pattern between these input sources
141 and the reference objects this class was created with. This pattern
142 informs of the shift and rotation needed to align the input source
143 objects into the frame of the references.
147 source_array : `numpy.ndarray`, (N, 3)
148 An array of spherical x,y,z coordinates and a magnitude in units
149 of objects having a lower value for sorting. The array should be
152 Number of sources to create a pattern from. Not all objects may be
153 checked if n_match criteria is before looping through all n_check
156 Number of objects to use in constructing a pattern to match.
158 Number of found patterns that must agree on their shift and
159 rotation before exiting. Set this value to 1 to recover the
160 expected behavior of Optimistic Pattern Matcher B.
161 max_n_patters : `int`
162 Number of patterns to create from the input source objects to
163 attempt to match into the reference objects.
165 Maximum allowed shift to match patterns in arcseconds.
166 max_rotation : `float`
167 Maximum allowed rotation between patterns in degrees.
169 Maximum distance in arcseconds allowed between candidate spokes in
170 the source and reference objects. Also sets that maximum distance
171 in the intermediate verify, pattern shift/rotation agreement, and
173 pattern_skip_array : `int`
174 Patterns we would like to skip. This could be due to the pattern
175 being matched on a previous iteration that we now consider invalid.
176 This assumes the ordering of the source objects is the same
177 between different runs of the matcher which, assuming no object
178 has been inserted or the magnitudes have changed, it should be.
182 output_struct : `lsst.pipe.base.Struct`
183 Result struct with components
185 - ``matches`` : (N, 2) array of matched ids for pairs. Empty list if no
186 match found (`numpy.ndarray`, (N, 2) or `list`)
187 - ``distances_rad`` : Radian distances between the matched objects.
188 Empty list if no match found (`numpy.ndarray`, (N,))
189 - ``pattern_idx``: Index of matched pattern. None if no match found
191 - ``shift`` : Magnitude for the shift between the source and reference
192 objects in arcseconds. None if no match found (`float`).
196 sorted_source_array = source_array[source_array[:, -1].argsort(), :3]
197 n_source = len(sorted_source_array)
200 output_match_struct = pipeBase.Struct(
208 self.
loglog.warning(
"Source object array is empty. Unable to match. Exiting matcher.")
221 max_cos_shift = np.cos(np.radians(max_shift / 3600.))
222 max_cos_rot_sq = np.cos(np.radians(max_rotation)) ** 2
223 max_dist_rad = np.radians(max_dist / 3600.)
227 for pattern_idx
in range(np.min((max_n_patterns,
228 n_source - n_match))):
232 if pattern_skip_array
is not None and \
233 np.any(pattern_skip_array == pattern_idx):
235 "Skipping previously matched bad pattern %i...",
239 pattern = sorted_source_array[
240 pattern_idx: np.min((pattern_idx + n_check, n_source)), :3]
245 construct_return_struct = \
247 pattern, n_match, max_cos_shift, max_cos_rot_sq,
251 if construct_return_struct.ref_candidates
is None or \
252 construct_return_struct.shift_rot_matrix
is None or \
253 construct_return_struct.cos_shift
is None or \
254 construct_return_struct.sin_rot
is None:
258 ref_candidates = construct_return_struct.ref_candidates
259 shift_rot_matrix = construct_return_struct.shift_rot_matrix
260 cos_shift = construct_return_struct.cos_shift
261 sin_rot = construct_return_struct.sin_rot
265 if len(ref_candidates) < n_match:
271 tmp_rot_vect_list = []
272 for test_vect
in test_vectors:
273 tmp_rot_vect_list.append(np.dot(shift_rot_matrix, test_vect))
280 tmp_rot_vect_list.append(pattern_idx)
281 rot_vect_list.append(tmp_rot_vect_list)
290 match_struct = self.
_final_verify_final_verify(source_array[:, :3],
294 if match_struct.match_ids
is None or \
295 match_struct.distances_rad
is None or \
296 match_struct.max_dist_rad
is None:
300 shift = np.degrees(np.arccos(cos_shift)) * 3600.
302 self.
loglog.debug(
"Succeeded after %i patterns.", pattern_idx)
303 self.
loglog.debug(
"\tShift %.4f arcsec", shift)
304 self.
loglog.debug(
"\tRotation: %.4f deg",
305 np.degrees(np.arcsin(sin_rot)))
308 output_match_struct.match_ids = \
309 match_struct.match_ids
310 output_match_struct.distances_rad = \
311 match_struct.distances_rad
312 output_match_struct.pattern_idx = pattern_idx
313 output_match_struct.shift = shift
314 output_match_struct.max_dist_rad = match_struct.max_dist_rad
315 return output_match_struct
317 self.
loglog.debug(
"Failed after %i patterns.", pattern_idx + 1)
318 return output_match_struct
320 def _compute_test_vectors(self, source_array):
321 """Compute spherical 3 vectors at the edges of the x, y, z extent
322 of the input source catalog.
326 source_array : `numpy.ndarray`, (N, 3)
327 array of 3 vectors representing positions on the unit
332 test_vectors : `numpy.ndarray`, (N, 3)
333 Array of vectors representing the maximum extents in x, y, z
334 of the input source array. These are used with the rotations
335 the code finds to test for agreement from different patterns
336 when the code is running in pessimistic mode.
340 if np.any(np.logical_not(np.isfinite(source_array))):
341 self.
loglog.warning(
"Input source objects contain non-finite values. "
342 "This could end badly.")
343 center_vect = np.nanmean(source_array, axis=0)
347 xbtm_vect = np.array([np.min(source_array[:, 0]), center_vect[1],
348 center_vect[2]], dtype=np.float64)
349 xtop_vect = np.array([np.max(source_array[:, 0]), center_vect[1],
350 center_vect[2]], dtype=np.float64)
351 xbtm_vect /= np.sqrt(np.dot(xbtm_vect, xbtm_vect))
352 xtop_vect /= np.sqrt(np.dot(xtop_vect, xtop_vect))
354 ybtm_vect = np.array([center_vect[0], np.min(source_array[:, 1]),
355 center_vect[2]], dtype=np.float64)
356 ytop_vect = np.array([center_vect[0], np.max(source_array[:, 1]),
357 center_vect[2]], dtype=np.float64)
358 ybtm_vect /= np.sqrt(np.dot(ybtm_vect, ybtm_vect))
359 ytop_vect /= np.sqrt(np.dot(ytop_vect, ytop_vect))
361 zbtm_vect = np.array([center_vect[0], center_vect[1],
362 np.min(source_array[:, 2])], dtype=np.float64)
363 ztop_vect = np.array([center_vect[0], center_vect[1],
364 np.max(source_array[:, 2])], dtype=np.float64)
365 zbtm_vect /= np.sqrt(np.dot(zbtm_vect, zbtm_vect))
366 ztop_vect /= np.sqrt(np.dot(ztop_vect, ztop_vect))
369 return np.array([xbtm_vect, xtop_vect, ybtm_vect, ytop_vect,
370 zbtm_vect, ztop_vect])
372 def _construct_pattern_and_shift_rot_matrix(self, src_pattern_array,
373 n_match, max_cos_theta_shift,
374 max_cos_rot_sq, max_dist_rad):
375 """Test an input source pattern against the reference catalog.
377 Returns the candidate matched patterns and their
378 implied rotation matrices or None.
382 src_pattern_array : `numpy.ndarray`, (N, 3)
383 Sub selection of source 3 vectors to create a pattern from
385 Number of points to attempt to create a pattern from. Must be
386 >= len(src_pattern_array)
387 max_cos_theta_shift : `float`
388 Maximum shift allowed between two patterns' centers.
389 max_cos_rot_sq : `float`
390 Maximum rotation between two patterns that have been shifted
391 to have their centers on top of each other.
392 max_dist_rad : `float`
393 Maximum delta distance allowed between the source and reference
394 pair distances to consider the reference pair a candidate for
395 the source pair. Also sets the tolerance between the opening
396 angles of the spokes when compared to the reference.
400 output_matched_pattern : `lsst.pipe.base.Struct`
401 Result struct with components:
403 - ``ref_candidates`` : ids of the matched pattern in the internal
404 reference_array object (`list` of `int`).
405 - ``src_candidates`` : Pattern ids of the sources matched
407 - ``shift_rot_matrix_src_to_ref`` : 3x3 matrix specifying the full
408 shift and rotation between the reference and source objects.
409 Rotates source into reference frame. `None` if match is not
410 found. (`numpy.ndarray`, (3, 3))
411 - ``shift_rot_matrix_ref_to_src`` : 3x3 matrix specifying the full
412 shift and rotation of the reference and source objects. Rotates
413 reference into source frame. None if match is not found
414 (`numpy.ndarray`, (3, 3)).
415 - ``cos_shift`` : Magnitude of the shift found between the two
416 patten centers. `None` if match is not found (`float`).
417 - ``sin_rot`` : float value of the rotation to align the already
418 shifted source pattern to the reference pattern. `None` if no match
426 output_matched_pattern = pipeBase.Struct(
429 shift_rot_matrix=
None,
435 src_delta_array = np.empty((len(src_pattern_array) - 1, 3))
436 src_delta_array[:, 0] = (src_pattern_array[1:, 0]
437 - src_pattern_array[0, 0])
438 src_delta_array[:, 1] = (src_pattern_array[1:, 1]
439 - src_pattern_array[0, 1])
440 src_delta_array[:, 2] = (src_pattern_array[1:, 2]
441 - src_pattern_array[0, 2])
442 src_dist_array = np.sqrt(src_delta_array[:, 0]**2
443 + src_delta_array[:, 1]**2
444 + src_delta_array[:, 2]**2)
450 src_dist_array[0], self.
_dist_array_dist_array, max_dist_rad)
453 for ref_dist_idx
in ref_dist_index_array:
457 tmp_ref_pair_list = self.
_id_array_id_array[ref_dist_idx]
458 for pair_idx, ref_id
in enumerate(tmp_ref_pair_list):
459 src_candidates = [0, 1]
461 shift_rot_matrix =
None
468 cos_shift = np.dot(src_pattern_array[0], ref_center)
469 if cos_shift < max_cos_theta_shift:
473 ref_candidates.append(ref_id)
476 ref_candidates.append(
477 tmp_ref_pair_list[1])
481 ref_candidates.append(
482 tmp_ref_pair_list[0])
492 src_pattern_array[0], ref_center, src_delta_array[0],
493 ref_delta, cos_shift, max_cos_rot_sq)
494 if test_rot_struct.cos_rot_sq
is None or \
495 test_rot_struct.proj_ref_ctr_delta
is None or \
496 test_rot_struct.shift_matrix
is None:
500 cos_rot_sq = test_rot_struct.cos_rot_sq
501 proj_ref_ctr_delta = test_rot_struct.proj_ref_ctr_delta
502 shift_matrix = test_rot_struct.shift_matrix
509 tmp_ref_dist_array = np.sqrt(
512 ** 2).sum(axis=1)).astype(
"float32")
513 tmp_sorted_args = np.argsort(tmp_ref_dist_array)
514 tmp_ref_id_array = tmp_ref_id_array[tmp_sorted_args]
515 tmp_ref_dist_array = tmp_ref_dist_array[tmp_sorted_args]
520 src_pattern_array[0], src_delta_array, src_dist_array,
522 tmp_ref_dist_array, tmp_ref_id_array, max_dist_rad,
527 if len(pattern_spoke_struct.ref_spoke_list) < n_match - 2
or \
528 len(pattern_spoke_struct.src_spoke_list) < n_match - 2:
532 ref_candidates.extend(pattern_spoke_struct.ref_spoke_list)
533 src_candidates.extend(pattern_spoke_struct.src_spoke_list)
539 cos_rot_sq, shift_matrix, src_delta_array[0],
543 if shift_rot_struct.sin_rot
is None or \
544 shift_rot_struct.shift_rot_matrix
is None:
548 sin_rot = shift_rot_struct.sin_rot
549 shift_rot_matrix = shift_rot_struct.shift_rot_matrix
558 src_pattern_array[src_candidates],
560 shift_rot_matrix, max_dist_rad)
562 if fit_shift_rot_matrix
is not None:
564 output_matched_pattern.ref_candidates = ref_candidates
565 output_matched_pattern.src_candidates = src_candidates
566 output_matched_pattern.shift_rot_matrix = \
568 output_matched_pattern.cos_shift = cos_shift
569 output_matched_pattern.sin_rot = sin_rot
570 return output_matched_pattern
572 return output_matched_pattern
574 def _find_candidate_reference_pairs(self, src_dist, ref_dist_array,
576 """Wrap numpy.searchsorted to find the range of reference spokes
577 within a spoke distance tolerance of our source spoke.
579 Returns an array sorted from the smallest absolute delta distance
580 between source and reference spoke length. This sorting increases the
581 speed for the pattern search greatly.
586 float value of the distance we would like to search for in
587 the reference array in radians.
588 ref_dist_array : `numpy.ndarray`, (N,)
589 sorted array of distances in radians.
590 max_dist_rad : `float`
591 maximum plus/minus search to find in the reference array in
596 tmp_diff_array : `numpy.ndarray`, (N,)
597 indices lookup into the input ref_dist_array sorted by the
598 difference in value to the src_dist from absolute value
603 start_idx = np.searchsorted(ref_dist_array, src_dist - max_dist_rad)
604 end_idx = np.searchsorted(ref_dist_array, src_dist + max_dist_rad,
608 if start_idx == end_idx:
609 return np.array([], dtype=
"int")
614 if end_idx > ref_dist_array.shape[0]:
615 end_idx = ref_dist_array.shape[0]
620 tmp_diff_array = np.fabs(ref_dist_array[start_idx:end_idx] - src_dist)
621 return tmp_diff_array.argsort() + start_idx
623 def _test_rotation(self, src_center, ref_center, src_delta, ref_delta,
624 cos_shift, max_cos_rot_sq):
625 """ Test if the rotation implied between the source
626 pattern and reference pattern is within tolerance. To test this
627 we need to create the first part of our spherical rotation matrix
628 which we also return for use later.
632 src_center : `numpy.ndarray`, (N, 3)
634 ref_center : `numpy.ndarray`, (N, 3)
635 3 vector defining the center of the candidate reference pinwheel
637 src_delta : `numpy.ndarray`, (N, 3)
638 3 vector delta between the source pattern center and the end of
640 ref_delta : `numpy.ndarray`, (N, 3)
641 3 vector delta of the candidate matched reference pair
643 Cosine of the angle between the source and reference candidate
645 max_cos_rot_sq : `float`
646 candidate reference pair after shifting the centers on top of each
647 other. The function will return None if the rotation implied is
648 greater than max_cos_rot_sq.
652 result : `lsst.pipe.base.Struct`
653 Result struct with components:
655 - ``cos_rot_sq`` : magnitude of the rotation needed to align the
656 two patterns after their centers are shifted on top of each
657 other. `None` if rotation test fails (`float`).
658 - ``shift_matrix`` : 3x3 rotation matrix describing the shift needed to
659 align the source and candidate reference center. `None` if rotation
660 test fails (`numpy.ndarray`, (N, 3)).
666 elif cos_shift < -1.0:
668 sin_shift = np.sqrt(1 - cos_shift ** 2)
674 rot_axis = np.cross(src_center, ref_center)
675 rot_axis /= sin_shift
677 rot_axis, cos_shift, sin_shift)
679 shift_matrix = np.identity(3)
683 rot_src_delta = np.dot(shift_matrix, src_delta)
684 proj_src_delta = (rot_src_delta
685 - np.dot(rot_src_delta, ref_center) * ref_center)
686 proj_ref_delta = (ref_delta
687 - np.dot(ref_delta, ref_center) * ref_center)
688 cos_rot_sq = (np.dot(proj_src_delta, proj_ref_delta) ** 2
689 / (np.dot(proj_src_delta, proj_src_delta)
690 * np.dot(proj_ref_delta, proj_ref_delta)))
692 if cos_rot_sq < max_cos_rot_sq:
693 return pipeBase.Struct(
695 proj_ref_ctr_delta=
None,
699 return pipeBase.Struct(
700 cos_rot_sq=cos_rot_sq,
701 proj_ref_ctr_delta=proj_ref_delta,
702 shift_matrix=shift_matrix,)
704 def _create_spherical_rotation_matrix(self, rot_axis, cos_rotation,
706 """Construct a generalized 3D rotation matrix about a given
711 rot_axis : `numpy.ndarray`, (3,)
712 3 vector defining the axis to rotate about.
713 cos_rotation : `float`
714 cosine of the rotation angle.
715 sin_rotation : `float`
716 sine of the rotation angle.
720 shift_matrix : `numpy.ndarray`, (3, 3)
721 3x3 spherical, rotation matrix.
724 rot_cross_matrix = np.array(
725 [[0., -rot_axis[2], rot_axis[1]],
726 [rot_axis[2], 0., -rot_axis[0]],
727 [-rot_axis[1], rot_axis[0], 0.]], dtype=np.float64)
728 shift_matrix = (cos_rotation*np.identity(3)
729 + sin_rotion*rot_cross_matrix
730 + (1. - cos_rotation)*np.outer(rot_axis, rot_axis))
734 def _create_pattern_spokes(self, src_ctr, src_delta_array, src_dist_array,
735 ref_ctr, ref_ctr_id, proj_ref_ctr_delta,
736 ref_dist_array, ref_id_array, max_dist_rad,
738 """ Create the individual spokes that make up the pattern now that the
739 shift and rotation are within tolerance.
741 If we can't create a valid pattern we exit early.
745 src_ctr : `numpy.ndarray`, (3,)
746 3 vector of the source pinwheel center
747 src_delta_array : `numpy.ndarray`, (N, 3)
748 Array of 3 vector deltas between the source center and the pairs
749 that make up the remaining spokes of the pinwheel
750 src_dist_array : `numpy.ndarray`, (N, 3)
751 Array of the distances of each src_delta in the pinwheel
752 ref_ctr : `numpy.ndarray`, (3,)
753 3 vector of the candidate reference center
755 id of the ref_ctr in the master reference array
756 proj_ref_ctr_delta : `numpy.ndarray`, (3,)
757 Plane projected 3 vector formed from the center point of the
758 candidate pin-wheel and the second point in the pattern to create
759 the first spoke pair. This is the candidate pair that was matched
760 in the main _construct_pattern_and_shift_rot_matrix loop
761 ref_dist_array : `numpy.ndarray`, (N,)
762 Array of vector distances for each of the reference pairs
763 ref_id_array : `numpy.ndarray`, (N,)
764 Array of id lookups into the master reference array that our
765 center id object is paired with.
766 max_dist_rad : `float`
767 Maximum search distance
769 Number of source deltas that must be matched into the reference
770 deltas in order to consider this a successful pattern match.
774 output_spokes : `lsst.pipe.base.Struct`
775 Result struct with components:
777 - ``ref_spoke_list`` : list of ints specifying ids into the master
778 reference array (`list` of `int`).
779 - ``src_spoke_list`` : list of ints specifying indices into the
780 current source pattern that is being tested (`list` of `int`).
783 output_spokes = pipeBase.Struct(
795 proj_src_ctr_delta = (src_delta_array[0]
796 - np.dot(src_delta_array[0], src_ctr) * src_ctr)
797 proj_src_ctr_dist_sq = np.dot(proj_src_ctr_delta, proj_src_ctr_delta)
800 proj_ref_ctr_dist_sq = np.dot(proj_ref_ctr_delta, proj_ref_ctr_delta)
803 for src_idx
in range(1, len(src_dist_array)):
804 if n_fail > len(src_dist_array) - (n_match - 1):
809 src_sin_tol = (max_dist_rad
810 / (src_dist_array[src_idx] + max_dist_rad))
817 if src_sin_tol > max_sin_tol:
818 src_sin_tol = max_sin_tol
823 src_delta_array[src_idx]
824 - np.dot(src_delta_array[src_idx], src_ctr) * src_ctr)
825 geom_dist_src = np.sqrt(
826 np.dot(proj_src_delta, proj_src_delta)
827 * proj_src_ctr_dist_sq)
830 cos_theta_src = (np.dot(proj_src_delta, proj_src_ctr_delta)
832 cross_src = (np.cross(proj_src_delta, proj_src_ctr_delta)
834 sin_theta_src = np.dot(cross_src, src_ctr)
839 src_dist_array[src_idx], ref_dist_array, max_dist_rad)
848 proj_ref_ctr_dist_sq,
859 ref_spoke_list.append(ref_id)
860 src_spoke_list.append(src_idx + 1)
864 if len(ref_spoke_list) >= n_match - 2:
866 output_spokes.ref_spoke_list = ref_spoke_list
867 output_spokes.src_spoke_list = src_spoke_list
872 def _test_spoke(self, cos_theta_src, sin_theta_src, ref_ctr, ref_ctr_id,
873 proj_ref_ctr_delta, proj_ref_ctr_dist_sq,
874 ref_dist_idx_array, ref_id_array, src_sin_tol):
875 """Test the opening angle between the first spoke of our pattern
876 for the source object against the reference object.
878 This method makes heavy use of the small angle approximation to perform
883 cos_theta_src : `float`
884 Cosine of the angle between the current candidate source spoke and
886 sin_theta_src : `float`
887 Sine of the angle between the current candidate source spoke and
889 ref_ctr : `numpy.ndarray`, (3,)
890 3 vector of the candidate reference center
892 id lookup of the ref_ctr into the master reference array
893 proj_ref_ctr_delta : `float`
894 Plane projected first spoke in the reference pattern using the
895 pattern center as normal.
896 proj_ref_ctr_dist_sq : `float`
897 Squared length of the projected vector.
898 ref_dist_idx_array : `numpy.ndarray`, (N,)
899 Indices sorted by the delta distance between the source
900 spoke we are trying to test and the candidate reference
902 ref_id_array : `numpy.ndarray`, (N,)
903 Array of id lookups into the master reference array that our
904 center id object is paired with.
905 src_sin_tol : `float`
906 Sine of tolerance allowed between source and reference spoke
912 If we can not find a candidate spoke we return `None` else we
913 return an int id into the master reference array.
917 for ref_dist_idx
in ref_dist_idx_array:
919 ref_delta = (self.
_reference_array_reference_array[ref_id_array[ref_dist_idx]]
923 proj_ref_delta = ref_delta - np.dot(ref_delta, ref_ctr) * ref_ctr
924 geom_dist_ref = np.sqrt(proj_ref_ctr_dist_sq
925 * np.dot(proj_ref_delta, proj_ref_delta))
926 cos_theta_ref = (np.dot(proj_ref_delta, proj_ref_ctr_delta)
931 if cos_theta_ref ** 2 < (1 - src_sin_tol ** 2):
932 cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2
933 / (1 - cos_theta_ref ** 2))
935 cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2
940 if cos_sq_comparison > src_sin_tol ** 2:
946 cross_ref = (np.cross(proj_ref_delta, proj_ref_ctr_delta)
948 sin_theta_ref = np.dot(cross_ref, ref_ctr)
952 if abs(cos_theta_src) < src_sin_tol:
953 sin_comparison = (sin_theta_src - sin_theta_ref) / src_sin_tol
956 (sin_theta_src - sin_theta_ref) / cos_theta_ref
957 if abs(sin_comparison) > src_sin_tol:
961 return ref_id_array[ref_dist_idx]
965 def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta,
967 """ Create the final part of our spherical rotation matrix.
972 cosine of the rotation needed to align our source and reference
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
980 ref_ctr : `numpy.ndarray`, (3,)
981 3 vector on the unit-sphere representing the center of our
983 ref_delta : `numpy.ndarray`, (3,)
984 3 vector delta made by the first pair of the reference pattern.
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)).
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)
1004 ref_ctr, cos_rot, sin_rot)
1006 shift_rot_matrix = np.dot(rot_matrix, shift_matrix)
1008 return pipeBase.Struct(
1010 shift_rot_matrix=shift_rot_matrix,)
1012 def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix,
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
1022 src_pattern : `numpy.ndarray`, (N,3)
1023 Array of 3 vectors representing the points that make up our source
1025 ref_pattern : `numpy.ndarray`, (N,3)
1026 Array of 3 vectors representing our candidate reference pinwheel
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.
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.
1041 if len(src_pattern) != len(ref_pattern):
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)))
1048 src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad):
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)
1063 src_pattern, ref_pattern, fit_shift_rot_matrix,
1065 return fit_shift_rot_matrix
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
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
1080 pattern_a : `numpy.ndarray`, (N,3)
1081 Array of 3 vectors representing the points that make up our source
1083 pattern_b : `numpy.ndarray`, (N,3)
1084 Array of 3 vectors representing our candidate reference pinwheel
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.
1096 True if all rotated source points are within max_dist_rad of
1097 the candidate references matches.
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
1111 This is useful for assuring the non unitary transforms do not contain
1112 too much distortion.
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
1127 dists = (test_pattern[:, 0] ** 2
1128 + test_pattern[:, 1] ** 2
1129 + test_pattern[:, 2] ** 2)
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
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
1152 Number of candidate rotations that agree for all of the rotated
1156 self.
loglog.debug(
"Comparing pattern %i to previous %i rotations...",
1157 rot_vects[-1][-1], len(rot_vects) - 1)
1160 for rot_idx
in range(max((len(rot_vects) - 1), 0)):
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):
1171 def _final_verify(self,
1176 """Match the all sources into the reference catalog using the shift/rot
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
1185 source_array : `numpy.ndarray` (N, 3)
1186 3-vector positions on the unit-sphere representing the sources to
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.
1194 Minimum number of matched objects required to consider the
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`).
1210 output_struct = pipeBase.Struct(
1217 match_sources_struct = self.
_match_sources_match_sources(source_array,
1219 cut_ids = match_sources_struct.match_ids[
1220 match_sources_struct.distances_rad < max_dist_rad]
1222 n_matched = len(cut_ids)
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
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],
1245 source_array, fit_shift_rot_matrix)
1250 match_sources_struct.distances_rad < max_dist_rad)
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
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
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.
1276 distances_rad : `numpy.ndarray`, (N,)
1277 Distances between pairs.
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.
1289 clipped_dists, _, clipped_max_dist = sigmaclip(
1295 if clipped_max_dist < 1e-16:
1296 clipped_max_dist = 1e-16
1297 n_matched_clipped = np.sum(distances_rad < clipped_max_dist)
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,
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
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.
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)).
1335 shifted_references = np.dot(
1336 np.linalg.inv(shift_rot_matrix),
1338 shifted_sources = np.dot(
1340 source_array.transpose()).transpose()
1342 ref_matches = np.empty((len(shifted_references), 2),
1344 src_matches = np.empty((len(shifted_sources), 2),
1347 ref_matches[:, 1] = np.arange(len(shifted_references),
1349 src_matches[:, 0] = np.arange(len(shifted_sources),
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_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'
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.
1383 handshake_mask_array : `numpy.ndarray`, (N,)
1384 Return the array positions where the two match catalogs agree.
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
def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad)
def _build_distances_and_angles(self)
def _clip_distances(self, distances_rad)
def _test_rotation(self, src_center, ref_center, src_delta, ref_delta, cos_shift, max_cos_rot_sq)
def _handshake_match(self, matches_src, matches_ref)
def _test_pattern_lengths(self, test_pattern, max_dist_rad)
def _create_pattern_spokes(self, src_ctr, src_delta_array, src_dist_array, ref_ctr, ref_ctr_id, proj_ref_ctr_delta, ref_dist_array, ref_id_array, max_dist_rad, n_match)
def _test_rotation_agreement(self, rot_vects, max_dist_rad)
def _create_spherical_rotation_matrix(self, rot_axis, cos_rotation, sin_rotion)
def _find_candidate_reference_pairs(self, src_dist, ref_dist_array, max_dist_rad)
def _construct_pattern_and_shift_rot_matrix(self, src_pattern_array, n_match, max_cos_theta_shift, max_cos_rot_sq, max_dist_rad)
def _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 _match_sources(self, source_array, shift_rot_matrix)
def __init__(self, reference_array, log)
def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta, ref_ctr, ref_delta)
def _compute_test_vectors(self, source_array)
def _intermediate_verify_comparison(self, pattern_a, pattern_b, shift_rot_matrix, max_dist_rad)
int check_spoke(double cos_theta_src, double sin_theta_src, ndarray::Array< double, 1, 1 > const &ref_ctr, ndarray::Array< double, 1, 1 > const &proj_ref_ctr_delta, double proj_ref_ctr_dist_sq, ndarray::Array< long int, 1, 1 > const &ref_dist_idx_array, ndarray::Array< uint16_t, 1, 1 > const &ref_id_array, ndarray::Array< double, 2, 1 > const &reference_array, double src_sin_tol)
Check the opening angle between the first spoke of our pattern for the source object against the refe...